Actual source code: ex11.c

  1: static char help[] = "Mesh Orientation Tutorial\n\n";

  3: #include <petscdmplex.h>
  4: #include <petscdmplextransform.h>

  6: typedef struct {
  7:   PetscBool genArr;        /* Generate all possible cell arrangements */
  8:   PetscBool refArr;        /* Refine all possible cell arrangements */
  9:   PetscBool printTable;    /* Print the CAyley table */
 10:   PetscInt  orntBounds[2]; /* Bounds for the orientation check */
 11:   PetscInt  numOrnt;       /* Number of specific orientations specified, or -1 for all orientations */
 12:   PetscInt  ornts[48];     /* Specific orientations if specified */
 13:   PetscInt  initOrnt;      /* Initial orientation for starting mesh */
 14: } AppCtx;

 16: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 17: {
 18:   PetscInt       n = 2;
 19:   PetscBool      flg;

 23:   options->genArr        = PETSC_FALSE;
 24:   options->refArr        = PETSC_FALSE;
 25:   options->printTable    = PETSC_FALSE;
 26:   options->orntBounds[0] = PETSC_MIN_INT;
 27:   options->orntBounds[1] = PETSC_MAX_INT;
 28:   options->numOrnt       = -1;
 29:   options->initOrnt      = 0;

 31:   PetscOptionsBegin(comm, "", "Mesh Orientation Tutorials Options", "DMPLEX");
 32:   PetscOptionsBool("-gen_arrangements", "Flag for generating all arrangements of the cell", "ex11.c", options->genArr, &options->genArr, NULL);
 33:   PetscOptionsBool("-ref_arrangements", "Flag for refining all arrangements of the cell", "ex11.c", options->refArr, &options->refArr, NULL);
 34:   PetscOptionsBool("-print_table", "Print the Cayley table", "ex11.c", options->printTable, &options->printTable, NULL);
 35:   PetscOptionsIntArray("-ornt_bounds", "Bounds for orientation checks", "ex11.c", options->orntBounds, &n, NULL);
 36:   n    = 48;
 37:   PetscOptionsIntArray("-ornts", "Specific orientations for checks", "ex11.c", options->ornts, &n, &flg);
 38:   if (flg) {
 39:     options->numOrnt = n;
 40:     PetscSortInt(n, options->ornts);
 41:   }
 42:   PetscOptionsInt("-init_ornt", "Initial orientation for starting mesh", "ex11.c", options->initOrnt, &options->initOrnt, NULL);
 43:   PetscOptionsEnd();
 44:   return 0;
 45: }

 47: static PetscBool ignoreOrnt(AppCtx *user, PetscInt o)
 48: {
 49:   PetscInt       loc;

 52:   if (user->numOrnt < 0) return PETSC_FALSE;
 53:   PetscFindInt(o, user->numOrnt, user->ornts, &loc);
 54:   if (loc < 0 || ierr)   return PETSC_TRUE;
 55:   return PETSC_FALSE;
 56: }

 58: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
 59: {
 61:   DMCreate(comm, dm);
 62:   DMSetType(*dm, DMPLEX);
 63:   DMSetFromOptions(*dm);
 64:   DMViewFromOptions(*dm, NULL, "-dm_view");
 65:   return 0;
 66: }

 68: static PetscErrorCode CheckCellVertices(DM dm, PetscInt cell, PetscInt o)
 69: {
 70:   DMPolytopeType  ct;
 71:   const PetscInt *arrVerts;
 72:   PetscInt       *closure = NULL;
 73:   PetscInt        Ncl, cl, Nv, vStart, vEnd, v;
 74:   MPI_Comm        comm;

 77:   PetscObjectGetComm((PetscObject) dm, &comm);
 78:   DMPlexGetCellType(dm, cell, &ct);
 79:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
 80:   DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);
 81:   for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) {
 82:     const PetscInt vertex = closure[cl];

 84:     if (vertex < vStart || vertex >= vEnd) continue;
 85:     closure[Nv++] = vertex;
 86:   }
 88:   arrVerts = DMPolytopeTypeGetVertexArrangment(ct, o);
 89:   for (v = 0; v < Nv; ++v) {
 91:   }
 92:   DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);
 93:   return 0;
 94: }

 96: /* Transform cell with group operation o */
 97: static PetscErrorCode ReorientCell(DM dm, PetscInt cell, PetscInt o, PetscBool swapCoords)
 98: {
 99:   DM              cdm;
100:   Vec             coordinates;
101:   PetscScalar    *coords, *ccoords = NULL;
102:   PetscInt       *closure = NULL;
103:   PetscInt        cdim, d, Nc, Ncl, cl, vStart, vEnd, Nv;

105:   /* Change vertex coordinates so that it plots as we expect */
106:   DMGetCoordinateDM(dm, &cdm);
107:   DMGetCoordinateDim(dm, &cdim);
108:   DMGetCoordinatesLocal(dm, &coordinates);
109:   DMPlexVecGetClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords);
110:   /* Reorient cone */
111:   DMPlexOrientPoint(dm, cell, o);
112:   /* Finish resetting coordinates */
113:   if (swapCoords) {
114:     DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
115:     VecGetArrayWrite(coordinates, &coords);
116:     DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);
117:     for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) {
118:       const PetscInt vertex = closure[cl];
119:       PetscScalar   *vcoords;

121:       if (vertex < vStart || vertex >= vEnd) continue;
122:       DMPlexPointLocalRef(cdm, vertex, coords, &vcoords);
123:       for (d = 0; d < cdim; ++d) vcoords[d] = ccoords[Nv*cdim + d];
124:       ++Nv;
125:     }
126:     DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);
127:     VecRestoreArrayWrite(coordinates, &coords);
128:   }
129:   DMPlexVecRestoreClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords);
130:   return 0;
131: }

