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