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 16d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 17d71ae5a4SJacob Faibussowitsch { 1802918637SMatthew G. Knepley PetscInt n = 2; 1902918637SMatthew G. Knepley PetscBool flg; 2002918637SMatthew G. Knepley 2102918637SMatthew G. Knepley PetscFunctionBeginUser; 2202918637SMatthew G. Knepley options->genArr = PETSC_FALSE; 2302918637SMatthew G. Knepley options->refArr = PETSC_FALSE; 2402918637SMatthew G. Knepley options->printTable = PETSC_FALSE; 25*1690c2aeSBarry Smith options->orntBounds[0] = PETSC_INT_MIN; 26*1690c2aeSBarry Smith options->orntBounds[1] = PETSC_INT_MAX; 2702918637SMatthew G. Knepley options->numOrnt = -1; 2802918637SMatthew G. Knepley options->initOrnt = 0; 2902918637SMatthew G. Knepley 30d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Mesh Orientation Tutorials Options", "DMPLEX"); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-gen_arrangements", "Flag for generating all arrangements of the cell", "ex11.c", options->genArr, &options->genArr, NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ref_arrangements", "Flag for refining all arrangements of the cell", "ex11.c", options->refArr, &options->refArr, NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-print_table", "Print the Cayley table", "ex11.c", options->printTable, &options->printTable, NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-ornt_bounds", "Bounds for orientation checks", "ex11.c", options->orntBounds, &n, NULL)); 3502918637SMatthew G. Knepley n = 48; 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-ornts", "Specific orientations for checks", "ex11.c", options->ornts, &n, &flg)); 3702918637SMatthew G. Knepley if (flg) { 3802918637SMatthew G. Knepley options->numOrnt = n; 399566063dSJacob Faibussowitsch PetscCall(PetscSortInt(n, options->ornts)); 4002918637SMatthew G. Knepley } 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-init_ornt", "Initial orientation for starting mesh", "ex11.c", options->initOrnt, &options->initOrnt, NULL)); 42d0609cedSBarry Smith PetscOptionsEnd(); 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4402918637SMatthew G. Knepley } 4502918637SMatthew G. Knepley 46d71ae5a4SJacob Faibussowitsch static PetscBool ignoreOrnt(AppCtx *user, PetscInt o) 47d71ae5a4SJacob Faibussowitsch { 4802918637SMatthew G. Knepley PetscInt loc; 4902918637SMatthew G. Knepley PetscErrorCode ierr; 5002918637SMatthew G. Knepley 5102918637SMatthew G. Knepley if (user->numOrnt < 0) return PETSC_FALSE; 5202918637SMatthew G. Knepley ierr = PetscFindInt(o, user->numOrnt, user->ornts, &loc); 5302918637SMatthew G. Knepley if (loc < 0 || ierr) return PETSC_TRUE; 5402918637SMatthew G. Knepley return PETSC_FALSE; 5502918637SMatthew G. Knepley } 5602918637SMatthew G. Knepley 57d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 58d71ae5a4SJacob Faibussowitsch { 5902918637SMatthew G. Knepley PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 619566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 629566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 639566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6502918637SMatthew G. Knepley } 6602918637SMatthew G. Knepley 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckCellVertices(DM dm, PetscInt cell, PetscInt o) 68d71ae5a4SJacob Faibussowitsch { 6902918637SMatthew G. Knepley DMPolytopeType ct; 7002918637SMatthew G. Knepley const PetscInt *arrVerts; 7102918637SMatthew G. Knepley PetscInt *closure = NULL; 7202918637SMatthew G. Knepley PetscInt Ncl, cl, Nv, vStart, vEnd, v; 7302918637SMatthew G. Knepley MPI_Comm comm; 7402918637SMatthew G. Knepley 7502918637SMatthew G. Knepley PetscFunctionBeginUser; 769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 779566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cell, &ct)); 789566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 799566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 8002918637SMatthew G. Knepley for (cl = 0, Nv = 0; cl < Ncl * 2; cl += 2) { 8102918637SMatthew G. Knepley const PetscInt vertex = closure[cl]; 8202918637SMatthew G. Knepley 8302918637SMatthew G. Knepley if (vertex < vStart || vertex >= vEnd) continue; 8402918637SMatthew G. Knepley closure[Nv++] = vertex; 8502918637SMatthew G. Knepley } 8663a3b9bcSJacob Faibussowitsch PetscCheck(Nv == DMPolytopeTypeGetNumVertices(ct), comm, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " has %" PetscInt_FMT " vertices != %" PetscInt_FMT " vertices in a %s", cell, Nv, DMPolytopeTypeGetNumVertices(ct), DMPolytopeTypes[ct]); 8785036b15SMatthew G. Knepley arrVerts = DMPolytopeTypeGetVertexArrangement(ct, o); 8802918637SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 8963a3b9bcSJacob Faibussowitsch PetscCheck(closure[v] == arrVerts[v] + vStart, comm, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " vertex[%" PetscInt_FMT "]: %" PetscInt_FMT " should be %" PetscInt_FMT " for arrangement %" PetscInt_FMT, cell, v, closure[v], arrVerts[v] + vStart, o); 9002918637SMatthew G. Knepley } 919566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9302918637SMatthew G. Knepley } 9402918637SMatthew G. Knepley 9502918637SMatthew G. Knepley /* Transform cell with group operation o */ 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode ReorientCell(DM dm, PetscInt cell, PetscInt o, PetscBool swapCoords) 97d71ae5a4SJacob Faibussowitsch { 9802918637SMatthew G. Knepley DM cdm; 9902918637SMatthew G. Knepley Vec coordinates; 10002918637SMatthew G. Knepley PetscScalar *coords, *ccoords = NULL; 10102918637SMatthew G. Knepley PetscInt *closure = NULL; 10202918637SMatthew G. Knepley PetscInt cdim, d, Nc, Ncl, cl, vStart, vEnd, Nv; 10302918637SMatthew G. Knepley 10402918637SMatthew G. Knepley PetscFunctionBegin; 10502918637SMatthew G. Knepley /* Change vertex coordinates so that it plots as we expect */ 1069566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 1089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1099566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords)); 11002918637SMatthew G. Knepley /* Reorient cone */ 1119566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(dm, cell, o)); 11202918637SMatthew G. Knepley /* Finish resetting coordinates */ 11302918637SMatthew G. Knepley if (swapCoords) { 1149566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(coordinates, &coords)); 1169566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 11702918637SMatthew G. Knepley for (cl = 0, Nv = 0; cl < Ncl * 2; cl += 2) { 11802918637SMatthew G. Knepley const PetscInt vertex = closure[cl]; 11902918637SMatthew G. Knepley PetscScalar *vcoords; 12002918637SMatthew G. Knepley 12102918637SMatthew G. Knepley if (vertex < vStart || vertex >= vEnd) continue; 1229566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(cdm, vertex, coords, &vcoords)); 12302918637SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = ccoords[Nv * cdim + d]; 12402918637SMatthew G. Knepley ++Nv; 12502918637SMatthew G. Knepley } 1269566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(coordinates, &coords)); 12802918637SMatthew G. Knepley } 1299566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords)); 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13102918637SMatthew G. Knepley } 13202918637SMatthew G. Knepley 13385036b15SMatthew G. Knepley static PetscErrorCode GenerateArrangements(DM dm, AppCtx *user) 134d71ae5a4SJacob Faibussowitsch { 13502918637SMatthew G. Knepley DM odm; 13602918637SMatthew G. Knepley DMPolytopeType ct; 13702918637SMatthew G. Knepley PetscInt No, o; 13802918637SMatthew G. Knepley const char *name; 13902918637SMatthew G. Knepley 14002918637SMatthew G. Knepley PetscFunctionBeginUser; 1413ba16761SJacob Faibussowitsch if (!user->genArr) PetscFunctionReturn(PETSC_SUCCESS); 1429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 1439566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, 0, &ct)); 14485036b15SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangements(ct) / 2; 14502918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 14602918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 1479566063dSJacob Faibussowitsch PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &odm)); 1489566063dSJacob Faibussowitsch PetscCall(ReorientCell(odm, 0, o, PETSC_TRUE)); 14963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%s orientation %" PetscInt_FMT "\n", name, o)); 1509566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(odm, NULL, "-gen_dm_view")); 1519566063dSJacob Faibussowitsch PetscCall(CheckCellVertices(odm, 0, o)); 1529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm)); 15302918637SMatthew G. Knepley } 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15502918637SMatthew G. Knepley } 15602918637SMatthew G. Knepley 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user) 158d71ae5a4SJacob Faibussowitsch { 15902918637SMatthew G. Knepley DM dm1, dm2; 16002918637SMatthew G. Knepley DMPolytopeType ct; 16102918637SMatthew G. Knepley const PetscInt *refcone, *cone; 16202918637SMatthew G. Knepley PetscInt No, o1, o2, o3, o4; 16302918637SMatthew G. Knepley PetscBool equal; 16402918637SMatthew G. Knepley const char *name; 16502918637SMatthew G. Knepley 16602918637SMatthew G. Knepley PetscFunctionBeginUser; 1673ba16761SJacob Faibussowitsch if (!user->genArr) PetscFunctionReturn(PETSC_SUCCESS); 1689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 1699566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, 0, &ct)); 1709566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, 0, &refcone)); 17185036b15SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangements(ct) / 2; 1729566063dSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cayley Table for %s\n", DMPolytopeTypes[ct])); 17302918637SMatthew G. Knepley for (o1 = PetscMax(-No, user->orntBounds[0]); o1 < PetscMin(No, user->orntBounds[1]); ++o1) { 17402918637SMatthew G. Knepley for (o2 = PetscMax(-No, user->orntBounds[0]); o2 < PetscMin(No, user->orntBounds[1]); ++o2) { 1759566063dSJacob Faibussowitsch PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm1)); 1769566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(dm1, 0, o2)); 1779566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(dm1, 0)); 1789566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(dm1, 0, o1)); 1799566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(dm1, 0)); 1805f80ce2aSJacob Faibussowitsch o3 = DMPolytopeTypeComposeOrientation(ct, o1, o2); 18102918637SMatthew G. Knepley /* First verification */ 1829566063dSJacob Faibussowitsch PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm2)); 1839566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(dm2, 0, o3)); 1849566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(dm2, 0)); 1859566063dSJacob Faibussowitsch PetscCall(DMPlexEqual(dm1, dm2, &equal)); 18602918637SMatthew G. Knepley if (!equal) { 1879566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm1, NULL, "-error_dm_view")); 1889566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm2, NULL, "-error_dm_view")); 18963a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " != %" PetscInt_FMT, DMPolytopeTypes[ct], o1, o2, o3); 19002918637SMatthew G. Knepley } 19102918637SMatthew G. Knepley /* Second verification */ 1929566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm1, 0, &cone)); 1939566063dSJacob Faibussowitsch PetscCall(DMPolytopeGetOrientation(ct, refcone, cone, &o4)); 19463a3b9bcSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ", ", o4)); 19563a3b9bcSJacob Faibussowitsch PetscCheck(o3 == o4, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " = %" PetscInt_FMT " != %" PetscInt_FMT, DMPolytopeTypes[ct], o1, o2, o3, o4); 1969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm1)); 1979566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm2)); 19802918637SMatthew G. Knepley } 1999566063dSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 20002918637SMatthew G. Knepley } 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20202918637SMatthew G. Knepley } 20302918637SMatthew G. Knepley 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode VerifyInverse(DM dm, AppCtx *user) 205d71ae5a4SJacob Faibussowitsch { 20602918637SMatthew G. Knepley DM dm1, dm2; 20702918637SMatthew G. Knepley DMPolytopeType ct; 20802918637SMatthew G. Knepley const PetscInt *refcone, *cone; 20902918637SMatthew G. Knepley PetscInt No, o, oi, o2; 21002918637SMatthew G. Knepley PetscBool equal; 21102918637SMatthew G. Knepley const char *name; 21202918637SMatthew G. Knepley 21302918637SMatthew G. Knepley PetscFunctionBeginUser; 2143ba16761SJacob Faibussowitsch if (!user->genArr) PetscFunctionReturn(PETSC_SUCCESS); 2159566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 2169566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, 0, &ct)); 2179566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, 0, &refcone)); 21885036b15SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangements(ct) / 2; 2199566063dSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Inverse table for %s\n", DMPolytopeTypes[ct])); 22002918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 22102918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 2225f80ce2aSJacob Faibussowitsch oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o); 2239566063dSJacob Faibussowitsch PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm1)); 2249566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(dm1, 0, o)); 2259566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(dm1, 0)); 2269566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(dm1, 0, oi)); 2279566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(dm1, 0)); 22802918637SMatthew G. Knepley /* First verification */ 2299566063dSJacob Faibussowitsch PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &dm2)); 2309566063dSJacob Faibussowitsch PetscCall(DMPlexEqual(dm1, dm2, &equal)); 23102918637SMatthew G. Knepley if (!equal) { 2329566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm1, NULL, "-error_dm_view")); 2339566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm2, NULL, "-error_dm_view")); 23463a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " != 0", DMPolytopeTypes[ct], o, oi); 23502918637SMatthew G. Knepley } 23602918637SMatthew G. Knepley /* Second verification */ 2379566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm1, 0, &cone)); 2389566063dSJacob Faibussowitsch PetscCall(DMPolytopeGetOrientation(ct, refcone, cone, &o2)); 23963a3b9bcSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ", ", oi)); 24063a3b9bcSJacob Faibussowitsch PetscCheck(o2 == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %" PetscInt_FMT " * %" PetscInt_FMT " = %" PetscInt_FMT " != 0", DMPolytopeTypes[ct], o, oi, o2); 2419566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm1)); 2429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm2)); 24302918637SMatthew G. Knepley } 2449566063dSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24602918637SMatthew G. Knepley } 24702918637SMatthew G. Knepley 24802918637SMatthew G. Knepley /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */ 249d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user) 250d71ae5a4SJacob Faibussowitsch { 25102918637SMatthew G. Knepley DMPlexTransform tr, otr; 25202918637SMatthew G. Knepley DMPolytopeType ct; 25302918637SMatthew G. Knepley DMPolytopeType *rct; 25402918637SMatthew G. Knepley const PetscInt *cone, *ornt, *ocone, *oornt; 25502918637SMatthew G. Knepley PetscInt *rsize, *rcone, *rornt; 25602918637SMatthew G. Knepley PetscInt Nct, n, oi, debug = 0; 25702918637SMatthew G. Knepley 25802918637SMatthew G. Knepley PetscFunctionBeginUser; 2599566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 2609566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetDM(tr, dm)); 2619566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetFromOptions(tr)); 2629566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(tr)); 26302918637SMatthew G. Knepley 2649566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)odm), &otr)); 2659566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetDM(otr, odm)); 2669566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetFromOptions(otr)); 2679566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(otr)); 26802918637SMatthew G. Knepley 2699566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, p, &ct)); 2709566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 2719566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, p, &ornt)); 2729566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(odm, p, &ocone)); 2739566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(odm, p, &oornt)); 2745f80ce2aSJacob Faibussowitsch oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o); 27563a3b9bcSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Orientation %" PetscInt_FMT "\n", oi)); 27602918637SMatthew G. Knepley 2779566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt)); 27802918637SMatthew G. Knepley for (n = 0; n < Nct; ++n) { 27902918637SMatthew G. Knepley DMPolytopeType ctNew = rct[n]; 28002918637SMatthew G. Knepley PetscInt r, ro; 28102918637SMatthew G. Knepley 2829566063dSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Checking type %s\n", DMPolytopeTypes[ctNew])); 28302918637SMatthew G. Knepley for (r = 0; r < rsize[n]; ++r) { 28402918637SMatthew G. Knepley const PetscInt *qcone, *qornt, *oqcone, *oqornt; 28502918637SMatthew G. Knepley PetscInt pNew, opNew, oo, pr, fo; 28602918637SMatthew G. Knepley PetscBool restore = PETSC_TRUE; 28702918637SMatthew G. Knepley 2889566063dSJacob Faibussowitsch PetscCall(DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew)); 2899566063dSJacob Faibussowitsch PetscCall(DMPlexTransformGetCone(tr, pNew, &qcone, &qornt)); 29002918637SMatthew G. Knepley if (debug) { 29102918637SMatthew G. Knepley PetscInt c; 29202918637SMatthew G. Knepley 29363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Checking replica %" PetscInt_FMT " (%" PetscInt_FMT ")\n Original Cone", r, pNew)); 29463a3b9bcSJacob Faibussowitsch for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT " (%" PetscInt_FMT ")", qcone[c], qornt[c])); 2959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 29602918637SMatthew G. Knepley } 29702918637SMatthew G. Knepley for (ro = 0; ro < rsize[n]; ++ro) { 29802918637SMatthew G. Knepley PetscBool found; 29902918637SMatthew G. Knepley 3009566063dSJacob Faibussowitsch PetscCall(DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew)); 3019566063dSJacob Faibussowitsch PetscCall(DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt)); 3029566063dSJacob Faibussowitsch PetscCall(DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found)); 30302918637SMatthew G. Knepley if (found) break; 3049566063dSJacob Faibussowitsch PetscCall(DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt)); 30502918637SMatthew G. Knepley } 30602918637SMatthew G. Knepley if (debug) { 30702918637SMatthew G. Knepley PetscInt c; 30802918637SMatthew G. Knepley 30963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Checking transform replica %" PetscInt_FMT " (%" PetscInt_FMT ") (%" PetscInt_FMT ")\n Transform Cone", ro, opNew, o)); 31063a3b9bcSJacob Faibussowitsch for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT " (%" PetscInt_FMT ")", oqcone[c], oqornt[c])); 3119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 31263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Matched %" PetscInt_FMT "\n", oo)); 31302918637SMatthew G. Knepley } 31402918637SMatthew G. Knepley if (ro == rsize[n]) { 31502918637SMatthew G. Knepley /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */ 31602918637SMatthew G. Knepley if (ct == DM_POLYTOPE_TETRAHEDRON) { 31702918637SMatthew G. Knepley /* The segment in a tetrahedron does not map into itself under the group action */ 3189371c9d4SSatish Balay if (ctNew == DM_POLYTOPE_SEGMENT) { 3199371c9d4SSatish Balay restore = PETSC_FALSE; 3209371c9d4SSatish Balay ro = r; 3219371c9d4SSatish Balay oo = 0; 3229371c9d4SSatish Balay } 32302918637SMatthew G. Knepley /* The last four interior faces do not map into themselves under the group action */ 3249371c9d4SSatish Balay if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) { 3259371c9d4SSatish Balay restore = PETSC_FALSE; 3269371c9d4SSatish Balay ro = r; 3279371c9d4SSatish Balay oo = 0; 3289371c9d4SSatish Balay } 32902918637SMatthew G. Knepley /* The last four interior faces do not map into themselves under the group action */ 3309371c9d4SSatish Balay if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) { 3319371c9d4SSatish Balay restore = PETSC_FALSE; 3329371c9d4SSatish Balay ro = r; 3339371c9d4SSatish Balay oo = 0; 3349371c9d4SSatish Balay } 33502918637SMatthew G. Knepley } 33663a3b9bcSJacob Faibussowitsch PetscCheck(!restore, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to find matching %s %" PetscInt_FMT " orientation for cell orientation %" PetscInt_FMT, DMPolytopeTypes[ctNew], r, o); 33702918637SMatthew G. Knepley } 33863a3b9bcSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ", %" PetscInt_FMT ", ", ro, oo)); 3399566063dSJacob Faibussowitsch PetscCall(DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo)); 34002918637SMatthew G. Knepley if (!user->printTable) { 34163a3b9bcSJacob Faibussowitsch PetscCheck(pr == ro, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong replica %" PetscInt_FMT " != %" PetscInt_FMT, pr, ro); 34263a3b9bcSJacob Faibussowitsch PetscCheck(fo == oo, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong orientation %" PetscInt_FMT " != %" PetscInt_FMT, fo, oo); 34302918637SMatthew G. Knepley } 3449566063dSJacob Faibussowitsch PetscCall(DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt)); 3459566063dSJacob Faibussowitsch if (restore) PetscCall(DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt)); 34602918637SMatthew G. Knepley } 3479566063dSJacob Faibussowitsch if (user->printTable) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 34802918637SMatthew G. Knepley } 3499566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 3509566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&otr)); 3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35202918637SMatthew G. Knepley } 35302918637SMatthew G. Knepley 35485036b15SMatthew G. Knepley static PetscErrorCode RefineArrangements(DM dm, AppCtx *user) 355d71ae5a4SJacob Faibussowitsch { 35602918637SMatthew G. Knepley DM odm, rdm; 35702918637SMatthew G. Knepley DMPolytopeType ct; 35802918637SMatthew G. Knepley PetscInt No, o; 35902918637SMatthew G. Knepley const char *name; 36002918637SMatthew G. Knepley 36102918637SMatthew G. Knepley PetscFunctionBeginUser; 3623ba16761SJacob Faibussowitsch if (!user->refArr) PetscFunctionReturn(PETSC_SUCCESS); 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 3649566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, 0, &ct)); 36585036b15SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangements(ct) / 2; 36602918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 36702918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 3689566063dSJacob Faibussowitsch PetscCall(CreateMesh(PetscObjectComm((PetscObject)dm), user, &odm)); 3699566063dSJacob Faibussowitsch if (user->initOrnt) PetscCall(ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE)); 3709566063dSJacob Faibussowitsch PetscCall(ReorientCell(odm, 0, o, PETSC_TRUE)); 3719566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(odm, NULL, "-orig_dm_view")); 3729566063dSJacob Faibussowitsch PetscCall(DMRefine(odm, MPI_COMM_NULL, &rdm)); 3739566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(rdm)); 37463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%s orientation %" PetscInt_FMT "\n", name, o)); 3759566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(rdm, NULL, "-ref_dm_view")); 3769566063dSJacob Faibussowitsch PetscCall(CheckSubcells(dm, odm, 0, o, user)); 3779566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm)); 3789566063dSJacob Faibussowitsch PetscCall(DMDestroy(&rdm)); 37902918637SMatthew G. Knepley } 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38102918637SMatthew G. Knepley } 38202918637SMatthew G. Knepley 383d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 384d71ae5a4SJacob Faibussowitsch { 38502918637SMatthew G. Knepley DM dm; 38602918637SMatthew G. Knepley AppCtx user; 38702918637SMatthew G. Knepley 388327415f7SBarry Smith PetscFunctionBeginUser; 3899566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3909566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 3919566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 39202918637SMatthew G. Knepley if (user.initOrnt) { 3939566063dSJacob Faibussowitsch PetscCall(ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE)); 3949566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-ornt_dm_view")); 39502918637SMatthew G. Knepley } 39685036b15SMatthew G. Knepley PetscCall(GenerateArrangements(dm, &user)); 3979566063dSJacob Faibussowitsch PetscCall(VerifyCayleyTable(dm, &user)); 3989566063dSJacob Faibussowitsch PetscCall(VerifyInverse(dm, &user)); 39985036b15SMatthew G. Knepley PetscCall(RefineArrangements(dm, &user)); 4009566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4019566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 402b122ec5aSJacob Faibussowitsch return 0; 40302918637SMatthew G. Knepley } 40402918637SMatthew G. Knepley 40502918637SMatthew G. Knepley /*TEST 40602918637SMatthew G. Knepley 407a01d01faSMatthew G. Knepley testset: 408a01d01faSMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_reference_cell_domain -gen_arrangements \ 409a01d01faSMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_tikzscale 0.5 410a01d01faSMatthew G. Knepley 41102918637SMatthew G. Knepley test: 41202918637SMatthew G. Knepley suffix: segment 413a01d01faSMatthew G. Knepley args: -dm_plex_cell segment \ 414a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 41502918637SMatthew G. Knepley 41602918637SMatthew G. Knepley test: 41702918637SMatthew G. Knepley suffix: triangle 418a01d01faSMatthew G. Knepley args: -dm_plex_cell triangle \ 419a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 42002918637SMatthew G. Knepley 42102918637SMatthew G. Knepley test: 42202918637SMatthew G. Knepley suffix: quadrilateral 423a01d01faSMatthew G. Knepley args: -dm_plex_cell quadrilateral \ 424a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 42502918637SMatthew G. Knepley 42602918637SMatthew G. Knepley test: 42702918637SMatthew G. Knepley suffix: tensor_segment 428a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quad \ 429a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 43002918637SMatthew G. Knepley 43102918637SMatthew G. Knepley test: 43202918637SMatthew G. Knepley suffix: tetrahedron 433a01d01faSMatthew G. Knepley args: -dm_plex_cell tetrahedron \ 434a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 43502918637SMatthew G. Knepley 43602918637SMatthew G. Knepley test: 43702918637SMatthew G. Knepley suffix: hexahedron 438a01d01faSMatthew G. Knepley args: -dm_plex_cell hexahedron \ 439a01d01faSMatthew 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 44002918637SMatthew G. Knepley 44102918637SMatthew G. Knepley test: 44202918637SMatthew G. Knepley suffix: triangular_prism 443a01d01faSMatthew G. Knepley args: -dm_plex_cell triangular_prism \ 444a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 44502918637SMatthew G. Knepley 44602918637SMatthew G. Knepley test: 44702918637SMatthew G. Knepley suffix: tensor_triangular_prism 448a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_triangular_prism \ 449a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 45002918637SMatthew G. Knepley 45102918637SMatthew G. Knepley test: 45202918637SMatthew G. Knepley suffix: tensor_quadrilateral_prism 453a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quadrilateral_prism \ 454a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 45502918637SMatthew G. Knepley 45602918637SMatthew G. Knepley test: 45702918637SMatthew G. Knepley suffix: pyramid 458a01d01faSMatthew G. Knepley args: -dm_plex_cell pyramid \ 459a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 460a01d01faSMatthew G. Knepley 461a01d01faSMatthew G. Knepley testset: 462a01d01faSMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_reference_cell_domain -ref_arrangements -dm_plex_check_all \ 463a01d01faSMatthew G. Knepley -ref_dm_view ::ascii_latex -dm_plex_view_tikzscale 1.0 46402918637SMatthew G. Knepley 46502918637SMatthew G. Knepley test: 46602918637SMatthew G. Knepley suffix: ref_segment 467a01d01faSMatthew G. Knepley args: -dm_plex_cell segment \ 468a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 46902918637SMatthew G. Knepley 47002918637SMatthew G. Knepley test: 47102918637SMatthew G. Knepley suffix: ref_triangle 472a01d01faSMatthew G. Knepley args: -dm_plex_cell triangle \ 473a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 47402918637SMatthew G. Knepley 47502918637SMatthew G. Knepley test: 47602918637SMatthew G. Knepley suffix: ref_quadrilateral 477a01d01faSMatthew G. Knepley args: -dm_plex_cell quadrilateral \ 478a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 47902918637SMatthew G. Knepley 48002918637SMatthew G. Knepley test: 48102918637SMatthew G. Knepley suffix: ref_tensor_segment 482a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quad \ 483a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 48402918637SMatthew G. Knepley 48502918637SMatthew G. Knepley test: 48602918637SMatthew G. Knepley suffix: ref_tetrahedron 487a01d01faSMatthew G. Knepley args: -dm_plex_cell tetrahedron \ 488a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 48902918637SMatthew G. Knepley 49002918637SMatthew G. Knepley test: 49102918637SMatthew G. Knepley suffix: ref_hexahedron 492a01d01faSMatthew G. Knepley args: -dm_plex_cell hexahedron \ 493a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 49402918637SMatthew G. Knepley 49502918637SMatthew G. Knepley test: 49602918637SMatthew G. Knepley suffix: ref_triangular_prism 497a01d01faSMatthew G. Knepley args: -dm_plex_cell triangular_prism \ 498a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 49902918637SMatthew G. Knepley 50002918637SMatthew G. Knepley test: 50102918637SMatthew G. Knepley suffix: ref_tensor_triangular_prism 502a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_triangular_prism \ 503a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 50402918637SMatthew G. Knepley 50502918637SMatthew G. Knepley test: 50602918637SMatthew G. Knepley suffix: ref_tensor_quadrilateral_prism 507a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quadrilateral_prism \ 508a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 50902918637SMatthew G. Knepley 51002918637SMatthew G. Knepley # ToBox should recreate the coordinate space since the cell shape changes 51102918637SMatthew G. Knepley testset: 51202918637SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all 51302918637SMatthew G. Knepley 51402918637SMatthew G. Knepley test: 51502918637SMatthew G. Knepley suffix: tobox_triangle 51602918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \ 51702918637SMatthew 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 51802918637SMatthew G. Knepley 51902918637SMatthew G. Knepley test: 52002918637SMatthew G. Knepley suffix: tobox_tensor_segment 52102918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \ 52202918637SMatthew 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 52302918637SMatthew G. Knepley 52402918637SMatthew G. Knepley test: 52502918637SMatthew G. Knepley suffix: tobox_tetrahedron 52602918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \ 52702918637SMatthew 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 52802918637SMatthew G. Knepley 52902918637SMatthew G. Knepley test: 53002918637SMatthew G. Knepley suffix: tobox_triangular_prism 53102918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 53202918637SMatthew 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 53302918637SMatthew G. Knepley 53402918637SMatthew G. Knepley test: 53502918637SMatthew G. Knepley suffix: tobox_tensor_triangular_prism 53602918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \ 53702918637SMatthew 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 53802918637SMatthew G. Knepley 53902918637SMatthew G. Knepley test: 54002918637SMatthew G. Knepley suffix: tobox_tensor_quadrilateral_prism 54102918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \ 54202918637SMatthew 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 54302918637SMatthew G. Knepley 54402918637SMatthew G. Knepley testset: 54502918637SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all 54602918637SMatthew G. Knepley 54702918637SMatthew G. Knepley test: 54802918637SMatthew G. Knepley suffix: alfeld_triangle 54902918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \ 55002918637SMatthew 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 55102918637SMatthew G. Knepley 55202918637SMatthew G. Knepley test: 55302918637SMatthew G. Knepley suffix: alfeld_tetrahedron 55402918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \ 55502918637SMatthew 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 55602918637SMatthew G. Knepley 55702918637SMatthew G. Knepley testset: 55802918637SMatthew G. Knepley args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all 55902918637SMatthew G. Knepley 56002918637SMatthew G. Knepley # This splits edge 1 of the triangle, and reflects about the added edge 56102918637SMatthew G. Knepley test: 56202918637SMatthew G. Knepley suffix: sbr_triangle_0 56302918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \ 56402918637SMatthew 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 56502918637SMatthew G. Knepley 56602918637SMatthew G. Knepley # This splits edge 0 of the triangle, and reflects about the added edge 56702918637SMatthew G. Knepley test: 56802918637SMatthew G. Knepley suffix: sbr_triangle_1 56902918637SMatthew 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 \ 57002918637SMatthew 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 57102918637SMatthew G. Knepley 57202918637SMatthew G. Knepley # This splits edge 2 of the triangle, and reflects about the added edge 57302918637SMatthew G. Knepley test: 57402918637SMatthew G. Knepley suffix: sbr_triangle_2 57502918637SMatthew 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 \ 57602918637SMatthew 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 57702918637SMatthew G. Knepley 57802918637SMatthew G. Knepley # This splits edges 1 and 2 of the triangle 57902918637SMatthew G. Knepley test: 58002918637SMatthew G. Knepley suffix: sbr_triangle_3 58102918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 0 \ 58202918637SMatthew 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 58302918637SMatthew G. Knepley 58402918637SMatthew G. Knepley # This splits edges 1 and 0 of the triangle 58502918637SMatthew G. Knepley test: 58602918637SMatthew G. Knepley suffix: sbr_triangle_4 58702918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 0 \ 58802918637SMatthew 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 58902918637SMatthew G. Knepley 59002918637SMatthew G. Knepley # This splits edges 0 and 1 of the triangle 59102918637SMatthew G. Knepley test: 59202918637SMatthew G. Knepley suffix: sbr_triangle_5 59302918637SMatthew 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 \ 59402918637SMatthew 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 59502918637SMatthew G. Knepley 59602918637SMatthew G. Knepley # This splits edges 0 and 2 of the triangle 59702918637SMatthew G. Knepley test: 59802918637SMatthew G. Knepley suffix: sbr_triangle_6 59902918637SMatthew 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 \ 60002918637SMatthew 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 60102918637SMatthew G. Knepley 60202918637SMatthew G. Knepley # This splits edges 2 and 0 of the triangle 60302918637SMatthew G. Knepley test: 60402918637SMatthew G. Knepley suffix: sbr_triangle_7 60502918637SMatthew 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 \ 60602918637SMatthew 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 60702918637SMatthew G. Knepley 60802918637SMatthew G. Knepley # This splits edges 2 and 1 of the triangle 60902918637SMatthew G. Knepley test: 61002918637SMatthew G. Knepley suffix: sbr_triangle_8 61102918637SMatthew 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 \ 61202918637SMatthew 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 61302918637SMatthew G. Knepley 61402918637SMatthew G. Knepley testset: 61502918637SMatthew G. Knepley args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all 61602918637SMatthew G. Knepley 61702918637SMatthew G. Knepley test: 61802918637SMatthew G. Knepley suffix: bl_tensor_segment 61902918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \ 62002918637SMatthew 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 62102918637SMatthew G. Knepley 62202918637SMatthew G. Knepley # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check 62302918637SMatthew G. Knepley test: 62402918637SMatthew G. Knepley suffix: bl_tensor_triangular_prism 62502918637SMatthew G. Knepley requires: TODO 62602918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \ 62702918637SMatthew 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 62802918637SMatthew G. Knepley 62902918637SMatthew G. Knepley TEST*/ 630