133: static PetscErrorCode GenerateArrangments(DM dm, AppCtx *user)
134: {
135:   DM             odm;
136:   DMPolytopeType ct;
137:   PetscInt       No, o;
138:   const char    *name;

141:   if (!user->genArr) return 0;
142:   PetscObjectGetName((PetscObject) dm, &name);
143:   DMPlexGetCellType(dm, 0, &ct);
144:   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
145:   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
146:     if (ignoreOrnt(user, o)) continue;
147:     CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);
148:     ReorientCell(odm, 0, o, PETSC_TRUE);
149:     PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);
150:     DMViewFromOptions(odm, NULL, "-gen_dm_view");
151:     CheckCellVertices(odm, 0, o);
152:     DMDestroy(&odm);
153:   }
154:   return 0;
155: }

157: static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user)
158: {
159:   DM              dm1, dm2;
160:   DMPolytopeType  ct;
161:   const PetscInt *refcone, *cone;
162:   PetscInt        No, o1, o2, o3, o4;
163:   PetscBool       equal;
164:   const char     *name;

167:   if (!user->genArr) return 0;
168:   PetscObjectGetName((PetscObject) dm, &name);
169:   DMPlexGetCellType(dm, 0, &ct);
170:   DMPlexGetCone(dm, 0, &refcone);
171:   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
172:   if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "Cayley Table for %s\n", DMPolytopeTypes[ct]);
173:   for (o1 = PetscMax(-No, user->orntBounds[0]); o1 < PetscMin(No, user->orntBounds[1]); ++o1) {
174:     for (o2 = PetscMax(-No, user->orntBounds[0]); o2 < PetscMin(No, user->orntBounds[1]); ++o2) {
175:       CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1);
176:       DMPlexOrientPoint(dm1, 0, o2);
177:       DMPlexCheckFaces(dm1, 0);
178:       DMPlexOrientPoint(dm1, 0, o1);
179:       DMPlexCheckFaces(dm1, 0);
180:       o3   = DMPolytopeTypeComposeOrientation(ct, o1, o2);
181:       /* First verification */
182:       CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2);
183:       DMPlexOrientPoint(dm2, 0, o3);
184:       DMPlexCheckFaces(dm2, 0);
185:       DMPlexEqual(dm1, dm2, &equal);
186:       if (!equal) {
187:         DMViewFromOptions(dm1, NULL, "-error_dm_view");
188:         DMViewFromOptions(dm2, NULL, "-error_dm_view");
189:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D != %D", DMPolytopeTypes[ct], o1, o2, o3);
190:       }
191:       /* Second verification */
192:       DMPlexGetCone(dm1, 0, &cone);
193:       DMPolytopeGetOrientation(ct, refcone, cone, &o4);
194:       if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "%D, ", o4);
196:       DMDestroy(&dm1);
197:       DMDestroy(&dm2);
198:     }
199:     if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "\n");
200:   }
201:   return 0;
202: }

