xref: /petsc/src/dm/impls/plex/tutorials/ex11.c (revision a01d01fa3fdd0823377c6ff796306a39b7cfc3fd)
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 
1602918637SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
1702918637SMatthew G. Knepley {
1802918637SMatthew G. Knepley   PetscInt       n = 2;
1902918637SMatthew G. Knepley   PetscBool      flg;
2002918637SMatthew G. Knepley   PetscErrorCode ierr;
2102918637SMatthew G. Knepley 
2202918637SMatthew G. Knepley   PetscFunctionBeginUser;
2302918637SMatthew G. Knepley   options->genArr        = PETSC_FALSE;
2402918637SMatthew G. Knepley   options->refArr        = PETSC_FALSE;
2502918637SMatthew G. Knepley   options->printTable    = PETSC_FALSE;
2602918637SMatthew G. Knepley   options->orntBounds[0] = PETSC_MIN_INT;
2702918637SMatthew G. Knepley   options->orntBounds[1] = PETSC_MAX_INT;
2802918637SMatthew G. Knepley   options->numOrnt       = -1;
2902918637SMatthew G. Knepley   options->initOrnt      = 0;
3002918637SMatthew G. Knepley 
3102918637SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Mesh Orientation Tutorials Options", "DMPLEX");CHKERRQ(ierr);
3202918637SMatthew G. Knepley   ierr = PetscOptionsBool("-gen_arrangements", "Flag for generating all arrangements of the cell", "ex11.c", options->genArr, &options->genArr, NULL);CHKERRQ(ierr);
3302918637SMatthew G. Knepley   ierr = PetscOptionsBool("-ref_arrangements", "Flag for refining all arrangements of the cell", "ex11.c", options->refArr, &options->refArr, NULL);CHKERRQ(ierr);
3402918637SMatthew G. Knepley   ierr = PetscOptionsBool("-print_table", "Print the Cayley table", "ex11.c", options->printTable, &options->printTable, NULL);CHKERRQ(ierr);
3502918637SMatthew G. Knepley   ierr = PetscOptionsIntArray("-ornt_bounds", "Bounds for orientation checks", "ex11.c", options->orntBounds, &n, NULL);CHKERRQ(ierr);
3602918637SMatthew G. Knepley   n    = 48;
3702918637SMatthew G. Knepley   ierr = PetscOptionsIntArray("-ornts", "Specific orientations for checks", "ex11.c", options->ornts, &n, &flg);CHKERRQ(ierr);
3802918637SMatthew G. Knepley   if (flg) {
3902918637SMatthew G. Knepley     options->numOrnt = n;
4002918637SMatthew G. Knepley     ierr = PetscSortInt(n, options->ornts);CHKERRQ(ierr);
4102918637SMatthew G. Knepley   }
4202918637SMatthew G. Knepley   ierr = PetscOptionsInt("-init_ornt", "Initial orientation for starting mesh", "ex11.c", options->initOrnt, &options->initOrnt, NULL);CHKERRQ(ierr);
4302918637SMatthew G. Knepley   ierr = PetscOptionsEnd();
4402918637SMatthew G. Knepley   PetscFunctionReturn(0);
4502918637SMatthew G. Knepley }
4602918637SMatthew G. Knepley 
4702918637SMatthew G. Knepley static PetscBool ignoreOrnt(AppCtx *user, PetscInt o)
4802918637SMatthew G. Knepley {
4902918637SMatthew G. Knepley   PetscInt       loc;
5002918637SMatthew G. Knepley   PetscErrorCode ierr;
5102918637SMatthew G. Knepley 
5202918637SMatthew G. Knepley   if (user->numOrnt < 0) return PETSC_FALSE;
5302918637SMatthew G. Knepley   ierr = PetscFindInt(o, user->numOrnt, user->ornts, &loc);
5402918637SMatthew G. Knepley   if (loc < 0 || ierr)   return PETSC_TRUE;
5502918637SMatthew G. Knepley   return PETSC_FALSE;
5602918637SMatthew G. Knepley }
5702918637SMatthew G. Knepley 
5802918637SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
5902918637SMatthew G. Knepley {
6002918637SMatthew G. Knepley   PetscErrorCode ierr;
6102918637SMatthew G. Knepley 
6202918637SMatthew G. Knepley   PetscFunctionBeginUser;
6302918637SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
6402918637SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
6502918637SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
6602918637SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
6702918637SMatthew G. Knepley   PetscFunctionReturn(0);
6802918637SMatthew G. Knepley }
6902918637SMatthew G. Knepley 
7002918637SMatthew G. Knepley static PetscErrorCode CheckCellVertices(DM dm, PetscInt cell, PetscInt o)
7102918637SMatthew G. Knepley {
7202918637SMatthew G. Knepley   DMPolytopeType  ct;
7302918637SMatthew G. Knepley   const PetscInt *arrVerts;
7402918637SMatthew G. Knepley   PetscInt       *closure = NULL;
7502918637SMatthew G. Knepley   PetscInt        Ncl, cl, Nv, vStart, vEnd, v;
7602918637SMatthew G. Knepley   MPI_Comm        comm;
7702918637SMatthew G. Knepley   PetscErrorCode  ierr;
7802918637SMatthew G. Knepley 
7902918637SMatthew G. Knepley   PetscFunctionBeginUser;
8002918637SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
8102918637SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr);
8202918637SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
8302918637SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
8402918637SMatthew G. Knepley   for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) {
8502918637SMatthew G. Knepley     const PetscInt vertex = closure[cl];
8602918637SMatthew G. Knepley 
8702918637SMatthew G. Knepley     if (vertex < vStart || vertex >= vEnd) continue;
8802918637SMatthew G. Knepley     closure[Nv++] = vertex;
8902918637SMatthew G. Knepley   }
902c71b3e2SJacob Faibussowitsch   PetscCheckFalse(Nv != DMPolytopeTypeGetNumVertices(ct),comm, PETSC_ERR_ARG_WRONG, "Cell %D has %D vertices != %D vertices in a %s", cell, Nv, DMPolytopeTypeGetNumVertices(ct), DMPolytopeTypes[ct]);
9102918637SMatthew G. Knepley   arrVerts = DMPolytopeTypeGetVertexArrangment(ct, o);
9202918637SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
932c71b3e2SJacob Faibussowitsch     PetscCheckFalse(closure[v] != arrVerts[v]+vStart,comm, PETSC_ERR_ARG_WRONG, "Cell %D vertex[%D]: %D should be %D for arrangement %D", cell, v, closure[v], arrVerts[v]+vStart, o);
9402918637SMatthew G. Knepley   }
9502918637SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
9602918637SMatthew G. Knepley   PetscFunctionReturn(0);
9702918637SMatthew G. Knepley }
9802918637SMatthew G. Knepley 
9902918637SMatthew G. Knepley /* Transform cell with group operation o */
10002918637SMatthew G. Knepley static PetscErrorCode ReorientCell(DM dm, PetscInt cell, PetscInt o, PetscBool swapCoords)
10102918637SMatthew G. Knepley {
10202918637SMatthew G. Knepley   DM              cdm;
10302918637SMatthew G. Knepley   Vec             coordinates;
10402918637SMatthew G. Knepley   PetscScalar    *coords, *ccoords = NULL;
10502918637SMatthew G. Knepley   PetscInt       *closure = NULL;
10602918637SMatthew G. Knepley   PetscInt        cdim, d, Nc, Ncl, cl, vStart, vEnd, Nv;
10702918637SMatthew G. Knepley   PetscErrorCode  ierr;
10802918637SMatthew G. Knepley 
10902918637SMatthew G. Knepley   PetscFunctionBegin;
11002918637SMatthew G. Knepley   /* Change vertex coordinates so that it plots as we expect */
11102918637SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
11202918637SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
11302918637SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
11402918637SMatthew G. Knepley   ierr = DMPlexVecGetClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords);CHKERRQ(ierr);
11502918637SMatthew G. Knepley   /* Reorient cone */
11602918637SMatthew G. Knepley   ierr = DMPlexOrientPoint(dm, cell, o);CHKERRQ(ierr);
11702918637SMatthew G. Knepley   /* Finish resetting coordinates */
11802918637SMatthew G. Knepley   if (swapCoords) {
11902918637SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12002918637SMatthew G. Knepley     ierr = VecGetArrayWrite(coordinates, &coords);CHKERRQ(ierr);
12102918637SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
12202918637SMatthew G. Knepley     for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) {
12302918637SMatthew G. Knepley       const PetscInt vertex = closure[cl];
12402918637SMatthew G. Knepley       PetscScalar   *vcoords;
12502918637SMatthew G. Knepley 
12602918637SMatthew G. Knepley       if (vertex < vStart || vertex >= vEnd) continue;
12702918637SMatthew G. Knepley       ierr = DMPlexPointLocalRef(cdm, vertex, coords, &vcoords);CHKERRQ(ierr);
12802918637SMatthew G. Knepley       for (d = 0; d < cdim; ++d) vcoords[d] = ccoords[Nv*cdim + d];
12902918637SMatthew G. Knepley       ++Nv;
13002918637SMatthew G. Knepley     }
13102918637SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
13202918637SMatthew G. Knepley     ierr = VecRestoreArrayWrite(coordinates, &coords);CHKERRQ(ierr);
13302918637SMatthew G. Knepley   }
13402918637SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords);CHKERRQ(ierr);
13502918637SMatthew G. Knepley   PetscFunctionReturn(0);
13602918637SMatthew G. Knepley }
13702918637SMatthew G. Knepley 
13802918637SMatthew G. Knepley static PetscErrorCode GenerateArrangments(DM dm, AppCtx *user)
13902918637SMatthew G. Knepley {
14002918637SMatthew G. Knepley   DM             odm;
14102918637SMatthew G. Knepley   DMPolytopeType ct;
14202918637SMatthew G. Knepley   PetscInt       No, o;
14302918637SMatthew G. Knepley   const char    *name;
14402918637SMatthew G. Knepley   PetscErrorCode ierr;
14502918637SMatthew G. Knepley 
14602918637SMatthew G. Knepley   PetscFunctionBeginUser;
14702918637SMatthew G. Knepley   if (!user->genArr) PetscFunctionReturn(0);
14802918637SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr);
14902918637SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr);
15002918637SMatthew G. Knepley   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
15102918637SMatthew G. Knepley   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
15202918637SMatthew G. Knepley     if (ignoreOrnt(user, o)) continue;
15302918637SMatthew G. Knepley     ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);CHKERRQ(ierr);
15402918637SMatthew G. Knepley     ierr = ReorientCell(odm, 0, o, PETSC_TRUE);CHKERRQ(ierr);
15502918637SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);CHKERRQ(ierr);
15602918637SMatthew G. Knepley     ierr = DMViewFromOptions(odm, NULL, "-gen_dm_view");CHKERRQ(ierr);
15702918637SMatthew G. Knepley     ierr = CheckCellVertices(odm, 0, o);CHKERRQ(ierr);
15802918637SMatthew G. Knepley     ierr = DMDestroy(&odm);CHKERRQ(ierr);
15902918637SMatthew G. Knepley   }
16002918637SMatthew G. Knepley   PetscFunctionReturn(0);
16102918637SMatthew G. Knepley }
16202918637SMatthew G. Knepley 
16302918637SMatthew G. Knepley static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user)
16402918637SMatthew G. Knepley {
16502918637SMatthew G. Knepley   DM              dm1, dm2;
16602918637SMatthew G. Knepley   DMPolytopeType  ct;
16702918637SMatthew G. Knepley   const PetscInt *refcone, *cone;
16802918637SMatthew G. Knepley   PetscInt        No, o1, o2, o3, o4;
16902918637SMatthew G. Knepley   PetscBool       equal;
17002918637SMatthew G. Knepley   const char     *name;
17102918637SMatthew G. Knepley   PetscErrorCode  ierr;
17202918637SMatthew G. Knepley 
17302918637SMatthew G. Knepley   PetscFunctionBeginUser;
17402918637SMatthew G. Knepley   if (!user->genArr) PetscFunctionReturn(0);
17502918637SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr);
17602918637SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr);
17702918637SMatthew G. Knepley   ierr = DMPlexGetCone(dm, 0, &refcone);CHKERRQ(ierr);
17802918637SMatthew G. Knepley   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
17902918637SMatthew G. Knepley   if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "Cayley Table for %s\n", DMPolytopeTypes[ct]);CHKERRQ(ierr);}
18002918637SMatthew G. Knepley   for (o1 = PetscMax(-No, user->orntBounds[0]); o1 < PetscMin(No, user->orntBounds[1]); ++o1) {
18102918637SMatthew G. Knepley     for (o2 = PetscMax(-No, user->orntBounds[0]); o2 < PetscMin(No, user->orntBounds[1]); ++o2) {
18202918637SMatthew G. Knepley       ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1);CHKERRQ(ierr);
18302918637SMatthew G. Knepley       ierr = DMPlexOrientPoint(dm1, 0, o2);CHKERRQ(ierr);
18402918637SMatthew G. Knepley       ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr);
18502918637SMatthew G. Knepley       ierr = DMPlexOrientPoint(dm1, 0, o1);CHKERRQ(ierr);
18602918637SMatthew G. Knepley       ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr);
18702918637SMatthew G. Knepley       o3   = DMPolytopeTypeComposeOrientation(ct, o1, o2);CHKERRQ(ierr);
18802918637SMatthew G. Knepley       /* First verification */
18902918637SMatthew G. Knepley       ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2);CHKERRQ(ierr);
19002918637SMatthew G. Knepley       ierr = DMPlexOrientPoint(dm2, 0, o3);CHKERRQ(ierr);
19102918637SMatthew G. Knepley       ierr = DMPlexCheckFaces(dm2, 0);CHKERRQ(ierr);
19202918637SMatthew G. Knepley       ierr = DMPlexEqual(dm1, dm2, &equal);CHKERRQ(ierr);
19302918637SMatthew G. Knepley       if (!equal) {
19402918637SMatthew G. Knepley         ierr = DMViewFromOptions(dm1, NULL, "-error_dm_view");CHKERRQ(ierr);
19502918637SMatthew G. Knepley         ierr = DMViewFromOptions(dm2, NULL, "-error_dm_view");CHKERRQ(ierr);
19698921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D != %D", DMPolytopeTypes[ct], o1, o2, o3);
19702918637SMatthew G. Knepley       }
19802918637SMatthew G. Knepley       /* Second verification */
19902918637SMatthew G. Knepley       ierr = DMPlexGetCone(dm1, 0, &cone);CHKERRQ(ierr);
20002918637SMatthew G. Knepley       ierr = DMPolytopeGetOrientation(ct, refcone, cone, &o4);CHKERRQ(ierr);
20102918637SMatthew G. Knepley       if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "%D, ", o4);CHKERRQ(ierr);}
2022c71b3e2SJacob Faibussowitsch       PetscCheckFalse(o3 != o4,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D = %D != %D", DMPolytopeTypes[ct], o1, o2, o3, o4);
20302918637SMatthew G. Knepley       ierr = DMDestroy(&dm1);CHKERRQ(ierr);
20402918637SMatthew G. Knepley       ierr = DMDestroy(&dm2);CHKERRQ(ierr);
20502918637SMatthew G. Knepley     }
20602918637SMatthew G. Knepley     if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);}
20702918637SMatthew G. Knepley   }
20802918637SMatthew G. Knepley   PetscFunctionReturn(0);
20902918637SMatthew G. Knepley }
21002918637SMatthew G. Knepley 
21102918637SMatthew G. Knepley static PetscErrorCode VerifyInverse(DM dm, AppCtx *user)
21202918637SMatthew G. Knepley {
21302918637SMatthew G. Knepley   DM              dm1, dm2;
21402918637SMatthew G. Knepley   DMPolytopeType  ct;
21502918637SMatthew G. Knepley   const PetscInt *refcone, *cone;
21602918637SMatthew G. Knepley   PetscInt        No, o, oi, o2;
21702918637SMatthew G. Knepley   PetscBool       equal;
21802918637SMatthew G. Knepley   const char     *name;
21902918637SMatthew G. Knepley   PetscErrorCode  ierr;
22002918637SMatthew G. Knepley 
22102918637SMatthew G. Knepley   PetscFunctionBeginUser;
22202918637SMatthew G. Knepley   if (!user->genArr) PetscFunctionReturn(0);
22302918637SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr);
22402918637SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr);
22502918637SMatthew G. Knepley   ierr = DMPlexGetCone(dm, 0, &refcone);CHKERRQ(ierr);
22602918637SMatthew G. Knepley   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
22702918637SMatthew G. Knepley   if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "Inverse table for %s\n", DMPolytopeTypes[ct]);CHKERRQ(ierr);}
22802918637SMatthew G. Knepley   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
22902918637SMatthew G. Knepley     if (ignoreOrnt(user, o)) continue;
23002918637SMatthew G. Knepley     oi   = DMPolytopeTypeComposeOrientationInv(ct, 0, o);CHKERRQ(ierr);
23102918637SMatthew G. Knepley     ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1);CHKERRQ(ierr);
23202918637SMatthew G. Knepley     ierr = DMPlexOrientPoint(dm1, 0, o);CHKERRQ(ierr);
23302918637SMatthew G. Knepley     ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr);
23402918637SMatthew G. Knepley     ierr = DMPlexOrientPoint(dm1, 0, oi);CHKERRQ(ierr);
23502918637SMatthew G. Knepley     ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr);
23602918637SMatthew G. Knepley     /* First verification */
23702918637SMatthew G. Knepley     ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2);CHKERRQ(ierr);
23802918637SMatthew G. Knepley     ierr = DMPlexEqual(dm1, dm2, &equal);CHKERRQ(ierr);
23902918637SMatthew G. Knepley     if (!equal) {
24002918637SMatthew G. Knepley       ierr = DMViewFromOptions(dm1, NULL, "-error_dm_view");CHKERRQ(ierr);
24102918637SMatthew G. Knepley       ierr = DMViewFromOptions(dm2, NULL, "-error_dm_view");CHKERRQ(ierr);
24298921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D != 0", DMPolytopeTypes[ct], o, oi);
24302918637SMatthew G. Knepley     }
24402918637SMatthew G. Knepley     /* Second verification */
24502918637SMatthew G. Knepley     ierr = DMPlexGetCone(dm1, 0, &cone);CHKERRQ(ierr);
24602918637SMatthew G. Knepley     ierr = DMPolytopeGetOrientation(ct, refcone, cone, &o2);CHKERRQ(ierr);
24702918637SMatthew G. Knepley     if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "%D, ", oi);CHKERRQ(ierr);}
2482c71b3e2SJacob Faibussowitsch     PetscCheckFalse(o2 != 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D = %D != 0", DMPolytopeTypes[ct], o, oi, o2);
24902918637SMatthew G. Knepley     ierr = DMDestroy(&dm1);CHKERRQ(ierr);
25002918637SMatthew G. Knepley     ierr = DMDestroy(&dm2);CHKERRQ(ierr);
25102918637SMatthew G. Knepley   }
25202918637SMatthew G. Knepley   if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);}
25302918637SMatthew G. Knepley   PetscFunctionReturn(0);
25402918637SMatthew G. Knepley }
25502918637SMatthew G. Knepley 
25602918637SMatthew G. Knepley /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */
25702918637SMatthew G. Knepley static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user)
25802918637SMatthew G. Knepley {
25902918637SMatthew G. Knepley   DMPlexTransform tr, otr;
26002918637SMatthew G. Knepley   DMPolytopeType  ct;
26102918637SMatthew G. Knepley   DMPolytopeType *rct;
26202918637SMatthew G. Knepley   const PetscInt *cone, *ornt, *ocone, *oornt;
26302918637SMatthew G. Knepley   PetscInt       *rsize, *rcone, *rornt;
26402918637SMatthew G. Knepley   PetscInt        Nct, n, oi, debug = 0;
26502918637SMatthew G. Knepley   PetscErrorCode  ierr;
26602918637SMatthew G. Knepley 
26702918637SMatthew G. Knepley   PetscFunctionBeginUser;
26802918637SMatthew G. Knepley   ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr);CHKERRQ(ierr);
26902918637SMatthew G. Knepley   ierr = DMPlexTransformSetDM(tr, dm);CHKERRQ(ierr);
27002918637SMatthew G. Knepley   ierr = DMPlexTransformSetFromOptions(tr);CHKERRQ(ierr);
27102918637SMatthew G. Knepley   ierr = DMPlexTransformSetUp(tr);CHKERRQ(ierr);
27202918637SMatthew G. Knepley 
27302918637SMatthew G. Knepley   ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) odm), &otr);CHKERRQ(ierr);
27402918637SMatthew G. Knepley   ierr = DMPlexTransformSetDM(otr, odm);CHKERRQ(ierr);
27502918637SMatthew G. Knepley   ierr = DMPlexTransformSetFromOptions(otr);CHKERRQ(ierr);
27602918637SMatthew G. Knepley   ierr = DMPlexTransformSetUp(otr);CHKERRQ(ierr);
27702918637SMatthew G. Knepley 
27802918637SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr);
27902918637SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
28002918637SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, p, &ornt);CHKERRQ(ierr);
28102918637SMatthew G. Knepley   ierr = DMPlexGetCone(odm, p, &ocone);CHKERRQ(ierr);
28202918637SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(odm, p, &oornt);CHKERRQ(ierr);
28302918637SMatthew G. Knepley   oi   = DMPolytopeTypeComposeOrientationInv(ct, 0, o);CHKERRQ(ierr);
28402918637SMatthew G. Knepley   if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "Orientation %D\n", oi);CHKERRQ(ierr);}
28502918637SMatthew G. Knepley 
28602918637SMatthew G. Knepley   ierr = DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt);CHKERRQ(ierr);
28702918637SMatthew G. Knepley   for (n = 0; n < Nct; ++n) {
28802918637SMatthew G. Knepley     DMPolytopeType ctNew = rct[n];
28902918637SMatthew G. Knepley     PetscInt       r, ro;
29002918637SMatthew G. Knepley 
29102918637SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  Checking type %s\n", DMPolytopeTypes[ctNew]);CHKERRQ(ierr);}
29202918637SMatthew G. Knepley     for (r = 0; r < rsize[n]; ++r) {
29302918637SMatthew G. Knepley       const PetscInt *qcone, *qornt, *oqcone, *oqornt;
29402918637SMatthew G. Knepley       PetscInt        pNew, opNew, oo, pr, fo;
29502918637SMatthew G. Knepley       PetscBool       restore = PETSC_TRUE;
29602918637SMatthew G. Knepley 
29702918637SMatthew G. Knepley       ierr = DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew);CHKERRQ(ierr);
29802918637SMatthew G. Knepley       ierr = DMPlexTransformGetCone(tr, pNew, &qcone, &qornt);CHKERRQ(ierr);
29902918637SMatthew G. Knepley       if (debug) {
30002918637SMatthew G. Knepley         PetscInt c;
30102918637SMatthew G. Knepley 
30202918637SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "    Checking replica %D (%D)\n      Original Cone", r, pNew);CHKERRQ(ierr);
30302918637SMatthew G. Knepley         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %D (%D)", qcone[c], qornt[c]);CHKERRQ(ierr);}
30402918637SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
30502918637SMatthew G. Knepley       }
30602918637SMatthew G. Knepley       for (ro = 0; ro < rsize[n]; ++ro) {
30702918637SMatthew G. Knepley         PetscBool found;
30802918637SMatthew G. Knepley 
30902918637SMatthew G. Knepley         ierr = DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew);CHKERRQ(ierr);
31002918637SMatthew G. Knepley         ierr = DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt);CHKERRQ(ierr);
31102918637SMatthew G. Knepley         ierr = DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found);CHKERRQ(ierr);
31202918637SMatthew G. Knepley         if (found) break;
31302918637SMatthew G. Knepley         ierr = DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);CHKERRQ(ierr);
31402918637SMatthew G. Knepley       }
31502918637SMatthew G. Knepley       if (debug) {
31602918637SMatthew G. Knepley         PetscInt c;
31702918637SMatthew G. Knepley 
31802918637SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "    Checking transform replica %D (%D) (%D)\n      Transform Cone", ro, opNew, o);CHKERRQ(ierr);
31902918637SMatthew G. Knepley         for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %D (%D)", oqcone[c], oqornt[c]);CHKERRQ(ierr);}
32002918637SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
32102918637SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "    Matched %D\n", oo);CHKERRQ(ierr);
32202918637SMatthew G. Knepley       }
32302918637SMatthew G. Knepley       if (ro == rsize[n]) {
32402918637SMatthew G. Knepley         /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */
32502918637SMatthew G. Knepley         if (ct == DM_POLYTOPE_TETRAHEDRON) {
32602918637SMatthew G. Knepley           /* The segment in a tetrahedron does not map into itself under the group action */
32702918637SMatthew G. Knepley           if (ctNew == DM_POLYTOPE_SEGMENT) {restore = PETSC_FALSE; ro = r; oo = 0;}
32802918637SMatthew G. Knepley           /* The last four interior faces do not map into themselves under the group action */
32902918637SMatthew G. Knepley           if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) {restore = PETSC_FALSE; ro = r; oo = 0;}
33002918637SMatthew G. Knepley           /* The last four interior faces do not map into themselves under the group action */
33102918637SMatthew G. Knepley           if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) {restore = PETSC_FALSE; ro = r; oo = 0;}
33202918637SMatthew G. Knepley         }
3332c71b3e2SJacob Faibussowitsch         PetscCheckFalse(restore,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to find matching %s %D orientation for cell orientation %D", DMPolytopeTypes[ctNew], r, o);
33402918637SMatthew G. Knepley       }
33502918637SMatthew G. Knepley       if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "%D, %D, ", ro, oo);CHKERRQ(ierr);}
33602918637SMatthew G. Knepley       ierr = DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo);CHKERRQ(ierr);
33702918637SMatthew G. Knepley       if (!user->printTable) {
3382c71b3e2SJacob Faibussowitsch         PetscCheckFalse(pr != ro,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong replica %D != %D", pr, ro);
3392c71b3e2SJacob Faibussowitsch         PetscCheckFalse(fo != oo,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong orientation %D != %D", fo, oo);
34002918637SMatthew G. Knepley       }
34102918637SMatthew G. Knepley       ierr = DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt);CHKERRQ(ierr);
34202918637SMatthew G. Knepley       if (restore) {ierr = DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);CHKERRQ(ierr);}
34302918637SMatthew G. Knepley     }
34402918637SMatthew G. Knepley     if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);}
34502918637SMatthew G. Knepley   }
34602918637SMatthew G. Knepley   ierr = DMPlexTransformDestroy(&tr);CHKERRQ(ierr);
34702918637SMatthew G. Knepley   ierr = DMPlexTransformDestroy(&otr);CHKERRQ(ierr);
34802918637SMatthew G. Knepley   PetscFunctionReturn(0);
34902918637SMatthew G. Knepley }
35002918637SMatthew G. Knepley 
35102918637SMatthew G. Knepley static PetscErrorCode RefineArrangments(DM dm, AppCtx *user)
35202918637SMatthew G. Knepley {
35302918637SMatthew G. Knepley   DM             odm, rdm;
35402918637SMatthew G. Knepley   DMPolytopeType ct;
35502918637SMatthew G. Knepley   PetscInt       No, o;
35602918637SMatthew G. Knepley   const char    *name;
35702918637SMatthew G. Knepley   PetscErrorCode ierr;
35802918637SMatthew G. Knepley 
35902918637SMatthew G. Knepley   PetscFunctionBeginUser;
36002918637SMatthew G. Knepley   if (!user->refArr) PetscFunctionReturn(0);
36102918637SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr);
36202918637SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr);
36302918637SMatthew G. Knepley   No   = DMPolytopeTypeGetNumArrangments(ct)/2;
36402918637SMatthew G. Knepley   for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) {
36502918637SMatthew G. Knepley     if (ignoreOrnt(user, o)) continue;
36602918637SMatthew G. Knepley     ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);CHKERRQ(ierr);
36702918637SMatthew G. Knepley     if (user->initOrnt) {ierr = ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE);CHKERRQ(ierr);}
36802918637SMatthew G. Knepley     ierr = ReorientCell(odm, 0, o, PETSC_TRUE);CHKERRQ(ierr);
36902918637SMatthew G. Knepley     ierr = DMViewFromOptions(odm, NULL, "-orig_dm_view");CHKERRQ(ierr);
37002918637SMatthew G. Knepley     ierr = DMRefine(odm, MPI_COMM_NULL, &rdm);CHKERRQ(ierr);
37102918637SMatthew G. Knepley     ierr = DMSetFromOptions(rdm);CHKERRQ(ierr);
37202918637SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);CHKERRQ(ierr);
37302918637SMatthew G. Knepley     ierr = DMViewFromOptions(rdm, NULL, "-ref_dm_view");CHKERRQ(ierr);
37402918637SMatthew G. Knepley     ierr = CheckSubcells(dm, odm, 0, o, user);CHKERRQ(ierr);
37502918637SMatthew G. Knepley     ierr = DMDestroy(&odm);CHKERRQ(ierr);
37602918637SMatthew G. Knepley     ierr = DMDestroy(&rdm);CHKERRQ(ierr);
37702918637SMatthew G. Knepley   }
37802918637SMatthew G. Knepley   PetscFunctionReturn(0);
37902918637SMatthew G. Knepley }
38002918637SMatthew G. Knepley 
38102918637SMatthew G. Knepley int main(int argc, char **argv)
38202918637SMatthew G. Knepley {
38302918637SMatthew G. Knepley   DM             dm;
38402918637SMatthew G. Knepley   AppCtx         user;
38502918637SMatthew G. Knepley   PetscErrorCode ierr;
38602918637SMatthew G. Knepley 
38702918637SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
38802918637SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
38902918637SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
39002918637SMatthew G. Knepley   if (user.initOrnt) {
39102918637SMatthew G. Knepley     ierr = ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE);CHKERRQ(ierr);
39202918637SMatthew G. Knepley     ierr = DMViewFromOptions(dm, NULL, "-ornt_dm_view");CHKERRQ(ierr);
39302918637SMatthew G. Knepley   }
39402918637SMatthew G. Knepley   ierr = GenerateArrangments(dm, &user);CHKERRQ(ierr);
39502918637SMatthew G. Knepley   ierr = VerifyCayleyTable(dm, &user);CHKERRQ(ierr);
39602918637SMatthew G. Knepley   ierr = VerifyInverse(dm, &user);CHKERRQ(ierr);
39702918637SMatthew G. Knepley   ierr = RefineArrangments(dm, &user);CHKERRQ(ierr);
39802918637SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
39902918637SMatthew G. Knepley   ierr = PetscFinalize();
40002918637SMatthew G. Knepley   return ierr;
40102918637SMatthew G. Knepley }
40202918637SMatthew G. Knepley 
40302918637SMatthew G. Knepley /*TEST
40402918637SMatthew G. Knepley 
405*a01d01faSMatthew G. Knepley   testset:
406*a01d01faSMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_reference_cell_domain -gen_arrangements \
407*a01d01faSMatthew G. Knepley           -gen_dm_view ::ascii_latex -dm_plex_view_tikzscale 0.5
408*a01d01faSMatthew G. Knepley 
40902918637SMatthew G. Knepley     test:
41002918637SMatthew G. Knepley       suffix: segment
411*a01d01faSMatthew G. Knepley       args: -dm_plex_cell segment \
412*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0
41302918637SMatthew G. Knepley 
41402918637SMatthew G. Knepley     test:
41502918637SMatthew G. Knepley       suffix: triangle
416*a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangle \
417*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
41802918637SMatthew G. Knepley 
41902918637SMatthew G. Knepley     test:
42002918637SMatthew G. Knepley       suffix: quadrilateral
421*a01d01faSMatthew G. Knepley       args: -dm_plex_cell quadrilateral \
422*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
42302918637SMatthew G. Knepley 
42402918637SMatthew G. Knepley     test:
42502918637SMatthew G. Knepley       suffix: tensor_segment
426*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quad \
427*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
42802918637SMatthew G. Knepley 
42902918637SMatthew G. Knepley     test:
43002918637SMatthew G. Knepley       suffix: tetrahedron
431*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tetrahedron \
432*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
43302918637SMatthew G. Knepley 
43402918637SMatthew G. Knepley     test:
43502918637SMatthew G. Knepley       suffix: hexahedron
436*a01d01faSMatthew G. Knepley       args: -dm_plex_cell hexahedron \
437*a01d01faSMatthew 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
43802918637SMatthew G. Knepley 
43902918637SMatthew G. Knepley     test:
44002918637SMatthew G. Knepley       suffix: triangular_prism
441*a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangular_prism \
442*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
44302918637SMatthew G. Knepley 
44402918637SMatthew G. Knepley     test:
44502918637SMatthew G. Knepley       suffix: tensor_triangular_prism
446*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_triangular_prism \
447*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
44802918637SMatthew G. Knepley 
44902918637SMatthew G. Knepley     test:
45002918637SMatthew G. Knepley       suffix: tensor_quadrilateral_prism
451*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quadrilateral_prism \
452*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
45302918637SMatthew G. Knepley 
45402918637SMatthew G. Knepley     test:
45502918637SMatthew G. Knepley       suffix: pyramid
456*a01d01faSMatthew G. Knepley       args: -dm_plex_cell pyramid \
457*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
458*a01d01faSMatthew G. Knepley 
459*a01d01faSMatthew G. Knepley   testset:
460*a01d01faSMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_reference_cell_domain -ref_arrangements -dm_plex_check_all \
461*a01d01faSMatthew G. Knepley           -ref_dm_view ::ascii_latex -dm_plex_view_tikzscale 1.0
46202918637SMatthew G. Knepley 
46302918637SMatthew G. Knepley     test:
46402918637SMatthew G. Knepley       suffix: ref_segment
465*a01d01faSMatthew G. Knepley       args: -dm_plex_cell segment \
466*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0
46702918637SMatthew G. Knepley 
46802918637SMatthew G. Knepley     test:
46902918637SMatthew G. Knepley       suffix: ref_triangle
470*a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangle \
471*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
47202918637SMatthew G. Knepley 
47302918637SMatthew G. Knepley     test:
47402918637SMatthew G. Knepley       suffix: ref_quadrilateral
475*a01d01faSMatthew G. Knepley       args: -dm_plex_cell quadrilateral \
476*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
47702918637SMatthew G. Knepley 
47802918637SMatthew G. Knepley     test:
47902918637SMatthew G. Knepley       suffix: ref_tensor_segment
480*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quad \
481*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0
48202918637SMatthew G. Knepley 
48302918637SMatthew G. Knepley     test:
48402918637SMatthew G. Knepley       suffix: ref_tetrahedron
485*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tetrahedron \
486*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
48702918637SMatthew G. Knepley 
48802918637SMatthew G. Knepley     test:
48902918637SMatthew G. Knepley       suffix: ref_hexahedron
490*a01d01faSMatthew G. Knepley       args: -dm_plex_cell hexahedron \
491*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
49202918637SMatthew G. Knepley 
49302918637SMatthew G. Knepley     test:
49402918637SMatthew G. Knepley       suffix: ref_triangular_prism
495*a01d01faSMatthew G. Knepley       args: -dm_plex_cell triangular_prism \
496*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
49702918637SMatthew G. Knepley 
49802918637SMatthew G. Knepley     test:
49902918637SMatthew G. Knepley       suffix: ref_tensor_triangular_prism
500*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_triangular_prism \
501*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
50202918637SMatthew G. Knepley 
50302918637SMatthew G. Knepley     test:
50402918637SMatthew G. Knepley       suffix: ref_tensor_quadrilateral_prism
505*a01d01faSMatthew G. Knepley       args: -dm_plex_cell tensor_quadrilateral_prism \
506*a01d01faSMatthew G. Knepley             -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0
50702918637SMatthew G. Knepley 
50802918637SMatthew G. Knepley   # ToBox should recreate the coordinate space since the cell shape changes
50902918637SMatthew G. Knepley   testset:
51002918637SMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all
51102918637SMatthew G. Knepley 
51202918637SMatthew G. Knepley     test:
51302918637SMatthew G. Knepley       suffix: tobox_triangle
51402918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
51502918637SMatthew 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
51602918637SMatthew G. Knepley 
51702918637SMatthew G. Knepley     test:
51802918637SMatthew G. Knepley       suffix: tobox_tensor_segment
51902918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
52002918637SMatthew 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
52102918637SMatthew G. Knepley 
52202918637SMatthew G. Knepley     test:
52302918637SMatthew G. Knepley       suffix: tobox_tetrahedron
52402918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
52502918637SMatthew 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
52602918637SMatthew G. Knepley 
52702918637SMatthew G. Knepley     test:
52802918637SMatthew G. Knepley       suffix: tobox_triangular_prism
52902918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
53002918637SMatthew 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
53102918637SMatthew G. Knepley 
53202918637SMatthew G. Knepley     test:
53302918637SMatthew G. Knepley       suffix: tobox_tensor_triangular_prism
53402918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
53502918637SMatthew 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
53602918637SMatthew G. Knepley 
53702918637SMatthew G. Knepley     test:
53802918637SMatthew G. Knepley       suffix: tobox_tensor_quadrilateral_prism
53902918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \
54002918637SMatthew 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
54102918637SMatthew G. Knepley 
54202918637SMatthew G. Knepley   testset:
54302918637SMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all
54402918637SMatthew G. Knepley 
54502918637SMatthew G. Knepley     test:
54602918637SMatthew G. Knepley       suffix: alfeld_triangle
54702918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \
54802918637SMatthew 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
54902918637SMatthew G. Knepley 
55002918637SMatthew G. Knepley     test:
55102918637SMatthew G. Knepley       suffix: alfeld_tetrahedron
55202918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \
55302918637SMatthew 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
55402918637SMatthew G. Knepley 
55502918637SMatthew G. Knepley   testset:
55602918637SMatthew G. Knepley     args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all
55702918637SMatthew G. Knepley 
55802918637SMatthew G. Knepley     # This splits edge 1 of the triangle, and reflects about the added edge
55902918637SMatthew G. Knepley     test:
56002918637SMatthew G. Knepley       suffix: sbr_triangle_0
56102918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \
56202918637SMatthew 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
56302918637SMatthew G. Knepley 
56402918637SMatthew G. Knepley     # This splits edge 0 of the triangle, and reflects about the added edge
56502918637SMatthew G. Knepley     test:
56602918637SMatthew G. Knepley       suffix: sbr_triangle_1
56702918637SMatthew 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 \
56802918637SMatthew 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
56902918637SMatthew G. Knepley 
57002918637SMatthew G. Knepley     # This splits edge 2 of the triangle, and reflects about the added edge
57102918637SMatthew G. Knepley     test:
57202918637SMatthew G. Knepley       suffix: sbr_triangle_2
57302918637SMatthew 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 \
57402918637SMatthew 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
57502918637SMatthew G. Knepley 
57602918637SMatthew G. Knepley     # This splits edges 1 and 2 of the triangle
57702918637SMatthew G. Knepley     test:
57802918637SMatthew G. Knepley       suffix: sbr_triangle_3
57902918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 0 \
58002918637SMatthew 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
58102918637SMatthew G. Knepley 
58202918637SMatthew G. Knepley     # This splits edges 1 and 0 of the triangle
58302918637SMatthew G. Knepley     test:
58402918637SMatthew G. Knepley       suffix: sbr_triangle_4
58502918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 0 \
58602918637SMatthew 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
58702918637SMatthew G. Knepley 
58802918637SMatthew G. Knepley     # This splits edges 0 and 1 of the triangle
58902918637SMatthew G. Knepley     test:
59002918637SMatthew G. Knepley       suffix: sbr_triangle_5
59102918637SMatthew 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 \
59202918637SMatthew 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
59302918637SMatthew G. Knepley 
59402918637SMatthew G. Knepley     # This splits edges 0 and 2 of the triangle
59502918637SMatthew G. Knepley     test:
59602918637SMatthew G. Knepley       suffix: sbr_triangle_6
59702918637SMatthew 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 \
59802918637SMatthew 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
59902918637SMatthew G. Knepley 
60002918637SMatthew G. Knepley     # This splits edges 2 and 0 of the triangle
60102918637SMatthew G. Knepley     test:
60202918637SMatthew G. Knepley       suffix: sbr_triangle_7
60302918637SMatthew 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 \
60402918637SMatthew 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
60502918637SMatthew G. Knepley 
60602918637SMatthew G. Knepley     # This splits edges 2 and 1 of the triangle
60702918637SMatthew G. Knepley     test:
60802918637SMatthew G. Knepley       suffix: sbr_triangle_8
60902918637SMatthew 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 \
61002918637SMatthew 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
61102918637SMatthew G. Knepley 
61202918637SMatthew G. Knepley   testset:
61302918637SMatthew G. Knepley     args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all
61402918637SMatthew G. Knepley 
61502918637SMatthew G. Knepley     test:
61602918637SMatthew G. Knepley       suffix: bl_tensor_segment
61702918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \
61802918637SMatthew 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
61902918637SMatthew G. Knepley 
62002918637SMatthew G. Knepley     # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check
62102918637SMatthew G. Knepley     test:
62202918637SMatthew G. Knepley       suffix: bl_tensor_triangular_prism
62302918637SMatthew G. Knepley       requires: TODO
62402918637SMatthew G. Knepley       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \
62502918637SMatthew 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
62602918637SMatthew G. Knepley 
62702918637SMatthew G. Knepley TEST*/
628