xref: /petsc/src/dm/impls/plex/tests/ex4.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "Tests for refinement of meshes created by hand\n\n";
2 
3 #include <petscdmplex.h>
4 
5 typedef struct {
6   PetscInt  debug;         /* The debugging level */
7   PetscInt  dim;           /* The topological mesh dimension */
8   PetscBool cellHybrid;    /* Use a hybrid mesh */
9   PetscBool cellSimplex;   /* Use simplices or hexes */
10   PetscBool testPartition; /* Use a fixed partitioning for testing */
11   PetscInt  testNum;       /* The particular mesh to test */
12   PetscBool uninterpolate; /* Uninterpolate the mesh at the end */
13   PetscBool reinterpolate; /* Reinterpolate the mesh at the end */
14 } AppCtx;
15 
16 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
17   PetscFunctionBegin;
18   options->debug         = 0;
19   options->dim           = 2;
20   options->cellHybrid    = PETSC_TRUE;
21   options->cellSimplex   = PETSC_TRUE;
22   options->testPartition = PETSC_TRUE;
23   options->testNum       = 0;
24   options->uninterpolate = PETSC_FALSE;
25   options->reinterpolate = PETSC_FALSE;
26 
27   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
28   PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex4.c", options->debug, &options->debug, NULL, 0));
29   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL, 1, 3));
30   PetscCall(PetscOptionsBool("-cell_hybrid", "Use a hybrid mesh", "ex4.c", options->cellHybrid, &options->cellHybrid, NULL));
31   PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex4.c", options->cellSimplex, &options->cellSimplex, NULL));
32   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex4.c", options->testPartition, &options->testPartition, NULL));
33   PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex4.c", options->testNum, &options->testNum, NULL, 0));
34   PetscCall(PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex4.c", options->uninterpolate, &options->uninterpolate, NULL));
35   PetscCall(PetscOptionsBool("-reinterpolate", "Reinterpolate the mesh at the end", "ex4.c", options->reinterpolate, &options->reinterpolate, NULL));
36   PetscOptionsEnd();
37   PetscFunctionReturn(0);
38 }
39 
40 /* Two segments
41 
42   2-------0-------3-------1-------4
43 
44 become
45 
46   4---0---7---1---5---2---8---3---6
47 
48 */
49 PetscErrorCode CreateSimplex_1D(MPI_Comm comm, DM *dm) {
50   PetscInt    depth = 1;
51   PetscMPIInt rank;
52 
53   PetscFunctionBegin;
54   PetscCallMPI(MPI_Comm_rank(comm, &rank));
55   if (rank == 0) {
56     PetscInt    numPoints[2]         = {3, 2};
57     PetscInt    coneSize[5]          = {2, 2, 0, 0, 0};
58     PetscInt    cones[4]             = {2, 3, 3, 4};
59     PetscInt    coneOrientations[16] = {0, 0, 0, 0};
60     PetscScalar vertexCoords[3]      = {-1.0, 0.0, 1.0};
61 
62     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
63   } else {
64     PetscInt numPoints[2] = {0, 0};
65 
66     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
67   }
68   PetscFunctionReturn(0);
69 }
70 
71 /* Two triangles
72         4
73       / | \
74      8  |  10
75     /   |   \
76    2  0 7  1 5
77     \   |   /
78      6  |  9
79       \ | /
80         3
81 
82 Becomes
83            10
84           / | \
85         21  |  26
86         /   |   \
87       14 2 20 4  16
88       /|\   |   /|\
89     22 | 28 | 32 | 25
90     /  |  \ | /  | 6\
91    8  29 3 13  7 31  11
92     \0 |  / | \  |  /
93     17 | 27 | 30 | 24
94       \|/   |   \|/
95       12 1 19 5  15
96         \   |   /
97         18  |  23
98           \ | /
99             9
100 */
101 PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *dm) {
102   PetscInt    depth = 2;
103   PetscMPIInt rank;
104 
105   PetscFunctionBegin;
106   PetscCallMPI(MPI_Comm_rank(comm, &rank));
107   if (rank == 0) {
108     PetscInt    numPoints[3]         = {4, 5, 2};
109     PetscInt    coneSize[11]         = {3, 3, 0, 0, 0, 0, 2, 2, 2, 2, 2};
110     PetscInt    cones[16]            = {6, 7, 8, 7, 9, 10, 2, 3, 3, 4, 4, 2, 3, 5, 5, 4};
111     PetscInt    coneOrientations[16] = {0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
112     PetscScalar vertexCoords[8]      = {-0.5, 0.0, 0.0, -0.5, 0.0, 0.5, 0.5, 0.0};
113 
114     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
115   } else {
116     PetscInt numPoints[3] = {0, 0, 0};
117 
118     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
119   }
120   PetscFunctionReturn(0);
121 }
122 
123 /* Two triangles separated by a zero-volume cell with 4 vertices/2 edges
124         5--16--8
125       / |      | \
126     11  |      |  12
127     /   |      |   \
128    3  0 10  2 14 1  6
129     \   |      |   /
130      9  |      |  13
131       \ |      | /
132         4--15--7
133 */
134 PetscErrorCode CreateSimplexHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) {
135   DM          idm, hdm = NULL;
136   DMLabel     faultLabel, hybridLabel;
137   PetscInt    p;
138   PetscMPIInt rank;
139 
140   PetscFunctionBegin;
141   PetscCallMPI(MPI_Comm_rank(comm, &rank));
142   if (rank == 0) {
143     switch (testNum) {
144     case 0: {
145       PetscInt    numPoints[2]        = {4, 2};
146       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
147       PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
148       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
149       PetscScalar vertexCoords[8]     = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, 1.0, 0.5};
150       PetscInt    faultPoints[2]      = {3, 4};
151 
152       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
153       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
154     } break;
155     case 1: {
156       PetscInt    numPoints[2]         = {5, 4};
157       PetscInt    coneSize[9]          = {3, 3, 3, 3, 0, 0, 0, 0, 0};
158       PetscInt    cones[12]            = {4, 5, 6, 6, 7, 4, 6, 5, 8, 6, 8, 7};
159       PetscInt    coneOrientations[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
160       PetscScalar vertexCoords[10]     = {-1.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0};
161       PetscInt    faultPoints[3]       = {5, 6, 7};
162 
163       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
164       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
165     } break;
166     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
167     }
168     PetscCall(DMPlexInterpolate(*dm, &idm));
169     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
170     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
171     PetscCall(DMSetFromOptions(idm));
172     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
173     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
174     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm));
175     PetscCall(DMLabelDestroy(&hybridLabel));
176     PetscCall(DMDestroy(&idm));
177     PetscCall(DMDestroy(dm));
178     *dm = hdm;
179   } else {
180     PetscInt numPoints[2] = {0, 0};
181 
182     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
183     PetscCall(DMPlexInterpolate(*dm, &idm));
184     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
185     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
186     PetscCall(DMSetFromOptions(idm));
187     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
188     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
189     PetscCall(DMDestroy(&idm));
190     PetscCall(DMDestroy(dm));
191     *dm = hdm;
192   }
193   PetscFunctionReturn(0);
194 }
195 
196 /* Two quadrilaterals
197 
198   5----10-----4----14-----7
199   |           |           |
200   |           |           |
201   |           |           |
202  11     0     9     1     13
203   |           |           |
204   |           |           |
205   |           |           |
206   2-----8-----3----12-----6
207 */
208 PetscErrorCode CreateTensorProduct_2D(MPI_Comm comm, PetscInt testNum, DM *dm) {
209   PetscInt    depth = 2;
210   PetscMPIInt rank;
211 
212   PetscFunctionBegin;
213   PetscCallMPI(MPI_Comm_rank(comm, &rank));
214   if (rank == 0) {
215     PetscInt    numPoints[3]         = {6, 7, 2};
216     PetscInt    coneSize[15]         = {4, 4, 0, 0, 0, 0, 0, 0, 2, 2, 2, 2, 2, 2, 2};
217     PetscInt    cones[22]            = {8, 9, 10, 11, 12, 13, 14, 9, 2, 3, 3, 4, 4, 5, 5, 2, 3, 6, 6, 7, 7, 4};
218     PetscInt    coneOrientations[22] = {0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
219     PetscScalar vertexCoords[12]     = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, -1.0, 0.5, 1.0, -0.5, 1.0, 0.5};
220 
221     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
222   } else {
223     PetscInt numPoints[3] = {0, 0, 0};
224 
225     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
226   }
227   PetscFunctionReturn(0);
228 }
229 
230 PetscErrorCode CreateTensorProductHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) {
231   DM          idm, hdm = NULL;
232   DMLabel     faultLabel, hybridLabel;
233   PetscInt    p;
234   PetscMPIInt rank;
235 
236   PetscFunctionBegin;
237   PetscCallMPI(MPI_Comm_rank(comm, &rank));
238   if (rank == 0) {
239     PetscInt numPoints[2] = {6, 2};
240     PetscInt coneSize[8]  = {4, 4, 0, 0, 0, 0, 0, 0};
241     PetscInt cones[8]     = {
242           2, 3, 4, 5, 3, 6, 7, 4,
243     };
244     PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
245     PetscScalar vertexCoords[12]    = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, -1.0, 0.5, 1.0, -0.5, 1.0, 0.5};
246     PetscInt    faultPoints[2]      = {3, 4};
247 
248     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
249     for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
250     PetscCall(DMPlexInterpolate(*dm, &idm));
251     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
252     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
253     PetscCall(DMSetFromOptions(idm));
254     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
255     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
256     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm));
257     PetscCall(DMLabelDestroy(&hybridLabel));
258   } else {
259     PetscInt numPoints[3] = {0, 0, 0};
260 
261     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
262     PetscCall(DMPlexInterpolate(*dm, &idm));
263     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
264     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
265     PetscCall(DMSetFromOptions(idm));
266     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
267     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
268   }
269   PetscCall(DMDestroy(&idm));
270   PetscCall(DMDestroy(dm));
271   *dm = hdm;
272   PetscFunctionReturn(0);
273 }
274 
275 /* Two tetrahedrons
276 
277  cell   5          5______    cell
278  0    / | \        |\      \     1
279     17  |  18      | 18 13  21
280     /8 19 10\     19  \      \
281    2-14-|----4     |   4--22--6
282     \ 9 | 7 /      |10 /      /
283     16  |  15      | 15  12 20
284       \ | /        |/      /
285         3          3------
286 */
287 PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscInt testNum, DM *dm) {
288   DM          idm;
289   PetscInt    depth = 3;
290   PetscMPIInt rank;
291 
292   PetscFunctionBegin;
293   PetscCallMPI(MPI_Comm_rank(comm, &rank));
294   if (rank == 0) {
295     switch (testNum) {
296     case 0: {
297       PetscInt    numPoints[4]         = {5, 9, 7, 2};
298       PetscInt    coneSize[23]         = {4, 4, 0, 0, 0, 0, 0, 3, 3, 3, 3, 3, 3, 3, 2, 2, 2, 2, 2, 2, 2, 2, 2};
299       PetscInt    cones[47]            = {7, 8, 9, 10, 10, 11, 12, 13, 14, 15, 16, 17, 18, 14, 16, 19, 17, 15, 18, 19, 20, 21, 19, 15, 22, 20, 18, 21, 22, 2, 4, 4, 3, 3, 2, 2, 5, 5, 4, 3, 5, 3, 6, 6, 5, 4, 6};
300       PetscInt    coneOrientations[47] = {0, 0, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, -1, -1, 0, -1, -1, -1, -1, 0, 0, -1, -1, 0, -1, -1, -1, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
301       PetscScalar vertexCoords[15]     = {0.0, 0.0, -0.5, 0.0, -0.5, 0.0, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5};
302 
303       PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
304     } break;
305     case 1: {
306       PetscInt    numPoints[2]        = {5, 2};
307       PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
308       PetscInt    cones[8]            = {4, 3, 5, 2, 5, 3, 4, 6};
309       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
310       PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0};
311 
312       depth = 1;
313       PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
314       PetscCall(DMPlexInterpolate(*dm, &idm));
315       PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
316       PetscCall(DMDestroy(dm));
317       *dm = idm;
318     } break;
319     case 2: {
320       PetscInt    numPoints[2]        = {4, 1};
321       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
322       PetscInt    cones[4]            = {2, 3, 4, 1};
323       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
324       PetscScalar vertexCoords[12]    = {0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
325 
326       depth = 1;
327       PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
328       PetscCall(DMPlexInterpolate(*dm, &idm));
329       PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
330       PetscCall(DMDestroy(dm));
331       *dm = idm;
332     } break;
333     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
334     }
335   } else {
336     PetscInt numPoints[4] = {0, 0, 0, 0};
337 
338     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
339     switch (testNum) {
340     case 1:
341       PetscCall(DMPlexInterpolate(*dm, &idm));
342       PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
343       PetscCall(DMDestroy(dm));
344       *dm = idm;
345       break;
346     }
347   }
348   PetscFunctionReturn(0);
349 }
350 
351 /* Two tetrahedrons separated by a zero-volume cell with 6 vertices
352 
353  cell   6 ___33___10______    cell
354  0    / | \        |\      \     1
355     21  |  23      | 29     27
356     /12 24 14\    30  \      \
357    3-20-|----5--32-|---9--26--7
358     \ 13| 11/      |18 /      /
359     19  |  22      | 28     25
360       \ | /        |/      /
361         4----31----8------
362          cell 2
363 */
364 PetscErrorCode CreateSimplexHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) {
365   DM          idm, hdm = NULL;
366   DMLabel     faultLabel, hybridLabel;
367   PetscInt    p;
368   PetscMPIInt rank;
369 
370   PetscFunctionBegin;
371   PetscCallMPI(MPI_Comm_rank(comm, &rank));
372   if (rank == 0) {
373     switch (testNum) {
374     case 0: {
375       PetscInt    numPoints[2]        = {5, 2};
376       PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
377       PetscInt    cones[8]            = {4, 3, 5, 2, 5, 3, 4, 6};
378       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
379       PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0};
380       PetscInt    faultPoints[3]      = {3, 4, 5};
381 
382       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
383       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
384     } break;
385     case 1: {
386       /* Tets 0,3,5 and 1,2,4 */
387       PetscInt    numPoints[2]         = {9, 6};
388       PetscInt    coneSize[15]         = {4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0};
389       PetscInt    cones[24]            = {7, 8, 9, 6, 11, 13, 9, 14, 10, 11, 13, 9, 9, 10, 11, 7, 9, 14, 13, 12, 7, 8, 11, 9};
390       PetscInt    coneOrientations[24] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
391       PetscScalar vertexCoords[27]     = {-2.0, -1.0, 0.0, -2.0, 0.0, 0.0, -2.0, 0.0, 1.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 2.0, -1.0, 0.0, 2.0, 0.0, 0.0, 2.0, 0.0, 1.0};
392       PetscInt    faultPoints[3]       = {9, 10, 11};
393 
394       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
395       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
396     } break;
397     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
398     }
399     PetscCall(DMPlexInterpolate(*dm, &idm));
400     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
401     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
402     PetscCall(DMSetFromOptions(idm));
403     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
404     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
405     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm));
406     PetscCall(DMLabelDestroy(&hybridLabel));
407     PetscCall(DMDestroy(&idm));
408     PetscCall(DMDestroy(dm));
409     *dm = hdm;
410   } else {
411     PetscInt numPoints[4] = {0, 0, 0, 0};
412 
413     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
414     PetscCall(DMPlexInterpolate(*dm, &idm));
415     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
416     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
417     PetscCall(DMSetFromOptions(idm));
418     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
419     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
420     PetscCall(DMDestroy(&idm));
421     PetscCall(DMDestroy(dm));
422     *dm = hdm;
423   }
424   PetscFunctionReturn(0);
425 }
426 
427 PetscErrorCode CreateTensorProduct_3D(MPI_Comm comm, PetscInt testNum, DM *dm) {
428   DM          idm;
429   PetscMPIInt rank;
430 
431   PetscFunctionBegin;
432   PetscCallMPI(MPI_Comm_rank(comm, &rank));
433   if (rank == 0) {
434     switch (testNum) {
435     case 0: {
436       PetscInt    numPoints[2]         = {12, 2};
437       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
438       PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8};
439       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
440       PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5, -1.0, 0.5, -0.5, 0.0, 0.5, -0.5, 0.0, -0.5, -0.5, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5};
441 
442       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
443     } break;
444     case 1: {
445       PetscInt    numPoints[2]        = {8, 1};
446       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
447       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
448       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
449       PetscScalar vertexCoords[24]    = {-1.0, -1.0, -1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, -1.0, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0};
450 
451       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
452     } break;
453     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
454     }
455   } else {
456     PetscInt numPoints[4] = {0, 0, 0, 0};
457 
458     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
459   }
460   PetscCall(DMPlexInterpolate(*dm, &idm));
461   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
462   PetscCall(DMDestroy(dm));
463   *dm = idm;
464   PetscFunctionReturn(0);
465 }
466 
467 PetscErrorCode CreateTensorProductHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) {
468   DM          idm, hdm = NULL;
469   DMLabel     faultLabel;
470   PetscInt    p;
471   PetscMPIInt rank;
472 
473   PetscFunctionBegin;
474   PetscCallMPI(MPI_Comm_rank(comm, &rank));
475   if (rank == 0) {
476     switch (testNum) {
477     case 0: {
478       PetscInt    numPoints[2]         = {12, 2};
479       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
480       PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8};
481       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
482       PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5, -1.0, 0.5, -0.5, 0.0, 0.5, -0.5, 0.0, -0.5, -0.5, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5};
483       PetscInt    faultPoints[4]       = {2, 3, 5, 6};
484 
485       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
486       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
487     } break;
488     case 1: {
489       PetscInt numPoints[2] = {30, 7};
490       PetscInt coneSize[37] = {8, 8, 8, 8, 8, 8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
491       PetscInt cones[56] = {8, 21, 20, 7, 13, 12, 23, 24, 14, 15, 10, 9, 13, 8, 21, 24, 15, 16, 11, 10, 24, 21, 22, 25, 30, 29, 28, 21, 35, 24, 33, 34, 24, 21, 30, 35, 25, 36, 31, 22, 27, 20, 21, 28, 32, 33, 24, 23, 15, 24, 13, 14, 19, 18, 17, 26};
492       PetscInt coneOrientations[56] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
493       PetscScalar vertexCoords[90]  = {-2.0, -2.0, -2.0, -2.0, -1.0, -2.0, -3.0, 0.0, -2.0, -2.0, 1.0,  -2.0, -2.0, 2.0, -2.0, -2.0, -2.0, 0.0,  -2.0, -1.0, 0.0, -3.0, 0.0, 0.0, -2.0, 1.0, 0.0, -2.0, 2.0, 0.0,
494                                        -2.0, -1.0, 2.0,  -3.0, 0.0,  2.0,  -2.0, 1.0, 2.0,  0.0,  -2.0, -2.0, 0.0,  0.0, -2.0, 0.0,  2.0,  -2.0, 0.0,  -2.0, 0.0, 0.0,  0.0, 0.0, 0.0,  2.0, 0.0, 0.0,  0.0, 2.0,
495                                        2.0,  -2.0, -2.0, 2.0,  -1.0, -2.0, 3.0,  0.0, -2.0, 2.0,  1.0,  -2.0, 2.0,  2.0, -2.0, 2.0,  -2.0, 0.0,  2.0,  -1.0, 0.0, 3.0,  0.0, 0.0, 2.0,  1.0, 0.0, 2.0,  2.0, 0.0};
496       PetscInt    faultPoints[6]    = {20, 21, 22, 23, 24, 25};
497 
498       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
499       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
500     } break;
501     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
502     }
503     PetscCall(DMPlexInterpolate(*dm, &idm));
504     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
505     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
506     PetscCall(DMSetFromOptions(idm));
507     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
508     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
509     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, NULL, NULL, NULL, &hdm));
510     PetscCall(DMDestroy(&idm));
511     PetscCall(DMDestroy(dm));
512     *dm = hdm;
513   } else {
514     PetscInt numPoints[4] = {0, 0, 0, 0};
515 
516     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
517     PetscCall(DMPlexInterpolate(*dm, &idm));
518     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
519     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
520     PetscCall(DMSetFromOptions(idm));
521     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
522     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
523     PetscCall(DMDestroy(&idm));
524     PetscCall(DMDestroy(dm));
525     *dm = hdm;
526   }
527   PetscFunctionReturn(0);
528 }
529 
530 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
531   PetscInt    dim         = user->dim;
532   PetscBool   cellHybrid  = user->cellHybrid;
533   PetscBool   cellSimplex = user->cellSimplex;
534   PetscMPIInt rank, size;
535 
536   PetscFunctionBegin;
537   PetscCallMPI(MPI_Comm_rank(comm, &rank));
538   PetscCallMPI(MPI_Comm_size(comm, &size));
539   PetscCall(DMCreate(comm, dm));
540   PetscCall(DMSetType(*dm, DMPLEX));
541   PetscCall(DMSetDimension(*dm, dim));
542   switch (dim) {
543   case 1:
544     PetscCheck(!cellHybrid, comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim);
545     PetscCall(CreateSimplex_1D(comm, dm));
546     break;
547   case 2:
548     if (cellSimplex) {
549       if (cellHybrid) {
550         PetscCall(CreateSimplexHybrid_2D(comm, user->testNum, dm));
551       } else {
552         PetscCall(CreateSimplex_2D(comm, dm));
553       }
554     } else {
555       if (cellHybrid) {
556         PetscCall(CreateTensorProductHybrid_2D(comm, user->testNum, dm));
557       } else {
558         PetscCall(CreateTensorProduct_2D(comm, user->testNum, dm));
559       }
560     }
561     break;
562   case 3:
563     if (cellSimplex) {
564       if (cellHybrid) {
565         PetscCall(CreateSimplexHybrid_3D(comm, user->testNum, dm));
566       } else {
567         PetscCall(CreateSimplex_3D(comm, user->testNum, dm));
568       }
569     } else {
570       if (cellHybrid) {
571         PetscCall(CreateTensorProductHybrid_3D(comm, user->testNum, dm));
572       } else {
573         PetscCall(CreateTensorProduct_3D(comm, user->testNum, dm));
574       }
575     }
576     break;
577   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
578   }
579   if (user->testPartition && size > 1) {
580     PetscPartitioner part;
581     PetscInt        *sizes  = NULL;
582     PetscInt        *points = NULL;
583 
584     if (rank == 0) {
585       if (dim == 2 && cellSimplex && !cellHybrid && size == 2) {
586         switch (user->testNum) {
587         case 0: {
588           PetscInt triSizes_p2[2]  = {1, 1};
589           PetscInt triPoints_p2[2] = {0, 1};
590 
591           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
592           PetscCall(PetscArraycpy(sizes, triSizes_p2, 2));
593           PetscCall(PetscArraycpy(points, triPoints_p2, 2));
594           break;
595         }
596         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum);
597         }
598       } else if (dim == 2 && cellSimplex && cellHybrid && size == 2) {
599         switch (user->testNum) {
600         case 0: {
601           PetscInt triSizes_p2[2]  = {1, 2};
602           PetscInt triPoints_p2[3] = {0, 1, 2};
603 
604           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
605           PetscCall(PetscArraycpy(sizes, triSizes_p2, 2));
606           PetscCall(PetscArraycpy(points, triPoints_p2, 3));
607           break;
608         }
609         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular hybrid mesh on 2 procs", user->testNum);
610         }
611       } else if (dim == 2 && !cellSimplex && !cellHybrid && size == 2) {
612         switch (user->testNum) {
613         case 0: {
614           PetscInt quadSizes_p2[2]  = {1, 1};
615           PetscInt quadPoints_p2[2] = {0, 1};
616 
617           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
618           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
619           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));
620           break;
621         }
622         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral mesh on 2 procs", user->testNum);
623         }
624       } else if (dim == 2 && !cellSimplex && cellHybrid && size == 2) {
625         switch (user->testNum) {
626         case 0: {
627           PetscInt quadSizes_p2[2]  = {1, 2};
628           PetscInt quadPoints_p2[3] = {0, 1, 2};
629 
630           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
631           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
632           PetscCall(PetscArraycpy(points, quadPoints_p2, 3));
633           break;
634         }
635         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral hybrid mesh on 2 procs", user->testNum);
636         }
637       } else if (dim == 3 && cellSimplex && !cellHybrid && size == 2) {
638         switch (user->testNum) {
639         case 0: {
640           PetscInt tetSizes_p2[2]  = {1, 1};
641           PetscInt tetPoints_p2[2] = {0, 1};
642 
643           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
644           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
645           PetscCall(PetscArraycpy(points, tetPoints_p2, 2));
646           break;
647         }
648         case 1: {
649           PetscInt tetSizes_p2[2]  = {1, 1};
650           PetscInt tetPoints_p2[2] = {0, 1};
651 
652           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
653           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
654           PetscCall(PetscArraycpy(points, tetPoints_p2, 2));
655           break;
656         }
657         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for tetrahedral mesh on 2 procs", user->testNum);
658         }
659       } else if (dim == 3 && cellSimplex && cellHybrid && size == 2) {
660         switch (user->testNum) {
661         case 0: {
662           PetscInt tetSizes_p2[2]  = {1, 2};
663           PetscInt tetPoints_p2[3] = {0, 1, 2};
664 
665           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
666           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
667           PetscCall(PetscArraycpy(points, tetPoints_p2, 3));
668           break;
669         }
670         case 1: {
671           PetscInt tetSizes_p2[2]  = {3, 4};
672           PetscInt tetPoints_p2[7] = {0, 3, 5, 1, 2, 4, 6};
673 
674           PetscCall(PetscMalloc2(2, &sizes, 7, &points));
675           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
676           PetscCall(PetscArraycpy(points, tetPoints_p2, 7));
677           break;
678         }
679         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for tetrahedral hybrid mesh on 2 procs", user->testNum);
680         }
681       } else if (dim == 3 && !cellSimplex && !cellHybrid && size == 2) {
682         switch (user->testNum) {
683         case 0: {
684           PetscInt hexSizes_p2[2]  = {1, 1};
685           PetscInt hexPoints_p2[2] = {0, 1};
686 
687           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
688           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
689           PetscCall(PetscArraycpy(points, hexPoints_p2, 2));
690           break;
691         }
692         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral mesh on 2 procs", user->testNum);
693         }
694       } else if (dim == 3 && !cellSimplex && cellHybrid && size == 2) {
695         switch (user->testNum) {
696         case 0: {
697           PetscInt hexSizes_p2[2]  = {1, 1};
698           PetscInt hexPoints_p2[2] = {0, 1};
699 
700           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
701           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
702           PetscCall(PetscArraycpy(points, hexPoints_p2, 2));
703           break;
704         }
705         case 1: {
706           PetscInt hexSizes_p2[2]  = {5, 4};
707           PetscInt hexPoints_p2[9] = {3, 4, 5, 7, 8, 0, 1, 2, 6};
708 
709           PetscCall(PetscMalloc2(2, &sizes, 9, &points));
710           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
711           PetscCall(PetscArraycpy(points, hexPoints_p2, 9));
712           break;
713         }
714         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral hybrid mesh on 2 procs", user->testNum);
715         }
716       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
717     }
718     PetscCall(DMPlexGetPartitioner(*dm, &part));
719     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
720     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
721     PetscCall(PetscFree2(sizes, points));
722   } else {
723     PetscPartitioner part;
724 
725     PetscCall(DMPlexGetPartitioner(*dm, &part));
726     PetscCall(PetscPartitionerSetFromOptions(part));
727   }
728   {
729     DM pdm = NULL;
730 
731     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
732     if (pdm) {
733       PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
734       PetscCall(DMDestroy(dm));
735       *dm = pdm;
736     }
737   }
738   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
739   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
740   PetscCall(DMSetFromOptions(*dm));
741   if (user->uninterpolate || user->reinterpolate) {
742     DM udm = NULL;
743 
744     PetscCall(DMPlexUninterpolate(*dm, &udm));
745     PetscCall(DMPlexCopyCoordinates(*dm, udm));
746     PetscCall(DMDestroy(dm));
747     *dm = udm;
748   }
749   if (user->reinterpolate) {
750     DM idm = NULL;
751 
752     PetscCall(DMPlexInterpolate(*dm, &idm));
753     PetscCall(DMPlexCopyCoordinates(*dm, idm));
754     PetscCall(DMDestroy(dm));
755     *dm = idm;
756   }
757   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
758   PetscCall(PetscObjectSetName((PetscObject)*dm, "Hybrid Mesh"));
759   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
760   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "hyb_"));
761   PetscCall(DMSetFromOptions(*dm));
762   PetscFunctionReturn(0);
763 }
764 
765 int main(int argc, char **argv) {
766   DM     dm;
767   AppCtx user; /* user-defined work context */
768 
769   PetscFunctionBeginUser;
770   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
771   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
772   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
773   PetscCall(DMDestroy(&dm));
774   PetscCall(PetscFinalize());
775   return 0;
776 }
777 
778 /*TEST
779 
780   # 1D Simplex 29-31
781   testset:
782     args: -dim 1 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
783     test:
784       suffix: 29
785     test:
786       suffix: 30
787       args: -dm_refine 1
788     test:
789       suffix: 31
790       args: -dm_refine 5
791 
792   # 2D Simplex 0-3
793   testset:
794     args: -dim 2 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
795     test:
796       suffix: 0
797     test:
798       suffix: 1
799       args: -dm_refine 1
800     test:
801       suffix: 2
802       nsize: 2
803     test:
804       suffix: 3
805       nsize: 2
806       args: -dm_refine 1
807     test:
808       suffix: 32
809       args: -dm_refine 1 -uninterpolate
810     test:
811       suffix: 33
812       nsize: 2
813       args: -dm_refine 1 -uninterpolate
814     test:
815       suffix: 34
816       nsize: 2
817       args: -dm_refine 3 -uninterpolate
818 
819   # 2D Hybrid Simplex 4-7
820   testset:
821     args: -dim 2 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
822     test:
823       suffix: 4
824     test:
825       suffix: 5
826       args: -dm_refine 1
827     test:
828       suffix: 6
829       nsize: 2
830     test:
831       suffix: 7
832       nsize: 2
833       args: -dm_refine 1
834     test:
835       suffix: 24
836       args: -test_num 1 -dm_refine 1
837 
838   # 2D Quad 12-13
839   testset:
840     args: -dim 2 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
841     test:
842       suffix: 12
843       args: -dm_refine 1
844     test:
845       suffix: 13
846       nsize: 2
847       args: -dm_refine 1
848 
849   # 2D Hybrid Quad 27-28
850   testset:
851     args: -dim 2 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
852     test:
853       suffix: 27
854       args: -dm_refine 1
855     test:
856       suffix: 28
857       nsize: 2
858       args: -dm_refine 1
859 
860   # 3D Simplex 8-11
861   testset:
862     args: -dim 3 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
863     test:
864       suffix: 8
865       args: -dm_refine 1
866     test:
867       suffix: 9
868       nsize: 2
869       args: -dm_refine 1
870     test:
871       suffix: 10
872       args: -test_num 1 -dm_refine 1
873     test:
874       suffix: 11
875       nsize: 2
876       args: -test_num 1 -dm_refine 1
877     test:
878       suffix: 25
879       args: -test_num 2 -dm_refine 2
880 
881   # 3D Hybrid Simplex 16-19
882   testset:
883     args: -dim 3 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
884     test:
885       suffix: 16
886       args: -dm_refine 1
887     test:
888       suffix: 17
889       nsize: 2
890       args: -dm_refine 1
891     test:
892       suffix: 18
893       args: -test_num 1 -dm_refine 1
894     test:
895       suffix: 19
896       nsize: 2
897       args: -test_num 1 -dm_refine 1
898 
899   # 3D Hex 14-15
900   testset:
901     args: -dim 3 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
902     test:
903       suffix: 14
904       args: -dm_refine 1
905     test:
906       suffix: 15
907       nsize: 2
908      args: -dm_refine 1
909     test:
910       suffix: 26
911       args: -test_num 1 -dm_refine 2
912 
913   # 3D Hybrid Hex 20-23
914   testset:
915     args: -dim 3 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
916     test:
917       suffix: 20
918       args: -dm_refine 1
919     test:
920       suffix: 21
921       nsize: 2
922       args: -dm_refine 1
923     test:
924       suffix: 22
925       args: -test_num 1 -dm_refine 1
926     test:
927       suffix: 23
928       nsize: 2
929       args: -test_num 1 -dm_refine 1
930 
931   # Hybrid interpolation
932   #   TODO Setup new tests (like -reinterpolate) that interpolate hybrid cells
933   testset:
934     nsize: 2
935     args: -test_partition 0 -petscpartitioner_type simple -dm_view -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
936     test:
937       suffix: hybint_2d_0
938       args: -dim 2 -dm_refine 2
939     test:
940       suffix: hybint_2d_1
941       args: -dim 2 -dm_refine 2 -test_num 1
942     test:
943       suffix: hybint_3d_0
944       args: -dim 3 -dm_refine 1
945     test:
946       suffix: hybint_3d_1
947       args: -dim 3 -dm_refine 1 -test_num 1
948 
949 TEST*/
950