xref: /petsc/src/dm/impls/plex/tutorials/ex11.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
102918637SMatthew G. Knepley static char help[] = "Mesh Orientation Tutorial\n\n";
202918637SMatthew G. Knepley 
302918637SMatthew G. Knepley #include <petscdmplex.h>
402918637SMatthew G. Knepley #include <petscdmplextransform.h>
502918637SMatthew G. Knepley 
602918637SMatthew G. Knepley typedef struct {
702918637SMatthew G. Knepley   PetscBool genArr;        /* Generate all possible cell arrangements */
802918637SMatthew G. Knepley   PetscBool refArr;        /* Refine all possible cell arrangements */
902918637SMatthew G. Knepley   PetscBool printTable;    /* Print the CAyley table */
1002918637SMatthew G. Knepley   PetscInt  orntBounds[2]; /* Bounds for the orientation check */
1102918637SMatthew G. Knepley   PetscInt  numOrnt;       /* Number of specific orientations specified, or -1 for all orientations */
1202918637SMatthew G. Knepley   PetscInt  ornts[48];     /* Specific orientations if specified */
1302918637SMatthew G. Knepley   PetscInt  initOrnt;      /* Initial orientation for starting mesh */
1402918637SMatthew G. Knepley } AppCtx;
1502918637SMatthew G. Knepley 
16*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
1702918637SMatthew G. Knepley   PetscInt  n = 2;
1802918637SMatthew G. Knepley   PetscBool flg;
1902918637SMatthew G. Knepley 
2002918637SMatthew G. Knepley   PetscFunctionBeginUser;
2102918637SMatthew G. Knepley   options->genArr        = PETSC_FALSE;
2202918637SMatthew G. Knepley   options->refArr        = PETSC_FALSE;
2302918637SMatthew G. Knepley   options->printTable    = PETSC_FALSE;
2402918637SMatthew G. Knepley   options->orntBounds[0] = PETSC_MIN_INT;
2502918637SMatthew G. Knepley   options->orntBounds[1] = PETSC_MAX_INT;
2602918637SMatthew G. Knepley   options->numOrnt       = -1;
2702918637SMatthew G. Knepley   options->initOrnt      = 0;
2802918637SMatthew G. Knepley 
29d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Mesh Orientation Tutorials Options", "DMPLEX");
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-gen_arrangements", "Flag for generating all arrangements of the cell", "ex11.c", options->genArr, &options->genArr, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ref_arrangements", "Flag for refining all arrangements of the cell", "ex11.c", options->refArr, &options->refArr, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-print_table", "Print the Cayley table", "ex11.c", options->printTable, &options->printTable, NULL));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-ornt_bounds", "Bounds for orientation checks", "ex11.c", options->orntBounds, &n, NULL));
3402918637SMatthew G. Knepley   n = 48;
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-ornts", "Specific orientations for checks", "ex11.c", options->ornts, &n, &flg));
3602918637SMatthew G. Knepley   if (flg) {
3702918637SMatthew G. Knepley     options->numOrnt = n;
389566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(n, options->ornts));
3902918637SMatthew G. Knepley   }
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-init_ornt", "Initial orientation for starting mesh", "ex11.c", options->initOrnt, &options->initOrnt, NULL));
41d0609cedSBarry Smith   PetscOptionsEnd();
4202918637SMatthew G. Knepley   PetscFunctionReturn(0);
4302918637SMatthew G. Knepley }
4402918637SMatthew G. Knepley 
45*9371c9d4SSatish Balay static PetscBool ignoreOrnt(AppCtx *user, PetscInt o) {
4602918637SMatthew G. Knepley   PetscInt       loc;
4702918637SMatthew G. Knepley   PetscErrorCode ierr;
4802918637SMatthew G. Knepley 
4902918637SMatthew G. Knepley   if (user->numOrnt < 0) return PETSC_FALSE;
5002918637SMatthew G. Knepley   ierr = PetscFindInt(o, user->numOrnt, user->ornts, &loc);
5102918637SMatthew G. Knepley   if (loc < 0 || ierr) return PETSC_TRUE;
5202918637SMatthew G. Knepley   return PETSC_FALSE;
5302918637SMatthew G. Knepley }
5402918637SMatthew G. Knepley 
55*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
5602918637SMatthew G. Knepley   PetscFunctionBeginUser;
579566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
589566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
599566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
609566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
6102918637SMatthew G. Knepley   PetscFunctionReturn(0);
6202918637SMatthew G. Knepley }
6302918637SMatthew G. Knepley 
64*9371c9d4SSatish Balay static PetscErrorCode CheckCellVertices(DM dm, PetscInt cell, PetscInt o) {
6502918637SMatthew G. Knepley   DMPolytopeType  ct;
6602918637SMatthew G. Knepley   const PetscInt *arrVerts;
6702918637SMatthew G. Knepley   PetscInt       *closure = NULL;
6802918637SMatthew G. Knepley   PetscInt        Ncl, cl, Nv, vStart, vEnd, v;
6902918637SMatthew G. Knepley   MPI_Comm        comm;
7002918637SMatthew G. Knepley 
7102918637SMatthew G. Knepley   PetscFunctionBeginUser;
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
739566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cell, &ct));
749566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
759566063dSJacob Faibussowitsch   PetscCall(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure));
7602918637SMatthew G. Knepley   for (cl = 0, Nv = 0; cl < Ncl * 2; cl += 2) {
7702918637SMatthew G. Knepley     const PetscInt vertex = closure[cl];
7802918637SMatthew G. Knepley 
7902918637SMatthew G. Knepley     if (vertex < vStart || vertex >= vEnd) continue;
8002918637SMatthew G. Knepley     closure[Nv++] = vertex;
8102918637SMatthew G. Knepley   }
8263a3b9bcSJacob Faibussowitsch   PetscCheck(Nv == DMPolytopeTypeGetNumVertices(ct), comm, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " has %" PetscInt_FMT " vertices != %" PetscInt_FMT " vertices in a %s", cell, Nv, DMPolytopeTypeGetNumVertices(ct), DMPolytopeTypes[ct]);
8302918637SMatthew G. Knepley   arrVerts = DMPolytopeTypeGetVertexArrangment(ct, o);
8402918637SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
8563a3b9bcSJacob Faibussowitsch     PetscCheck(closure[v] == arrVerts[v] + vStart, comm, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " vertex[%" PetscInt_FMT "]: %" PetscInt_FMT " should be %" PetscInt_FMT " for arrangement %" PetscInt_FMT, cell, v, closure[v], arrVerts[v] + vStart, o);
8602918637SMatthew G. Knepley   }
879566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure));
8802918637SMatthew G. Knepley   PetscFunctionReturn(0);
8902918637SMatthew G. Knepley }
9002918637SMatthew G. Knepley 
9102918637SMatthew G. Knepley /* Transform cell with group operation o */
92*9371c9d4SSatish Balay static PetscErrorCode ReorientCell(DM dm, PetscInt cell, PetscInt o, PetscBool swapCoords) {
9302918637SMatthew G. Knepley   DM           cdm;
9402918637SMatthew G. Knepley   Vec          coordinates;
9502918637SMatthew G. Knepley   PetscScalar *coords, *ccoords = NULL;
9602918637SMatthew G. Knepley   PetscInt    *closure = NULL;
9702918637SMatthew G. Knepley   PetscInt     cdim, d, Nc, Ncl, cl, vStart, vEnd, Nv;
9802918637SMatthew G. Knepley 
9902918637SMatthew G. Knepley   PetscFunctionBegin;
10002918637SMatthew G. Knepley   /* Change vertex coordinates so that it plots as we expect */
1019566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
1029566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
1039566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1049566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords));
10502918637SMatthew G. Knepley   /* Reorient cone */
1069566063dSJacob Faibussowitsch   PetscCall(DMPlexOrientPoint(dm, cell, o));
10702918637SMatthew G. Knepley   /* Finish resetting coordinates */
10802918637SMatthew G. Knepley   if (swapCoords) {
1099566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1109566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(coordinates, &coords));
1119566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure));
11202918637SMatthew G. Knepley     for (cl = 0, Nv = 0; cl < Ncl * 2; cl += 2) {
11302918637SMatthew G. Knepley       const PetscInt vertex = closure[cl];
11402918637SMatthew G. Knepley       PetscScalar   *vcoords;
11502918637SMatthew G. Knepley 
11602918637SMatthew G. Knepley       if (vertex < vStart || vertex >= vEnd) continue;
1179566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(cdm, vertex, coords, &vcoords));
11802918637SMatthew G. Knepley       for (d = 0; d < cdim; ++d) vcoords[d] = ccoords[Nv * cdim + d];
11902918637SMatthew G. Knepley       ++Nv;
12002918637SMatthew G. Knepley     }
1219566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure));
1229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(coordinates, &coords));
12302918637SMatthew G. Knepley   }
1249566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords));
12502918637SMatthew G. Knepley   PetscFunctionReturn(0);
12602918637SMatthew G. Knepley }
12702918637SMatthew G. Knepley 
128*9371c9d4SSatish Balay static PetscErrorCode GenerateArrangments(DM dm, AppCtx *user) {
12902918637SMatthew G. Knepley   DM             odm;
13002918637SMatthew G. Knepley   DMPolytopeType ct;
13102918637SMatthew G. Knepley   PetscInt       No, o;
13202918637SMatthew G. Knepley   const char    *name;
13302918637SMatthew G. Knepley 
13402918637SMatthew G. Knepley   PetscFunctionBeginUser;
13502918637SMatthew G. Knepley   if (!user->genArr) PetscFunctionReturn(0);
1369566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1379566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, 0, &ct));
13802918637SMatthew G. Knepley   No = DMPolytopeTypeGetNumArrangments(ct) / 2;
13902918637SMatthew G. Knepley   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
14002918637SMatthew G. Knepley     if (ignoreOrnt(user, o)) continue;
1419566063dSJacob Faibussowitsch     PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &odm));
1429566063dSJacob Faibussowitsch     PetscCall(ReorientCell(odm, 0, o, PETSC_TRUE));
14363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%s orientation %" PetscInt_FMT "\n", name, o));
1449566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(odm, NULL, "-gen_dm_view"));
1459566063dSJacob Faibussowitsch     PetscCall(CheckCellVertices(odm, 0, o));
1469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&odm));
14702918637SMatthew G. Knepley   }
14802918637SMatthew G. Knepley   PetscFunctionReturn(0);
14902918637SMatthew G. Knepley }
15002918637SMatthew G. Knepley 
151*9371c9d4SSatish Balay static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user) {
15202918637SMatthew G. Knepley   DM              dm1, dm2;
15302918637SMatthew G. Knepley   DMPolytopeType  ct;
15402918637SMatthew G. Knepley   const PetscInt *refcone, *cone;
15502918637SMatthew G. Knepley   PetscInt        No, o1, o2, o3, o4;
15602918637SMatthew G. Knepley   PetscBool       equal;
15702918637SMatthew G. Knepley   const char     *name;
15802918637SMatthew G. Knepley 
15902918637SMatthew G. Knepley   PetscFunctionBeginUser;
16002918637SMatthew G. Knepley   if (!user->genArr) PetscFunctionReturn(0);
1619566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, 0, &ct));
1639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, 0, &refcone));
16402918637SMatthew G. Knepley   No = DMPolytopeTypeGetNumArrangments(ct) / 2;
1659566063dSJacob Faibussowitsch   if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cayley Table for %s\n", DMPolytopeTypes[ct]));
16602918637SMatthew G. Knepley   for (o1 = PetscMax(-No, user->orntBounds[0]); o1 < PetscMin(No, user->orntBounds[1]); ++o1) {
16702918637SMatthew G. Knepley     for (o2 = PetscMax(-No, user->orntBounds[0]); o2 < PetscMin(No, user->orntBounds[1]); ++o2) {
1689566063dSJacob Faibussowitsch       PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm1));
1699566063dSJacob Faibussowitsch       PetscCall(DMPlexOrientPoint(dm1, 0, o2));
1709566063dSJacob Faibussowitsch       PetscCall(DMPlexCheckFaces(dm1, 0));
1719566063dSJacob Faibussowitsch       PetscCall(DMPlexOrientPoint(dm1, 0, o1));
1729566063dSJacob Faibussowitsch       PetscCall(DMPlexCheckFaces(dm1, 0));
1735f80ce2aSJacob Faibussowitsch       o3 = DMPolytopeTypeComposeOrientation(ct, o1, o2);
17402918637SMatthew G. Knepley       /* First verification */
1759566063dSJacob Faibussowitsch       PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm2));
1769566063dSJacob Faibussowitsch       PetscCall(DMPlexOrientPoint(dm2, 0, o3));
1779566063dSJacob Faibussowitsch       PetscCall(DMPlexCheckFaces(dm2, 0));
1789566063dSJacob Faibussowitsch       PetscCall(DMPlexEqual(dm1, dm2, &equal));
17902918637SMatthew G. Knepley       if (!equal) {
1809566063dSJacob Faibussowitsch         PetscCall(DMViewFromOptions(dm1, NULL, "-error_dm_view"));
1819566063dSJacob Faibussowitsch         PetscCall(DMViewFromOptions(dm2, NULL, "-error_dm_view"));
18263a3b9bcSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " != %" PetscInt_FMT, DMPolytopeTypes[ct], o1, o2, o3);
18302918637SMatthew G. Knepley       }
18402918637SMatthew G. Knepley       /* Second verification */
1859566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm1, 0, &cone));
1869566063dSJacob Faibussowitsch       PetscCall(DMPolytopeGetOrientation(ct, refcone, cone, &o4));
18763a3b9bcSJacob Faibussowitsch       if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ", ", o4));
18863a3b9bcSJacob Faibussowitsch       PetscCheck(o3 == o4, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " = %" PetscInt_FMT " != %" PetscInt_FMT, DMPolytopeTypes[ct], o1, o2, o3, o4);
1899566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm1));
1909566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm2));
19102918637SMatthew G. Knepley     }
1929566063dSJacob Faibussowitsch     if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
19302918637SMatthew G. Knepley   }
19402918637SMatthew G. Knepley   PetscFunctionReturn(0);
19502918637SMatthew G. Knepley }
19602918637SMatthew G. Knepley 
197*9371c9d4SSatish Balay static PetscErrorCode VerifyInverse(DM dm, AppCtx *user) {
19802918637SMatthew G. Knepley   DM              dm1, dm2;
19902918637SMatthew G. Knepley   DMPolytopeType  ct;
20002918637SMatthew G. Knepley   const PetscInt *refcone, *cone;
20102918637SMatthew G. Knepley   PetscInt        No, o, oi, o2;
20202918637SMatthew G. Knepley   PetscBool       equal;
20302918637SMatthew G. Knepley   const char     *name;
20402918637SMatthew G. Knepley 
20502918637SMatthew G. Knepley   PetscFunctionBeginUser;
20602918637SMatthew G. Knepley   if (!user->genArr) PetscFunctionReturn(0);
2079566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
2089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, 0, &ct));
2099566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, 0, &refcone));
21002918637SMatthew G. Knepley   No = DMPolytopeTypeGetNumArrangments(ct) / 2;
2119566063dSJacob Faibussowitsch   if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Inverse table for %s\n", DMPolytopeTypes[ct]));
21202918637SMatthew G. Knepley   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
21302918637SMatthew G. Knepley     if (ignoreOrnt(user, o)) continue;
2145f80ce2aSJacob Faibussowitsch     oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o);
2159566063dSJacob Faibussowitsch     PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm1));
2169566063dSJacob Faibussowitsch     PetscCall(DMPlexOrientPoint(dm1, 0, o));
2179566063dSJacob Faibussowitsch     PetscCall(DMPlexCheckFaces(dm1, 0));
2189566063dSJacob Faibussowitsch     PetscCall(DMPlexOrientPoint(dm1, 0, oi));
2199566063dSJacob Faibussowitsch     PetscCall(DMPlexCheckFaces(dm1, 0));
22002918637SMatthew G. Knepley     /* First verification */
2219566063dSJacob Faibussowitsch     PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm2));
2229566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dm1, dm2, &equal));
22302918637SMatthew G. Knepley     if (!equal) {
2249566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dm1, NULL, "-error_dm_view"));
2259566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dm2, NULL, "-error_dm_view"));
22663a3b9bcSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " != 0", DMPolytopeTypes[ct], o, oi);
22702918637SMatthew G. Knepley     }
22802918637SMatthew G. Knepley     /* Second verification */
2299566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm1, 0, &cone));
2309566063dSJacob Faibussowitsch     PetscCall(DMPolytopeGetOrientation(ct, refcone, cone, &o2));
23163a3b9bcSJacob Faibussowitsch     if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ", ", oi));
23263a3b9bcSJacob Faibussowitsch     PetscCheck(o2 == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " = %" PetscInt_FMT " != 0", DMPolytopeTypes[ct], o, oi, o2);
2339566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm1));
2349566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm2));
23502918637SMatthew G. Knepley   }
2369566063dSJacob Faibussowitsch   if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
23702918637SMatthew G. Knepley   PetscFunctionReturn(0);
23802918637SMatthew G. Knepley }
23902918637SMatthew G. Knepley 
24002918637SMatthew G. Knepley /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */
241*9371c9d4SSatish Balay static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user) {
24202918637SMatthew G. Knepley   DMPlexTransform tr, otr;
24302918637SMatthew G. Knepley   DMPolytopeType  ct;
24402918637SMatthew G. Knepley   DMPolytopeType *rct;
24502918637SMatthew G. Knepley   const PetscInt *cone, *ornt, *ocone, *oornt;
24602918637SMatthew G. Knepley   PetscInt       *rsize, *rcone, *rornt;
24702918637SMatthew G. Knepley   PetscInt        Nct, n, oi, debug = 0;
24802918637SMatthew G. Knepley 
24902918637SMatthew G. Knepley   PetscFunctionBeginUser;
2509566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
2519566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetDM(tr, dm));
2529566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetFromOptions(tr));
2539566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetUp(tr));
25402918637SMatthew G. Knepley 
2559566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)odm), &otr));
2569566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetDM(otr, odm));
2579566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetFromOptions(otr));
2589566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetUp(otr));
25902918637SMatthew G. Knepley 
2609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, p, &ct));
2619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, p, &cone));
2629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
2639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(odm, p, &ocone));
2649566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(odm, p, &oornt));
2655f80ce2aSJacob Faibussowitsch   oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o);
26663a3b9bcSJacob Faibussowitsch   if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Orientation %" PetscInt_FMT "\n", oi));
26702918637SMatthew G. Knepley 
2689566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
26902918637SMatthew G. Knepley   for (n = 0; n < Nct; ++n) {
27002918637SMatthew G. Knepley     DMPolytopeType ctNew = rct[n];
27102918637SMatthew G. Knepley     PetscInt       r, ro;
27202918637SMatthew G. Knepley 
2739566063dSJacob Faibussowitsch     if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Checking type %s\n", DMPolytopeTypes[ctNew]));
27402918637SMatthew G. Knepley     for (r = 0; r < rsize[n]; ++r) {
27502918637SMatthew G. Knepley       const PetscInt *qcone, *qornt, *oqcone, *oqornt;
27602918637SMatthew G. Knepley       PetscInt        pNew, opNew, oo, pr, fo;
27702918637SMatthew G. Knepley       PetscBool       restore = PETSC_TRUE;
27802918637SMatthew G. Knepley 
2799566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew));
2809566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformGetCone(tr, pNew, &qcone, &qornt));
28102918637SMatthew G. Knepley       if (debug) {
28202918637SMatthew G. Knepley         PetscInt c;
28302918637SMatthew G. Knepley 
28463a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    Checking replica %" PetscInt_FMT " (%" PetscInt_FMT ")\n      Original Cone", r, pNew));
28563a3b9bcSJacob Faibussowitsch         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT " (%" PetscInt_FMT ")", qcone[c], qornt[c]));
2869566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
28702918637SMatthew G. Knepley       }
28802918637SMatthew G. Knepley       for (ro = 0; ro < rsize[n]; ++ro) {
28902918637SMatthew G. Knepley         PetscBool found;
29002918637SMatthew G. Knepley 
2919566063dSJacob Faibussowitsch         PetscCall(DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew));
2929566063dSJacob Faibussowitsch         PetscCall(DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt));
2939566063dSJacob Faibussowitsch         PetscCall(DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found));
29402918637SMatthew G. Knepley         if (found) break;
2959566063dSJacob Faibussowitsch         PetscCall(DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt));
29602918637SMatthew G. Knepley       }
29702918637SMatthew G. Knepley       if (debug) {
29802918637SMatthew G. Knepley         PetscInt c;
29902918637SMatthew G. Knepley 
30063a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    Checking transform replica %" PetscInt_FMT " (%" PetscInt_FMT ") (%" PetscInt_FMT ")\n      Transform Cone", ro, opNew, o));
30163a3b9bcSJacob Faibussowitsch         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT " (%" PetscInt_FMT ")", oqcone[c], oqornt[c]));
3029566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
30363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    Matched %" PetscInt_FMT "\n", oo));
30402918637SMatthew G. Knepley       }
30502918637SMatthew G. Knepley       if (ro == rsize[n]) {
30602918637SMatthew G. Knepley         /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */
30702918637SMatthew G. Knepley         if (ct == DM_POLYTOPE_TETRAHEDRON) {
30802918637SMatthew G. Knepley           /* The segment in a tetrahedron does not map into itself under the group action */
309*9371c9d4SSatish Balay           if (ctNew == DM_POLYTOPE_SEGMENT) {
310*9371c9d4SSatish Balay             restore = PETSC_FALSE;
311*9371c9d4SSatish Balay             ro      = r;
312*9371c9d4SSatish Balay             oo      = 0;
313*9371c9d4SSatish Balay           }
31402918637SMatthew G. Knepley           /* The last four interior faces do not map into themselves under the group action */
315*9371c9d4SSatish Balay           if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) {
316*9371c9d4SSatish Balay             restore = PETSC_FALSE;
317*9371c9d4SSatish Balay             ro      = r;
318*9371c9d4SSatish Balay             oo      = 0;
319*9371c9d4SSatish Balay           }
32002918637SMatthew G. Knepley           /* The last four interior faces do not map into themselves under the group action */
321*9371c9d4SSatish Balay           if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) {
322*9371c9d4SSatish Balay             restore = PETSC_FALSE;
323*9371c9d4SSatish Balay             ro      = r;
324*9371c9d4SSatish Balay             oo      = 0;
325*9371c9d4SSatish Balay           }
32602918637SMatthew G. Knepley         }
32763a3b9bcSJacob Faibussowitsch         PetscCheck(!restore, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to find matching %s %" PetscInt_FMT " orientation for cell orientation %" PetscInt_FMT, DMPolytopeTypes[ctNew], r, o);
32802918637SMatthew G. Knepley       }
32963a3b9bcSJacob Faibussowitsch       if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ", %" PetscInt_FMT ", ", ro, oo));
3309566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo));
33102918637SMatthew G. Knepley       if (!user->printTable) {
33263a3b9bcSJacob Faibussowitsch         PetscCheck(pr == ro, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong replica %" PetscInt_FMT " != %" PetscInt_FMT, pr, ro);
33363a3b9bcSJacob Faibussowitsch         PetscCheck(fo == oo, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong orientation %" PetscInt_FMT " != %" PetscInt_FMT, fo, oo);
33402918637SMatthew G. Knepley       }
3359566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt));
3369566063dSJacob Faibussowitsch       if (restore) PetscCall(DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt));
33702918637SMatthew G. Knepley     }
3389566063dSJacob Faibussowitsch     if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
33902918637SMatthew G. Knepley   }
3409566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&tr));
3419566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&otr));
34202918637SMatthew G. Knepley   PetscFunctionReturn(0);
34302918637SMatthew G. Knepley }
34402918637SMatthew G. Knepley 
345*9371c9d4SSatish Balay static PetscErrorCode RefineArrangments(DM dm, AppCtx *user) {
34602918637SMatthew G. Knepley   DM             odm, rdm;
34702918637SMatthew G. Knepley   DMPolytopeType ct;
34802918637SMatthew G. Knepley   PetscInt       No, o;
34902918637SMatthew G. Knepley   const char    *name;
35002918637SMatthew G. Knepley 
35102918637SMatthew G. Knepley   PetscFunctionBeginUser;
35202918637SMatthew G. Knepley   if (!user->refArr) PetscFunctionReturn(0);
3539566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
3549566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, 0, &ct));
35502918637SMatthew G. Knepley   No = DMPolytopeTypeGetNumArrangments(ct) / 2;
35602918637SMatthew G. Knepley   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
35702918637SMatthew G. Knepley     if (ignoreOrnt(user, o)) continue;
3589566063dSJacob Faibussowitsch     PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &odm));
3599566063dSJacob Faibussowitsch     if (user->initOrnt) PetscCall(ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE));
3609566063dSJacob Faibussowitsch     PetscCall(ReorientCell(odm, 0, o, PETSC_TRUE));
3619566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(odm, NULL, "-orig_dm_view"));
3629566063dSJacob Faibussowitsch     PetscCall(DMRefine(odm, MPI_COMM_NULL, &rdm));
3639566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(rdm));
36463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%s orientation %" PetscInt_FMT "\n", name, o));
3659566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(rdm, NULL, "-ref_dm_view"));
3669566063dSJacob Faibussowitsch     PetscCall(CheckSubcells(dm, odm, 0, o, user));
3679566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&odm));
3689566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&rdm));
36902918637SMatthew G. Knepley   }
37002918637SMatthew G. Knepley   PetscFunctionReturn(0);
37102918637SMatthew G. Knepley }
37202918637SMatthew G. Knepley 
373*9371c9d4SSatish Balay int main(int argc, char **argv) {
37402918637SMatthew G. Knepley   DM     dm;
37502918637SMatthew G. Knepley   AppCtx user;
37602918637SMatthew G. Knepley 
377327415f7SBarry Smith   PetscFunctionBeginUser;
3789566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3799566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3809566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
38102918637SMatthew G. Knepley   if (user.initOrnt) {
3829566063dSJacob Faibussowitsch     PetscCall(ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE));
3839566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm, NULL, "-ornt_dm_view"));
38402918637SMatthew G. Knepley   }
3859566063dSJacob Faibussowitsch   PetscCall(GenerateArrangments(dm, &user));
3869566063dSJacob Faibussowitsch   PetscCall(VerifyCayleyTable(dm, &user));
3879566063dSJacob Faibussowitsch   PetscCall(VerifyInverse(dm, &user));
3889566063dSJacob Faibussowitsch   PetscCall(RefineArrangments(dm, &user));
3899566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3909566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
391b122ec5aSJacob Faibussowitsch   return 0;
39202918637SMatthew G. Knepley }
39302918637SMatthew G. Knepley 
39402918637SMatthew G. Knepley /*TEST
39502918637SMatthew G. Knepley 
396a01d01faSMatthew G. Knepley   testset:
397a01d01faSMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_reference_cell_domain -gen_arrangements \
398a01d01faSMatthew G. Knepley           -gen_dm_view ::ascii_latex -dm_plex_view_tikzscale 0.5
399a01d01faSMatthew G. Knepley 
40002918637SMatthew G. Knepley     test:
40102918637SMatthew G. Knepley       suffix: segment
402a01d01faSMatthew G. Knepley       args: -dm_plex_cell segment \
403a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0
40402918637SMatthew G. Knepley 
40502918637SMatthew G. Knepley     test:
40602918637SMatthew G. Knepley       suffix: triangle
407a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangle \
408a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
40902918637SMatthew G. Knepley 
41002918637SMatthew G. Knepley     test:
41102918637SMatthew G. Knepley       suffix: quadrilateral
412a01d01faSMatthew G. Knepley       args: -dm_plex_cell quadrilateral \
413a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
41402918637SMatthew G. Knepley 
41502918637SMatthew G. Knepley     test:
41602918637SMatthew G. Knepley       suffix: tensor_segment
417a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quad \
418a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
41902918637SMatthew G. Knepley 
42002918637SMatthew G. Knepley     test:
42102918637SMatthew G. Knepley       suffix: tetrahedron
422a01d01faSMatthew G. Knepley       args: -dm_plex_cell tetrahedron \
423a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
42402918637SMatthew G. Knepley 
42502918637SMatthew G. Knepley     test:
42602918637SMatthew G. Knepley       suffix: hexahedron
427a01d01faSMatthew G. Knepley       args: -dm_plex_cell hexahedron \
428a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.3
42902918637SMatthew G. Knepley 
43002918637SMatthew G. Knepley     test:
43102918637SMatthew G. Knepley       suffix: triangular_prism
432a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangular_prism \
433a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
43402918637SMatthew G. Knepley 
43502918637SMatthew G. Knepley     test:
43602918637SMatthew G. Knepley       suffix: tensor_triangular_prism
437a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_triangular_prism \
438a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
43902918637SMatthew G. Knepley 
44002918637SMatthew G. Knepley     test:
44102918637SMatthew G. Knepley       suffix: tensor_quadrilateral_prism
442a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quadrilateral_prism \
443a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
44402918637SMatthew G. Knepley 
44502918637SMatthew G. Knepley     test:
44602918637SMatthew G. Knepley       suffix: pyramid
447a01d01faSMatthew G. Knepley       args: -dm_plex_cell pyramid \
448a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
449a01d01faSMatthew G. Knepley 
450a01d01faSMatthew G. Knepley   testset:
451a01d01faSMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_reference_cell_domain -ref_arrangements -dm_plex_check_all \
452a01d01faSMatthew G. Knepley           -ref_dm_view ::ascii_latex -dm_plex_view_tikzscale 1.0
45302918637SMatthew G. Knepley 
45402918637SMatthew G. Knepley     test:
45502918637SMatthew G. Knepley       suffix: ref_segment
456a01d01faSMatthew G. Knepley       args: -dm_plex_cell segment \
457a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0
45802918637SMatthew G. Knepley 
45902918637SMatthew G. Knepley     test:
46002918637SMatthew G. Knepley       suffix: ref_triangle
461a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangle \
462a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
46302918637SMatthew G. Knepley 
46402918637SMatthew G. Knepley     test:
46502918637SMatthew G. Knepley       suffix: ref_quadrilateral
466a01d01faSMatthew G. Knepley       args: -dm_plex_cell quadrilateral \
467a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
46802918637SMatthew G. Knepley 
46902918637SMatthew G. Knepley     test:
47002918637SMatthew G. Knepley       suffix: ref_tensor_segment
471a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quad \
472a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
47302918637SMatthew G. Knepley 
47402918637SMatthew G. Knepley     test:
47502918637SMatthew G. Knepley       suffix: ref_tetrahedron
476a01d01faSMatthew G. Knepley       args: -dm_plex_cell tetrahedron \
477a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
47802918637SMatthew G. Knepley 
47902918637SMatthew G. Knepley     test:
48002918637SMatthew G. Knepley       suffix: ref_hexahedron
481a01d01faSMatthew G. Knepley       args: -dm_plex_cell hexahedron \
482a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
48302918637SMatthew G. Knepley 
48402918637SMatthew G. Knepley     test:
48502918637SMatthew G. Knepley       suffix: ref_triangular_prism
486a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangular_prism \
487a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
48802918637SMatthew G. Knepley 
48902918637SMatthew G. Knepley     test:
49002918637SMatthew G. Knepley       suffix: ref_tensor_triangular_prism
491a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_triangular_prism \
492a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
49302918637SMatthew G. Knepley 
49402918637SMatthew G. Knepley     test:
49502918637SMatthew G. Knepley       suffix: ref_tensor_quadrilateral_prism
496a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quadrilateral_prism \
497a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
49802918637SMatthew G. Knepley 
49902918637SMatthew G. Knepley   # ToBox should recreate the coordinate space since the cell shape changes
50002918637SMatthew G. Knepley   testset:
50102918637SMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all
50202918637SMatthew G. Knepley 
50302918637SMatthew G. Knepley     test:
50402918637SMatthew G. Knepley       suffix: tobox_triangle
50502918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
50602918637SMatthew G. Knepley             -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
50702918637SMatthew G. Knepley 
50802918637SMatthew G. Knepley     test:
50902918637SMatthew G. Knepley       suffix: tobox_tensor_segment
51002918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
51102918637SMatthew G. Knepley             -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
51202918637SMatthew G. Knepley 
51302918637SMatthew G. Knepley     test:
51402918637SMatthew G. Knepley       suffix: tobox_tetrahedron
51502918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
51602918637SMatthew G. Knepley             -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
51702918637SMatthew G. Knepley 
51802918637SMatthew G. Knepley     test:
51902918637SMatthew G. Knepley       suffix: tobox_triangular_prism
52002918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
52102918637SMatthew G. Knepley             -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
52202918637SMatthew G. Knepley 
52302918637SMatthew G. Knepley     test:
52402918637SMatthew G. Knepley       suffix: tobox_tensor_triangular_prism
52502918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
52602918637SMatthew G. Knepley             -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
52702918637SMatthew G. Knepley 
52802918637SMatthew G. Knepley     test:
52902918637SMatthew G. Knepley       suffix: tobox_tensor_quadrilateral_prism
53002918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \
53102918637SMatthew G. Knepley             -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
53202918637SMatthew G. Knepley 
53302918637SMatthew G. Knepley   testset:
53402918637SMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all
53502918637SMatthew G. Knepley 
53602918637SMatthew G. Knepley     test:
53702918637SMatthew G. Knepley       suffix: alfeld_triangle
53802918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
53902918637SMatthew G. Knepley             -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
54002918637SMatthew G. Knepley 
54102918637SMatthew G. Knepley     test:
54202918637SMatthew G. Knepley       suffix: alfeld_tetrahedron
54302918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
54402918637SMatthew G. Knepley             -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
54502918637SMatthew G. Knepley 
54602918637SMatthew G. Knepley   testset:
54702918637SMatthew G. Knepley     args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all
54802918637SMatthew G. Knepley 
54902918637SMatthew G. Knepley     # This splits edge 1 of the triangle, and reflects about the added edge
55002918637SMatthew G. Knepley     test:
55102918637SMatthew G. Knepley       suffix: sbr_triangle_0
55202918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \
55302918637SMatthew G. Knepley             -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
55402918637SMatthew G. Knepley 
55502918637SMatthew G. Knepley     # This splits edge 0 of the triangle, and reflects about the added edge
55602918637SMatthew G. Knepley     test:
55702918637SMatthew G. Knepley       suffix: sbr_triangle_1
55802918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -init_ornt 1 -ornts -3,0 \
55902918637SMatthew G. Knepley             -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
56002918637SMatthew G. Knepley 
56102918637SMatthew G. Knepley     # This splits edge 2 of the triangle, and reflects about the added edge
56202918637SMatthew G. Knepley     test:
56302918637SMatthew G. Knepley       suffix: sbr_triangle_2
56402918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -init_ornt 2 -ornts -1,0 \
56502918637SMatthew G. Knepley             -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
56602918637SMatthew G. Knepley 
56702918637SMatthew G. Knepley     # This splits edges 1 and 2 of the triangle
56802918637SMatthew G. Knepley     test:
56902918637SMatthew G. Knepley       suffix: sbr_triangle_3
57002918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 0 \
57102918637SMatthew G. Knepley             -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
57202918637SMatthew G. Knepley 
57302918637SMatthew G. Knepley     # This splits edges 1 and 0 of the triangle
57402918637SMatthew G. Knepley     test:
57502918637SMatthew G. Knepley       suffix: sbr_triangle_4
57602918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 0 \
57702918637SMatthew G. Knepley             -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
57802918637SMatthew G. Knepley 
57902918637SMatthew G. Knepley     # This splits edges 0 and 1 of the triangle
58002918637SMatthew G. Knepley     test:
58102918637SMatthew G. Knepley       suffix: sbr_triangle_5
58202918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -init_ornt 1 -ornts 0 \
58302918637SMatthew G. Knepley             -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
58402918637SMatthew G. Knepley 
58502918637SMatthew G. Knepley     # This splits edges 0 and 2 of the triangle
58602918637SMatthew G. Knepley     test:
58702918637SMatthew G. Knepley       suffix: sbr_triangle_6
58802918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -init_ornt 1 -ornts 0 \
58902918637SMatthew G. Knepley             -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
59002918637SMatthew G. Knepley 
59102918637SMatthew G. Knepley     # This splits edges 2 and 0 of the triangle
59202918637SMatthew G. Knepley     test:
59302918637SMatthew G. Knepley       suffix: sbr_triangle_7
59402918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -init_ornt 2 -ornts 0 \
59502918637SMatthew G. Knepley             -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
59602918637SMatthew G. Knepley 
59702918637SMatthew G. Knepley     # This splits edges 2 and 1 of the triangle
59802918637SMatthew G. Knepley     test:
59902918637SMatthew G. Knepley       suffix: sbr_triangle_8
60002918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -init_ornt 2 -ornts 0 \
60102918637SMatthew G. Knepley             -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
60202918637SMatthew G. Knepley 
60302918637SMatthew G. Knepley   testset:
60402918637SMatthew G. Knepley     args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all
60502918637SMatthew G. Knepley 
60602918637SMatthew G. Knepley     test:
60702918637SMatthew G. Knepley       suffix: bl_tensor_segment
60802918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
60902918637SMatthew G. Knepley             -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
61002918637SMatthew G. Knepley 
61102918637SMatthew G. Knepley     # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check
61202918637SMatthew G. Knepley     test:
61302918637SMatthew G. Knepley       suffix: bl_tensor_triangular_prism
61402918637SMatthew G. Knepley       requires: TODO
61502918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
61602918637SMatthew G. Knepley             -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
61702918637SMatthew G. Knepley 
61802918637SMatthew G. Knepley TEST*/
619