xref: /petsc/src/dm/impls/plex/tests/ex12.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown static char help[] = "Partition a mesh in parallel, perhaps with overlap\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscsf.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /* Sample usage:
7c4762a1bSJed Brown 
8c4762a1bSJed Brown Load a file in serial and distribute it on 24 processes:
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   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
11c4762a1bSJed Brown 
12c4762a1bSJed Brown Load a file in serial and distribute it on 24 processes using a custom partitioner:
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   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
15c4762a1bSJed Brown 
16c4762a1bSJed Brown Load a file in serial, distribute it, and then redistribute it on 24 processes using two different partitioners:
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   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
19c4762a1bSJed Brown 
20c4762a1bSJed Brown 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:
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   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
23c4762a1bSJed Brown 
24c4762a1bSJed Brown */
25c4762a1bSJed Brown 
26c4762a1bSJed Brown enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_REDISTRIBUTE};
27c4762a1bSJed Brown 
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown   /* Domain and mesh definition */
30c4762a1bSJed Brown   PetscInt  overlap;                      /* The cell overlap to use during partitioning */
31c4762a1bSJed Brown   PetscBool testPartition;                /* Use a fixed partitioning for testing */
32c4762a1bSJed Brown   PetscBool testRedundant;                /* Use a redundant partitioning for testing */
33c4762a1bSJed Brown   PetscBool loadBalance;                  /* Load balance via a second distribute step */
34c4762a1bSJed Brown   PetscBool partitionBalance;             /* Balance shared point partition */
35c4762a1bSJed Brown   PetscLogStage stages[4];
36c4762a1bSJed Brown } AppCtx;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   PetscFunctionBegin;
41c4762a1bSJed Brown   options->overlap          = 0;
42c4762a1bSJed Brown   options->testPartition    = PETSC_FALSE;
43c4762a1bSJed Brown   options->testRedundant    = PETSC_FALSE;
44c4762a1bSJed Brown   options->loadBalance      = PETSC_FALSE;
45c4762a1bSJed Brown   options->partitionBalance = PETSC_FALSE;
46c4762a1bSJed Brown 
47*d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL,0));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL));
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL));
53*d0609cedSBarry Smith   PetscOptionsEnd();
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshLoad",         &options->stages[STAGE_LOAD]));
569566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshDistribute",   &options->stages[STAGE_DISTRIBUTE]));
579566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRefine",       &options->stages[STAGE_REFINE]));
589566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]));
59c4762a1bSJed Brown   PetscFunctionReturn(0);
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
63c4762a1bSJed Brown {
64c4762a1bSJed Brown   DM             pdm             = NULL;
65c4762a1bSJed Brown   PetscInt       triSizes_n2[2]  = {4, 4};
66c4762a1bSJed Brown   PetscInt       triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7};
67c4762a1bSJed Brown   PetscInt       triSizes_n3[3]  = {3, 2, 3};
68c4762a1bSJed Brown   PetscInt       triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4};
69c4762a1bSJed Brown   PetscInt       triSizes_n4[4]  = {2, 2, 2, 2};
70c4762a1bSJed Brown   PetscInt       triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6};
71c4762a1bSJed Brown   PetscInt       triSizes_n8[8]  = {1, 1, 1, 1, 1, 1, 1, 1};
72c4762a1bSJed Brown   PetscInt       triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7};
73c4762a1bSJed Brown   PetscInt       quadSizes[2]    = {2, 2};
74c4762a1bSJed Brown   PetscInt       quadPoints[4]   = {2, 3, 0, 1};
75c4762a1bSJed Brown   PetscInt       overlap         = user->overlap >= 0 ? user->overlap : 0;
7630602db0SMatthew G. Knepley   PetscInt       dim;
7730602db0SMatthew G. Knepley   PetscBool      simplex;
78c4762a1bSJed Brown   PetscMPIInt    rank, size;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   PetscFunctionBegin;
819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
839566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD]));
849566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
859566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
869566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
879566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
889566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
899566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
909566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &dim));
919566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(*dm, &simplex));
929566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
93c4762a1bSJed Brown   if (!user->testRedundant) {
94c4762a1bSJed Brown     PetscPartitioner part;
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
979566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
989566063dSJacob Faibussowitsch     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
99c4762a1bSJed Brown     if (user->testPartition) {
100c4762a1bSJed Brown       const PetscInt *sizes = NULL;
101c4762a1bSJed Brown       const PetscInt *points = NULL;
102c4762a1bSJed Brown 
103dd400576SPatrick Sanan       if (rank == 0) {
10430602db0SMatthew G. Knepley         if (dim == 2 && simplex && size == 2) {
105c4762a1bSJed Brown           sizes = triSizes_n2; points = triPoints_n2;
10630602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 3) {
107c4762a1bSJed Brown           sizes = triSizes_n3; points = triPoints_n3;
10830602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 4) {
109c4762a1bSJed Brown           sizes = triSizes_n4; points = triPoints_n4;
11030602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 8) {
111c4762a1bSJed Brown           sizes = triSizes_n8; points = triPoints_n8;
11230602db0SMatthew G. Knepley         } else if (dim == 2 && !simplex && size == 2) {
113c4762a1bSJed Brown           sizes = quadSizes; points = quadPoints;
114c4762a1bSJed Brown         }
115c4762a1bSJed Brown       }
1169566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1179566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
118c4762a1bSJed Brown     }
1199566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
120c4762a1bSJed Brown   } else {
121c4762a1bSJed Brown     PetscSF sf;
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetRedundantDM(*dm, &sf, &pdm));
124c4762a1bSJed Brown     if (sf) {
125c4762a1bSJed Brown       DM test;
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch       PetscCall(DMPlexCreate(comm,&test));
1289566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh"));
1299566063dSJacob Faibussowitsch       PetscCall(DMPlexMigrate(*dm, sf, test));
1309566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view"));
1319566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&test));
132c4762a1bSJed Brown     }
1339566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown   if (pdm) {
1369566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
137c4762a1bSJed Brown     *dm  = pdm;
138c4762a1bSJed Brown   }
1399566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1409566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
141c4762a1bSJed Brown   if (user->loadBalance) {
142c4762a1bSJed Brown     PetscPartitioner part;
143c4762a1bSJed Brown 
1449566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-prelb_dm_view"));
1459566063dSJacob Faibussowitsch     PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_"));
1469566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]));
1479566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
1489566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) part, "lb_"));
1499566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
150c4762a1bSJed Brown     if (user->testPartition) {
151c4762a1bSJed Brown       PetscInt         reSizes_n2[2]  = {2, 2};
152c4762a1bSJed Brown       PetscInt         rePoints_n2[4] = {2, 3, 0, 1};
153c4762a1bSJed Brown       if (rank) {rePoints_n2[0] = 1; rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;}
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1569566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2));
157c4762a1bSJed Brown     }
1589566063dSJacob Faibussowitsch     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
1599566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
160c4762a1bSJed Brown     if (pdm) {
1619566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
162c4762a1bSJed Brown       *dm  = pdm;
163c4762a1bSJed Brown     }
1649566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
165c4762a1bSJed Brown   }
1669566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE]));
1679566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1689566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
169c4762a1bSJed Brown   PetscFunctionReturn(0);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown int main(int argc, char **argv)
173c4762a1bSJed Brown {
174c4762a1bSJed Brown   DM             dm;
175c4762a1bSJed Brown   AppCtx         user; /* user-defined work context */
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1789566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1799566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
182b122ec5aSJacob Faibussowitsch   return 0;
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /*TEST
186c4762a1bSJed Brown   # Parallel, no overlap tests 0-2
187c4762a1bSJed Brown   test:
188c4762a1bSJed Brown     suffix: 0
189c4762a1bSJed Brown     requires: triangle
19030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex
191c4762a1bSJed Brown   test:
192c4762a1bSJed Brown     suffix: 1
193c4762a1bSJed Brown     requires: triangle
194c4762a1bSJed Brown     nsize: 3
19530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
196c4762a1bSJed Brown   test:
197c4762a1bSJed Brown     suffix: 2
198c4762a1bSJed Brown     requires: triangle
199c4762a1bSJed Brown     nsize: 8
20030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
201c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
202c4762a1bSJed Brown   test:
203c4762a1bSJed Brown     suffix: 3
204c4762a1bSJed Brown     requires: triangle
205c4762a1bSJed Brown     nsize: 3
20630602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
207c4762a1bSJed Brown   test:
208c4762a1bSJed Brown     suffix: 4
209c4762a1bSJed Brown     requires: triangle
210c4762a1bSJed Brown     nsize: 8
21130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
212c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
213c4762a1bSJed Brown   test:
214c4762a1bSJed Brown     suffix: 5
215c4762a1bSJed Brown     requires: triangle
216c4762a1bSJed Brown     nsize: 8
21730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail
218c4762a1bSJed Brown   # Parallel load balancing, test 6-7
219c4762a1bSJed Brown   test:
220c4762a1bSJed Brown     suffix: 6
221c4762a1bSJed Brown     requires: triangle
222c4762a1bSJed Brown     nsize: 2
22330602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
224c4762a1bSJed Brown   test:
225c4762a1bSJed Brown     suffix: 7
226c4762a1bSJed Brown     requires: triangle
227c4762a1bSJed Brown     nsize: 2
22830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail
229c4762a1bSJed Brown   # Parallel redundant copying, test 8
230c4762a1bSJed Brown   test:
231c4762a1bSJed Brown     suffix: 8
232c4762a1bSJed Brown     requires: triangle
233c4762a1bSJed Brown     nsize: 2
23430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail
235c4762a1bSJed Brown   test:
236c4762a1bSJed Brown     suffix: lb_0
237c4762a1bSJed Brown     requires: parmetis
238c4762a1bSJed Brown     nsize: 4
23930602db0SMatthew G. Knepley     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
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   # Same tests as above, but with balancing of the shared point partition
242c4762a1bSJed Brown   test:
243c4762a1bSJed Brown     suffix: 9
244c4762a1bSJed Brown     requires: triangle
24530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance
246c4762a1bSJed Brown   test:
247c4762a1bSJed Brown     suffix: 10
248c4762a1bSJed Brown     requires: triangle
249c4762a1bSJed Brown     nsize: 3
25030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
251c4762a1bSJed Brown   test:
252c4762a1bSJed Brown     suffix: 11
253c4762a1bSJed Brown     requires: triangle
254c4762a1bSJed Brown     nsize: 8
25530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
256c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
257c4762a1bSJed Brown   test:
258c4762a1bSJed Brown     suffix: 12
259c4762a1bSJed Brown     requires: triangle
260c4762a1bSJed Brown     nsize: 3
26130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
262c4762a1bSJed Brown   test:
263c4762a1bSJed Brown     suffix: 13
264c4762a1bSJed Brown     requires: triangle
265c4762a1bSJed Brown     nsize: 8
26630602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
267c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
268c4762a1bSJed Brown   test:
269c4762a1bSJed Brown     suffix: 14
270c4762a1bSJed Brown     requires: triangle
271c4762a1bSJed Brown     nsize: 8
27230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance
273c4762a1bSJed Brown   # Parallel load balancing, test 6-7
274c4762a1bSJed Brown   test:
275c4762a1bSJed Brown     suffix: 15
276c4762a1bSJed Brown     requires: triangle
277c4762a1bSJed Brown     nsize: 2
27830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
279c4762a1bSJed Brown   test:
280c4762a1bSJed Brown     suffix: 16
281c4762a1bSJed Brown     requires: triangle
282c4762a1bSJed Brown     nsize: 2
28330602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance
284c4762a1bSJed Brown   # Parallel redundant copying, test 8
285c4762a1bSJed Brown   test:
286c4762a1bSJed Brown     suffix: 17
287c4762a1bSJed Brown     requires: triangle
288c4762a1bSJed Brown     nsize: 2
28930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance
290c4762a1bSJed Brown   test:
291c4762a1bSJed Brown     suffix: lb_1
292c4762a1bSJed Brown     requires: parmetis
293c4762a1bSJed Brown     nsize: 4
29430602db0SMatthew G. Knepley     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
295c4762a1bSJed Brown TEST*/
296