xref: /petsc/src/dm/impls/plex/tests/ex12.c (revision 1b37a2a7cc4a4fb30c3e967db1c694c0a1013f51)
1 static char help[] = "Partition a mesh in parallel, perhaps with overlap\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscsf.h>
5 
6 /* Sample usage:
7 
8 Load a file in serial and distribute it on 24 processes:
9 
10   make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/squaremotor-30.exo -orig_dm_view -dm_view" NP=24
11 
12 Load a file in serial, distribute it, and then redistribute it on 24 processes using two different partitioners:
13 
14   make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/squaremotor-30.exo -petscpartitioner_type simple -load_balance -lb_petscpartitioner_type parmetis -orig_dm_view -dm_view" NP=24
15 
16 Load a file in serial, distribute it randomly, refine it in parallel, and then redistribute it on 24 processes using two different partitioners, and view to VTK:
17 
18   make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/squaremotor-30.exo -petscpartitioner_type shell -petscpartitioner_shell_random -dm_refine 1 -load_balance -lb_petscpartitioner_type parmetis -prelb_dm_view vtk:$PWD/prelb.vtk -dm_view vtk:$PWD/balance.vtk -dm_partition_view" NP=24
19 
20 */
21 
22 enum {
23   STAGE_LOAD,
24   STAGE_DISTRIBUTE,
25   STAGE_REFINE,
26   STAGE_REDISTRIBUTE
27 };
28 
29 typedef struct {
30   /* Domain and mesh definition */
31   PetscInt      overlap;          /* The cell overlap to use during partitioning */
32   PetscBool     testPartition;    /* Use a fixed partitioning for testing */
33   PetscBool     testRedundant;    /* Use a redundant partitioning for testing */
34   PetscBool     loadBalance;      /* Load balance via a second distribute step */
35   PetscBool     partitionBalance; /* Balance shared point partition */
36   PetscLogStage stages[4];
37 } AppCtx;
38 
39 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
40 {
41   PetscFunctionBegin;
42   options->overlap          = 0;
43   options->testPartition    = PETSC_FALSE;
44   options->testRedundant    = PETSC_FALSE;
45   options->loadBalance      = PETSC_FALSE;
46   options->partitionBalance = PETSC_FALSE;
47 
48   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
49   PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL, 0));
50   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL));
51   PetscCall(PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL));
52   PetscCall(PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL));
53   PetscCall(PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL));
54   PetscOptionsEnd();
55 
56   PetscCall(PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]));
57   PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
58   PetscCall(PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]));
59   PetscCall(PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]));
60   PetscFunctionReturn(PETSC_SUCCESS);
61 }
62 
63 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
64 {
65   DM          pdm             = NULL;
66   PetscInt    triSizes_n2[2]  = {4, 4};
67   PetscInt    triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7};
68   PetscInt    triSizes_n3[3]  = {3, 2, 3};
69   PetscInt    triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4};
70   PetscInt    triSizes_n4[4]  = {2, 2, 2, 2};
71   PetscInt    triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6};
72   PetscInt    triSizes_n8[8]  = {1, 1, 1, 1, 1, 1, 1, 1};
73   PetscInt    triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7};
74   PetscInt    quadSizes[2]    = {2, 2};
75   PetscInt    quadPoints[4]   = {2, 3, 0, 1};
76   PetscInt    overlap         = user->overlap >= 0 ? user->overlap : 0;
77   PetscInt    dim;
78   PetscBool   simplex;
79   PetscMPIInt rank, size;
80 
81   PetscFunctionBegin;
82   PetscCallMPI(MPI_Comm_rank(comm, &rank));
83   PetscCallMPI(MPI_Comm_size(comm, &size));
84   PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD]));
85   PetscCall(DMCreate(comm, dm));
86   PetscCall(DMSetType(*dm, DMPLEX));
87   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
88   PetscCall(DMSetFromOptions(*dm));
89   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
90   PetscCall(PetscLogStagePop());
91   PetscCall(DMGetDimension(*dm, &dim));
92   PetscCall(DMPlexIsSimplex(*dm, &simplex));
93   PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
94   if (!user->testRedundant) {
95     PetscPartitioner part;
96 
97     PetscCall(DMPlexGetPartitioner(*dm, &part));
98     PetscCall(PetscPartitionerSetFromOptions(part));
99     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
100     if (user->testPartition) {
101       const PetscInt *sizes  = NULL;
102       const PetscInt *points = NULL;
103 
104       if (rank == 0) {
105         if (dim == 2 && simplex && size == 2) {
106           sizes  = triSizes_n2;
107           points = triPoints_n2;
108         } else if (dim == 2 && simplex && size == 3) {
109           sizes  = triSizes_n3;
110           points = triPoints_n3;
111         } else if (dim == 2 && simplex && size == 4) {
112           sizes  = triSizes_n4;
113           points = triPoints_n4;
114         } else if (dim == 2 && simplex && size == 8) {
115           sizes  = triSizes_n8;
116           points = triPoints_n8;
117         } else if (dim == 2 && !simplex && size == 2) {
118           sizes  = quadSizes;
119           points = quadPoints;
120         }
121       }
122       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
123       PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
124     }
125     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
126   } else {
127     PetscSF sf;
128 
129     PetscCall(DMPlexGetRedundantDM(*dm, &sf, &pdm));
130     if (sf) {
131       DM test;
132 
133       PetscCall(DMPlexCreate(comm, &test));
134       PetscCall(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh"));
135       PetscCall(DMPlexMigrate(*dm, sf, test));
136       PetscCall(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view"));
137       PetscCall(DMDestroy(&test));
138     }
139     PetscCall(PetscSFDestroy(&sf));
140   }
141   if (pdm) {
142     PetscCall(DMDestroy(dm));
143     *dm = pdm;
144   }
145   PetscCall(PetscLogStagePop());
146   PetscCall(DMSetFromOptions(*dm));
147   if (user->loadBalance) {
148     PetscPartitioner part;
149 
150     PetscCall(DMViewFromOptions(*dm, NULL, "-prelb_dm_view"));
151     PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_"));
152     PetscCall(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]));
153     PetscCall(DMPlexGetPartitioner(*dm, &part));
154     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "lb_"));
155     PetscCall(PetscPartitionerSetFromOptions(part));
156     if (user->testPartition) {
157       PetscInt reSizes_n2[2]  = {2, 2};
158       PetscInt rePoints_n2[4] = {2, 3, 0, 1};
159       if (rank) {
160         rePoints_n2[0] = 1;
161         rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;
162       }
163 
164       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
165       PetscCall(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2));
166     }
167     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
168     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
169     if (pdm) {
170       PetscCall(DMDestroy(dm));
171       *dm = pdm;
172     }
173     PetscCall(PetscLogStagePop());
174   }
175   PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE]));
176   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
177   PetscCall(PetscLogStagePop());
178   PetscFunctionReturn(PETSC_SUCCESS);
179 }
180 
181 int main(int argc, char **argv)
182 {
183   DM     dm;
184   AppCtx user; /* user-defined work context */
185 
186   PetscFunctionBeginUser;
187   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
188   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
189   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
190   PetscCall(DMDestroy(&dm));
191   PetscCall(PetscFinalize());
192   return 0;
193 }
194 
195 /*TEST
196   # Parallel, no overlap tests 0-2
197   test:
198     suffix: 0
199     requires: triangle
200     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex
201   test:
202     suffix: 1
203     requires: triangle
204     nsize: 3
205     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
206   test:
207     suffix: 2
208     requires: triangle
209     nsize: 8
210     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
211   # Parallel, level-1 overlap tests 3-4
212   test:
213     suffix: 3
214     requires: triangle
215     nsize: 3
216     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
217   test:
218     suffix: 4
219     requires: triangle
220     nsize: 8
221     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
222   # Parallel, level-2 overlap test 5
223   test:
224     suffix: 5
225     requires: triangle
226     nsize: 8
227     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail
228   # Parallel load balancing, test 6-7
229   test:
230     suffix: 6
231     requires: triangle
232     nsize: 2
233     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
234   test:
235     suffix: 7
236     requires: triangle
237     nsize: 2
238     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail
239   # Parallel redundant copying, test 8
240   test:
241     suffix: 8
242     requires: triangle
243     nsize: 2
244     args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail
245   test:
246     suffix: lb_0
247     requires: parmetis
248     nsize: 4
249     args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type shell -petscpartitioner_shell_random -lb_petscpartitioner_type parmetis -load_balance -lb_petscpartitioner_view -prelb_dm_view ::load_balance -dm_view ::load_balance
250 
251   # Same tests as above, but with balancing of the shared point partition
252   test:
253     suffix: 9
254     requires: triangle
255     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance
256   test:
257     suffix: 10
258     requires: triangle
259     nsize: 3
260     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
261   test:
262     suffix: 11
263     requires: triangle
264     nsize: 8
265     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
266   # Parallel, level-1 overlap tests 3-4
267   test:
268     suffix: 12
269     requires: triangle
270     nsize: 3
271     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
272   test:
273     suffix: 13
274     requires: triangle
275     nsize: 8
276     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
277   # Parallel, level-2 overlap test 5
278   test:
279     suffix: 14
280     requires: triangle
281     nsize: 8
282     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance
283   # Parallel load balancing, test 6-7
284   test:
285     suffix: 15
286     requires: triangle
287     nsize: 2
288     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
289   test:
290     suffix: 16
291     requires: triangle
292     nsize: 2
293     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance
294   # Parallel redundant copying, test 8
295   test:
296     suffix: 17
297     requires: triangle
298     nsize: 2
299     args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance
300   test:
301     suffix: lb_1
302     requires: parmetis
303     nsize: 4
304     args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type shell -petscpartitioner_shell_random -lb_petscpartitioner_type parmetis -load_balance -lb_petscpartitioner_view -partition_balance -prelb_dm_view ::load_balance -dm_view ::load_balance
305 TEST*/
306