xref: /petsc/src/dm/impls/plex/tests/ex12.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode ierr;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBegin;
43c4762a1bSJed Brown   options->overlap          = 0;
44c4762a1bSJed Brown   options->testPartition    = PETSC_FALSE;
45c4762a1bSJed Brown   options->testRedundant    = PETSC_FALSE;
46c4762a1bSJed Brown   options->loadBalance      = PETSC_FALSE;
47c4762a1bSJed Brown   options->partitionBalance = PETSC_FALSE;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL,0));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL));
551e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
56c4762a1bSJed Brown 
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshLoad",         &options->stages[STAGE_LOAD]));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshDistribute",   &options->stages[STAGE_DISTRIBUTE]));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshRefine",       &options->stages[STAGE_REFINE]));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]));
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   DM             pdm             = NULL;
67c4762a1bSJed Brown   PetscInt       triSizes_n2[2]  = {4, 4};
68c4762a1bSJed Brown   PetscInt       triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7};
69c4762a1bSJed Brown   PetscInt       triSizes_n3[3]  = {3, 2, 3};
70c4762a1bSJed Brown   PetscInt       triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4};
71c4762a1bSJed Brown   PetscInt       triSizes_n4[4]  = {2, 2, 2, 2};
72c4762a1bSJed Brown   PetscInt       triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6};
73c4762a1bSJed Brown   PetscInt       triSizes_n8[8]  = {1, 1, 1, 1, 1, 1, 1, 1};
74c4762a1bSJed Brown   PetscInt       triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7};
75c4762a1bSJed Brown   PetscInt       quadSizes[2]    = {2, 2};
76c4762a1bSJed Brown   PetscInt       quadPoints[4]   = {2, 3, 0, 1};
77c4762a1bSJed Brown   PetscInt       overlap         = user->overlap >= 0 ? user->overlap : 0;
7830602db0SMatthew G. Knepley   PetscInt       dim;
7930602db0SMatthew G. Knepley   PetscBool      simplex;
80c4762a1bSJed Brown   PetscMPIInt    rank, size;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBegin;
835f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
845f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_LOAD]));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(*dm, &dim));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(*dm, &simplex));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
95c4762a1bSJed Brown   if (!user->testRedundant) {
96c4762a1bSJed Brown     PetscPartitioner part;
97c4762a1bSJed Brown 
985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPartitioner(*dm, &part));
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerSetFromOptions(part));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
101c4762a1bSJed Brown     if (user->testPartition) {
102c4762a1bSJed Brown       const PetscInt *sizes = NULL;
103c4762a1bSJed Brown       const PetscInt *points = NULL;
104c4762a1bSJed Brown 
105dd400576SPatrick Sanan       if (rank == 0) {
10630602db0SMatthew G. Knepley         if (dim == 2 && simplex && size == 2) {
107c4762a1bSJed Brown           sizes = triSizes_n2; points = triPoints_n2;
10830602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 3) {
109c4762a1bSJed Brown           sizes = triSizes_n3; points = triPoints_n3;
11030602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 4) {
111c4762a1bSJed Brown           sizes = triSizes_n4; points = triPoints_n4;
11230602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 8) {
113c4762a1bSJed Brown           sizes = triSizes_n8; points = triPoints_n8;
11430602db0SMatthew G. Knepley         } else if (dim == 2 && !simplex && size == 2) {
115c4762a1bSJed Brown           sizes = quadSizes; points = quadPoints;
116c4762a1bSJed Brown         }
117c4762a1bSJed Brown       }
1185f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPartitionerShellSetPartition(part, size, sizes, points));
120c4762a1bSJed Brown     }
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistribute(*dm, overlap, NULL, &pdm));
122c4762a1bSJed Brown   } else {
123c4762a1bSJed Brown     PetscSF sf;
124c4762a1bSJed Brown 
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetRedundantDM(*dm, &sf, &pdm));
126c4762a1bSJed Brown     if (sf) {
127c4762a1bSJed Brown       DM test;
128c4762a1bSJed Brown 
1295f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCreate(comm,&test));
1305f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh"));
1315f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexMigrate(*dm, sf, test));
1325f80ce2aSJacob Faibussowitsch       CHKERRQ(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view"));
1335f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&test));
134c4762a1bSJed Brown     }
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFDestroy(&sf));
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown   if (pdm) {
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(dm));
139c4762a1bSJed Brown     *dm  = pdm;
140c4762a1bSJed Brown   }
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
143c4762a1bSJed Brown   if (user->loadBalance) {
144c4762a1bSJed Brown     PetscPartitioner part;
145c4762a1bSJed Brown 
1465f80ce2aSJacob Faibussowitsch     CHKERRQ(DMViewFromOptions(*dm, NULL, "-prelb_dm_view"));
1475f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexSetOptionsPrefix(*dm, "lb_"));
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]));
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPartitioner(*dm, &part));
1505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) part, "lb_"));
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerSetFromOptions(part));
152c4762a1bSJed Brown     if (user->testPartition) {
153c4762a1bSJed Brown       PetscInt         reSizes_n2[2]  = {2, 2};
154c4762a1bSJed Brown       PetscInt         rePoints_n2[4] = {2, 3, 0, 1};
155c4762a1bSJed Brown       if (rank) {rePoints_n2[0] = 1; rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;}
156c4762a1bSJed Brown 
1575f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1585f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2));
159c4762a1bSJed Brown     }
1605f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
1615f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistribute(*dm, overlap, NULL, &pdm));
162c4762a1bSJed Brown     if (pdm) {
1635f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(dm));
164c4762a1bSJed Brown       *dm  = pdm;
165c4762a1bSJed Brown     }
1665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
167c4762a1bSJed Brown   }
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_REFINE]));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
171c4762a1bSJed Brown   PetscFunctionReturn(0);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown int main(int argc, char **argv)
175c4762a1bSJed Brown {
176c4762a1bSJed Brown   DM             dm;
177c4762a1bSJed Brown   AppCtx         user; /* user-defined work context */
178c4762a1bSJed Brown 
179*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
183*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
184*b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*TEST
188c4762a1bSJed Brown   # Parallel, no overlap tests 0-2
189c4762a1bSJed Brown   test:
190c4762a1bSJed Brown     suffix: 0
191c4762a1bSJed Brown     requires: triangle
19230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex
193c4762a1bSJed Brown   test:
194c4762a1bSJed Brown     suffix: 1
195c4762a1bSJed Brown     requires: triangle
196c4762a1bSJed Brown     nsize: 3
19730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
198c4762a1bSJed Brown   test:
199c4762a1bSJed Brown     suffix: 2
200c4762a1bSJed Brown     requires: triangle
201c4762a1bSJed Brown     nsize: 8
20230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
203c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
204c4762a1bSJed Brown   test:
205c4762a1bSJed Brown     suffix: 3
206c4762a1bSJed Brown     requires: triangle
207c4762a1bSJed Brown     nsize: 3
20830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
209c4762a1bSJed Brown   test:
210c4762a1bSJed Brown     suffix: 4
211c4762a1bSJed Brown     requires: triangle
212c4762a1bSJed Brown     nsize: 8
21330602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
214c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
215c4762a1bSJed Brown   test:
216c4762a1bSJed Brown     suffix: 5
217c4762a1bSJed Brown     requires: triangle
218c4762a1bSJed Brown     nsize: 8
21930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail
220c4762a1bSJed Brown   # Parallel load balancing, test 6-7
221c4762a1bSJed Brown   test:
222c4762a1bSJed Brown     suffix: 6
223c4762a1bSJed Brown     requires: triangle
224c4762a1bSJed Brown     nsize: 2
22530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
226c4762a1bSJed Brown   test:
227c4762a1bSJed Brown     suffix: 7
228c4762a1bSJed Brown     requires: triangle
229c4762a1bSJed Brown     nsize: 2
23030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail
231c4762a1bSJed Brown   # Parallel redundant copying, test 8
232c4762a1bSJed Brown   test:
233c4762a1bSJed Brown     suffix: 8
234c4762a1bSJed Brown     requires: triangle
235c4762a1bSJed Brown     nsize: 2
23630602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail
237c4762a1bSJed Brown   test:
238c4762a1bSJed Brown     suffix: lb_0
239c4762a1bSJed Brown     requires: parmetis
240c4762a1bSJed Brown     nsize: 4
24130602db0SMatthew 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
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   # Same tests as above, but with balancing of the shared point partition
244c4762a1bSJed Brown   test:
245c4762a1bSJed Brown     suffix: 9
246c4762a1bSJed Brown     requires: triangle
24730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance
248c4762a1bSJed Brown   test:
249c4762a1bSJed Brown     suffix: 10
250c4762a1bSJed Brown     requires: triangle
251c4762a1bSJed Brown     nsize: 3
25230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
253c4762a1bSJed Brown   test:
254c4762a1bSJed Brown     suffix: 11
255c4762a1bSJed Brown     requires: triangle
256c4762a1bSJed Brown     nsize: 8
25730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
258c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
259c4762a1bSJed Brown   test:
260c4762a1bSJed Brown     suffix: 12
261c4762a1bSJed Brown     requires: triangle
262c4762a1bSJed Brown     nsize: 3
26330602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
264c4762a1bSJed Brown   test:
265c4762a1bSJed Brown     suffix: 13
266c4762a1bSJed Brown     requires: triangle
267c4762a1bSJed Brown     nsize: 8
26830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
269c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
270c4762a1bSJed Brown   test:
271c4762a1bSJed Brown     suffix: 14
272c4762a1bSJed Brown     requires: triangle
273c4762a1bSJed Brown     nsize: 8
27430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance
275c4762a1bSJed Brown   # Parallel load balancing, test 6-7
276c4762a1bSJed Brown   test:
277c4762a1bSJed Brown     suffix: 15
278c4762a1bSJed Brown     requires: triangle
279c4762a1bSJed Brown     nsize: 2
28030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
281c4762a1bSJed Brown   test:
282c4762a1bSJed Brown     suffix: 16
283c4762a1bSJed Brown     requires: triangle
284c4762a1bSJed Brown     nsize: 2
28530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance
286c4762a1bSJed Brown   # Parallel redundant copying, test 8
287c4762a1bSJed Brown   test:
288c4762a1bSJed Brown     suffix: 17
289c4762a1bSJed Brown     requires: triangle
290c4762a1bSJed Brown     nsize: 2
29130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance
292c4762a1bSJed Brown   test:
293c4762a1bSJed Brown     suffix: lb_1
294c4762a1bSJed Brown     requires: parmetis
295c4762a1bSJed Brown     nsize: 4
29630602db0SMatthew 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
297c4762a1bSJed Brown TEST*/
298