204: static PetscErrorCode VerifyInverse(DM dm, AppCtx *user)
205: {
206:   DM              dm1, dm2;
207:   DMPolytopeType  ct;
208:   const PetscInt *refcone, *cone;
209:   PetscInt        No, o, oi, o2;
210:   PetscBool       equal;
211:   const char     *name;

214:   if (!user->genArr) return 0;
215:   PetscObjectGetName((PetscObject) dm, &name);
216:   DMPlexGetCellType(dm, 0, &ct);
217:   DMPlexGetCone(dm, 0, &refcone);
218:   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
219:   if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "Inverse table for %s\n", DMPolytopeTypes[ct]);
220:   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
221:     if (ignoreOrnt(user, o)) continue;
222:     oi   = DMPolytopeTypeComposeOrientationInv(ct, 0, o);
223:     CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1);
224:     DMPlexOrientPoint(dm1, 0, o);
225:     DMPlexCheckFaces(dm1, 0);
226:     DMPlexOrientPoint(dm1, 0, oi);
227:     DMPlexCheckFaces(dm1, 0);
228:     /* First verification */
229:     CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2);
230:     DMPlexEqual(dm1, dm2, &equal);
231:     if (!equal) {
232:       DMViewFromOptions(dm1, NULL, "-error_dm_view");
233:       DMViewFromOptions(dm2, NULL, "-error_dm_view");
234:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D != 0", DMPolytopeTypes[ct], o, oi);
235:     }
236:     /* Second verification */
237:     DMPlexGetCone(dm1, 0, &cone);
238:     DMPolytopeGetOrientation(ct, refcone, cone, &o2);
239:     if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "%D, ", oi);
241:     DMDestroy(&dm1);
242:     DMDestroy(&dm2);
243:   }
244:   if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "\n");
245:   return 0;
246: }

248: /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */
249: static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user)
250: {
251:   DMPlexTransform tr, otr;
252:   DMPolytopeType  ct;
253:   DMPolytopeType *rct;
254:   const PetscInt *cone, *ornt, *ocone, *oornt;
255:   PetscInt       *rsize, *rcone, *rornt;
256:   PetscInt        Nct, n, oi, debug = 0;

259:   DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr);
260:   DMPlexTransformSetDM(tr, dm);
261:   DMPlexTransformSetFromOptions(tr);
262:   DMPlexTransformSetUp(tr);

264:   DMPlexTransformCreate(PetscObjectComm((PetscObject) odm), &otr);
265:   DMPlexTransformSetDM(otr, odm);
266:   DMPlexTransformSetFromOptions(otr);
267:   DMPlexTransformSetUp(otr);

269:   DMPlexGetCellType(dm, p, &ct);
270:   DMPlexGetCone(dm, p, &cone);
271:   DMPlexGetConeOrientation(dm, p, &ornt);
272:   DMPlexGetCone(odm, p, &ocone);
273:   DMPlexGetConeOrientation(odm, p, &oornt);
274:   oi   = DMPolytopeTypeComposeOrientationInv(ct, 0, o);
275:   if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "Orientation %D\n", oi);

