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