xref: /petsc/src/dm/impls/plex/tests/ex24.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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);
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare_is", "Compare ISs and PetscSections?", FILENAME, options->compare_is, &options->compare_is, NULL));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare_dm", "Compare DMs?", FILENAME, options->compare_dm, &options->compare_dm, NULL));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(options->partitioning,MATPARTITIONINGPARMETIS,sizeof(options->partitioning)));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-partitioning","The mat partitioning type to test","None",options->partitioning, options->partitioning,sizeof(options->partitioning),NULL));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-repartition", "Partition again after the first partition?", FILENAME, repartition, &repartition, NULL));
36c4762a1bSJed Brown   if (repartition) {
37*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrncpy(options->repartitioning,MATPARTITIONINGPARMETIS,64));
38*5f80ce2aSJacob 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   }
42*5f80ce2aSJacob 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;
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
63*5f80ce2aSJacob 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   PetscErrorCode ierr;
81c4762a1bSJed Brown   PetscMPIInt    size;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
84c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
85*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(comm, &user));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &user, &dm1));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &user, &dm2));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   if (user.tpw) {
91*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(comm, &tpws));
92*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(tpws, 0, size));
93c4762a1bSJed Brown     for (i=0;i<size;i++) {
94c4762a1bSJed Brown       PetscInt tdof = i%2 ? 2*i -1 : i+2;
95*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(tpws, i, tdof));
96c4762a1bSJed Brown     }
97c4762a1bSJed Brown     if (size > 1) { /* test zero tpw entry */
98*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(tpws, 0, 0));
99c4762a1bSJed Brown     }
100*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(tpws));
101c4762a1bSJed Brown   }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* partition dm1 using PETSCPARTITIONERPARMETIS */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dm1, &part1));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part1,"p1_"));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part1, user.partitioning));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part1));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s1));
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part1, dm1, tpws, s1, &is1));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* partition dm2 using PETSCPARTITIONERMATPARTITIONING with MATPARTITIONINGPARMETIS */
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dm2, &part2));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part2,"p2_"));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningSetType(mp, user.partitioning));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part2));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s2));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part2, dm2, tpws, s2, &is2));
122c4762a1bSJed Brown 
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is1g, NULL, "-seq_is1_view"));
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is2g, NULL, "-seq_is2_view"));
127c4762a1bSJed Brown   /* compare the two ISs */
128c4762a1bSJed Brown   if (user.compare_is) {
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqualUnsorted(is1g, is2g, &flg));
130*5f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "ISs are not equal with type %s with size %d.\n",user.partitioning,size));
131c4762a1bSJed Brown   }
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1g));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2g));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* compare the two PetscSections */
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s1, NULL, "-seq_s1_view"));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s2, NULL, "-seq_s2_view"));
138c4762a1bSJed Brown   if (user.compare_is) {
139*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCompare(s1, s2, &flg));
140*5f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "PetscSections are not equal with %s with size %d.\n",user.partitioning,size));
141c4762a1bSJed Brown   }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* distribute both DMs */
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dm1, 0, NULL, &dmdist1));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dm2, 0, NULL, &dmdist2));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* cleanup */
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&tpws));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s1));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s2));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm1));
156*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm2));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* if distributed DMs are NULL (sequential case), then quit */
159c4762a1bSJed Brown   if (!dmdist1 && !dmdist2) return ierr;
160c4762a1bSJed Brown 
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dmdist1, NULL, "-dm_dist1_view"));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dmdist2, NULL, "-dm_dist2_view"));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* compare the two distributed DMs */
165c4762a1bSJed Brown   if (user.compare_dm) {
166*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexEqual(dmdist1, dmdist2, &flg));
167*5f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Distributed DMs are not equal %s with size %d.\n",user.partitioning,size));
168c4762a1bSJed Brown   }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* if repartitioning is disabled, then quit */
171c4762a1bSJed Brown   if (user.repartitioning[0] == '\0') return ierr;
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   if (user.tpw) {
174*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(comm, &tpws));
175*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(tpws, 0, size));
176c4762a1bSJed Brown     for (i=0;i<size;i++) {
177c4762a1bSJed Brown       PetscInt tdof = i%2 ? i+1 : size - i;
178*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(tpws, i, tdof));
179c4762a1bSJed Brown     }
180*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(tpws));
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* repartition distributed DM dmdist1 */
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dmdist1, &part1));
186*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part1,"dp1_"));
187*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part1, user.repartitioning));
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part1));
189*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s1));
190*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part1, dmdist1, tpws, s1, &is1));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* repartition distributed DM dmdist2 */
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dmdist2, &part2));
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part2,"dp2_"));
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING));
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp));
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningSetType(mp, user.repartitioning));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part2));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &s2));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDMPlexPartition(part2, dmdist2, tpws, s2, &is2));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* compare the two ISs */
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g));
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is1g, NULL, "-dist_is1_view"));
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISViewFromOptions(is2g, NULL, "-dist_is2_view"));
208c4762a1bSJed Brown   if (user.compare_is) {
209*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqualUnsorted(is1g, is2g, &flg));
210*5f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Distributed ISs are not equal, with %s with size %d.\n",user.repartitioning,size));
211c4762a1bSJed Brown   }
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1g));
213*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2g));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* compare the two PetscSections */
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s1, NULL, "-dist_s1_view"));
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s2, NULL, "-dist_s2_view"));
218c4762a1bSJed Brown   if (user.compare_is) {
219*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCompare(s1, s2, &flg));
220*5f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Distributed PetscSections are not equal, with %s with size %d.\n",user.repartitioning,size));
221c4762a1bSJed Brown   }
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* redistribute both distributed DMs */
224*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dmdist1, 0, NULL, &dm1));
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ScotchResetRandomSeed());
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dmdist2, 0, NULL, &dm2));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* compare the two distributed DMs */
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsInterpolated(dm1, &interp));
23130602db0SMatthew G. Knepley   if (interp == DMPLEX_INTERPOLATED_NONE) {
232*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexEqual(dm1, dm2, &flg));
233*5f80ce2aSJacob Faibussowitsch     if (!flg) CHKERRQ(PetscPrintf(comm, "Redistributed DMs are not equal, with %s with size %d.\n",user.repartitioning,size));
234c4762a1bSJed Brown   }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* cleanup */
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&tpws));
238*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s1));
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s2));
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm1));
243*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm2));
244*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmdist1));
245*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmdist2));
246c4762a1bSJed Brown   ierr = PetscFinalize();
247c4762a1bSJed Brown   return ierr;
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown /*TEST
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   test:
253c4762a1bSJed Brown     # partition sequential mesh loaded from Exodus file
254c4762a1bSJed Brown     suffix: 0
255c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
256c4762a1bSJed Brown     requires: chaco parmetis ptscotch exodusii
25730602db0SMatthew G. Knepley     args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
258c4762a1bSJed Brown     args: -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}}
259c4762a1bSJed Brown   test:
260c4762a1bSJed Brown     # repartition mesh already partitioned naively by MED loader
261c4762a1bSJed Brown     suffix: 1
262c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
263c4762a1bSJed Brown     TODO: MED
264c4762a1bSJed Brown     requires: parmetis ptscotch med
26530602db0SMatthew G. Knepley     args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med
266c4762a1bSJed Brown     args: -repartition 0 -partitioning {{parmetis ptscotch}}
267c4762a1bSJed Brown   test:
268c4762a1bSJed Brown     # partition mesh generated by ctetgen using scotch, then repartition with scotch, diff view
269c4762a1bSJed Brown     suffix: 3
270c4762a1bSJed Brown     nsize: 4
271c4762a1bSJed Brown     requires: ptscotch ctetgen
27230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,2 -partitioning ptscotch -repartitioning ptscotch
273c4762a1bSJed Brown     args: -p1_petscpartitioner_view -p2_petscpartitioner_view -dp1_petscpartitioner_view -dp2_petscpartitioner_view -tpweight {{0 1}}
274c4762a1bSJed Brown   test:
275c4762a1bSJed Brown     # partition mesh generated by ctetgen using partitioners supported both by MatPartitioning and PetscPartitioner
276c4762a1bSJed Brown     suffix: 4
277c4762a1bSJed Brown     nsize: {{1 2 3 4 8}}
278c4762a1bSJed Brown     requires: chaco parmetis ptscotch ctetgen
27930602db0SMatthew 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}}
280c4762a1bSJed Brown 
281c4762a1bSJed Brown TEST*/
282