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