xref: /petsc/src/dm/impls/plex/tests/ex12.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
269371c9d4SSatish Balay enum {
279371c9d4SSatish Balay   STAGE_LOAD,
289371c9d4SSatish Balay   STAGE_DISTRIBUTE,
299371c9d4SSatish Balay   STAGE_REFINE,
309371c9d4SSatish Balay   STAGE_REDISTRIBUTE
319371c9d4SSatish Balay };
32c4762a1bSJed Brown 
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown   /* Domain and mesh definition */
35c4762a1bSJed Brown   PetscInt      overlap;          /* The cell overlap to use during partitioning */
36c4762a1bSJed Brown   PetscBool     testPartition;    /* Use a fixed partitioning for testing */
37c4762a1bSJed Brown   PetscBool     testRedundant;    /* Use a redundant partitioning for testing */
38c4762a1bSJed Brown   PetscBool     loadBalance;      /* Load balance via a second distribute step */
39c4762a1bSJed Brown   PetscBool     partitionBalance; /* Balance shared point partition */
40c4762a1bSJed Brown   PetscLogStage stages[4];
41c4762a1bSJed Brown } AppCtx;
42c4762a1bSJed Brown 
43d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown   PetscFunctionBegin;
46c4762a1bSJed Brown   options->overlap          = 0;
47c4762a1bSJed Brown   options->testPartition    = PETSC_FALSE;
48c4762a1bSJed Brown   options->testRedundant    = PETSC_FALSE;
49c4762a1bSJed Brown   options->loadBalance      = PETSC_FALSE;
50c4762a1bSJed Brown   options->partitionBalance = PETSC_FALSE;
51c4762a1bSJed Brown 
52d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL, 0));
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL));
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL));
58d0609cedSBarry Smith   PetscOptionsEnd();
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]));
619566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
629566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]));
639566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]));
64*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
68d71ae5a4SJacob Faibussowitsch {
69c4762a1bSJed Brown   DM          pdm             = NULL;
70c4762a1bSJed Brown   PetscInt    triSizes_n2[2]  = {4, 4};
71c4762a1bSJed Brown   PetscInt    triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7};
72c4762a1bSJed Brown   PetscInt    triSizes_n3[3]  = {3, 2, 3};
73c4762a1bSJed Brown   PetscInt    triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4};
74c4762a1bSJed Brown   PetscInt    triSizes_n4[4]  = {2, 2, 2, 2};
75c4762a1bSJed Brown   PetscInt    triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6};
76c4762a1bSJed Brown   PetscInt    triSizes_n8[8]  = {1, 1, 1, 1, 1, 1, 1, 1};
77c4762a1bSJed Brown   PetscInt    triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7};
78c4762a1bSJed Brown   PetscInt    quadSizes[2]    = {2, 2};
79c4762a1bSJed Brown   PetscInt    quadPoints[4]   = {2, 3, 0, 1};
80c4762a1bSJed Brown   PetscInt    overlap         = user->overlap >= 0 ? user->overlap : 0;
8130602db0SMatthew G. Knepley   PetscInt    dim;
8230602db0SMatthew G. Knepley   PetscBool   simplex;
83c4762a1bSJed Brown   PetscMPIInt rank, size;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   PetscFunctionBegin;
869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
889566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD]));
899566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
909566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
919566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
929566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
939566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
949566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
959566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &dim));
969566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(*dm, &simplex));
979566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
98c4762a1bSJed Brown   if (!user->testRedundant) {
99c4762a1bSJed Brown     PetscPartitioner part;
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
1029566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
1039566063dSJacob Faibussowitsch     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
104c4762a1bSJed Brown     if (user->testPartition) {
105c4762a1bSJed Brown       const PetscInt *sizes  = NULL;
106c4762a1bSJed Brown       const PetscInt *points = NULL;
107c4762a1bSJed Brown 
108dd400576SPatrick Sanan       if (rank == 0) {
10930602db0SMatthew G. Knepley         if (dim == 2 && simplex && size == 2) {
1109371c9d4SSatish Balay           sizes  = triSizes_n2;
1119371c9d4SSatish Balay           points = triPoints_n2;
11230602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 3) {
1139371c9d4SSatish Balay           sizes  = triSizes_n3;
1149371c9d4SSatish Balay           points = triPoints_n3;
11530602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 4) {
1169371c9d4SSatish Balay           sizes  = triSizes_n4;
1179371c9d4SSatish Balay           points = triPoints_n4;
11830602db0SMatthew G. Knepley         } else if (dim == 2 && simplex && size == 8) {
1199371c9d4SSatish Balay           sizes  = triSizes_n8;
1209371c9d4SSatish Balay           points = triPoints_n8;
12130602db0SMatthew G. Knepley         } else if (dim == 2 && !simplex && size == 2) {
1229371c9d4SSatish Balay           sizes  = quadSizes;
1239371c9d4SSatish Balay           points = quadPoints;
124c4762a1bSJed Brown         }
125c4762a1bSJed Brown       }
1269566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1279566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
128c4762a1bSJed Brown     }
1299566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
130c4762a1bSJed Brown   } else {
131c4762a1bSJed Brown     PetscSF sf;
132c4762a1bSJed Brown 
1339566063dSJacob Faibussowitsch     PetscCall(DMPlexGetRedundantDM(*dm, &sf, &pdm));
134c4762a1bSJed Brown     if (sf) {
135c4762a1bSJed Brown       DM test;
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch       PetscCall(DMPlexCreate(comm, &test));
1389566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh"));
1399566063dSJacob Faibussowitsch       PetscCall(DMPlexMigrate(*dm, sf, test));
1409566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view"));
1419566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&test));
142c4762a1bSJed Brown     }
1439566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown   if (pdm) {
1469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
147c4762a1bSJed Brown     *dm = pdm;
148c4762a1bSJed Brown   }
1499566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1509566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
151c4762a1bSJed Brown   if (user->loadBalance) {
152c4762a1bSJed Brown     PetscPartitioner part;
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-prelb_dm_view"));
1559566063dSJacob Faibussowitsch     PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_"));
1569566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]));
1579566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
1589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "lb_"));
1599566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
160c4762a1bSJed Brown     if (user->testPartition) {
161c4762a1bSJed Brown       PetscInt reSizes_n2[2]  = {2, 2};
162c4762a1bSJed Brown       PetscInt rePoints_n2[4] = {2, 3, 0, 1};
1639371c9d4SSatish Balay       if (rank) {
1649371c9d4SSatish Balay         rePoints_n2[0] = 1;
1659371c9d4SSatish Balay         rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;
1669371c9d4SSatish Balay       }
167c4762a1bSJed Brown 
1689566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
1699566063dSJacob Faibussowitsch       PetscCall(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2));
170c4762a1bSJed Brown     }
1719566063dSJacob Faibussowitsch     PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance));
1729566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
173c4762a1bSJed Brown     if (pdm) {
1749566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
175c4762a1bSJed Brown       *dm = pdm;
176c4762a1bSJed Brown     }
1779566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
178c4762a1bSJed Brown   }
1799566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE]));
1809566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1819566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
182*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
186d71ae5a4SJacob Faibussowitsch {
187c4762a1bSJed Brown   DM     dm;
188c4762a1bSJed Brown   AppCtx user; /* user-defined work context */
189c4762a1bSJed Brown 
190327415f7SBarry Smith   PetscFunctionBeginUser;
1919566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1929566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1939566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1949566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
196b122ec5aSJacob Faibussowitsch   return 0;
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /*TEST
200c4762a1bSJed Brown   # Parallel, no overlap tests 0-2
201c4762a1bSJed Brown   test:
202c4762a1bSJed Brown     suffix: 0
203c4762a1bSJed Brown     requires: triangle
20430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex
205c4762a1bSJed Brown   test:
206c4762a1bSJed Brown     suffix: 1
207c4762a1bSJed Brown     requires: triangle
208c4762a1bSJed Brown     nsize: 3
20930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
210c4762a1bSJed Brown   test:
211c4762a1bSJed Brown     suffix: 2
212c4762a1bSJed Brown     requires: triangle
213c4762a1bSJed Brown     nsize: 8
21430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail
215c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
216c4762a1bSJed Brown   test:
217c4762a1bSJed Brown     suffix: 3
218c4762a1bSJed Brown     requires: triangle
219c4762a1bSJed Brown     nsize: 3
22030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
221c4762a1bSJed Brown   test:
222c4762a1bSJed Brown     suffix: 4
223c4762a1bSJed Brown     requires: triangle
224c4762a1bSJed Brown     nsize: 8
22530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
226c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
227c4762a1bSJed Brown   test:
228c4762a1bSJed Brown     suffix: 5
229c4762a1bSJed Brown     requires: triangle
230c4762a1bSJed Brown     nsize: 8
23130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail
232c4762a1bSJed Brown   # Parallel load balancing, test 6-7
233c4762a1bSJed Brown   test:
234c4762a1bSJed Brown     suffix: 6
235c4762a1bSJed Brown     requires: triangle
236c4762a1bSJed Brown     nsize: 2
23730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
238c4762a1bSJed Brown   test:
239c4762a1bSJed Brown     suffix: 7
240c4762a1bSJed Brown     requires: triangle
241c4762a1bSJed Brown     nsize: 2
24230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail
243c4762a1bSJed Brown   # Parallel redundant copying, test 8
244c4762a1bSJed Brown   test:
245c4762a1bSJed Brown     suffix: 8
246c4762a1bSJed Brown     requires: triangle
247c4762a1bSJed Brown     nsize: 2
24830602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail
249c4762a1bSJed Brown   test:
250c4762a1bSJed Brown     suffix: lb_0
251c4762a1bSJed Brown     requires: parmetis
252c4762a1bSJed Brown     nsize: 4
25330602db0SMatthew 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
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   # Same tests as above, but with balancing of the shared point partition
256c4762a1bSJed Brown   test:
257c4762a1bSJed Brown     suffix: 9
258c4762a1bSJed Brown     requires: triangle
25930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance
260c4762a1bSJed Brown   test:
261c4762a1bSJed Brown     suffix: 10
262c4762a1bSJed Brown     requires: triangle
263c4762a1bSJed Brown     nsize: 3
26430602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
265c4762a1bSJed Brown   test:
266c4762a1bSJed Brown     suffix: 11
267c4762a1bSJed Brown     requires: triangle
268c4762a1bSJed Brown     nsize: 8
26930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance
270c4762a1bSJed Brown   # Parallel, level-1 overlap tests 3-4
271c4762a1bSJed Brown   test:
272c4762a1bSJed Brown     suffix: 12
273c4762a1bSJed Brown     requires: triangle
274c4762a1bSJed Brown     nsize: 3
27530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
276c4762a1bSJed Brown   test:
277c4762a1bSJed Brown     suffix: 13
278c4762a1bSJed Brown     requires: triangle
279c4762a1bSJed Brown     nsize: 8
28030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
281c4762a1bSJed Brown   # Parallel, level-2 overlap test 5
282c4762a1bSJed Brown   test:
283c4762a1bSJed Brown     suffix: 14
284c4762a1bSJed Brown     requires: triangle
285c4762a1bSJed Brown     nsize: 8
28630602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance
287c4762a1bSJed Brown   # Parallel load balancing, test 6-7
288c4762a1bSJed Brown   test:
289c4762a1bSJed Brown     suffix: 15
290c4762a1bSJed Brown     requires: triangle
291c4762a1bSJed Brown     nsize: 2
29230602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance
293c4762a1bSJed Brown   test:
294c4762a1bSJed Brown     suffix: 16
295c4762a1bSJed Brown     requires: triangle
296c4762a1bSJed Brown     nsize: 2
29730602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance
298c4762a1bSJed Brown   # Parallel redundant copying, test 8
299c4762a1bSJed Brown   test:
300c4762a1bSJed Brown     suffix: 17
301c4762a1bSJed Brown     requires: triangle
302c4762a1bSJed Brown     nsize: 2
30330602db0SMatthew G. Knepley     args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance
304c4762a1bSJed Brown   test:
305c4762a1bSJed Brown     suffix: lb_1
306c4762a1bSJed Brown     requires: parmetis
307c4762a1bSJed Brown     nsize: 4
30830602db0SMatthew 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
309c4762a1bSJed Brown TEST*/
310