277:   DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt);
278:   for (n = 0; n < Nct; ++n) {
279:     DMPolytopeType ctNew = rct[n];
280:     PetscInt       r, ro;

282:     if (debug) PetscPrintf(PETSC_COMM_SELF, "  Checking type %s\n", DMPolytopeTypes[ctNew]);
283:     for (r = 0; r < rsize[n]; ++r) {
284:       const PetscInt *qcone, *qornt, *oqcone, *oqornt;
285:       PetscInt        pNew, opNew, oo, pr, fo;
286:       PetscBool       restore = PETSC_TRUE;

288:       DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew);
289:       DMPlexTransformGetCone(tr, pNew, &qcone, &qornt);
290:       if (debug) {
291:         PetscInt c;

293:         PetscPrintf(PETSC_COMM_SELF, "    Checking replica %D (%D)\n      Original Cone", r, pNew);
294:         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) PetscPrintf(PETSC_COMM_SELF, " %D (%D)", qcone[c], qornt[c]);
295:         PetscPrintf(PETSC_COMM_SELF, "\n");
296:       }
297:       for (ro = 0; ro < rsize[n]; ++ro) {
298:         PetscBool found;

300:         DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew);
301:         DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt);
302:         DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found);
303:         if (found) break;
304:         DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);
305:       }
306:       if (debug) {
307:         PetscInt c;

309:         PetscPrintf(PETSC_COMM_SELF, "    Checking transform replica %D (%D) (%D)\n      Transform Cone", ro, opNew, o);
310:         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) PetscPrintf(PETSC_COMM_SELF, " %D (%D)", oqcone[c], oqornt[c]);
311:         PetscPrintf(PETSC_COMM_SELF, "\n");
312:         PetscPrintf(PETSC_COMM_SELF, "    Matched %D\n", oo);
313:       }
314:       if (ro == rsize[n]) {
315:         /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */
316:         if (ct == DM_POLYTOPE_TETRAHEDRON) {
317:           /* The segment in a tetrahedron does not map into itself under the group action */
318:           if (ctNew == DM_POLYTOPE_SEGMENT) {restore = PETSC_FALSE; ro = r; oo = 0;}
319:           /* The last four interior faces do not map into themselves under the group action */
320:           if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) {restore = PETSC_FALSE; ro = r; oo = 0;}
321:           /* The last four interior faces do not map into themselves under the group action */
322:           if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) {restore = PETSC_FALSE; ro = r; oo = 0;}
323:         }
325:       }
326:       if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "%D, %D, ", ro, oo);
327:       DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo);
328:       if (!user->printTable) {
331:       }
332:       DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt);
333:       if (restore) DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);
334:     }
335:     if (user->printTable) PetscPrintf(PETSC_COMM_SELF, "\n");
336:   }
337:   DMPlexTransformDestroy(&tr);
338:   DMPlexTransformDestroy(&otr);
339:   return 0;
340: }

342: static PetscErrorCode RefineArrangments(DM dm, AppCtx *user)
343: {
344:   DM             odm, rdm;
345:   DMPolytopeType ct;
346:   PetscInt       No, o;
347:   const char    *name;

350:   if (!user->refArr) return 0;
351:   PetscObjectGetName((PetscObject) dm, &name);
352:   DMPlexGetCellType(dm, 0, &ct);
353:   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
354:   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
355:     if (ignoreOrnt(user, o)) continue;
356:     CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);
357:     if (user->initOrnt) ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE);
358:     ReorientCell(odm, 0, o, PETSC_TRUE);
359:     DMViewFromOptions(odm, NULL, "-orig_dm_view");
360:     DMRefine(odm, MPI_COMM_NULL, &rdm);
361:     DMSetFromOptions(rdm);
362:     PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);
363:     DMViewFromOptions(rdm, NULL, "-ref_dm_view");
364:     CheckSubcells(dm, odm, 0, o, user);
365:     DMDestroy(&odm);
366:     DMDestroy(&rdm);
367:   }
368:   return 0;
369: }

371: int main(int argc, char **argv)
372: {
373:   DM             dm;
374:   AppCtx         user;

376:   PetscInitialize(&argc, &argv, NULL, help);
377:   ProcessOptions(PETSC_COMM_WORLD, &user);
378:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
379:   if (user.initOrnt) {
380:     ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE);
381:     DMViewFromOptions(dm, NULL, "-ornt_dm_view");
382:   }
383:   GenerateArrangments(dm, &user);
384:   VerifyCayleyTable(dm, &user);
385:   VerifyInverse(dm, &user);
386:   RefineArrangments(dm, &user);
387:   DMDestroy(&dm);
388:   PetscFinalize();
389:   return 0;
390: }

