xref: /petsc/src/dm/impls/plex/tests/ex71.c (revision 90a64276df3bb793126cd17906ff016f564b4b70)
1*90a64276Sksagiyam static char help[] = "Tests for submesh creation\n\n";
2*90a64276Sksagiyam 
3*90a64276Sksagiyam #include <petscdmplex.h>
4*90a64276Sksagiyam #include <petscsf.h>
5*90a64276Sksagiyam 
6*90a64276Sksagiyam /* Submesh of a 2 x 2 mesh using 3 processes.
7*90a64276Sksagiyam 
8*90a64276Sksagiyam Local numbering on each rank:
9*90a64276Sksagiyam 
10*90a64276Sksagiyam      (6)(16)-(7)(17)-(8)       5--14---6--15---7      (10)(20)(11)(21)(12)
11*90a64276Sksagiyam       |       |       |        |       |       |        |       |       |
12*90a64276Sksagiyam     (18) (1)(19) (2)(20)      16   0  17   1  18      (22) (2)(23) (3)(24)
13*90a64276Sksagiyam       |       |       |        |       |       |        |       |       |
14*90a64276Sksagiyam      (5)(15)(11)(22)(12)       4--13-(11)(22)(12)      (9)(19)--6--14---7
15*90a64276Sksagiyam       |       |       |        |       |       |        |       |       |
16*90a64276Sksagiyam      14   0 (23) (3)(24)     (20) (2)(23) (3)(24)     (18) (1) 15   0  16
17*90a64276Sksagiyam       |       |       |        |       |       |        |       |       |
18*90a64276Sksagiyam       4--13--(9)(21)(10)      (8)(19)-(9)(21)(10)      (8)(17)--4--13---5
19*90a64276Sksagiyam 
20*90a64276Sksagiyam            mesh_0                   mesh_1                   mesh_2
21*90a64276Sksagiyam 
22*90a64276Sksagiyam where () represents ghost points. We extract the left 2 cells.
23*90a64276Sksagiyam With sanitize_submesh = PETSC_FALSE, we get:
24*90a64276Sksagiyam 
25*90a64276Sksagiyam      (4)(11)-(5)               3---9---4
26*90a64276Sksagiyam       |       |                |       |
27*90a64276Sksagiyam     (12) (1)(13)              10   0  11
28*90a64276Sksagiyam       |       |                |       |
29*90a64276Sksagiyam      (3)(10) (7)               2---8---7
30*90a64276Sksagiyam       |       |                |       |
31*90a64276Sksagiyam       9   0 (14)             (13) (1) 14
32*90a64276Sksagiyam       |       |                |       |
33*90a64276Sksagiyam       2---8--(6)              (5)(12)--6
34*90a64276Sksagiyam 
35*90a64276Sksagiyam On the other hand, with sanitize_submesh = PETSC_TRUE, we get:
36*90a64276Sksagiyam 
37*90a64276Sksagiyam      (4)(11)-(5)               3---9---4
38*90a64276Sksagiyam       |       |                |       |
39*90a64276Sksagiyam     (12) (1)(13)              10   0  11
40*90a64276Sksagiyam       |       |                |       |
41*90a64276Sksagiyam      (3)(10) (7)               2---8---7
42*90a64276Sksagiyam       |       |                |       |
43*90a64276Sksagiyam       9   0  14              (13) (1)(14)
44*90a64276Sksagiyam       |       |                |       |
45*90a64276Sksagiyam       2---8---6               (5)(12)-(6)
46*90a64276Sksagiyam 
47*90a64276Sksagiyam         submesh_0                submesh_1               submesh_2
48*90a64276Sksagiyam 
49*90a64276Sksagiyam as points 15 and 4 of mesh_2 are in the closure of a submesh cell owned by rank 0 (point 0 of submesh_0),
50*90a64276Sksagiyam and not in the closure of any submesh cell owned by rank 1.
51*90a64276Sksagiyam 
52*90a64276Sksagiyam */
53*90a64276Sksagiyam 
54*90a64276Sksagiyam typedef struct {
55*90a64276Sksagiyam   PetscBool ignoreLabelHalo; /* Ignore filter values in the halo. */
56*90a64276Sksagiyam   PetscBool sanitizeSubmesh; /* Sanitize submesh. */
57*90a64276Sksagiyam } AppCtx;
58*90a64276Sksagiyam 
59*90a64276Sksagiyam PetscErrorCode ProcessOptions(AppCtx *options)
60*90a64276Sksagiyam {
61*90a64276Sksagiyam   PetscFunctionBegin;
62*90a64276Sksagiyam   options->ignoreLabelHalo = PETSC_FALSE;
63*90a64276Sksagiyam   options->sanitizeSubmesh = PETSC_FALSE;
64*90a64276Sksagiyam 
65*90a64276Sksagiyam   PetscOptionsBegin(PETSC_COMM_SELF, "", "Filtering Problem Options", "DMPLEX");
66*90a64276Sksagiyam   PetscCall(PetscOptionsBool("-ignore_label_halo", "Ignore filter values in the halo", "ex80.c", options->ignoreLabelHalo, &options->ignoreLabelHalo, NULL));
67*90a64276Sksagiyam   PetscCall(PetscOptionsBool("-sanitize_submesh", "Sanitize submesh", "ex80.c", options->sanitizeSubmesh, &options->sanitizeSubmesh, NULL));
68*90a64276Sksagiyam   PetscOptionsEnd();
69*90a64276Sksagiyam   PetscFunctionReturn(PETSC_SUCCESS);
70*90a64276Sksagiyam }
71*90a64276Sksagiyam int main(int argc, char **argv)
72*90a64276Sksagiyam {
73*90a64276Sksagiyam   DM             dm, subdm;
74*90a64276Sksagiyam   PetscSF        ownershipTransferSF;
75*90a64276Sksagiyam   DMLabel        filter;
76*90a64276Sksagiyam   const PetscInt filterValue = 1;
77*90a64276Sksagiyam   MPI_Comm       comm;
78*90a64276Sksagiyam   PetscMPIInt    size, rank;
79*90a64276Sksagiyam   AppCtx         user;
80*90a64276Sksagiyam 
81*90a64276Sksagiyam   PetscFunctionBeginUser;
82*90a64276Sksagiyam   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
83*90a64276Sksagiyam   PetscCall(ProcessOptions(&user));
84*90a64276Sksagiyam   comm = PETSC_COMM_WORLD;
85*90a64276Sksagiyam   PetscCallMPI(MPI_Comm_size(comm, &size));
86*90a64276Sksagiyam   if (size != 3) {
87*90a64276Sksagiyam     PetscCall(PetscPrintf(comm, "This example is specifically designed for comm size == 3.\n"));
88*90a64276Sksagiyam     PetscCall(PetscFinalize());
89*90a64276Sksagiyam     return 0;
90*90a64276Sksagiyam   }
91*90a64276Sksagiyam   PetscCallMPI(MPI_Comm_rank(comm, &rank));
92*90a64276Sksagiyam   {
93*90a64276Sksagiyam     DM             pdm;
94*90a64276Sksagiyam     const PetscInt faces[2] = {2, 2};
95*90a64276Sksagiyam     PetscInt       overlap  = 1;
96*90a64276Sksagiyam 
97*90a64276Sksagiyam     PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_FALSE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm));
98*90a64276Sksagiyam     {
99*90a64276Sksagiyam       PetscPartitioner part;
100*90a64276Sksagiyam       PetscInt        *sizes  = NULL;
101*90a64276Sksagiyam       PetscInt        *points = NULL;
102*90a64276Sksagiyam 
103*90a64276Sksagiyam       if (rank == 0) {
104*90a64276Sksagiyam         PetscInt sizes1[3]  = {1, 2, 1};
105*90a64276Sksagiyam         PetscInt points1[4] = {0, 2, 3, 1};
106*90a64276Sksagiyam 
107*90a64276Sksagiyam         PetscCall(PetscMalloc2(3, &sizes, 4, &points));
108*90a64276Sksagiyam         PetscCall(PetscArraycpy(sizes, sizes1, 3));
109*90a64276Sksagiyam         PetscCall(PetscArraycpy(points, points1, 4));
110*90a64276Sksagiyam       }
111*90a64276Sksagiyam       PetscCall(DMPlexGetPartitioner(dm, &part));
112*90a64276Sksagiyam       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
113*90a64276Sksagiyam       PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
114*90a64276Sksagiyam       PetscCall(PetscFree2(sizes, points));
115*90a64276Sksagiyam     }
116*90a64276Sksagiyam     PetscCall(DMSetAdjacency(dm, -1, PETSC_FALSE, PETSC_TRUE));
117*90a64276Sksagiyam     PetscCall(DMPlexDistribute(dm, overlap, NULL, &pdm));
118*90a64276Sksagiyam     if (pdm) {
119*90a64276Sksagiyam       PetscCall(DMDestroy(&dm));
120*90a64276Sksagiyam       dm = pdm;
121*90a64276Sksagiyam     }
122*90a64276Sksagiyam   }
123*90a64276Sksagiyam   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "filter", &filter));
124*90a64276Sksagiyam   switch (rank) {
125*90a64276Sksagiyam   case 0:
126*90a64276Sksagiyam     PetscCall(DMLabelSetValue(filter, 0, filterValue));
127*90a64276Sksagiyam     PetscCall(DMLabelSetValue(filter, 1, filterValue));
128*90a64276Sksagiyam     break;
129*90a64276Sksagiyam   case 1:
130*90a64276Sksagiyam     PetscCall(DMLabelSetValue(filter, 0, filterValue));
131*90a64276Sksagiyam     PetscCall(DMLabelSetValue(filter, 2, filterValue));
132*90a64276Sksagiyam     break;
133*90a64276Sksagiyam   case 2:
134*90a64276Sksagiyam     break;
135*90a64276Sksagiyam   }
136*90a64276Sksagiyam   PetscCall(PetscObjectSetName((PetscObject)dm, "Example_DM"));
137*90a64276Sksagiyam   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
138*90a64276Sksagiyam   PetscCall(DMPlexFilter(dm, filter, filterValue, user.ignoreLabelHalo, user.sanitizeSubmesh, &ownershipTransferSF, &subdm));
139*90a64276Sksagiyam   PetscCall(DMLabelDestroy(&filter));
140*90a64276Sksagiyam   PetscCall(DMDestroy(&dm));
141*90a64276Sksagiyam   PetscCall(PetscObjectSetName((PetscObject)subdm, "Example_SubDM"));
142*90a64276Sksagiyam   PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view"));
143*90a64276Sksagiyam   PetscCall(DMDestroy(&subdm));
144*90a64276Sksagiyam   PetscCall(PetscObjectSetName((PetscObject)ownershipTransferSF, "Example_Ownership_Transfer_SF"));
145*90a64276Sksagiyam   PetscCall(PetscSFView(ownershipTransferSF, PETSC_VIEWER_STDOUT_WORLD));
146*90a64276Sksagiyam   PetscCall(PetscSFDestroy(&ownershipTransferSF));
147*90a64276Sksagiyam   PetscCall(PetscFinalize());
148*90a64276Sksagiyam   return 0;
149*90a64276Sksagiyam }
150*90a64276Sksagiyam 
151*90a64276Sksagiyam /*TEST
152*90a64276Sksagiyam 
153*90a64276Sksagiyam   testset:
154*90a64276Sksagiyam     nsize: 3
155*90a64276Sksagiyam     args: -dm_view ascii::ascii_info_detail
156*90a64276Sksagiyam 
157*90a64276Sksagiyam     test:
158*90a64276Sksagiyam       suffix: 0
159*90a64276Sksagiyam       args:
160*90a64276Sksagiyam 
161*90a64276Sksagiyam     test:
162*90a64276Sksagiyam       suffix: 1
163*90a64276Sksagiyam       args: -sanitize_submesh
164*90a64276Sksagiyam 
165*90a64276Sksagiyam     test:
166*90a64276Sksagiyam       suffix: 2
167*90a64276Sksagiyam       args: -ignore_label_halo
168*90a64276Sksagiyam 
169*90a64276Sksagiyam TEST*/
170