xref: /petsc/src/dm/impls/plex/tests/ex24.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[]     = "Test that MatPartitioning and PetscPartitioner interfaces are equivalent when using PETSCPARTITIONERMATPARTITIONING\n\n";
2c4762a1bSJed Brown static char FILENAME[] = "ex24.c";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown #include <petscviewerhdf5.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #if defined(PETSC_HAVE_PTSCOTCH)
8c4762a1bSJed Brown EXTERN_C_BEGIN
9c4762a1bSJed Brown   #include <ptscotch.h>
10c4762a1bSJed Brown EXTERN_C_END
11c4762a1bSJed Brown #endif
12c4762a1bSJed Brown 
13c4762a1bSJed Brown typedef struct {
14c4762a1bSJed Brown   PetscBool compare_is; /* Compare ISs and PetscSections */
15c4762a1bSJed Brown   PetscBool compare_dm; /* Compare DM */
16c4762a1bSJed Brown   PetscBool tpw;        /* Use target partition weights */
17c4762a1bSJed Brown   char      partitioning[64];
18c4762a1bSJed Brown   char      repartitioning[64];
19c4762a1bSJed Brown } AppCtx;
20c4762a1bSJed Brown 
ProcessOptions(MPI_Comm comm,AppCtx * options)21d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown   PetscBool repartition = PETSC_TRUE;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   PetscFunctionBegin;
26c4762a1bSJed Brown   options->compare_is = PETSC_FALSE;
27c4762a1bSJed Brown   options->compare_dm = PETSC_FALSE;
2830602db0SMatthew G. Knepley 
29d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_is", "Compare ISs and PetscSections?", FILENAME, options->compare_is, &options->compare_is, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_dm", "Compare DMs?", FILENAME, options->compare_dm, &options->compare_dm, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(options->partitioning, MATPARTITIONINGPARMETIS, sizeof(options->partitioning)));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-partitioning", "The mat partitioning type to test", "None", options->partitioning, options->partitioning, sizeof(options->partitioning), NULL));
349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-repartition", "Partition again after the first partition?", FILENAME, repartition, &repartition, NULL));
35c4762a1bSJed Brown   if (repartition) {
369566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(options->repartitioning, MATPARTITIONINGPARMETIS, 64));
379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-repartitioning", "The mat partitioning type to test (second partitioning)", "None", options->repartitioning, options->repartitioning, sizeof(options->repartitioning), NULL));
38c4762a1bSJed Brown   } else {
39c4762a1bSJed Brown     options->repartitioning[0] = '\0';
40c4762a1bSJed Brown   }
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tpweight", "Use target partition weights", FILENAME, options->tpw, &options->tpw, NULL));
42d0609cedSBarry Smith   PetscOptionsEnd();
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
ScotchResetRandomSeed()46d71ae5a4SJacob Faibussowitsch static PetscErrorCode ScotchResetRandomSeed()
47d71ae5a4SJacob Faibussowitsch {
48362febeeSStefano Zampini   PetscFunctionBegin;
49c4762a1bSJed Brown #if defined(PETSC_HAVE_PTSCOTCH)
50c4762a1bSJed Brown   SCOTCH_randomReset();
51c4762a1bSJed Brown #endif
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)55d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
599566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
609566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
619566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
629566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
main(int argc,char ** argv)66d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
67d71ae5a4SJacob Faibussowitsch {
68c4762a1bSJed Brown   MPI_Comm               comm;
69c4762a1bSJed Brown   DM                     dm1, dm2, dmdist1, dmdist2;
7030602db0SMatthew G. Knepley   DMPlexInterpolatedFlag interp;
71c4762a1bSJed Brown   MatPartitioning        mp;
72c4762a1bSJed Brown   PetscPartitioner       part1, part2;
73c4762a1bSJed Brown   AppCtx                 user;
74c4762a1bSJed Brown   IS                     is1 = NULL, is2 = NULL;
75c4762a1bSJed Brown   IS                     is1g, is2g;
76c4762a1bSJed Brown   PetscSection           s1 = NULL, s2 = NULL, tpws = NULL;
77c4762a1bSJed Brown   PetscInt               i;
78c4762a1bSJed Brown   PetscBool              flg;
79c4762a1bSJed Brown   PetscMPIInt            size;
80c4762a1bSJed Brown 
81327415f7SBarry Smith   PetscFunctionBeginUser;
829566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
83c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
859566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
869566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &user, &dm1));
879566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &user, &dm2));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   if (user.tpw) {
90fd599e24SMatthew G. Knepley     PetscBool isscotch;
91fd599e24SMatthew G. Knepley 
929566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(comm, &tpws));
939566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(tpws, 0, size));
94c4762a1bSJed Brown     for (i = 0; i < size; i++) {
95c4762a1bSJed Brown       PetscInt tdof = i % 2 ? 2 * i - 1 : i + 2;
969566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(tpws, i, tdof));
97c4762a1bSJed Brown     }
98fd599e24SMatthew G. Knepley     // PTScotch cannot have a zero partition weight
99fd599e24SMatthew G. Knepley     PetscCall(PetscStrncmp(user.partitioning, PETSCPARTITIONERPTSCOTCH, sizeof(PETSCPARTITIONERPTSCOTCH), &isscotch));
100fd599e24SMatthew G. Knepley     if (!isscotch && size > 1) { /* test zero tpw entry */
1019566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(tpws, 0, 0));
102c4762a1bSJed Brown     }
1039566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(tpws));
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* partition dm1 using PETSCPARTITIONERPARMETIS */
1079566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
1089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dm1, &part1));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part1, "p1_"));
1109566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetType(part1, user.partitioning));
1119566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part1));
1129566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &s1));
1139566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDMPlexPartition(part1, dm1, tpws, s1, &is1));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* partition dm2 using PETSCPARTITIONERMATPARTITIONING with MATPARTITIONINGPARMETIS */
1169566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
1179566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dm2, &part2));
1189566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part2, "p2_"));
1199566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING));
1209566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp));
1219566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetType(mp, user.partitioning));
1229566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part2));
1239566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &s2));
1249566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDMPlexPartition(part2, dm2, tpws, s2, &is2));
125c4762a1bSJed Brown 
1269566063dSJacob Faibussowitsch   PetscCall(ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g));
1279566063dSJacob Faibussowitsch   PetscCall(ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g));
1289566063dSJacob Faibussowitsch   PetscCall(ISViewFromOptions(is1g, NULL, "-seq_is1_view"));
1299566063dSJacob Faibussowitsch   PetscCall(ISViewFromOptions(is2g, NULL, "-seq_is2_view"));
130c4762a1bSJed Brown   /* compare the two ISs */
131c4762a1bSJed Brown   if (user.compare_is) {
1329566063dSJacob Faibussowitsch     PetscCall(ISEqualUnsorted(is1g, is2g, &flg));
1339566063dSJacob Faibussowitsch     if (!flg) PetscCall(PetscPrintf(comm, "ISs are not equal with type %s with size %d.\n", user.partitioning, size));
134c4762a1bSJed Brown   }
1359566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1g));
1369566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2g));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* compare the two PetscSections */
1399566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(s1, NULL, "-seq_s1_view"));
1409566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(s2, NULL, "-seq_s2_view"));
141c4762a1bSJed Brown   if (user.compare_is) {
1429566063dSJacob Faibussowitsch     PetscCall(PetscSectionCompare(s1, s2, &flg));
1439566063dSJacob Faibussowitsch     if (!flg) PetscCall(PetscPrintf(comm, "PetscSections are not equal with %s with size %d.\n", user.partitioning, size));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* distribute both DMs */
1479566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
1489566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dm1, 0, NULL, &dmdist1));
1499566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
1509566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dm2, 0, NULL, &dmdist2));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* cleanup */
1539566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&tpws));
1549566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s1));
1559566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s2));
1569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
1579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
1589566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm1));
1599566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm2));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* if distributed DMs are NULL (sequential case), then quit */
162b122ec5aSJacob Faibussowitsch   if (!dmdist1 && !dmdist2) return 0;
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dmdist1, NULL, "-dm_dist1_view"));
1659566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dmdist2, NULL, "-dm_dist2_view"));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* compare the two distributed DMs */
168c4762a1bSJed Brown   if (user.compare_dm) {
1699566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dmdist1, dmdist2, &flg));
1709566063dSJacob Faibussowitsch     if (!flg) PetscCall(PetscPrintf(comm, "Distributed DMs are not equal %s with size %d.\n", user.partitioning, size));
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* if repartitioning is disabled, then quit */
174b122ec5aSJacob Faibussowitsch   if (user.repartitioning[0] == '\0') return 0;
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   if (user.tpw) {
1779566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(comm, &tpws));
1789566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(tpws, 0, size));
179c4762a1bSJed Brown     for (i = 0; i < size; i++) {
180c4762a1bSJed Brown       PetscInt tdof = i % 2 ? i + 1 : size - i;
1819566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(tpws, i, tdof));
182c4762a1bSJed Brown     }
1839566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(tpws));
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* repartition distributed DM dmdist1 */
1879566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
1889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dmdist1, &part1));
1899566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part1, "dp1_"));
1909566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetType(part1, user.repartitioning));
1919566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part1));
1929566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &s1));
1939566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDMPlexPartition(part1, dmdist1, tpws, s1, &is1));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* repartition distributed DM dmdist2 */
1969566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
1979566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dmdist2, &part2));
1989566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part2, "dp2_"));
1999566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING));
2009566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp));
2019566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetType(mp, user.repartitioning));
2029566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part2));
2039566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &s2));
2049566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDMPlexPartition(part2, dmdist2, tpws, s2, &is2));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* compare the two ISs */
2079566063dSJacob Faibussowitsch   PetscCall(ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g));
2089566063dSJacob Faibussowitsch   PetscCall(ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g));
2099566063dSJacob Faibussowitsch   PetscCall(ISViewFromOptions(is1g, NULL, "-dist_is1_view"));
2109566063dSJacob Faibussowitsch   PetscCall(ISViewFromOptions(is2g, NULL, "-dist_is2_view"));
211c4762a1bSJed Brown   if (user.compare_is) {
2129566063dSJacob Faibussowitsch     PetscCall(ISEqualUnsorted(is1g, is2g, &flg));
2139566063dSJacob Faibussowitsch     if (!flg) PetscCall(PetscPrintf(comm, "Distributed ISs are not equal, with %s with size %d.\n", user.repartitioning, size));
214c4762a1bSJed Brown   }
2159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1g));
2169566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2g));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* compare the two PetscSections */
2199566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(s1, NULL, "-dist_s1_view"));
2209566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(s2, NULL, "-dist_s2_view"));
221c4762a1bSJed Brown   if (user.compare_is) {
2229566063dSJacob Faibussowitsch     PetscCall(PetscSectionCompare(s1, s2, &flg));
2239566063dSJacob Faibussowitsch     if (!flg) PetscCall(PetscPrintf(comm, "Distributed PetscSections are not equal, with %s with size %d.\n", user.repartitioning, size));
224c4762a1bSJed Brown   }
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* redistribute both distributed DMs */
2279566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
2289566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dmdist1, 0, NULL, &dm1));
2299566063dSJacob Faibussowitsch   PetscCall(ScotchResetRandomSeed());
2309566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dmdist2, 0, NULL, &dm2));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* compare the two distributed DMs */
2339566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolated(dm1, &interp));
23430602db0SMatthew G. Knepley   if (interp == DMPLEX_INTERPOLATED_NONE) {
2359566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dm1, dm2, &flg));
2369566063dSJacob Faibussowitsch     if (!flg) PetscCall(PetscPrintf(comm, "Redistributed DMs are not equal, with %s with size %d.\n", user.repartitioning, size));
237c4762a1bSJed Brown   }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /* cleanup */
2409566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&tpws));
2419566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s1));
2429566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s2));
2439566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
2449566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
2459566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm1));
2469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm2));
2479566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmdist1));
2489566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmdist2));
2499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
250b122ec5aSJacob Faibussowitsch   return 0;
251c4762a1bSJed Brown }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown /*TEST
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   test:
256c4762a1bSJed Brown     # partition sequential mesh loaded from Exodus file
257c4762a1bSJed Brown     suffix: 0
258c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
259c4762a1bSJed Brown     requires: chaco parmetis ptscotch exodusii
26046ac1a18SMatthew G. Knepley     args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
261c4762a1bSJed Brown     args: -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}}
262*3886731fSPierre Jolivet     output_file: output/empty.out
263c4762a1bSJed Brown   test:
264c4762a1bSJed Brown     # partition mesh generated by ctetgen using scotch, then repartition with scotch, diff view
265c4762a1bSJed Brown     suffix: 3
266c4762a1bSJed Brown     nsize: 4
267c4762a1bSJed Brown     requires: ptscotch ctetgen
26830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,2 -partitioning ptscotch -repartitioning ptscotch
269c4762a1bSJed Brown     args: -p1_petscpartitioner_view -p2_petscpartitioner_view -dp1_petscpartitioner_view -dp2_petscpartitioner_view -tpweight {{0 1}}
270c4762a1bSJed Brown   test:
271c4762a1bSJed Brown     # partition mesh generated by ctetgen using partitioners supported both by MatPartitioning and PetscPartitioner
272c4762a1bSJed Brown     suffix: 4
273c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
274c4762a1bSJed Brown     requires: chaco parmetis ptscotch ctetgen
27530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces {{2,3,4 5,4,3 7,11,5}} -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}}
276*3886731fSPierre Jolivet     output_file: output/empty.out
277c4762a1bSJed Brown 
278c4762a1bSJed Brown TEST*/
279