392: /*TEST

394:   testset:
395:     args: -dm_coord_space 0 -dm_plex_reference_cell_domain -gen_arrangements \
396:           -gen_dm_view ::ascii_latex -dm_plex_view_tikzscale 0.5

398:     test:
399:       suffix: segment
400:       args: -dm_plex_cell segment \
401:             -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0

403:     test:
404:       suffix: triangle
405:       args: -dm_plex_cell triangle \
406:             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0

408:     test:
409:       suffix: quadrilateral
410:       args: -dm_plex_cell quadrilateral \
411:             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0

413:     test:
414:       suffix: tensor_segment
415:       args: -dm_plex_cell tensor_quad \
416:             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0

418:     test:
419:       suffix: tetrahedron
420:       args: -dm_plex_cell tetrahedron \
421:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

423:     test:
424:       suffix: hexahedron
425:       args: -dm_plex_cell hexahedron \
426:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.3

428:     test:
429:       suffix: triangular_prism
430:       args: -dm_plex_cell triangular_prism \
431:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

433:     test:
434:       suffix: tensor_triangular_prism
435:       args: -dm_plex_cell tensor_triangular_prism \
436:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

438:     test:
439:       suffix: tensor_quadrilateral_prism
440:       args: -dm_plex_cell tensor_quadrilateral_prism \
441:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

443:     test:
444:       suffix: pyramid
445:       args: -dm_plex_cell pyramid \
446:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

448:   testset:
449:     args: -dm_coord_space 0 -dm_plex_reference_cell_domain -ref_arrangements -dm_plex_check_all \
450:           -ref_dm_view ::ascii_latex -dm_plex_view_tikzscale 1.0

452:     test:
453:       suffix: ref_segment
454:       args: -dm_plex_cell segment \
455:             -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0

457:     test:
458:       suffix: ref_triangle
459:       args: -dm_plex_cell triangle \
460:             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0

462:     test:
463:       suffix: ref_quadrilateral
464:       args: -dm_plex_cell quadrilateral \
465:             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0

467:     test:
468:       suffix: ref_tensor_segment
469:       args: -dm_plex_cell tensor_quad \
470:             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0

472:     test:
473:       suffix: ref_tetrahedron
474:       args: -dm_plex_cell tetrahedron \
475:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

477:     test:
478:       suffix: ref_hexahedron
479:       args: -dm_plex_cell hexahedron \
480:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

482:     test:
483:       suffix: ref_triangular_prism
484:       args: -dm_plex_cell triangular_prism \
485:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

487:     test:
488:       suffix: ref_tensor_triangular_prism
489:       args: -dm_plex_cell tensor_triangular_prism \
490:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

492:     test:
493:       suffix: ref_tensor_quadrilateral_prism
494:       args: -dm_plex_cell tensor_quadrilateral_prism \
495:             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0

497:   # ToBox should recreate the coordinate space since the cell shape changes
498:   testset:
499:     args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all

501:     test:
502:       suffix: tobox_triangle
503:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
504:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

506:     test:
507:       suffix: tobox_tensor_segment
508:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
509:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

511:     test:
512:       suffix: tobox_tetrahedron
513:       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
514:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0

516:     test:
517:       suffix: tobox_triangular_prism
518:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
519:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0

521:     test:
522:       suffix: tobox_tensor_triangular_prism
523:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
524:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0

526:     test:
527:       suffix: tobox_tensor_quadrilateral_prism
528:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \
529:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0

531:   testset:
532:     args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all

534:     test:
535:       suffix: alfeld_triangle
536:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
537:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

539:     test:
540:       suffix: alfeld_tetrahedron
541:       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
542:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0

544:   testset:
545:     args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all

547:     # This splits edge 1 of the triangle, and reflects about the added edge
548:     test:
549:       suffix: sbr_triangle_0
550:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \
551:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

553:     # This splits edge 0 of the triangle, and reflects about the added edge
554:     test:
555:       suffix: sbr_triangle_1
556:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -init_ornt 1 -ornts -3,0 \
557:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

559:     # This splits edge 2 of the triangle, and reflects about the added edge
560:     test:
561:       suffix: sbr_triangle_2
562:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -init_ornt 2 -ornts -1,0 \
563:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

565:     # This splits edges 1 and 2 of the triangle
566:     test:
567:       suffix: sbr_triangle_3
568:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 0 \
569:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

571:     # This splits edges 1 and 0 of the triangle
572:     test:
573:       suffix: sbr_triangle_4
574:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 0 \
575:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

577:     # This splits edges 0 and 1 of the triangle
578:     test:
579:       suffix: sbr_triangle_5
580:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -init_ornt 1 -ornts 0 \
581:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

583:     # This splits edges 0 and 2 of the triangle
584:     test:
585:       suffix: sbr_triangle_6
586:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -init_ornt 1 -ornts 0 \
587:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

589:     # This splits edges 2 and 0 of the triangle
590:     test:
591:       suffix: sbr_triangle_7
592:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -init_ornt 2 -ornts 0 \
593:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

595:     # This splits edges 2 and 1 of the triangle
596:     test:
597:       suffix: sbr_triangle_8
598:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -init_ornt 2 -ornts 0 \
599:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

601:   testset:
602:     args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all

604:     test:
605:       suffix: bl_tensor_segment
606:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
607:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 1.0

609:     # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check
610:     test:
611:       suffix: bl_tensor_triangular_prism
612:       requires: TODO
613:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
614:             -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0

616: TEST*/