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: {

 63:   DMCreate(comm, dm);
 64:   DMSetType(*dm, DMPLEX);
 65:   DMSetFromOptions(*dm);
 66:   DMViewFromOptions(*dm, NULL, "-dm_view");
 67:   return(0);
 68: }

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

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

 87:     if (vertex < vStart || vertex >= vEnd) continue;
 88:     closure[Nv++] = vertex;
 89:   }
 90:   if (Nv != DMPolytopeTypeGetNumVertices(ct)) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Cell %D has %D vertices != %D vertices in a %s", cell, Nv, DMPolytopeTypeGetNumVertices(ct), DMPolytopeTypes[ct]);
 91:   arrVerts = DMPolytopeTypeGetVertexArrangment(ct, o);
 92:   for (v = 0; v < Nv; ++v) {
 93:     if (closure[v] != arrVerts[v]+vStart) SETERRQ5(comm, PETSC_ERR_ARG_WRONG, "Cell %D vertex[%D]: %D should be %D for arrangement %D", cell, v, closure[v], arrVerts[v]+vStart, o);
 94:   }
 95:   DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);
 96:   return(0);
 97: }

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

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

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

138: static PetscErrorCode GenerateArrangments(DM dm, AppCtx *user)
139: {
140:   DM             odm;
141:   DMPolytopeType ct;
142:   PetscInt       No, o;
143:   const char    *name;

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

163: static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user)
164: {
165:   DM              dm1, dm2;
166:   DMPolytopeType  ct;
167:   const PetscInt *refcone, *cone;
168:   PetscInt        No, o1, o2, o3, o4;
169:   PetscBool       equal;
170:   const char     *name;
171:   PetscErrorCode  ierr;

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

211: static PetscErrorCode VerifyInverse(DM dm, AppCtx *user)
212: {
213:   DM              dm1, dm2;
214:   DMPolytopeType  ct;
215:   const PetscInt *refcone, *cone;
216:   PetscInt        No, o, oi, o2;
217:   PetscBool       equal;
218:   const char     *name;
219:   PetscErrorCode  ierr;

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

256: /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */
257: static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user)
258: {
259:   DMPlexTransform tr, otr;
260:   DMPolytopeType  ct;
261:   DMPolytopeType *rct;
262:   const PetscInt *cone, *ornt, *ocone, *oornt;
263:   PetscInt       *rsize, *rcone, *rornt;
264:   PetscInt        Nct, n, oi, debug = 0;
265:   PetscErrorCode  ierr;

268:   DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr);
269:   DMPlexTransformSetDM(tr, dm);
270:   DMPlexTransformSetFromOptions(tr);
271:   DMPlexTransformSetUp(tr);

273:   DMPlexTransformCreate(PetscObjectComm((PetscObject) odm), &otr);
274:   DMPlexTransformSetDM(otr, odm);
275:   DMPlexTransformSetFromOptions(otr);
276:   DMPlexTransformSetUp(otr);

278:   DMPlexGetCellType(dm, p, &ct);
279:   DMPlexGetCone(dm, p, &cone);
280:   DMPlexGetConeOrientation(dm, p, &ornt);
281:   DMPlexGetCone(odm, p, &ocone);
282:   DMPlexGetConeOrientation(odm, p, &oornt);
283:   oi   = DMPolytopeTypeComposeOrientationInv(ct, 0, o);
284:   if (user->printTable) {PetscPrintf(PETSC_COMM_SELF, "Orientation %D\n", oi);}

286:   DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt);
287:   for (n = 0; n < Nct; ++n) {
288:     DMPolytopeType ctNew = rct[n];
289:     PetscInt       r, ro;

291:     if (debug) {PetscPrintf(PETSC_COMM_SELF, "  Checking type %s\n", DMPolytopeTypes[ctNew]);}
292:     for (r = 0; r < rsize[n]; ++r) {
293:       const PetscInt *qcone, *qornt, *oqcone, *oqornt;
294:       PetscInt        pNew, opNew, oo, pr, fo;
295:       PetscBool       restore = PETSC_TRUE;

297:       DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew);
298:       DMPlexTransformGetCone(tr, pNew, &qcone, &qornt);
299:       if (debug) {
300:         PetscInt c;

302:         PetscPrintf(PETSC_COMM_SELF, "    Checking replica %D (%D)\n      Original Cone", r, pNew);
303:         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) {PetscPrintf(PETSC_COMM_SELF, " %D (%D)", qcone[c], qornt[c]);}
304:         PetscPrintf(PETSC_COMM_SELF, "\n");
305:       }
306:       for (ro = 0; ro < rsize[n]; ++ro) {
307:         PetscBool found;

309:         DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew);
310:         DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt);
311:         DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found);
312:         if (found) break;
313:         DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);
314:       }
315:       if (debug) {
316:         PetscInt c;

318:         PetscPrintf(PETSC_COMM_SELF, "    Checking transform replica %D (%D) (%D)\n      Transform Cone", ro, opNew, o);
319:         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) {PetscPrintf(PETSC_COMM_SELF, " %D (%D)", oqcone[c], oqornt[c]);}
320:         PetscPrintf(PETSC_COMM_SELF, "\n");
321:         PetscPrintf(PETSC_COMM_SELF, "    Matched %D\n", oo);
322:       }
323:       if (ro == rsize[n]) {
324:         /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */
325:         if (ct == DM_POLYTOPE_TETRAHEDRON) {
326:           /* The segment in a tetrahedron does not map into itself under the group action */
327:           if (ctNew == DM_POLYTOPE_SEGMENT) {restore = PETSC_FALSE; ro = r; oo = 0;}
328:           /* The last four interior faces do not map into themselves under the group action */
329:           if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) {restore = PETSC_FALSE; ro = r; oo = 0;}
330:           /* The last four interior faces do not map into themselves under the group action */
331:           if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) {restore = PETSC_FALSE; ro = r; oo = 0;}
332:         }
333:         if (restore) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to find matching %s %D orientation for cell orientation %D", DMPolytopeTypes[ctNew], r, o);
334:       }
335:       if (user->printTable) {PetscPrintf(PETSC_COMM_SELF, "%D, %D, ", ro, oo);}
336:       DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo);
337:       if (!user->printTable) {
338:         if (pr != ro) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong replica %D != %D", pr, ro);
339:         if (fo != oo)  SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong orientation %D != %D", fo, oo);
340:       }
341:       DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt);
342:       if (restore) {DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);}
343:     }
344:     if (user->printTable) {PetscPrintf(PETSC_COMM_SELF, "\n");}
345:   }
346:   DMPlexTransformDestroy(&tr);
347:   DMPlexTransformDestroy(&otr);
348:   return(0);
349: }

