xref: /petsc/src/dm/impls/plex/tests/ex16.c (revision d3ef4daaf99d19e98efa749367c59a1afaa53b93)
1c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5*d3ef4daaSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
6c4762a1bSJed Brown {
7*d3ef4daaSMatthew G. Knepley   PetscFunctionBeginUser;
89566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
99566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
109566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
119566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
12c4762a1bSJed Brown   PetscFunctionReturn(0);
13c4762a1bSJed Brown }
14c4762a1bSJed Brown 
15*d3ef4daaSMatthew G. Knepley // Label half of the cells
16*d3ef4daaSMatthew G. Knepley static PetscErrorCode CreateHalfCellsLabel(DM dm, PetscBool lower, DMLabel *label)
17*d3ef4daaSMatthew G. Knepley {
18*d3ef4daaSMatthew G. Knepley   PetscInt cStart, cEnd, cStartSub, cEndSub;
19*d3ef4daaSMatthew G. Knepley 
20*d3ef4daaSMatthew G. Knepley   PetscFunctionBeginUser;
21*d3ef4daaSMatthew G. Knepley   PetscCall(DMCreateLabel(dm, "cells"));
22*d3ef4daaSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "cells", label));
23*d3ef4daaSMatthew G. Knepley   PetscCall(DMLabelClearStratum(*label, 1));
24*d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
25*d3ef4daaSMatthew G. Knepley   if (lower) {cStartSub = cStart; cEndSub = cEnd/2;}
26*d3ef4daaSMatthew G. Knepley   else       {cStartSub = cEnd/2; cEndSub = cEnd;}
27*d3ef4daaSMatthew G. Knepley   for (PetscInt c = cStartSub; c < cEndSub; ++c) PetscCall(DMLabelSetValue(*label, c, 1));
28*d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexLabelComplete(dm, *label));
29*d3ef4daaSMatthew G. Knepley   PetscFunctionReturn(0);
30*d3ef4daaSMatthew G. Knepley }
31*d3ef4daaSMatthew G. Knepley 
32*d3ef4daaSMatthew G. Knepley // Label everything on the right half of the domain
33*d3ef4daaSMatthew G. Knepley static PetscErrorCode CreateHalfDomainLabel(DM dm, PetscBool lower, DMLabel *label)
34*d3ef4daaSMatthew G. Knepley {
35*d3ef4daaSMatthew G. Knepley   PetscReal centroid[3];
36*d3ef4daaSMatthew G. Knepley   PetscInt  cStart, cEnd, cdim;
37*d3ef4daaSMatthew G. Knepley 
38*d3ef4daaSMatthew G. Knepley   PetscFunctionBeginUser;
39*d3ef4daaSMatthew G. Knepley   PetscCall(DMCreateLabel(dm, "cells"));
40*d3ef4daaSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "cells", label));
41*d3ef4daaSMatthew G. Knepley   PetscCall(DMLabelClearStratum(*label, 1));
42*d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
43*d3ef4daaSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
44*d3ef4daaSMatthew G. Knepley   for (PetscInt c = cStart; c < cEnd; ++c) {
45*d3ef4daaSMatthew G. Knepley     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
46*d3ef4daaSMatthew G. Knepley     if (lower) {if (centroid[0] < 0.5) PetscCall(DMLabelSetValue(*label, c, 1));}
47*d3ef4daaSMatthew G. Knepley     else       {if (centroid[0] > 0.5) PetscCall(DMLabelSetValue(*label, c, 1));}
48*d3ef4daaSMatthew G. Knepley   }
49*d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexLabelComplete(dm, *label));
50*d3ef4daaSMatthew G. Knepley   PetscFunctionReturn(0);
51*d3ef4daaSMatthew G. Knepley }
52*d3ef4daaSMatthew G. Knepley 
53*d3ef4daaSMatthew G. Knepley PetscErrorCode CreateSubmesh(DM dm, PetscBool domain, PetscBool lower, DM *subdm)
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   DMLabel label, map;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   PetscFunctionBegin;
58*d3ef4daaSMatthew G. Knepley   if (domain) PetscCall(CreateHalfDomainLabel(dm, lower, &label));
59*d3ef4daaSMatthew G. Knepley   else        PetscCall(CreateHalfCellsLabel(dm, lower, &label));
609566063dSJacob Faibussowitsch   PetscCall(DMPlexFilter(dm, label, 1, subdm));
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *subdm, "Submesh"));
62*d3ef4daaSMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *subdm, "sub_"));
639566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view"));
649566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSubpointMap(*subdm, &map));
65*d3ef4daaSMatthew G. Knepley   PetscCall(PetscObjectViewFromOptions((PetscObject) map, NULL, "-map_view"));
66c4762a1bSJed Brown   PetscFunctionReturn(0);
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown int main(int argc, char **argv)
70c4762a1bSJed Brown {
71c4762a1bSJed Brown   DM        dm, subdm;
72*d3ef4daaSMatthew G. Knepley   PetscBool domain = PETSC_FALSE;
73c4762a1bSJed Brown 
74327415f7SBarry Smith   PetscFunctionBeginUser;
759566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
76*d3ef4daaSMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-domain", &domain, NULL));
779566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
78*d3ef4daaSMatthew G. Knepley   PetscCall(CreateSubmesh(dm, domain, PETSC_TRUE, &subdm));
799566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(subdm));
809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&subdm));
81*d3ef4daaSMatthew G. Knepley   PetscCall(CreateSubmesh(dm, domain, PETSC_FALSE, &subdm));
829566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(subdm));
839566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&subdm));
849566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
86b122ec5aSJacob Faibussowitsch   return 0;
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown /*TEST
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   test:
92c4762a1bSJed Brown     suffix: 0
93c4762a1bSJed Brown     requires: triangle
94*d3ef4daaSMatthew G. Knepley     args: -dm_coord_space 0 -sub_dm_plex_check_all \
95*d3ef4daaSMatthew G. Knepley           -dm_view ascii::ascii_info_detail -sub_dm_view ascii::ascii_info_detail -map_view
96*d3ef4daaSMatthew G. Knepley 
97*d3ef4daaSMatthew G. Knepley   # These tests check that filtering is stable when boundary point ownership could change, so it needs 3 processes
98*d3ef4daaSMatthew G. Knepley   testset:
99*d3ef4daaSMatthew G. Knepley     nsize: 3
100*d3ef4daaSMatthew G. Knepley     requires: parmetis
101*d3ef4daaSMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 20,20 -petscpartitioner_type parmetis -dm_distribute_overlap 1 -sub_dm_distribute 0 \
102*d3ef4daaSMatthew G. Knepley           -sub_dm_plex_check_all -dm_view -sub_dm_view
103*d3ef4daaSMatthew G. Knepley 
104*d3ef4daaSMatthew G. Knepley     test:
105*d3ef4daaSMatthew G. Knepley       suffix: 1
106*d3ef4daaSMatthew G. Knepley 
107*d3ef4daaSMatthew G. Knepley     test:
108*d3ef4daaSMatthew G. Knepley       suffix: 2
109*d3ef4daaSMatthew G. Knepley       args: -domain
110c4762a1bSJed Brown 
111c4762a1bSJed Brown TEST*/
112