xref: /petsc/src/dm/impls/plex/tests/ex12.c (revision e0008cae9d8b220d3a191ff57f84cb339f1acef2)
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, distribute it, and then redistribute it on 24 processes using two different partitioners:
13c4762a1bSJed Brown 
14c4762a1bSJed 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
15c4762a1bSJed Brown 
16c4762a1bSJed 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:
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 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
19c4762a1bSJed Brown 
20c4762a1bSJed Brown */
21c4762a1bSJed Brown 
229371c9d4SSatish Balay enum {
239371c9d4SSatish Balay   STAGE_LOAD,
249371c9d4SSatish Balay   STAGE_DISTRIBUTE,
259371c9d4SSatish Balay   STAGE_REFINE,
269371c9d4SSatish Balay   STAGE_REDISTRIBUTE
279371c9d4SSatish Balay };
28c4762a1bSJed Brown 
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   /* Domain and mesh definition */
31c4762a1bSJed Brown   PetscInt      overlap;          /* The cell overlap to use during partitioning */
32c4762a1bSJed Brown   PetscBool     testPartition;    /* Use a fixed partitioning for testing */
33c4762a1bSJed Brown   PetscBool     testRedundant;    /* Use a redundant partitioning for testing */
34c4762a1bSJed Brown   PetscBool     loadBalance;      /* Load balance via a second distribute step */
35c4762a1bSJed Brown   PetscBool     partitionBalance; /* Balance shared point partition */
36c4762a1bSJed Brown   PetscLogStage stages[4];
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
39d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
40d71ae5a4SJacob Faibussowitsch {
41c4762a1bSJed Brown   PetscFunctionBegin;
42c4762a1bSJed Brown   options->overlap          = 0;
43c4762a1bSJed Brown   options->testPartition    = PETSC_FALSE;
44c4762a1bSJed Brown   options->testRedundant    = PETSC_FALSE;
45c4762a1bSJed Brown   options->loadBalance      = PETSC_FALSE;
46c4762a1bSJed Brown   options->partitionBalance = PETSC_FALSE;
47c4762a1bSJed Brown 
48d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL, 0));
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL));
539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL));
54d0609cedSBarry Smith   PetscOptionsEnd();
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]));
579566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
589566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]));
599566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]));
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
63d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
64d71ae5a4SJacob Faibussowitsch {
65c4762a1bSJed Brown   DM          pdm             = NULL;
66c4762a1bSJed Brown   PetscInt    triSizes_n2[2]  = {4, 4};
67c4762a1bSJed Brown   PetscInt    triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7};
68c4762a1bSJed Brown   PetscInt    triSizes_n3[3]  = {3, 2, 3};
69c4762a1bSJed Brown   PetscInt    triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4};
70c4762a1bSJed Brown   PetscInt    triSizes_n4[4]  = {2, 2, 2, 2};
71c4762a1bSJed Brown   PetscInt    triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6};
72c4762a1bSJed Brown   PetscInt    triSizes_n8[8]  = {1, 1, 1, 1, 1, 1, 1, 1};
73c4762a1bSJed Brown   PetscInt    triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7};
74c4762a1bSJed Brown   PetscInt    quadSizes[2]    = {2, 2};
75c4762a1bSJed Brown   PetscInt    quadPoints[4]   = {2, 3, 0, 1};
76c4762a1bSJed Brown   PetscInt    overlap         = user->overlap >= 0 ? user->overlap : 0;
7730602db0SMatthew G. Knepley   PetscInt    dim;
7830602db0SMatthew G. Knepley   PetscBool   simplex;
79c4762a1bSJed Brown   PetscMPIInt rank, size;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   PetscFunctionBegin;
829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
849566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD]));
859566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
869566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
879566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
889566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
899566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
909566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
919566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &dim));
929566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(*dm, &simplex));
939566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
94c4762a1bSJed Brown   if (!user->testRedundant) {
95c4762a1bSJed Brown     PetscPartitioner part;
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
989566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
999566063dSJacob Faibussowitsch     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
100c4762a1bSJed Brown     if (user->testPartition) {
101c4762a1bSJed Brown       const PetscInt *sizes  = NULL;
102c4762a1bSJed Brown       const PetscInt *points = NULL;
103c4762a1bSJed Brown 
104dd400576SPatrick Sanan       if (rank == 0) {
10530602db0SMatthew G. Knepley         if (dim == 2 && simplex && size == 2) {
1069371c9d4SSatish Balay           sizes  = triSizes_n2;
1079371c9d4SSatish Balay           points = triPoints_n2;
10830602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 3) {
1099371c9d4SSatish Balay           sizes  = triSizes_n3;
1109371c9d4SSatish Balay           points = triPoints_n3;
11130602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 4) {
1129371c9d4SSatish Balay           sizes  = triSizes_n4;
1139371c9d4SSatish Balay           points = triPoints_n4;
11430602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 8) {
1159371c9d4SSatish Balay           sizes  = triSizes_n8;
1169371c9d4SSatish Balay           points = triPoints_n8;
11730602db0SMatthew G. Knepley         } else if (dim == 2 && !simplex && size == 2) {
1189371c9d4SSatish Balay           sizes  = quadSizes;
1199371c9d4SSatish Balay           points = quadPoints;
120c4762a1bSJed Brown         }
121c4762a1bSJed Brown       }
1229566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1239566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
124c4762a1bSJed Brown     }
1259566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
126c4762a1bSJed Brown   } else {
127c4762a1bSJed Brown     PetscSF sf;
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch     PetscCall(DMPlexGetRedundantDM(*dm, &sf, &pdm));
130c4762a1bSJed Brown     if (sf) {
131c4762a1bSJed Brown       DM test;
132c4762a1bSJed Brown 
1339566063dSJacob Faibussowitsch       PetscCall(DMPlexCreate(comm, &test));
1349566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh"));
1359566063dSJacob Faibussowitsch       PetscCall(DMPlexMigrate(*dm, sf, test));
1369566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view"));
1379566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&test));
138c4762a1bSJed Brown     }
1399566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
140c4762a1bSJed Brown   }
141c4762a1bSJed Brown   if (pdm) {
1429566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
143c4762a1bSJed Brown     *dm = pdm;
144c4762a1bSJed Brown   }
1459566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1469566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
147c4762a1bSJed Brown   if (user->loadBalance) {
148c4762a1bSJed Brown     PetscPartitioner part;
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-prelb_dm_view"));
1519566063dSJacob Faibussowitsch     PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_"));
1529566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]));
1539566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
1549566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "lb_"));
1559566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
156c4762a1bSJed Brown     if (user->testPartition) {
157c4762a1bSJed Brown       PetscInt reSizes_n2[2]  = {2, 2};
158c4762a1bSJed Brown       PetscInt rePoints_n2[4] = {2, 3, 0, 1};
1599371c9d4SSatish Balay       if (rank) {
1609371c9d4SSatish Balay         rePoints_n2[0] = 1;
1619371c9d4SSatish Balay         rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;
1629371c9d4SSatish Balay       }
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1659566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2));
166c4762a1bSJed Brown     }
1679566063dSJacob Faibussowitsch     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
1689566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
169c4762a1bSJed Brown     if (pdm) {
1709566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
171c4762a1bSJed Brown       *dm = pdm;
172c4762a1bSJed Brown     }
1739566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
174c4762a1bSJed Brown   }
1759566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE]));
1769566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1779566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
182d71ae5a4SJacob Faibussowitsch {
183c4762a1bSJed Brown   DM     dm;
184c4762a1bSJed Brown   AppCtx user; /* user-defined work context */
185c4762a1bSJed Brown 
186327415f7SBarry Smith   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1889566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1899566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1909566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1919566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
192b122ec5aSJacob Faibussowitsch   return 0;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown /*TEST
196c4762a1bSJed Brown   # Parallel, no overlap tests 0-2
197c4762a1bSJed Brown   test:
198c4762a1bSJed Brown     suffix: 0
199c4762a1bSJed Brown     requires: triangle
20030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex
201*e0008caeSPierre Jolivet     output_file: output/empty.out
202c4762a1bSJed Brown   test:
203c4762a1bSJed Brown     suffix: 1
204c4762a1bSJed Brown     requires: triangle
205c4762a1bSJed Brown     nsize: 3
20630602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
207c4762a1bSJed Brown   test:
208c4762a1bSJed Brown     suffix: 2
209c4762a1bSJed Brown     requires: triangle
210c4762a1bSJed Brown     nsize: 8
21130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
212c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
213c4762a1bSJed Brown   test:
214c4762a1bSJed Brown     suffix: 3
215c4762a1bSJed Brown     requires: triangle
216c4762a1bSJed Brown     nsize: 3
21730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
218c4762a1bSJed Brown   test:
219c4762a1bSJed Brown     suffix: 4
220c4762a1bSJed Brown     requires: triangle
221c4762a1bSJed Brown     nsize: 8
22230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
223c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
224c4762a1bSJed Brown   test:
225c4762a1bSJed Brown     suffix: 5
226c4762a1bSJed Brown     requires: triangle
227c4762a1bSJed Brown     nsize: 8
22830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail
229c4762a1bSJed Brown   # Parallel load balancing, test 6-7
230c4762a1bSJed Brown   test:
231c4762a1bSJed Brown     suffix: 6
232c4762a1bSJed Brown     requires: triangle
233c4762a1bSJed Brown     nsize: 2
23430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
235c4762a1bSJed Brown   test:
236c4762a1bSJed Brown     suffix: 7
237c4762a1bSJed Brown     requires: triangle
238c4762a1bSJed Brown     nsize: 2
23930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail
240c4762a1bSJed Brown   # Parallel redundant copying, test 8
241c4762a1bSJed Brown   test:
242c4762a1bSJed Brown     suffix: 8
243c4762a1bSJed Brown     requires: triangle
244c4762a1bSJed Brown     nsize: 2
24530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail
246c4762a1bSJed Brown   test:
247c4762a1bSJed Brown     suffix: lb_0
248c4762a1bSJed Brown     requires: parmetis
249c4762a1bSJed Brown     nsize: 4
25030602db0SMatthew 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
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   # Same tests as above, but with balancing of the shared point partition
253c4762a1bSJed Brown   test:
254c4762a1bSJed Brown     suffix: 9
255c4762a1bSJed Brown     requires: triangle
25630602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance
257*e0008caeSPierre Jolivet     output_file: output/empty.out
258c4762a1bSJed Brown   test:
259c4762a1bSJed Brown     suffix: 10
260c4762a1bSJed Brown     requires: triangle
261c4762a1bSJed Brown     nsize: 3
26230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
263c4762a1bSJed Brown   test:
264c4762a1bSJed Brown     suffix: 11
265c4762a1bSJed Brown     requires: triangle
266c4762a1bSJed Brown     nsize: 8
26730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
268c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
269c4762a1bSJed Brown   test:
270c4762a1bSJed Brown     suffix: 12
271c4762a1bSJed Brown     requires: triangle
272c4762a1bSJed Brown     nsize: 3
27330602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
274c4762a1bSJed Brown   test:
275c4762a1bSJed Brown     suffix: 13
276c4762a1bSJed Brown     requires: triangle
277c4762a1bSJed Brown     nsize: 8
27830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
279c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
280c4762a1bSJed Brown   test:
281c4762a1bSJed Brown     suffix: 14
282c4762a1bSJed Brown     requires: triangle
283c4762a1bSJed Brown     nsize: 8
28430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance
285c4762a1bSJed Brown   # Parallel load balancing, test 6-7
286c4762a1bSJed Brown   test:
287c4762a1bSJed Brown     suffix: 15
288c4762a1bSJed Brown     requires: triangle
289c4762a1bSJed Brown     nsize: 2
29030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
291c4762a1bSJed Brown   test:
292c4762a1bSJed Brown     suffix: 16
293c4762a1bSJed Brown     requires: triangle
294c4762a1bSJed Brown     nsize: 2
29530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance
296c4762a1bSJed Brown   # Parallel redundant copying, test 8
297c4762a1bSJed Brown   test:
298c4762a1bSJed Brown     suffix: 17
299c4762a1bSJed Brown     requires: triangle
300c4762a1bSJed Brown     nsize: 2
30130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance
302c4762a1bSJed Brown   test:
303c4762a1bSJed Brown     suffix: lb_1
304c4762a1bSJed Brown     requires: parmetis
305c4762a1bSJed Brown     nsize: 4
30630602db0SMatthew 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
307c4762a1bSJed Brown TEST*/
308