351: static PetscErrorCode RefineArrangments(DM dm, AppCtx *user)
352: {
353:   DM             odm, rdm;
354:   DMPolytopeType ct;
355:   PetscInt       No, o;
356:   const char    *name;

360:   if (!user->refArr) return(0);
361:   PetscObjectGetName((PetscObject) dm, &name);
362:   DMPlexGetCellType(dm, 0, &ct);
363:   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
364:   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
365:     if (ignoreOrnt(user, o)) continue;
366:     CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);
367:     if (user->initOrnt) {ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE);}
368:     ReorientCell(odm, 0, o, PETSC_TRUE);
369:     DMViewFromOptions(odm, NULL, "-orig_dm_view");
370:     DMRefine(odm, MPI_COMM_NULL, &rdm);
371:     DMSetFromOptions(rdm);
372:     PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);
373:     DMViewFromOptions(rdm, NULL, "-ref_dm_view");
374:     CheckSubcells(dm, odm, 0, o, user);
375:     DMDestroy(&odm);
376:     DMDestroy(&rdm);
377:   }
378:   return(0);
379: }

381: int main(int argc, char **argv)
382: {
383:   DM             dm;
384:   AppCtx         user;

387:   PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
388:   ProcessOptions(PETSC_COMM_WORLD, &user);
389:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
390:   if (user.initOrnt) {
391:     ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE);
392:     DMViewFromOptions(dm, NULL, "-ornt_dm_view");
393:   }
394:   GenerateArrangments(dm, &user);
395:   VerifyCayleyTable(dm, &user);
396:   VerifyInverse(dm, &user);
397:   RefineArrangments(dm, &user);
398:   DMDestroy(&dm);
399:   PetscFinalize();
400:   return ierr;
401: }

