102918637SMatthew G. Knepley static char help[] = "Mesh Orientation Tutorial\n\n"; 202918637SMatthew G. Knepley 302918637SMatthew G. Knepley #include <petscdmplex.h> 402918637SMatthew G. Knepley #include <petscdmplextransform.h> 502918637SMatthew G. Knepley 602918637SMatthew G. Knepley typedef struct { 702918637SMatthew G. Knepley PetscBool genArr; /* Generate all possible cell arrangements */ 802918637SMatthew G. Knepley PetscBool refArr; /* Refine all possible cell arrangements */ 902918637SMatthew G. Knepley PetscBool printTable; /* Print the CAyley table */ 1002918637SMatthew G. Knepley PetscInt orntBounds[2]; /* Bounds for the orientation check */ 1102918637SMatthew G. Knepley PetscInt numOrnt; /* Number of specific orientations specified, or -1 for all orientations */ 1202918637SMatthew G. Knepley PetscInt ornts[48]; /* Specific orientations if specified */ 1302918637SMatthew G. Knepley PetscInt initOrnt; /* Initial orientation for starting mesh */ 1402918637SMatthew G. Knepley } AppCtx; 1502918637SMatthew G. Knepley 1602918637SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 1702918637SMatthew G. Knepley { 1802918637SMatthew G. Knepley PetscInt n = 2; 1902918637SMatthew G. Knepley PetscBool flg; 2002918637SMatthew G. Knepley PetscErrorCode ierr; 2102918637SMatthew G. Knepley 2202918637SMatthew G. Knepley PetscFunctionBeginUser; 2302918637SMatthew G. Knepley options->genArr = PETSC_FALSE; 2402918637SMatthew G. Knepley options->refArr = PETSC_FALSE; 2502918637SMatthew G. Knepley options->printTable = PETSC_FALSE; 2602918637SMatthew G. Knepley options->orntBounds[0] = PETSC_MIN_INT; 2702918637SMatthew G. Knepley options->orntBounds[1] = PETSC_MAX_INT; 2802918637SMatthew G. Knepley options->numOrnt = -1; 2902918637SMatthew G. Knepley options->initOrnt = 0; 3002918637SMatthew G. Knepley 3102918637SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Mesh Orientation Tutorials Options", "DMPLEX");CHKERRQ(ierr); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-gen_arrangements", "Flag for generating all arrangements of the cell", "ex11.c", options->genArr, &options->genArr, NULL)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-ref_arrangements", "Flag for refining all arrangements of the cell", "ex11.c", options->refArr, &options->refArr, NULL)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-print_table", "Print the Cayley table", "ex11.c", options->printTable, &options->printTable, NULL)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-ornt_bounds", "Bounds for orientation checks", "ex11.c", options->orntBounds, &n, NULL)); 3602918637SMatthew G. Knepley n = 48; 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-ornts", "Specific orientations for checks", "ex11.c", options->ornts, &n, &flg)); 3802918637SMatthew G. Knepley if (flg) { 3902918637SMatthew G. Knepley options->numOrnt = n; 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortInt(n, options->ornts)); 4102918637SMatthew G. Knepley } 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-init_ornt", "Initial orientation for starting mesh", "ex11.c", options->initOrnt, &options->initOrnt, NULL)); 4302918637SMatthew G. Knepley ierr = PetscOptionsEnd(); 4402918637SMatthew G. Knepley PetscFunctionReturn(0); 4502918637SMatthew G. Knepley } 4602918637SMatthew G. Knepley 4702918637SMatthew G. Knepley static PetscBool ignoreOrnt(AppCtx *user, PetscInt o) 4802918637SMatthew G. Knepley { 4902918637SMatthew G. Knepley PetscInt loc; 5002918637SMatthew G. Knepley PetscErrorCode ierr; 5102918637SMatthew G. Knepley 5202918637SMatthew G. Knepley if (user->numOrnt < 0) return PETSC_FALSE; 5302918637SMatthew G. Knepley ierr = PetscFindInt(o, user->numOrnt, user->ornts, &loc); 5402918637SMatthew G. Knepley if (loc < 0 || ierr) return PETSC_TRUE; 5502918637SMatthew G. Knepley return PETSC_FALSE; 5602918637SMatthew G. Knepley } 5702918637SMatthew G. Knepley 5802918637SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 5902918637SMatthew G. Knepley { 6002918637SMatthew G. Knepley PetscFunctionBeginUser; 615f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 6502918637SMatthew G. Knepley PetscFunctionReturn(0); 6602918637SMatthew G. Knepley } 6702918637SMatthew G. Knepley 6802918637SMatthew G. Knepley static PetscErrorCode CheckCellVertices(DM dm, PetscInt cell, PetscInt o) 6902918637SMatthew G. Knepley { 7002918637SMatthew G. Knepley DMPolytopeType ct; 7102918637SMatthew G. Knepley const PetscInt *arrVerts; 7202918637SMatthew G. Knepley PetscInt *closure = NULL; 7302918637SMatthew G. Knepley PetscInt Ncl, cl, Nv, vStart, vEnd, v; 7402918637SMatthew G. Knepley MPI_Comm comm; 7502918637SMatthew G. Knepley 7602918637SMatthew G. Knepley PetscFunctionBeginUser; 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cell, &ct)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 8102918637SMatthew G. Knepley for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) { 8202918637SMatthew G. Knepley const PetscInt vertex = closure[cl]; 8302918637SMatthew G. Knepley 8402918637SMatthew G. Knepley if (vertex < vStart || vertex >= vEnd) continue; 8502918637SMatthew G. Knepley closure[Nv++] = vertex; 8602918637SMatthew G. Knepley } 872c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nv != DMPolytopeTypeGetNumVertices(ct),comm, PETSC_ERR_ARG_WRONG, "Cell %D has %D vertices != %D vertices in a %s", cell, Nv, DMPolytopeTypeGetNumVertices(ct), DMPolytopeTypes[ct]); 8802918637SMatthew G. Knepley arrVerts = DMPolytopeTypeGetVertexArrangment(ct, o); 8902918637SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 902c71b3e2SJacob Faibussowitsch PetscCheckFalse(closure[v] != arrVerts[v]+vStart,comm, PETSC_ERR_ARG_WRONG, "Cell %D vertex[%D]: %D should be %D for arrangement %D", cell, v, closure[v], arrVerts[v]+vStart, o); 9102918637SMatthew G. Knepley } 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 9302918637SMatthew G. Knepley PetscFunctionReturn(0); 9402918637SMatthew G. Knepley } 9502918637SMatthew G. Knepley 9602918637SMatthew G. Knepley /* Transform cell with group operation o */ 9702918637SMatthew G. Knepley static PetscErrorCode ReorientCell(DM dm, PetscInt cell, PetscInt o, PetscBool swapCoords) 9802918637SMatthew G. Knepley { 9902918637SMatthew G. Knepley DM cdm; 10002918637SMatthew G. Knepley Vec coordinates; 10102918637SMatthew G. Knepley PetscScalar *coords, *ccoords = NULL; 10202918637SMatthew G. Knepley PetscInt *closure = NULL; 10302918637SMatthew G. Knepley PetscInt cdim, d, Nc, Ncl, cl, vStart, vEnd, Nv; 10402918637SMatthew G. Knepley 10502918637SMatthew G. Knepley PetscFunctionBegin; 10602918637SMatthew G. Knepley /* Change vertex coordinates so that it plots as we expect */ 1075f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &cdim)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords)); 11102918637SMatthew G. Knepley /* Reorient cone */ 1125f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexOrientPoint(dm, cell, o)); 11302918637SMatthew G. Knepley /* Finish resetting coordinates */ 11402918637SMatthew G. Knepley if (swapCoords) { 1155f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(coordinates, &coords)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 11802918637SMatthew G. Knepley for (cl = 0, Nv = 0; cl < Ncl*2; cl += 2) { 11902918637SMatthew G. Knepley const PetscInt vertex = closure[cl]; 12002918637SMatthew G. Knepley PetscScalar *vcoords; 12102918637SMatthew G. Knepley 12202918637SMatthew G. Knepley if (vertex < vStart || vertex >= vEnd) continue; 1235f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRef(cdm, vertex, coords, &vcoords)); 12402918637SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = ccoords[Nv*cdim + d]; 12502918637SMatthew G. Knepley ++Nv; 12602918637SMatthew G. Knepley } 1275f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &Ncl, &closure)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(coordinates, &coords)); 12902918637SMatthew G. Knepley } 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(cdm, NULL, coordinates, cell, &Nc, &ccoords)); 13102918637SMatthew G. Knepley PetscFunctionReturn(0); 13202918637SMatthew G. Knepley } 13302918637SMatthew G. Knepley 13402918637SMatthew G. Knepley static PetscErrorCode GenerateArrangments(DM dm, AppCtx *user) 13502918637SMatthew G. Knepley { 13602918637SMatthew G. Knepley DM odm; 13702918637SMatthew G. Knepley DMPolytopeType ct; 13802918637SMatthew G. Knepley PetscInt No, o; 13902918637SMatthew G. Knepley const char *name; 14002918637SMatthew G. Knepley 14102918637SMatthew G. Knepley PetscFunctionBeginUser; 14202918637SMatthew G. Knepley if (!user->genArr) PetscFunctionReturn(0); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) dm, &name)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, 0, &ct)); 14502918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 14602918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 14702918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 1485f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(ReorientCell(odm, 0, o, PETSC_TRUE)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(odm, NULL, "-gen_dm_view")); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(CheckCellVertices(odm, 0, o)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&odm)); 15402918637SMatthew G. Knepley } 15502918637SMatthew G. Knepley PetscFunctionReturn(0); 15602918637SMatthew G. Knepley } 15702918637SMatthew G. Knepley 15802918637SMatthew G. Knepley static PetscErrorCode VerifyCayleyTable(DM dm, AppCtx *user) 15902918637SMatthew G. Knepley { 16002918637SMatthew G. Knepley DM dm1, dm2; 16102918637SMatthew G. Knepley DMPolytopeType ct; 16202918637SMatthew G. Knepley const PetscInt *refcone, *cone; 16302918637SMatthew G. Knepley PetscInt No, o1, o2, o3, o4; 16402918637SMatthew G. Knepley PetscBool equal; 16502918637SMatthew G. Knepley const char *name; 16602918637SMatthew G. Knepley 16702918637SMatthew G. Knepley PetscFunctionBeginUser; 16802918637SMatthew G. Knepley if (!user->genArr) PetscFunctionReturn(0); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) dm, &name)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, 0, &ct)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, 0, &refcone)); 17202918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 1735f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Cayley Table for %s\n", DMPolytopeTypes[ct])); 17402918637SMatthew G. Knepley for (o1 = PetscMax(-No, user->orntBounds[0]); o1 < PetscMin(No, user->orntBounds[1]); ++o1) { 17502918637SMatthew G. Knepley for (o2 = PetscMax(-No, user->orntBounds[0]); o2 < PetscMin(No, user->orntBounds[1]); ++o2) { 1765f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexOrientPoint(dm1, 0, o2)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCheckFaces(dm1, 0)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexOrientPoint(dm1, 0, o1)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCheckFaces(dm1, 0)); 1815f80ce2aSJacob Faibussowitsch o3 = DMPolytopeTypeComposeOrientation(ct, o1, o2); 18202918637SMatthew G. Knepley /* First verification */ 1835f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexOrientPoint(dm2, 0, o3)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCheckFaces(dm2, 0)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexEqual(dm1, dm2, &equal)); 18702918637SMatthew G. Knepley if (!equal) { 1885f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm1, NULL, "-error_dm_view")); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm2, NULL, "-error_dm_view")); 19098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D != %D", DMPolytopeTypes[ct], o1, o2, o3); 19102918637SMatthew G. Knepley } 19202918637SMatthew G. Knepley /* Second verification */ 1935f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm1, 0, &cone)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMPolytopeGetOrientation(ct, refcone, cone, &o4)); 1955f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%D, ", o4)); 1962c71b3e2SJacob Faibussowitsch PetscCheckFalse(o3 != o4,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cayley table error for %s: %D * %D = %D != %D", DMPolytopeTypes[ct], o1, o2, o3, o4); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm1)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm2)); 19902918637SMatthew G. Knepley } 2005f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n")); 20102918637SMatthew G. Knepley } 20202918637SMatthew G. Knepley PetscFunctionReturn(0); 20302918637SMatthew G. Knepley } 20402918637SMatthew G. Knepley 20502918637SMatthew G. Knepley static PetscErrorCode VerifyInverse(DM dm, AppCtx *user) 20602918637SMatthew G. Knepley { 20702918637SMatthew G. Knepley DM dm1, dm2; 20802918637SMatthew G. Knepley DMPolytopeType ct; 20902918637SMatthew G. Knepley const PetscInt *refcone, *cone; 21002918637SMatthew G. Knepley PetscInt No, o, oi, o2; 21102918637SMatthew G. Knepley PetscBool equal; 21202918637SMatthew G. Knepley const char *name; 21302918637SMatthew G. Knepley 21402918637SMatthew G. Knepley PetscFunctionBeginUser; 21502918637SMatthew G. Knepley if (!user->genArr) PetscFunctionReturn(0); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) dm, &name)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, 0, &ct)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, 0, &refcone)); 21902918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 2205f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Inverse table for %s\n", DMPolytopeTypes[ct])); 22102918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 22202918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 2235f80ce2aSJacob Faibussowitsch oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm1)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexOrientPoint(dm1, 0, o)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCheckFaces(dm1, 0)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexOrientPoint(dm1, 0, oi)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCheckFaces(dm1, 0)); 22902918637SMatthew G. Knepley /* First verification */ 2305f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PetscObjectComm((PetscObject) dm), user, &dm2)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexEqual(dm1, dm2, &equal)); 23202918637SMatthew G. Knepley if (!equal) { 2335f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm1, NULL, "-error_dm_view")); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm2, NULL, "-error_dm_view")); 23598921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D != 0", DMPolytopeTypes[ct], o, oi); 23602918637SMatthew G. Knepley } 23702918637SMatthew G. Knepley /* Second verification */ 2385f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm1, 0, &cone)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPolytopeGetOrientation(ct, refcone, cone, &o2)); 2405f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%D, ", oi)); 2412c71b3e2SJacob Faibussowitsch PetscCheckFalse(o2 != 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inverse error for %s: %D * %D = %D != 0", DMPolytopeTypes[ct], o, oi, o2); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm1)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm2)); 24402918637SMatthew G. Knepley } 2455f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n")); 24602918637SMatthew G. Knepley PetscFunctionReturn(0); 24702918637SMatthew G. Knepley } 24802918637SMatthew G. Knepley 24902918637SMatthew G. Knepley /* Suppose that point p has the same arrangement as o from canonical, compare the subcells to canonical subcells */ 25002918637SMatthew G. Knepley static PetscErrorCode CheckSubcells(DM dm, DM odm, PetscInt p, PetscInt o, AppCtx *user) 25102918637SMatthew G. Knepley { 25202918637SMatthew G. Knepley DMPlexTransform tr, otr; 25302918637SMatthew G. Knepley DMPolytopeType ct; 25402918637SMatthew G. Knepley DMPolytopeType *rct; 25502918637SMatthew G. Knepley const PetscInt *cone, *ornt, *ocone, *oornt; 25602918637SMatthew G. Knepley PetscInt *rsize, *rcone, *rornt; 25702918637SMatthew G. Knepley PetscInt Nct, n, oi, debug = 0; 25802918637SMatthew G. Knepley 25902918637SMatthew G. Knepley PetscFunctionBeginUser; 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetDM(tr, dm)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetFromOptions(tr)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetUp(tr)); 26402918637SMatthew G. Knepley 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformCreate(PetscObjectComm((PetscObject) odm), &otr)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetDM(otr, odm)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetFromOptions(otr)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetUp(otr)); 26902918637SMatthew G. Knepley 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, p, &ct)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, p, &cone)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeOrientation(dm, p, &ornt)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(odm, p, &ocone)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeOrientation(odm, p, &oornt)); 2755f80ce2aSJacob Faibussowitsch oi = DMPolytopeTypeComposeOrientationInv(ct, 0, o); 2765f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Orientation %D\n", oi)); 27702918637SMatthew G. Knepley 2785f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt)); 27902918637SMatthew G. Knepley for (n = 0; n < Nct; ++n) { 28002918637SMatthew G. Knepley DMPolytopeType ctNew = rct[n]; 28102918637SMatthew G. Knepley PetscInt r, ro; 28202918637SMatthew G. Knepley 2835f80ce2aSJacob Faibussowitsch if (debug) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Checking type %s\n", DMPolytopeTypes[ctNew])); 28402918637SMatthew G. Knepley for (r = 0; r < rsize[n]; ++r) { 28502918637SMatthew G. Knepley const PetscInt *qcone, *qornt, *oqcone, *oqornt; 28602918637SMatthew G. Knepley PetscInt pNew, opNew, oo, pr, fo; 28702918637SMatthew G. Knepley PetscBool restore = PETSC_TRUE; 28802918637SMatthew G. Knepley 2895f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformGetTargetPoint(tr, ct, ctNew, p, r, &pNew)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformGetCone(tr, pNew, &qcone, &qornt)); 29102918637SMatthew G. Knepley if (debug) { 29202918637SMatthew G. Knepley PetscInt c; 29302918637SMatthew G. Knepley 2945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Checking replica %D (%D)\n Original Cone", r, pNew)); 2955f80ce2aSJacob Faibussowitsch for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " %D (%D)", qcone[c], qornt[c])); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n")); 29702918637SMatthew G. Knepley } 29802918637SMatthew G. Knepley for (ro = 0; ro < rsize[n]; ++ro) { 29902918637SMatthew G. Knepley PetscBool found; 30002918637SMatthew G. Knepley 3015f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformGetTargetPoint(otr, ct, ctNew, p, ro, &opNew)); 3025f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformGetConeOriented(otr, opNew, o, &oqcone, &oqornt)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(DMPolytopeMatchOrientation(ctNew, oqcone, qcone, &oo, &found)); 30402918637SMatthew G. Knepley if (found) break; 3055f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt)); 30602918637SMatthew G. Knepley } 30702918637SMatthew G. Knepley if (debug) { 30802918637SMatthew G. Knepley PetscInt c; 30902918637SMatthew G. Knepley 3105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Checking transform replica %D (%D) (%D)\n Transform Cone", ro, opNew, o)); 3115f80ce2aSJacob Faibussowitsch for (c = 0; c < DMPolytopeTypeGetConeSize(ctNew); ++c) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " %D (%D)", oqcone[c], oqornt[c])); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n")); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Matched %D\n", oo)); 31402918637SMatthew G. Knepley } 31502918637SMatthew G. Knepley if (ro == rsize[n]) { 31602918637SMatthew G. Knepley /* The tetrahedron has 3 pairs of opposing edges, and any pair can be connected by the interior segment */ 31702918637SMatthew G. Knepley if (ct == DM_POLYTOPE_TETRAHEDRON) { 31802918637SMatthew G. Knepley /* The segment in a tetrahedron does not map into itself under the group action */ 31902918637SMatthew G. Knepley if (ctNew == DM_POLYTOPE_SEGMENT) {restore = PETSC_FALSE; ro = r; oo = 0;} 32002918637SMatthew G. Knepley /* The last four interior faces do not map into themselves under the group action */ 32102918637SMatthew G. Knepley if (r > 3 && ctNew == DM_POLYTOPE_TRIANGLE) {restore = PETSC_FALSE; ro = r; oo = 0;} 32202918637SMatthew G. Knepley /* The last four interior faces do not map into themselves under the group action */ 32302918637SMatthew G. Knepley if (r > 3 && ctNew == DM_POLYTOPE_TETRAHEDRON) {restore = PETSC_FALSE; ro = r; oo = 0;} 32402918637SMatthew G. Knepley } 32528b400f6SJacob Faibussowitsch PetscCheck(!restore,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to find matching %s %D orientation for cell orientation %D", DMPolytopeTypes[ctNew], r, o); 32602918637SMatthew G. Knepley } 3275f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%D, %D, ", ro, oo)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformGetSubcellOrientation(tr, ct, p, oi, ctNew, r, 0, &pr, &fo)); 32902918637SMatthew G. Knepley if (!user->printTable) { 3302c71b3e2SJacob Faibussowitsch PetscCheckFalse(pr != ro,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong replica %D != %D", pr, ro); 3312c71b3e2SJacob Faibussowitsch PetscCheckFalse(fo != oo,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Choose wrong orientation %D != %D", fo, oo); 33202918637SMatthew G. Knepley } 3335f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformRestoreCone(tr, pNew, &qcone, &qornt)); 3345f80ce2aSJacob Faibussowitsch if (restore) CHKERRQ(DMPlexTransformRestoreCone(otr, pNew, &oqcone, &oqornt)); 33502918637SMatthew G. Knepley } 3365f80ce2aSJacob Faibussowitsch if (user->printTable) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n")); 33702918637SMatthew G. Knepley } 3385f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformDestroy(&tr)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformDestroy(&otr)); 34002918637SMatthew G. Knepley PetscFunctionReturn(0); 34102918637SMatthew G. Knepley } 34202918637SMatthew G. Knepley 34302918637SMatthew G. Knepley static PetscErrorCode RefineArrangments(DM dm, AppCtx *user) 34402918637SMatthew G. Knepley { 34502918637SMatthew G. Knepley DM odm, rdm; 34602918637SMatthew G. Knepley DMPolytopeType ct; 34702918637SMatthew G. Knepley PetscInt No, o; 34802918637SMatthew G. Knepley const char *name; 34902918637SMatthew G. Knepley 35002918637SMatthew G. Knepley PetscFunctionBeginUser; 35102918637SMatthew G. Knepley if (!user->refArr) PetscFunctionReturn(0); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) dm, &name)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, 0, &ct)); 35402918637SMatthew G. Knepley No = DMPolytopeTypeGetNumArrangments(ct)/2; 35502918637SMatthew G. Knepley for (o = PetscMax(-No, user->orntBounds[0]); o < PetscMin(No, user->orntBounds[1]); ++o) { 35602918637SMatthew G. Knepley if (ignoreOrnt(user, o)) continue; 3575f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PetscObjectComm((PetscObject) dm), user, &odm)); 3585f80ce2aSJacob Faibussowitsch if (user->initOrnt) CHKERRQ(ReorientCell(odm, 0, user->initOrnt, PETSC_FALSE)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(ReorientCell(odm, 0, o, PETSC_TRUE)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(odm, NULL, "-orig_dm_view")); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(odm, MPI_COMM_NULL, &rdm)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(rdm)); 3635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "%s orientation %D\n", name, o)); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(rdm, NULL, "-ref_dm_view")); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(CheckSubcells(dm, odm, 0, o, user)); 3665f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&odm)); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&rdm)); 36802918637SMatthew G. Knepley } 36902918637SMatthew G. Knepley PetscFunctionReturn(0); 37002918637SMatthew G. Knepley } 37102918637SMatthew G. Knepley 37202918637SMatthew G. Knepley int main(int argc, char **argv) 37302918637SMatthew G. Knepley { 37402918637SMatthew G. Knepley DM dm; 37502918637SMatthew G. Knepley AppCtx user; 37602918637SMatthew G. Knepley 377*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 3795f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 38002918637SMatthew G. Knepley if (user.initOrnt) { 3815f80ce2aSJacob Faibussowitsch CHKERRQ(ReorientCell(dm, 0, user.initOrnt, PETSC_FALSE)); 3825f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-ornt_dm_view")); 38302918637SMatthew G. Knepley } 3845f80ce2aSJacob Faibussowitsch CHKERRQ(GenerateArrangments(dm, &user)); 3855f80ce2aSJacob Faibussowitsch CHKERRQ(VerifyCayleyTable(dm, &user)); 3865f80ce2aSJacob Faibussowitsch CHKERRQ(VerifyInverse(dm, &user)); 3875f80ce2aSJacob Faibussowitsch CHKERRQ(RefineArrangments(dm, &user)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 389*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 390*b122ec5aSJacob Faibussowitsch return 0; 39102918637SMatthew G. Knepley } 39202918637SMatthew G. Knepley 39302918637SMatthew G. Knepley /*TEST 39402918637SMatthew G. Knepley 395a01d01faSMatthew G. Knepley testset: 396a01d01faSMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_reference_cell_domain -gen_arrangements \ 397a01d01faSMatthew G. Knepley -gen_dm_view ::ascii_latex -dm_plex_view_tikzscale 0.5 398a01d01faSMatthew G. Knepley 39902918637SMatthew G. Knepley test: 40002918637SMatthew G. Knepley suffix: segment 401a01d01faSMatthew G. Knepley args: -dm_plex_cell segment \ 402a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 40302918637SMatthew G. Knepley 40402918637SMatthew G. Knepley test: 40502918637SMatthew G. Knepley suffix: triangle 406a01d01faSMatthew G. Knepley args: -dm_plex_cell triangle \ 407a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 40802918637SMatthew G. Knepley 40902918637SMatthew G. Knepley test: 41002918637SMatthew G. Knepley suffix: quadrilateral 411a01d01faSMatthew G. Knepley args: -dm_plex_cell quadrilateral \ 412a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 41302918637SMatthew G. Knepley 41402918637SMatthew G. Knepley test: 41502918637SMatthew G. Knepley suffix: tensor_segment 416a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quad \ 417a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 41802918637SMatthew G. Knepley 41902918637SMatthew G. Knepley test: 42002918637SMatthew G. Knepley suffix: tetrahedron 421a01d01faSMatthew G. Knepley args: -dm_plex_cell tetrahedron \ 422a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 42302918637SMatthew G. Knepley 42402918637SMatthew G. Knepley test: 42502918637SMatthew G. Knepley suffix: hexahedron 426a01d01faSMatthew G. Knepley args: -dm_plex_cell hexahedron \ 427a01d01faSMatthew 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 42802918637SMatthew G. Knepley 42902918637SMatthew G. Knepley test: 43002918637SMatthew G. Knepley suffix: triangular_prism 431a01d01faSMatthew G. Knepley args: -dm_plex_cell triangular_prism \ 432a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 43302918637SMatthew G. Knepley 43402918637SMatthew G. Knepley test: 43502918637SMatthew G. Knepley suffix: tensor_triangular_prism 436a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_triangular_prism \ 437a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 43802918637SMatthew G. Knepley 43902918637SMatthew G. Knepley test: 44002918637SMatthew G. Knepley suffix: tensor_quadrilateral_prism 441a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quadrilateral_prism \ 442a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 44302918637SMatthew G. Knepley 44402918637SMatthew G. Knepley test: 44502918637SMatthew G. Knepley suffix: pyramid 446a01d01faSMatthew G. Knepley args: -dm_plex_cell pyramid \ 447a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 448a01d01faSMatthew G. Knepley 449a01d01faSMatthew G. Knepley testset: 450a01d01faSMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_reference_cell_domain -ref_arrangements -dm_plex_check_all \ 451a01d01faSMatthew G. Knepley -ref_dm_view ::ascii_latex -dm_plex_view_tikzscale 1.0 45202918637SMatthew G. Knepley 45302918637SMatthew G. Knepley test: 45402918637SMatthew G. Knepley suffix: ref_segment 455a01d01faSMatthew G. Knepley args: -dm_plex_cell segment \ 456a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0 -dm_plex_view_colors_depth 1,0 45702918637SMatthew G. Knepley 45802918637SMatthew G. Knepley test: 45902918637SMatthew G. Knepley suffix: ref_triangle 460a01d01faSMatthew G. Knepley args: -dm_plex_cell triangle \ 461a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 46202918637SMatthew G. Knepley 46302918637SMatthew G. Knepley test: 46402918637SMatthew G. Knepley suffix: ref_quadrilateral 465a01d01faSMatthew G. Knepley args: -dm_plex_cell quadrilateral \ 466a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 46702918637SMatthew G. Knepley 46802918637SMatthew G. Knepley test: 46902918637SMatthew G. Knepley suffix: ref_tensor_segment 470a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quad \ 471a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,1,0 -dm_plex_view_colors_depth 1,1,0 47202918637SMatthew G. Knepley 47302918637SMatthew G. Knepley test: 47402918637SMatthew G. Knepley suffix: ref_tetrahedron 475a01d01faSMatthew G. Knepley args: -dm_plex_cell tetrahedron \ 476a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 47702918637SMatthew G. Knepley 47802918637SMatthew G. Knepley test: 47902918637SMatthew G. Knepley suffix: ref_hexahedron 480a01d01faSMatthew G. Knepley args: -dm_plex_cell hexahedron \ 481a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 48202918637SMatthew G. Knepley 48302918637SMatthew G. Knepley test: 48402918637SMatthew G. Knepley suffix: ref_triangular_prism 485a01d01faSMatthew G. Knepley args: -dm_plex_cell triangular_prism \ 486a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 48702918637SMatthew G. Knepley 48802918637SMatthew G. Knepley test: 48902918637SMatthew G. Knepley suffix: ref_tensor_triangular_prism 490a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_triangular_prism \ 491a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 49202918637SMatthew G. Knepley 49302918637SMatthew G. Knepley test: 49402918637SMatthew G. Knepley suffix: ref_tensor_quadrilateral_prism 495a01d01faSMatthew G. Knepley args: -dm_plex_cell tensor_quadrilateral_prism \ 496a01d01faSMatthew G. Knepley -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 49702918637SMatthew G. Knepley 49802918637SMatthew G. Knepley # ToBox should recreate the coordinate space since the cell shape changes 49902918637SMatthew G. Knepley testset: 50002918637SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -ref_arrangements -dm_plex_check_all 50102918637SMatthew G. Knepley 50202918637SMatthew G. Knepley test: 50302918637SMatthew G. Knepley suffix: tobox_triangle 50402918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \ 50502918637SMatthew 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 50602918637SMatthew G. Knepley 50702918637SMatthew G. Knepley test: 50802918637SMatthew G. Knepley suffix: tobox_tensor_segment 50902918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \ 51002918637SMatthew 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 51102918637SMatthew G. Knepley 51202918637SMatthew G. Knepley test: 51302918637SMatthew G. Knepley suffix: tobox_tetrahedron 51402918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \ 51502918637SMatthew 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 51602918637SMatthew G. Knepley 51702918637SMatthew G. Knepley test: 51802918637SMatthew G. Knepley suffix: tobox_triangular_prism 51902918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 52002918637SMatthew 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 52102918637SMatthew G. Knepley 52202918637SMatthew G. Knepley test: 52302918637SMatthew G. Knepley suffix: tobox_tensor_triangular_prism 52402918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \ 52502918637SMatthew G. Knepley -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0 52602918637SMatthew G. Knepley 52702918637SMatthew G. Knepley test: 52802918637SMatthew G. Knepley suffix: tobox_tensor_quadrilateral_prism 52902918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism \ 53002918637SMatthew G. Knepley -ref_dm_view ::ascii_latex -dm_plex_view_numbers_depth 1,0,0,0 -dm_plex_view_colors_depth 1,0,0,0 -dm_plex_view_tikzscale 1.0 53102918637SMatthew G. Knepley 53202918637SMatthew G. Knepley testset: 53302918637SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_transform_type refine_alfeld -ref_arrangements -dm_plex_check_all 53402918637SMatthew G. Knepley 53502918637SMatthew G. Knepley test: 53602918637SMatthew G. Knepley suffix: alfeld_triangle 53702918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle \ 53802918637SMatthew 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 53902918637SMatthew G. Knepley 54002918637SMatthew G. Knepley test: 54102918637SMatthew G. Knepley suffix: alfeld_tetrahedron 54202918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron \ 54302918637SMatthew 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 54402918637SMatthew G. Knepley 54502918637SMatthew G. Knepley testset: 54602918637SMatthew G. Knepley args: -dm_plex_transform_type refine_sbr -ref_arrangements -dm_plex_check_all 54702918637SMatthew G. Knepley 54802918637SMatthew G. Knepley # This splits edge 1 of the triangle, and reflects about the added edge 54902918637SMatthew G. Knepley test: 55002918637SMatthew G. Knepley suffix: sbr_triangle_0 55102918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5 -ornts -2,0 \ 55202918637SMatthew 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 55302918637SMatthew G. Knepley 55402918637SMatthew G. Knepley # This splits edge 0 of the triangle, and reflects about the added edge 55502918637SMatthew G. Knepley test: 55602918637SMatthew G. Knepley suffix: sbr_triangle_1 55702918637SMatthew 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 \ 55802918637SMatthew 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 55902918637SMatthew G. Knepley 56002918637SMatthew G. Knepley # This splits edge 2 of the triangle, and reflects about the added edge 56102918637SMatthew G. Knepley test: 56202918637SMatthew G. Knepley suffix: sbr_triangle_2 56302918637SMatthew 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 \ 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 edges 1 and 2 of the triangle 56702918637SMatthew G. Knepley test: 56802918637SMatthew G. Knepley suffix: sbr_triangle_3 56902918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 5,6 -ornts 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 edges 1 and 0 of the triangle 57302918637SMatthew G. Knepley test: 57402918637SMatthew G. Knepley suffix: sbr_triangle_4 57502918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_transform_sbr_ref_cell 4,5 -ornts 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 0 and 1 of the triangle 57902918637SMatthew G. Knepley test: 58002918637SMatthew G. Knepley suffix: sbr_triangle_5 58102918637SMatthew 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 \ 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 0 and 2 of the triangle 58502918637SMatthew G. Knepley test: 58602918637SMatthew G. Knepley suffix: sbr_triangle_6 58702918637SMatthew 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 \ 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 2 and 0 of the triangle 59102918637SMatthew G. Knepley test: 59202918637SMatthew G. Knepley suffix: sbr_triangle_7 59302918637SMatthew 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 \ 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 2 and 1 of the triangle 59702918637SMatthew G. Knepley test: 59802918637SMatthew G. Knepley suffix: sbr_triangle_8 59902918637SMatthew 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 \ 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 testset: 60302918637SMatthew G. Knepley args: -dm_plex_transform_type refine_boundary_layer -dm_plex_transform_bl_splits 2 -ref_arrangements -dm_plex_check_all 60402918637SMatthew G. Knepley 60502918637SMatthew G. Knepley test: 60602918637SMatthew G. Knepley suffix: bl_tensor_segment 60702918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad \ 60802918637SMatthew 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 60902918637SMatthew G. Knepley 61002918637SMatthew G. Knepley # The subcell check is broken because at orientation 3, the internal triangles do not get properly permuted for the check 61102918637SMatthew G. Knepley test: 61202918637SMatthew G. Knepley suffix: bl_tensor_triangular_prism 61302918637SMatthew G. Knepley requires: TODO 61402918637SMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism \ 61502918637SMatthew 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 61602918637SMatthew G. Knepley 61702918637SMatthew G. Knepley TEST*/ 618