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