xref: /petsc/src/dm/impls/plex/tests/ex24.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
21c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22c4762a1bSJed Brown {
23c4762a1bSJed Brown   PetscBool      repartition = PETSC_TRUE;
24c4762a1bSJed Brown   PetscErrorCode ierr;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   PetscFunctionBegin;
27c4762a1bSJed Brown   options->compare_is = PETSC_FALSE;
28c4762a1bSJed Brown   options->compare_dm = PETSC_FALSE;
2930602db0SMatthew G. Knepley 
30c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare_is", "Compare ISs and PetscSections?", FILENAME, options->compare_is, &options->compare_is, NULL));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare_dm", "Compare DMs?", FILENAME, options->compare_dm, &options->compare_dm, NULL));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(options->partitioning,MATPARTITIONINGPARMETIS,sizeof(options->partitioning)));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-partitioning","The mat partitioning type to test","None",options->partitioning, options->partitioning,sizeof(options->partitioning),NULL));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-repartition", "Partition again after the first partition?", FILENAME, repartition, &repartition, NULL));
36c4762a1bSJed Brown   if (repartition) {
375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrncpy(options->repartitioning,MATPARTITIONINGPARMETIS,64));
385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsString("-repartitioning","The mat partitioning type to test (second partitioning)","None", options->repartitioning, options->repartitioning,sizeof(options->repartitioning),NULL));
39c4762a1bSJed Brown   } else {
40c4762a1bSJed Brown     options->repartitioning[0] = '\0';
41c4762a1bSJed Brown   }
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-tpweight", "Use target partition weights", FILENAME, options->tpw, &options->tpw, NULL));
431e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
44c4762a1bSJed Brown   PetscFunctionReturn(0);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown static PetscErrorCode ScotchResetRandomSeed()
48c4762a1bSJed Brown {
49362febeeSStefano Zampini   PetscFunctionBegin;
50c4762a1bSJed Brown #if defined(PETSC_HAVE_PTSCOTCH)
51c4762a1bSJed Brown   SCOTCH_randomReset();
52c4762a1bSJed Brown #endif
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
57c4762a1bSJed Brown {
58c4762a1bSJed Brown   PetscFunctionBegin;
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
64c4762a1bSJed Brown   PetscFunctionReturn(0);
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown int main(int argc, char **argv)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   MPI_Comm       comm;
70c4762a1bSJed Brown   DM             dm1, dm2, dmdist1, dmdist2;
7130602db0SMatthew G. Knepley   DMPlexInterpolatedFlag interp;
72c4762a1bSJed Brown   MatPartitioning mp;
73c4762a1bSJed Brown   PetscPartitioner part1, part2;
74c4762a1bSJed Brown   AppCtx         user;
75c4762a1bSJed Brown   IS             is1=NULL, is2=NULL;
76c4762a1bSJed Brown   IS             is1g, is2g;
77c4762a1bSJed Brown   PetscSection   s1=NULL, s2=NULL, tpws = NULL;
78c4762a1bSJed Brown   PetscInt       i;
79c4762a1bSJed Brown   PetscBool      flg;
80c4762a1bSJed Brown   PetscMPIInt    size;
81c4762a1bSJed Brown 
82*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
83c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
845f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(comm, &user));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &user, &dm1));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &user, &dm2));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   if (user.tpw) {
905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(comm, &tpws));
915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(tpws, 0, size));
92c4762a1bSJed Brown     for (i=0;i<size;i++) {
93c4762a1bSJed Brown       PetscInt tdof = i%2 ? 2*i -1 : i+2;
945f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(tpws, i, tdof));
95c4762a1bSJed Brown     }
96c4762a1bSJed Brown     if (size > 1) { /* test zero tpw entry */
975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(tpws, 0, 0));
98c4762a1bSJed Brown     }
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(tpws));
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* partition dm1 using PETSCPARTITIONERPARMETIS */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dm1, &part1));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part1,"p1_"));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part1, user.partitioning));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part1));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s1));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part1, dm1, tpws, s1, &is1));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* partition dm2 using PETSCPARTITIONERMATPARTITIONING with MATPARTITIONINGPARMETIS */
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dm2, &part2));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part2,"p2_"));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningSetType(mp, user.partitioning));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part2));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s2));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part2, dm2, tpws, s2, &is2));
121c4762a1bSJed Brown 
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is1g, NULL, "-seq_is1_view"));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is2g, NULL, "-seq_is2_view"));
126c4762a1bSJed Brown   /* compare the two ISs */
127c4762a1bSJed Brown   if (user.compare_is) {
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqualUnsorted(is1g, is2g, &flg));
1295f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "ISs are not equal with type %s with size %d.\n",user.partitioning,size));
130c4762a1bSJed Brown   }
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1g));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2g));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* compare the two PetscSections */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s1, NULL, "-seq_s1_view"));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s2, NULL, "-seq_s2_view"));
137c4762a1bSJed Brown   if (user.compare_is) {
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCompare(s1, s2, &flg));
1395f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "PetscSections are not equal with %s with size %d.\n",user.partitioning,size));
140c4762a1bSJed Brown   }
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* distribute both DMs */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dm1, 0, NULL, &dmdist1));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dm2, 0, NULL, &dmdist2));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* cleanup */
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&tpws));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s1));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s2));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm1));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm2));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* if distributed DMs are NULL (sequential case), then quit */
158*b122ec5aSJacob Faibussowitsch   if (!dmdist1 && !dmdist2) return 0;
159c4762a1bSJed Brown 
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dmdist1, NULL, "-dm_dist1_view"));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dmdist2, NULL, "-dm_dist2_view"));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* compare the two distributed DMs */
164c4762a1bSJed Brown   if (user.compare_dm) {
1655f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexEqual(dmdist1, dmdist2, &flg));
1665f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Distributed DMs are not equal %s with size %d.\n",user.partitioning,size));
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* if repartitioning is disabled, then quit */
170*b122ec5aSJacob Faibussowitsch   if (user.repartitioning[0] == '\0') return 0;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   if (user.tpw) {
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(comm, &tpws));
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(tpws, 0, size));
175c4762a1bSJed Brown     for (i=0;i<size;i++) {
176c4762a1bSJed Brown       PetscInt tdof = i%2 ? i+1 : size - i;
1775f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(tpws, i, tdof));
178c4762a1bSJed Brown     }
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(tpws));
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* repartition distributed DM dmdist1 */
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dmdist1, &part1));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part1,"dp1_"));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part1, user.repartitioning));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part1));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s1));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part1, dmdist1, tpws, s1, &is1));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* repartition distributed DM dmdist2 */
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dmdist2, &part2));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part2,"dp2_"));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningSetType(mp, user.repartitioning));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part2));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s2));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part2, dmdist2, tpws, s2, &is2));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* compare the two ISs */
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g));
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is1g, NULL, "-dist_is1_view"));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is2g, NULL, "-dist_is2_view"));
207c4762a1bSJed Brown   if (user.compare_is) {
2085f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqualUnsorted(is1g, is2g, &flg));
2095f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Distributed ISs are not equal, with %s with size %d.\n",user.repartitioning,size));
210c4762a1bSJed Brown   }
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1g));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2g));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* compare the two PetscSections */
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s1, NULL, "-dist_s1_view"));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s2, NULL, "-dist_s2_view"));
217c4762a1bSJed Brown   if (user.compare_is) {
2185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCompare(s1, s2, &flg));
2195f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Distributed PetscSections are not equal, with %s with size %d.\n",user.repartitioning,size));
220c4762a1bSJed Brown   }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /* redistribute both distributed DMs */
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dmdist1, 0, NULL, &dm1));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dmdist2, 0, NULL, &dm2));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* compare the two distributed DMs */
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsInterpolated(dm1, &interp));
23030602db0SMatthew G. Knepley   if (interp == DMPLEX_INTERPOLATED_NONE) {
2315f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexEqual(dm1, dm2, &flg));
2325f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Redistributed DMs are not equal, with %s with size %d.\n",user.repartitioning,size));
233c4762a1bSJed Brown   }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* cleanup */
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&tpws));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s1));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s2));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm1));
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm2));
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmdist1));
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmdist2));
245*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
246*b122ec5aSJacob Faibussowitsch   return 0;
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown /*TEST
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   test:
252c4762a1bSJed Brown     # partition sequential mesh loaded from Exodus file
253c4762a1bSJed Brown     suffix: 0
254c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
255c4762a1bSJed Brown     requires: chaco parmetis ptscotch exodusii
25630602db0SMatthew G. Knepley     args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
257c4762a1bSJed Brown     args: -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}}
258c4762a1bSJed Brown   test:
259c4762a1bSJed Brown     # repartition mesh already partitioned naively by MED loader
260c4762a1bSJed Brown     suffix: 1
261c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
262c4762a1bSJed Brown     TODO: MED
263c4762a1bSJed Brown     requires: parmetis ptscotch med
26430602db0SMatthew G. Knepley     args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med
265c4762a1bSJed Brown     args: -repartition 0 -partitioning {{parmetis ptscotch}}
266c4762a1bSJed Brown   test:
267c4762a1bSJed Brown     # partition mesh generated by ctetgen using scotch, then repartition with scotch, diff view
268c4762a1bSJed Brown     suffix: 3
269c4762a1bSJed Brown     nsize: 4
270c4762a1bSJed Brown     requires: ptscotch ctetgen
27130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,2 -partitioning ptscotch -repartitioning ptscotch
272c4762a1bSJed Brown     args: -p1_petscpartitioner_view -p2_petscpartitioner_view -dp1_petscpartitioner_view -dp2_petscpartitioner_view -tpweight {{0 1}}
273c4762a1bSJed Brown   test:
274c4762a1bSJed Brown     # partition mesh generated by ctetgen using partitioners supported both by MatPartitioning and PetscPartitioner
275c4762a1bSJed Brown     suffix: 4
276c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
277c4762a1bSJed Brown     requires: chaco parmetis ptscotch ctetgen
27830602db0SMatthew 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}}
279c4762a1bSJed Brown 
280c4762a1bSJed Brown TEST*/
281