403: /*TEST

405:   test:
406:     suffix: segment
407:     args: -dm_plex_reference_cell_domain -dm_plex_cell segment -gen_arrangements \
408:           -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 -dm_plex_view_tikzscale 0.5

410:   test:
411:     suffix: triangle
412:     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -gen_arrangements \
413:           -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 0.5

415:   test:
416:     suffix: quadrilateral
417:     args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -gen_arrangements \
418:           -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 0.5

420:   test:
421:     suffix: tensor_segment
422:     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -gen_arrangements \
423:           -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 0.5

425:   test:
426:     suffix: tetrahedron
427:     args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -gen_arrangements \
428:           -gen_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 0.5

430:   test:
431:     suffix: hexahedron
432:     args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -gen_arrangements \
433:           -gen_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 0.3

435:   test:
436:     suffix: triangular_prism
437:     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -gen_arrangements \
438:           -gen_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 0.5

440:   test:
441:     suffix: tensor_triangular_prism
442:     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -gen_arrangements \
443:           -gen_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 0.5

445:   test:
446:     suffix: tensor_quadrilateral_prism
447:     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -gen_arrangements \
448:           -gen_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 0.5

450:   test:
451:     suffix: pyramid
452:     args: -dm_plex_reference_cell_domain -dm_plex_cell pyramid -gen_arrangements \
453:           -gen_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 0.5

455:   test:
456:     suffix: ref_segment
457:     args: -dm_plex_reference_cell_domain -dm_plex_cell segment -ref_arrangements -dm_plex_check_all \
458:           -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 -dm_plex_view_tikzscale 1.0

460:   test:
461:     suffix: ref_triangle
462:     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -ref_arrangements -dm_plex_check_all \
463:           -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

465:   test:
466:     suffix: ref_quadrilateral
467:     args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -ref_arrangements -dm_plex_check_all \
468:           -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

470:   test:
471:     suffix: ref_tensor_segment
472:     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -ref_arrangements -dm_plex_check_all \
473:           -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

475:   test:
476:     suffix: ref_tetrahedron
477:     args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -ref_arrangements -dm_plex_check_all \
478:           -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

480:   test:
481:     suffix: ref_hexahedron
482:     args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -ref_arrangements -dm_plex_check_all \
483:           -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

485:   test:
486:     suffix: ref_triangular_prism
487:     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -ref_arrangements -dm_plex_check_all \
488:           -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

490:   test:
491:     suffix: ref_tensor_triangular_prism
492:     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -ref_arrangements -dm_plex_check_all \
493:           -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

495:   test:
496:     suffix: ref_tensor_quadrilateral_prism
497:     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -ref_arrangements -dm_plex_check_all \
498:           -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

500:   # ToBox should recreate the coordinate space since the cell shape changes
501:   testset:
502:     args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all

504:     test:
505:       suffix: tobox_triangle
506:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
507:             -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

509:     test:
510:       suffix: tobox_tensor_segment
511:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
512:             -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

514:     test:
515:       suffix: tobox_tetrahedron
516:       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
517:             -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

519:     test:
520:       suffix: tobox_triangular_prism
521:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
522:             -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

524:     test:
525:       suffix: tobox_tensor_triangular_prism
526:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
527:             -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

529:     test:
530:       suffix: tobox_tensor_quadrilateral_prism
531:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \
532:             -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

534:   testset:
535:     args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all

537:     test:
538:       suffix: alfeld_triangle
539:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
540:             -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

542:     test:
543:       suffix: alfeld_tetrahedron
544:       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
545:             -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

547:   testset:
548:     args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all

550:     # This splits edge 1 of the triangle, and reflects about the added edge
551:     test:
552:       suffix: sbr_triangle_0
553:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \
554:             -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

556:     # This splits edge 0 of the triangle, and reflects about the added edge
557:     test:
558:       suffix: sbr_triangle_1
559:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -init_ornt 1 -ornts -3,0 \
560:             -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

562:     # This splits edge 2 of the triangle, and reflects about the added edge
563:     test:
564:       suffix: sbr_triangle_2
565:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -init_ornt 2 -ornts -1,0 \
566:             -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

568:     # This splits edges 1 and 2 of the triangle
569:     test:
570:       suffix: sbr_triangle_3
571:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 0 \
572:             -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

574:     # This splits edges 1 and 0 of the triangle
575:     test:
576:       suffix: sbr_triangle_4
577:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 0 \
578:             -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

580:     # This splits edges 0 and 1 of the triangle
581:     test:
582:       suffix: sbr_triangle_5
583:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -init_ornt 1 -ornts 0 \
584:             -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

586:     # This splits edges 0 and 2 of the triangle
587:     test:
588:       suffix: sbr_triangle_6
589:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -init_ornt 1 -ornts 0 \
590:             -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

592:     # This splits edges 2 and 0 of the triangle
593:     test:
594:       suffix: sbr_triangle_7
595:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -init_ornt 2 -ornts 0 \
596:             -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

598:     # This splits edges 2 and 1 of the triangle
599:     test:
600:       suffix: sbr_triangle_8
601:       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -init_ornt 2 -ornts 0 \
602:             -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

604:   testset:
605:     args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all

607:     test:
608:       suffix: bl_tensor_segment
609:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
610:             -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

612:     # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check
613:     test:
614:       suffix: bl_tensor_triangular_prism
615:       requires: TODO
616:       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
617:             -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

619: TEST*/