1*02918637SMatthew G. Knepley static char help[] = "Mesh Orientation Tutorial\n\n"; 2*02918637SMatthew G. Knepley 3*02918637SMatthew G. Knepley #include <petscdmplex.h> 4*02918637SMatthew G. Knepley #include <petscdmplextransform.h> 5*02918637SMatthew G. Knepley 6*02918637SMatthew G. Knepley typedef struct { 7*02918637SMatthew G. Knepley PetscBool genArr; /* Generate all possible cell arrangements */ 8*02918637SMatthew G. Knepley PetscBool refArr; /* Refine all possible cell arrangements */ 9*02918637SMatthew G. Knepley PetscBool printTable; /* Print the CAyley table */ 10*02918637SMatthew G. Knepley PetscInt orntBounds[2]; /* Bounds for the orientation check */ 11*02918637SMatthew G. Knepley PetscInt numOrnt; /* Number of specific orientations specified, or -1 for all orientations */ 12*02918637SMatthew G. Knepley PetscInt ornts[48]; /* Specific orientations if specified */ 13*02918637SMatthew G. Knepley PetscInt initOrnt; /* Initial orientation for starting mesh */ 14*02918637SMatthew G. Knepley } AppCtx; 15*02918637SMatthew G. Knepley 16*02918637SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 17*02918637SMatthew G. Knepley { 18*02918637SMatthew G. Knepley PetscInt n = 2; 19*02918637SMatthew G. Knepley PetscBool flg; 20*02918637SMatthew G. Knepley PetscErrorCode ierr; 21*02918637SMatthew G. Knepley 22*02918637SMatthew G. Knepley PetscFunctionBeginUser; 23*02918637SMatthew G. Knepley options->genArr = PETSC_FALSE; 24*02918637SMatthew G. Knepley options->refArr = PETSC_FALSE; 25*02918637SMatthew G. Knepley options->printTable = PETSC_FALSE; 26*02918637SMatthew G. Knepley options->orntBounds[0] = PETSC_MIN_INT; 27*02918637SMatthew G. Knepley options->orntBounds[1] = PETSC_MAX_INT; 28*02918637SMatthew G. Knepley options->numOrnt = -1; 29*02918637SMatthew G. Knepley options->initOrnt = 0; 30*02918637SMatthew G. Knepley 31*02918637SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Mesh Orientation Tutorials Options", "DMPLEX");CHKERRQ(ierr); 32*02918637SMatthew G. Knepley ierr = PetscOptionsBool("-gen_arrangements", "Flag for generating all arrangements of the cell", "ex11.c", options->genArr, &options->genArr, NULL);CHKERRQ(ierr); 33*02918637SMatthew G. Knepley ierr = PetscOptionsBool("-ref_arrangements", "Flag for refining all arrangements of the cell", "ex11.c", options->refArr, &options->refArr, NULL);CHKERRQ(ierr); 34*02918637SMatthew G. Knepley ierr = PetscOptionsBool("-print_table", "Print the Cayley table", "ex11.c", options->printTable, &options->printTable, NULL);CHKERRQ(ierr); 35*02918637SMatthew G. Knepley ierr = PetscOptionsIntArray("-ornt_bounds", "Bounds for orientation checks", "ex11.c", options->orntBounds, &n, NULL);CHKERRQ(ierr); 36*02918637SMatthew G. Knepley n = 48; 37*02918637SMatthew G. Knepley ierr = PetscOptionsIntArray("-ornts", "Specific orientations for checks", "ex11.c", options->ornts, &n, &flg);CHKERRQ(ierr); 38*02918637SMatthew G. Knepley if (flg) { 39*02918637SMatthew G. Knepley options->numOrnt = n; 40*02918637SMatthew G. Knepley ierr = PetscSortInt(n, options->ornts);CHKERRQ(ierr); 41*02918637SMatthew G. Knepley } 42*02918637SMatthew G. Knepley ierr = PetscOptionsInt("-init_ornt", "Initial orientation for starting mesh", "ex11.c", options->initOrnt, &options->initOrnt, NULL);CHKERRQ(ierr); 43*02918637SMatthew G. Knepley ierr = PetscOptionsEnd(); 44*02918637SMatthew G. Knepley PetscFunctionReturn(0); 45*02918637SMatthew G. Knepley } 46*02918637SMatthew G. Knepley 47*02918637SMatthew G. Knepley static PetscBool ignoreOrnt(AppCtx *user, PetscInt o) 48*02918637SMatthew G. Knepley { 49*02918637SMatthew G. Knepley PetscInt loc; 50*02918637SMatthew G. Knepley PetscErrorCode ierr; 51*02918637SMatthew G. Knepley 52*02918637SMatthew G. Knepley if (user->numOrnt < 0) return PETSC_FALSE; 53*02918637SMatthew G. Knepley ierr = PetscFindInt(o, user->numOrnt, user->ornts, &loc); 54*02918637SMatthew G. Knepley if (loc < 0 || ierr) return PETSC_TRUE; 55*02918637SMatthew G. Knepley return PETSC_FALSE; 56*02918637SMatthew G. Knepley } 57*02918637SMatthew G. Knepley 58*02918637SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 59*02918637SMatthew G. Knepley { 60*02918637SMatthew G. Knepley PetscErrorCode ierr; 61*02918637SMatthew G. Knepley 62*02918637SMatthew G. Knepley PetscFunctionBeginUser; 63*02918637SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 64*02918637SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 65*02918637SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 66*02918637SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 67*02918637SMatthew G. Knepley PetscFunctionReturn(0); 68*02918637SMatthew G. Knepley } 69*02918637SMatthew G. Knepley 70*02918637SMatthew G. Knepley static PetscErrorCode CheckCellVertices(DM dm, PetscInt cell, PetscInt o) 71*02918637SMatthew G. Knepley { 72*02918637SMatthew G. Knepley DMPolytopeType ct; 73*02918637SMatthew G. Knepley const PetscInt *arrVerts; 74*02918637SMatthew G. Knepley PetscInt *closure = NULL; 75*02918637SMatthew G. Knepley PetscInt Ncl, cl, Nv, vStart, vEnd, v; 76*02918637SMatthew G. Knepley MPI_Comm comm; 77*02918637SMatthew G. Knepley PetscErrorCode ierr; 78*02918637SMatthew G. Knepley 79*02918637SMatthew G. Knepley PetscFunctionBeginUser; 80*02918637SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 81*02918637SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr); 82*02918637SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 83*02918637SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr); 84*02918637SMatthew G. Knepley for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) { 85*02918637SMatthew G. Knepley const PetscInt vertex = closure[cl]; 86*02918637SMatthew G. Knepley 87*02918637SMatthew G. Knepley if (vertex < vStart || vertex >= vEnd) continue; 88*02918637SMatthew G. Knepley closure[Nv++] = vertex; 89*02918637SMatthew G. Knepley } 90*02918637SMatthew G. Knepley if (Nv != DMPolytopeTypeGetNumVertices(ct)) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Cell %D has %D vertices != %D vertices in a %s", cell, Nv, DMPolytopeTypeGetNumVertices(ct), DMPolytopeTypes[ct]); 91*02918637SMatthew G. Knepley arrVerts = DMPolytopeTypeGetVertexArrangment(ct, o); 92*02918637SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 93*02918637SMatthew G. Knepley if (closure[v] != arrVerts[v]+vStart) SETERRQ5(comm, PETSC_ERR_ARG_WRONG, "Cell %D vertex[%D]: %D should be %D for arrangement %D", cell, v, closure[v], arrVerts[v]+vStart, o); 94*02918637SMatthew G. Knepley } 95*02918637SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr); 96*02918637SMatthew G. Knepley PetscFunctionReturn(0); 97*02918637SMatthew G. Knepley } 98*02918637SMatthew G. Knepley 99*02918637SMatthew G. Knepley /* Transform cell with group operation o */ 100*02918637SMatthew G. Knepley static PetscErrorCode ReorientCell(DM dm, PetscInt cell, PetscInt o, PetscBool swapCoords) 101*02918637SMatthew G. Knepley { 102*02918637SMatthew G. Knepley DM cdm; 103*02918637SMatthew G. Knepley Vec coordinates; 104*02918637SMatthew G. Knepley PetscScalar *coords, *ccoords = NULL; 105*02918637SMatthew G. Knepley PetscInt *closure = NULL; 106*02918637SMatthew G. Knepley PetscInt cdim, d, Nc, Ncl, cl, vStart, vEnd, Nv; 107*02918637SMatthew G. Knepley PetscErrorCode ierr; 108*02918637SMatthew G. Knepley 109*02918637SMatthew G. Knepley PetscFunctionBegin; 110*02918637SMatthew G. Knepley /* Change vertex coordinates so that it plots as we expect */ 111*02918637SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 112*02918637SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 113*02918637SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 114*02918637SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords);CHKERRQ(ierr); 115*02918637SMatthew G. Knepley /* Reorient cone */ 116*02918637SMatthew G. Knepley ierr = DMPlexOrientPoint(dm, cell, o);CHKERRQ(ierr); 117*02918637SMatthew G. Knepley /* Finish resetting coordinates */ 118*02918637SMatthew G. Knepley if (swapCoords) { 119*02918637SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 120*02918637SMatthew G. Knepley ierr = VecGetArrayWrite(coordinates, &coords);CHKERRQ(ierr); 121*02918637SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr); 122*02918637SMatthew G. Knepley for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) { 123*02918637SMatthew G. Knepley const PetscInt vertex = closure[cl]; 124*02918637SMatthew G. Knepley PetscScalar *vcoords; 125*02918637SMatthew G. Knepley 126*02918637SMatthew G. Knepley if (vertex < vStart || vertex >= vEnd) continue; 127*02918637SMatthew G. Knepley ierr = DMPlexPointLocalRef(cdm, vertex, coords, &vcoords);CHKERRQ(ierr); 128*02918637SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = ccoords[Nv*cdim + d]; 129*02918637SMatthew G. Knepley ++Nv; 130*02918637SMatthew G. Knepley } 131*02918637SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr); 132*02918637SMatthew G. Knepley ierr = VecRestoreArrayWrite(coordinates, &coords);CHKERRQ(ierr); 133*02918637SMatthew G. Knepley } 134*02918637SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords);CHKERRQ(ierr); 135*02918637SMatthew G. Knepley PetscFunctionReturn(0); 136*02918637SMatthew G. Knepley } 137*02918637SMatthew G. Knepley 138*02918637SMatthew G. Knepley static PetscErrorCode GenerateArrangments(DM dm, AppCtx *user) 139*02918637SMatthew G. Knepley { 140*02918637SMatthew G. Knepley DM odm; 141*02918637SMatthew G. Knepley DMPolytopeType ct; 142*02918637SMatthew G. Knepley PetscInt No, o; 143*02918637SMatthew G. Knepley const char *name; 144*02918637SMatthew G. Knepley PetscErrorCode ierr; 145*02918637SMatthew G. Knepley 146*02918637SMatthew G. Knepley PetscFunctionBeginUser; 147*02918637SMatthew G. Knepley if (!user->genArr) PetscFunctionReturn(0); 148*02918637SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 149*02918637SMatthew G. Knepley ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr); 150*02918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 151*02918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 152*02918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 153*02918637SMatthew G. Knepley ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);CHKERRQ(ierr); 154*02918637SMatthew G. Knepley ierr = ReorientCell(odm, 0, o, PETSC_TRUE);CHKERRQ(ierr); 155*02918637SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);CHKERRQ(ierr); 156*02918637SMatthew G. Knepley ierr = DMViewFromOptions(odm, NULL, "-gen_dm_view");CHKERRQ(ierr); 157*02918637SMatthew G. Knepley ierr = CheckCellVertices(odm, 0, o);CHKERRQ(ierr); 158*02918637SMatthew G. Knepley ierr = DMDestroy(&odm);CHKERRQ(ierr); 159*02918637SMatthew G. Knepley } 160*02918637SMatthew G. Knepley PetscFunctionReturn(0); 161*02918637SMatthew G. Knepley } 162*02918637SMatthew G. Knepley 163*02918637SMatthew G. Knepley static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user) 164*02918637SMatthew G. Knepley { 165*02918637SMatthew G. Knepley DM dm1, dm2; 166*02918637SMatthew G. Knepley DMPolytopeType ct; 167*02918637SMatthew G. Knepley const PetscInt *refcone, *cone; 168*02918637SMatthew G. Knepley PetscInt No, o1, o2, o3, o4; 169*02918637SMatthew G. Knepley PetscBool equal; 170*02918637SMatthew G. Knepley const char *name; 171*02918637SMatthew G. Knepley PetscErrorCode ierr; 172*02918637SMatthew G. Knepley 173*02918637SMatthew G. Knepley PetscFunctionBeginUser; 174*02918637SMatthew G. Knepley if (!user->genArr) PetscFunctionReturn(0); 175*02918637SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 176*02918637SMatthew G. Knepley ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr); 177*02918637SMatthew G. Knepley ierr = DMPlexGetCone(dm, 0, &refcone);CHKERRQ(ierr); 178*02918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 179*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "Cayley Table for %s\n", DMPolytopeTypes[ct]);CHKERRQ(ierr);} 180*02918637SMatthew G. Knepley for (o1 = PetscMax(-No, user->orntBounds[0]); o1 < PetscMin(No, user->orntBounds[1]); ++o1) { 181*02918637SMatthew G. Knepley for (o2 = PetscMax(-No, user->orntBounds[0]); o2 < PetscMin(No, user->orntBounds[1]); ++o2) { 182*02918637SMatthew G. Knepley ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1);CHKERRQ(ierr); 183*02918637SMatthew G. Knepley ierr = DMPlexOrientPoint(dm1, 0, o2);CHKERRQ(ierr); 184*02918637SMatthew G. Knepley ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr); 185*02918637SMatthew G. Knepley ierr = DMPlexOrientPoint(dm1, 0, o1);CHKERRQ(ierr); 186*02918637SMatthew G. Knepley ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr); 187*02918637SMatthew G. Knepley o3 = DMPolytopeTypeComposeOrientation(ct, o1, o2);CHKERRQ(ierr); 188*02918637SMatthew G. Knepley /* First verification */ 189*02918637SMatthew G. Knepley ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2);CHKERRQ(ierr); 190*02918637SMatthew G. Knepley ierr = DMPlexOrientPoint(dm2, 0, o3);CHKERRQ(ierr); 191*02918637SMatthew G. Knepley ierr = DMPlexCheckFaces(dm2, 0);CHKERRQ(ierr); 192*02918637SMatthew G. Knepley ierr = DMPlexEqual(dm1, dm2, &equal);CHKERRQ(ierr); 193*02918637SMatthew G. Knepley if (!equal) { 194*02918637SMatthew G. Knepley ierr = DMViewFromOptions(dm1, NULL, "-error_dm_view");CHKERRQ(ierr); 195*02918637SMatthew G. Knepley ierr = DMViewFromOptions(dm2, NULL, "-error_dm_view");CHKERRQ(ierr); 196*02918637SMatthew G. Knepley SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D != %D", DMPolytopeTypes[ct], o1, o2, o3); 197*02918637SMatthew G. Knepley } 198*02918637SMatthew G. Knepley /* Second verification */ 199*02918637SMatthew G. Knepley ierr = DMPlexGetCone(dm1, 0, &cone);CHKERRQ(ierr); 200*02918637SMatthew G. Knepley ierr = DMPolytopeGetOrientation(ct, refcone, cone, &o4);CHKERRQ(ierr); 201*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "%D, ", o4);CHKERRQ(ierr);} 202*02918637SMatthew G. Knepley if (o3 != o4) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D = %D != %D", DMPolytopeTypes[ct], o1, o2, o3, o4); 203*02918637SMatthew G. Knepley ierr = DMDestroy(&dm1);CHKERRQ(ierr); 204*02918637SMatthew G. Knepley ierr = DMDestroy(&dm2);CHKERRQ(ierr); 205*02918637SMatthew G. Knepley } 206*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);} 207*02918637SMatthew G. Knepley } 208*02918637SMatthew G. Knepley PetscFunctionReturn(0); 209*02918637SMatthew G. Knepley } 210*02918637SMatthew G. Knepley 211*02918637SMatthew G. Knepley static PetscErrorCode VerifyInverse(DM dm, AppCtx *user) 212*02918637SMatthew G. Knepley { 213*02918637SMatthew G. Knepley DM dm1, dm2; 214*02918637SMatthew G. Knepley DMPolytopeType ct; 215*02918637SMatthew G. Knepley const PetscInt *refcone, *cone; 216*02918637SMatthew G. Knepley PetscInt No, o, oi, o2; 217*02918637SMatthew G. Knepley PetscBool equal; 218*02918637SMatthew G. Knepley const char *name; 219*02918637SMatthew G. Knepley PetscErrorCode ierr; 220*02918637SMatthew G. Knepley 221*02918637SMatthew G. Knepley PetscFunctionBeginUser; 222*02918637SMatthew G. Knepley if (!user->genArr) PetscFunctionReturn(0); 223*02918637SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 224*02918637SMatthew G. Knepley ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr); 225*02918637SMatthew G. Knepley ierr = DMPlexGetCone(dm, 0, &refcone);CHKERRQ(ierr); 226*02918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 227*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "Inverse table for %s\n", DMPolytopeTypes[ct]);CHKERRQ(ierr);} 228*02918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 229*02918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 230*02918637SMatthew G. Knepley oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o);CHKERRQ(ierr); 231*02918637SMatthew G. Knepley ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1);CHKERRQ(ierr); 232*02918637SMatthew G. Knepley ierr = DMPlexOrientPoint(dm1, 0, o);CHKERRQ(ierr); 233*02918637SMatthew G. Knepley ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr); 234*02918637SMatthew G. Knepley ierr = DMPlexOrientPoint(dm1, 0, oi);CHKERRQ(ierr); 235*02918637SMatthew G. Knepley ierr = DMPlexCheckFaces(dm1, 0);CHKERRQ(ierr); 236*02918637SMatthew G. Knepley /* First verification */ 237*02918637SMatthew G. Knepley ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2);CHKERRQ(ierr); 238*02918637SMatthew G. Knepley ierr = DMPlexEqual(dm1, dm2, &equal);CHKERRQ(ierr); 239*02918637SMatthew G. Knepley if (!equal) { 240*02918637SMatthew G. Knepley ierr = DMViewFromOptions(dm1, NULL, "-error_dm_view");CHKERRQ(ierr); 241*02918637SMatthew G. Knepley ierr = DMViewFromOptions(dm2, NULL, "-error_dm_view");CHKERRQ(ierr); 242*02918637SMatthew G. Knepley SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D != 0", DMPolytopeTypes[ct], o, oi); 243*02918637SMatthew G. Knepley } 244*02918637SMatthew G. Knepley /* Second verification */ 245*02918637SMatthew G. Knepley ierr = DMPlexGetCone(dm1, 0, &cone);CHKERRQ(ierr); 246*02918637SMatthew G. Knepley ierr = DMPolytopeGetOrientation(ct, refcone, cone, &o2);CHKERRQ(ierr); 247*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "%D, ", oi);CHKERRQ(ierr);} 248*02918637SMatthew G. Knepley if (o2 != 0) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D = %D != 0", DMPolytopeTypes[ct], o, oi, o2); 249*02918637SMatthew G. Knepley ierr = DMDestroy(&dm1);CHKERRQ(ierr); 250*02918637SMatthew G. Knepley ierr = DMDestroy(&dm2);CHKERRQ(ierr); 251*02918637SMatthew G. Knepley } 252*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);} 253*02918637SMatthew G. Knepley PetscFunctionReturn(0); 254*02918637SMatthew G. Knepley } 255*02918637SMatthew G. Knepley 256*02918637SMatthew G. Knepley /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */ 257*02918637SMatthew G. Knepley static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user) 258*02918637SMatthew G. Knepley { 259*02918637SMatthew G. Knepley DMPlexTransform tr, otr; 260*02918637SMatthew G. Knepley DMPolytopeType ct; 261*02918637SMatthew G. Knepley DMPolytopeType *rct; 262*02918637SMatthew G. Knepley const PetscInt *cone, *ornt, *ocone, *oornt; 263*02918637SMatthew G. Knepley PetscInt *rsize, *rcone, *rornt; 264*02918637SMatthew G. Knepley PetscInt Nct, n, oi, debug = 0; 265*02918637SMatthew G. Knepley PetscErrorCode ierr; 266*02918637SMatthew G. Knepley 267*02918637SMatthew G. Knepley PetscFunctionBeginUser; 268*02918637SMatthew G. Knepley ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr);CHKERRQ(ierr); 269*02918637SMatthew G. Knepley ierr = DMPlexTransformSetDM(tr, dm);CHKERRQ(ierr); 270*02918637SMatthew G. Knepley ierr = DMPlexTransformSetFromOptions(tr);CHKERRQ(ierr); 271*02918637SMatthew G. Knepley ierr = DMPlexTransformSetUp(tr);CHKERRQ(ierr); 272*02918637SMatthew G. Knepley 273*02918637SMatthew G. Knepley ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) odm), &otr);CHKERRQ(ierr); 274*02918637SMatthew G. Knepley ierr = DMPlexTransformSetDM(otr, odm);CHKERRQ(ierr); 275*02918637SMatthew G. Knepley ierr = DMPlexTransformSetFromOptions(otr);CHKERRQ(ierr); 276*02918637SMatthew G. Knepley ierr = DMPlexTransformSetUp(otr);CHKERRQ(ierr); 277*02918637SMatthew G. Knepley 278*02918637SMatthew G. Knepley ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr); 279*02918637SMatthew G. Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 280*02918637SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, p, &ornt);CHKERRQ(ierr); 281*02918637SMatthew G. Knepley ierr = DMPlexGetCone(odm, p, &ocone);CHKERRQ(ierr); 282*02918637SMatthew G. Knepley ierr = DMPlexGetConeOrientation(odm, p, &oornt);CHKERRQ(ierr); 283*02918637SMatthew G. Knepley oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o);CHKERRQ(ierr); 284*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "Orientation %D\n", oi);CHKERRQ(ierr);} 285*02918637SMatthew G. Knepley 286*02918637SMatthew G. Knepley ierr = DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt);CHKERRQ(ierr); 287*02918637SMatthew G. Knepley for (n = 0; n < Nct; ++n) { 288*02918637SMatthew G. Knepley DMPolytopeType ctNew = rct[n]; 289*02918637SMatthew G. Knepley PetscInt r, ro; 290*02918637SMatthew G. Knepley 291*02918637SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " Checking type %s\n", DMPolytopeTypes[ctNew]);CHKERRQ(ierr);} 292*02918637SMatthew G. Knepley for (r = 0; r < rsize[n]; ++r) { 293*02918637SMatthew G. Knepley const PetscInt *qcone, *qornt, *oqcone, *oqornt; 294*02918637SMatthew G. Knepley PetscInt pNew, opNew, oo, pr, fo; 295*02918637SMatthew G. Knepley PetscBool restore = PETSC_TRUE; 296*02918637SMatthew G. Knepley 297*02918637SMatthew G. Knepley ierr = DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew);CHKERRQ(ierr); 298*02918637SMatthew G. Knepley ierr = DMPlexTransformGetCone(tr, pNew, &qcone, &qornt);CHKERRQ(ierr); 299*02918637SMatthew G. Knepley if (debug) { 300*02918637SMatthew G. Knepley PetscInt c; 301*02918637SMatthew G. Knepley 302*02918637SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Checking replica %D (%D)\n Original Cone", r, pNew);CHKERRQ(ierr); 303*02918637SMatthew G. Knepley for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %D (%D)", qcone[c], qornt[c]);CHKERRQ(ierr);} 304*02918637SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 305*02918637SMatthew G. Knepley } 306*02918637SMatthew G. Knepley for (ro = 0; ro < rsize[n]; ++ro) { 307*02918637SMatthew G. Knepley PetscBool found; 308*02918637SMatthew G. Knepley 309*02918637SMatthew G. Knepley ierr = DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew);CHKERRQ(ierr); 310*02918637SMatthew G. Knepley ierr = DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt);CHKERRQ(ierr); 311*02918637SMatthew G. Knepley ierr = DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found);CHKERRQ(ierr); 312*02918637SMatthew G. Knepley if (found) break; 313*02918637SMatthew G. Knepley ierr = DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);CHKERRQ(ierr); 314*02918637SMatthew G. Knepley } 315*02918637SMatthew G. Knepley if (debug) { 316*02918637SMatthew G. Knepley PetscInt c; 317*02918637SMatthew G. Knepley 318*02918637SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Checking transform replica %D (%D) (%D)\n Transform Cone", ro, opNew, o);CHKERRQ(ierr); 319*02918637SMatthew G. Knepley for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %D (%D)", oqcone[c], oqornt[c]);CHKERRQ(ierr);} 320*02918637SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 321*02918637SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Matched %D\n", oo);CHKERRQ(ierr); 322*02918637SMatthew G. Knepley } 323*02918637SMatthew G. Knepley if (ro == rsize[n]) { 324*02918637SMatthew G. Knepley /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */ 325*02918637SMatthew G. Knepley if (ct == DM_POLYTOPE_TETRAHEDRON) { 326*02918637SMatthew G. Knepley /* The segment in a tetrahedron does not map into itself under the group action */ 327*02918637SMatthew G. Knepley if (ctNew == DM_POLYTOPE_SEGMENT) {restore = PETSC_FALSE; ro = r; oo = 0;} 328*02918637SMatthew G. Knepley /* The last four interior faces do not map into themselves under the group action */ 329*02918637SMatthew G. Knepley if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) {restore = PETSC_FALSE; ro = r; oo = 0;} 330*02918637SMatthew G. Knepley /* The last four interior faces do not map into themselves under the group action */ 331*02918637SMatthew G. Knepley if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) {restore = PETSC_FALSE; ro = r; oo = 0;} 332*02918637SMatthew G. Knepley } 333*02918637SMatthew G. Knepley if (restore) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to find matching %s %D orientation for cell orientation %D", DMPolytopeTypes[ctNew], r, o); 334*02918637SMatthew G. Knepley } 335*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "%D, %D, ", ro, oo);CHKERRQ(ierr);} 336*02918637SMatthew G. Knepley ierr = DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo);CHKERRQ(ierr); 337*02918637SMatthew G. Knepley if (!user->printTable) { 338*02918637SMatthew G. Knepley if (pr != ro) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong replica %D != %D", pr, ro); 339*02918637SMatthew G. Knepley if (fo != oo) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong orientation %D != %D", fo, oo); 340*02918637SMatthew G. Knepley } 341*02918637SMatthew G. Knepley ierr = DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt);CHKERRQ(ierr); 342*02918637SMatthew G. Knepley if (restore) {ierr = DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt);CHKERRQ(ierr);} 343*02918637SMatthew G. Knepley } 344*02918637SMatthew G. Knepley if (user->printTable) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);} 345*02918637SMatthew G. Knepley } 346*02918637SMatthew G. Knepley ierr = DMPlexTransformDestroy(&tr);CHKERRQ(ierr); 347*02918637SMatthew G. Knepley ierr = DMPlexTransformDestroy(&otr);CHKERRQ(ierr); 348*02918637SMatthew G. Knepley PetscFunctionReturn(0); 349*02918637SMatthew G. Knepley } 350*02918637SMatthew G. Knepley 351*02918637SMatthew G. Knepley static PetscErrorCode RefineArrangments(DM dm, AppCtx *user) 352*02918637SMatthew G. Knepley { 353*02918637SMatthew G. Knepley DM odm, rdm; 354*02918637SMatthew G. Knepley DMPolytopeType ct; 355*02918637SMatthew G. Knepley PetscInt No, o; 356*02918637SMatthew G. Knepley const char *name; 357*02918637SMatthew G. Knepley PetscErrorCode ierr; 358*02918637SMatthew G. Knepley 359*02918637SMatthew G. Knepley PetscFunctionBeginUser; 360*02918637SMatthew G. Knepley if (!user->refArr) PetscFunctionReturn(0); 361*02918637SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 362*02918637SMatthew G. Knepley ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr); 363*02918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 364*02918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 365*02918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 366*02918637SMatthew G. Knepley ierr = CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm);CHKERRQ(ierr); 367*02918637SMatthew G. Knepley if (user->initOrnt) {ierr = ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE);CHKERRQ(ierr);} 368*02918637SMatthew G. Knepley ierr = ReorientCell(odm, 0, o, PETSC_TRUE);CHKERRQ(ierr); 369*02918637SMatthew G. Knepley ierr = DMViewFromOptions(odm, NULL, "-orig_dm_view");CHKERRQ(ierr); 370*02918637SMatthew G. Knepley ierr = DMRefine(odm, MPI_COMM_NULL, &rdm);CHKERRQ(ierr); 371*02918637SMatthew G. Knepley ierr = DMSetFromOptions(rdm);CHKERRQ(ierr); 372*02918637SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o);CHKERRQ(ierr); 373*02918637SMatthew G. Knepley ierr = DMViewFromOptions(rdm, NULL, "-ref_dm_view");CHKERRQ(ierr); 374*02918637SMatthew G. Knepley ierr = CheckSubcells(dm, odm, 0, o, user);CHKERRQ(ierr); 375*02918637SMatthew G. Knepley ierr = DMDestroy(&odm);CHKERRQ(ierr); 376*02918637SMatthew G. Knepley ierr = DMDestroy(&rdm);CHKERRQ(ierr); 377*02918637SMatthew G. Knepley } 378*02918637SMatthew G. Knepley PetscFunctionReturn(0); 379*02918637SMatthew G. Knepley } 380*02918637SMatthew G. Knepley 381*02918637SMatthew G. Knepley int main(int argc, char **argv) 382*02918637SMatthew G. Knepley { 383*02918637SMatthew G. Knepley DM dm; 384*02918637SMatthew G. Knepley AppCtx user; 385*02918637SMatthew G. Knepley PetscErrorCode ierr; 386*02918637SMatthew G. Knepley 387*02918637SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 388*02918637SMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 389*02918637SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 390*02918637SMatthew G. Knepley if (user.initOrnt) { 391*02918637SMatthew G. Knepley ierr = ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE);CHKERRQ(ierr); 392*02918637SMatthew G. Knepley ierr = DMViewFromOptions(dm, NULL, "-ornt_dm_view");CHKERRQ(ierr); 393*02918637SMatthew G. Knepley } 394*02918637SMatthew G. Knepley ierr = GenerateArrangments(dm, &user);CHKERRQ(ierr); 395*02918637SMatthew G. Knepley ierr = VerifyCayleyTable(dm, &user);CHKERRQ(ierr); 396*02918637SMatthew G. Knepley ierr = VerifyInverse(dm, &user);CHKERRQ(ierr); 397*02918637SMatthew G. Knepley ierr = RefineArrangments(dm, &user);CHKERRQ(ierr); 398*02918637SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 399*02918637SMatthew G. Knepley ierr = PetscFinalize(); 400*02918637SMatthew G. Knepley return ierr; 401*02918637SMatthew G. Knepley } 402*02918637SMatthew G. Knepley 403*02918637SMatthew G. Knepley /*TEST 404*02918637SMatthew G. Knepley 405*02918637SMatthew G. Knepley test: 406*02918637SMatthew G. Knepley suffix: segment 407*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell segment -gen_arrangements \ 408*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 -dm_plex_view_tikzscale 0.5 409*02918637SMatthew G. Knepley 410*02918637SMatthew G. Knepley test: 411*02918637SMatthew G. Knepley suffix: triangle 412*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -gen_arrangements \ 413*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 0.5 414*02918637SMatthew G. Knepley 415*02918637SMatthew G. Knepley test: 416*02918637SMatthew G. Knepley suffix: quadrilateral 417*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -gen_arrangements \ 418*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 0.5 419*02918637SMatthew G. Knepley 420*02918637SMatthew G. Knepley test: 421*02918637SMatthew G. Knepley suffix: tensor_segment 422*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -gen_arrangements \ 423*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 -dm_plex_view_tikzscale 0.5 424*02918637SMatthew G. Knepley 425*02918637SMatthew G. Knepley test: 426*02918637SMatthew G. Knepley suffix: tetrahedron 427*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -gen_arrangements \ 428*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.5 429*02918637SMatthew G. Knepley 430*02918637SMatthew G. Knepley test: 431*02918637SMatthew G. Knepley suffix: hexahedron 432*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -gen_arrangements \ 433*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.3 434*02918637SMatthew G. Knepley 435*02918637SMatthew G. Knepley test: 436*02918637SMatthew G. Knepley suffix: triangular_prism 437*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -gen_arrangements \ 438*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.5 439*02918637SMatthew G. Knepley 440*02918637SMatthew G. Knepley test: 441*02918637SMatthew G. Knepley suffix: tensor_triangular_prism 442*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -gen_arrangements \ 443*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.5 444*02918637SMatthew G. Knepley 445*02918637SMatthew G. Knepley test: 446*02918637SMatthew G. Knepley suffix: tensor_quadrilateral_prism 447*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -gen_arrangements \ 448*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.5 449*02918637SMatthew G. Knepley 450*02918637SMatthew G. Knepley test: 451*02918637SMatthew G. Knepley suffix: pyramid 452*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell pyramid -gen_arrangements \ 453*02918637SMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 0.5 454*02918637SMatthew G. Knepley 455*02918637SMatthew G. Knepley test: 456*02918637SMatthew G. Knepley suffix: ref_segment 457*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell segment -ref_arrangements -dm_plex_check_all \ 458*02918637SMatthew G. Knepley -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 -dm_plex_view_tikzscale 1.0 459*02918637SMatthew G. Knepley 460*02918637SMatthew G. Knepley test: 461*02918637SMatthew G. Knepley suffix: ref_triangle 462*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -ref_arrangements -dm_plex_check_all \ 463*02918637SMatthew 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 464*02918637SMatthew G. Knepley 465*02918637SMatthew G. Knepley test: 466*02918637SMatthew G. Knepley suffix: ref_quadrilateral 467*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -ref_arrangements -dm_plex_check_all \ 468*02918637SMatthew 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 469*02918637SMatthew G. Knepley 470*02918637SMatthew G. Knepley test: 471*02918637SMatthew G. Knepley suffix: ref_tensor_segment 472*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -ref_arrangements -dm_plex_check_all \ 473*02918637SMatthew 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 474*02918637SMatthew G. Knepley 475*02918637SMatthew G. Knepley test: 476*02918637SMatthew G. Knepley suffix: ref_tetrahedron 477*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -ref_arrangements -dm_plex_check_all \ 478*02918637SMatthew 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 479*02918637SMatthew G. Knepley 480*02918637SMatthew G. Knepley test: 481*02918637SMatthew G. Knepley suffix: ref_hexahedron 482*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -ref_arrangements -dm_plex_check_all \ 483*02918637SMatthew 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 484*02918637SMatthew G. Knepley 485*02918637SMatthew G. Knepley test: 486*02918637SMatthew G. Knepley suffix: ref_triangular_prism 487*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -ref_arrangements -dm_plex_check_all \ 488*02918637SMatthew 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 489*02918637SMatthew G. Knepley 490*02918637SMatthew G. Knepley test: 491*02918637SMatthew G. Knepley suffix: ref_tensor_triangular_prism 492*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -ref_arrangements -dm_plex_check_all \ 493*02918637SMatthew 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 494*02918637SMatthew G. Knepley 495*02918637SMatthew G. Knepley test: 496*02918637SMatthew G. Knepley suffix: ref_tensor_quadrilateral_prism 497*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -ref_arrangements -dm_plex_check_all \ 498*02918637SMatthew 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 499*02918637SMatthew G. Knepley 500*02918637SMatthew G. Knepley # ToBox should recreate the coordinate space since the cell shape changes 501*02918637SMatthew G. Knepley testset: 502*02918637SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all 503*02918637SMatthew G. Knepley 504*02918637SMatthew G. Knepley test: 505*02918637SMatthew G. Knepley suffix: tobox_triangle 506*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \ 507*02918637SMatthew 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 508*02918637SMatthew G. Knepley 509*02918637SMatthew G. Knepley test: 510*02918637SMatthew G. Knepley suffix: tobox_tensor_segment 511*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \ 512*02918637SMatthew 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 513*02918637SMatthew G. Knepley 514*02918637SMatthew G. Knepley test: 515*02918637SMatthew G. Knepley suffix: tobox_tetrahedron 516*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \ 517*02918637SMatthew 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 518*02918637SMatthew G. Knepley 519*02918637SMatthew G. Knepley test: 520*02918637SMatthew G. Knepley suffix: tobox_triangular_prism 521*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 522*02918637SMatthew 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 523*02918637SMatthew G. Knepley 524*02918637SMatthew G. Knepley test: 525*02918637SMatthew G. Knepley suffix: tobox_tensor_triangular_prism 526*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \ 527*02918637SMatthew 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 528*02918637SMatthew G. Knepley 529*02918637SMatthew G. Knepley test: 530*02918637SMatthew G. Knepley suffix: tobox_tensor_quadrilateral_prism 531*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \ 532*02918637SMatthew 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 533*02918637SMatthew G. Knepley 534*02918637SMatthew G. Knepley testset: 535*02918637SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all 536*02918637SMatthew G. Knepley 537*02918637SMatthew G. Knepley test: 538*02918637SMatthew G. Knepley suffix: alfeld_triangle 539*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \ 540*02918637SMatthew 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 541*02918637SMatthew G. Knepley 542*02918637SMatthew G. Knepley test: 543*02918637SMatthew G. Knepley suffix: alfeld_tetrahedron 544*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \ 545*02918637SMatthew 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 546*02918637SMatthew G. Knepley 547*02918637SMatthew G. Knepley testset: 548*02918637SMatthew G. Knepley args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all 549*02918637SMatthew G. Knepley 550*02918637SMatthew G. Knepley # This splits edge 1 of the triangle, and reflects about the added edge 551*02918637SMatthew G. Knepley test: 552*02918637SMatthew G. Knepley suffix: sbr_triangle_0 553*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \ 554*02918637SMatthew 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 555*02918637SMatthew G. Knepley 556*02918637SMatthew G. Knepley # This splits edge 0 of the triangle, and reflects about the added edge 557*02918637SMatthew G. Knepley test: 558*02918637SMatthew G. Knepley suffix: sbr_triangle_1 559*02918637SMatthew 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 \ 560*02918637SMatthew 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 561*02918637SMatthew G. Knepley 562*02918637SMatthew G. Knepley # This splits edge 2 of the triangle, and reflects about the added edge 563*02918637SMatthew G. Knepley test: 564*02918637SMatthew G. Knepley suffix: sbr_triangle_2 565*02918637SMatthew 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 \ 566*02918637SMatthew 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 567*02918637SMatthew G. Knepley 568*02918637SMatthew G. Knepley # This splits edges 1 and 2 of the triangle 569*02918637SMatthew G. Knepley test: 570*02918637SMatthew G. Knepley suffix: sbr_triangle_3 571*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 0 \ 572*02918637SMatthew 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 573*02918637SMatthew G. Knepley 574*02918637SMatthew G. Knepley # This splits edges 1 and 0 of the triangle 575*02918637SMatthew G. Knepley test: 576*02918637SMatthew G. Knepley suffix: sbr_triangle_4 577*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 0 \ 578*02918637SMatthew 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 579*02918637SMatthew G. Knepley 580*02918637SMatthew G. Knepley # This splits edges 0 and 1 of the triangle 581*02918637SMatthew G. Knepley test: 582*02918637SMatthew G. Knepley suffix: sbr_triangle_5 583*02918637SMatthew 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 \ 584*02918637SMatthew 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 585*02918637SMatthew G. Knepley 586*02918637SMatthew G. Knepley # This splits edges 0 and 2 of the triangle 587*02918637SMatthew G. Knepley test: 588*02918637SMatthew G. Knepley suffix: sbr_triangle_6 589*02918637SMatthew 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 \ 590*02918637SMatthew 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 591*02918637SMatthew G. Knepley 592*02918637SMatthew G. Knepley # This splits edges 2 and 0 of the triangle 593*02918637SMatthew G. Knepley test: 594*02918637SMatthew G. Knepley suffix: sbr_triangle_7 595*02918637SMatthew 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 \ 596*02918637SMatthew 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 597*02918637SMatthew G. Knepley 598*02918637SMatthew G. Knepley # This splits edges 2 and 1 of the triangle 599*02918637SMatthew G. Knepley test: 600*02918637SMatthew G. Knepley suffix: sbr_triangle_8 601*02918637SMatthew 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 \ 602*02918637SMatthew 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 603*02918637SMatthew G. Knepley 604*02918637SMatthew G. Knepley testset: 605*02918637SMatthew G. Knepley args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all 606*02918637SMatthew G. Knepley 607*02918637SMatthew G. Knepley test: 608*02918637SMatthew G. Knepley suffix: bl_tensor_segment 609*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \ 610*02918637SMatthew 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 611*02918637SMatthew G. Knepley 612*02918637SMatthew G. Knepley # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check 613*02918637SMatthew G. Knepley test: 614*02918637SMatthew G. Knepley suffix: bl_tensor_triangular_prism 615*02918637SMatthew G. Knepley requires: TODO 616*02918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \ 617*02918637SMatthew 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 618*02918637SMatthew G. Knepley 619*02918637SMatthew G. Knepley TEST*/ 620