154fcfd0cSMatthew G. Knepley static char help[] = "Tests for refinement of meshes created by hand\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown PetscInt debug; /* The debugging level */ 7c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 8c4762a1bSJed Brown PetscBool cellHybrid; /* Use a hybrid mesh */ 9c4762a1bSJed Brown PetscBool cellSimplex; /* Use simplices or hexes */ 10c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 11c4762a1bSJed Brown PetscInt testNum; /* The particular mesh to test */ 12c4762a1bSJed Brown PetscBool uninterpolate; /* Uninterpolate the mesh at the end */ 13c4762a1bSJed Brown PetscBool reinterpolate; /* Reinterpolate the mesh at the end */ 14c4762a1bSJed Brown } AppCtx; 15c4762a1bSJed Brown 16*9371c9d4SSatish Balay PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 17c4762a1bSJed Brown PetscFunctionBegin; 18c4762a1bSJed Brown options->debug = 0; 19c4762a1bSJed Brown options->dim = 2; 20c4762a1bSJed Brown options->cellHybrid = PETSC_TRUE; 21c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 22c4762a1bSJed Brown options->testPartition = PETSC_TRUE; 23c4762a1bSJed Brown options->testNum = 0; 24c4762a1bSJed Brown options->uninterpolate = PETSC_FALSE; 25c4762a1bSJed Brown options->reinterpolate = PETSC_FALSE; 26c4762a1bSJed Brown 27d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex4.c", options->debug, &options->debug, NULL, 0)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL, 1, 3)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cell_hybrid", "Use a hybrid mesh", "ex4.c", options->cellHybrid, &options->cellHybrid, NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex4.c", options->cellSimplex, &options->cellSimplex, NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex4.c", options->testPartition, &options->testPartition, NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex4.c", options->testNum, &options->testNum, NULL, 0)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex4.c", options->uninterpolate, &options->uninterpolate, NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-reinterpolate", "Reinterpolate the mesh at the end", "ex4.c", options->reinterpolate, &options->reinterpolate, NULL)); 36d0609cedSBarry Smith PetscOptionsEnd(); 37c4762a1bSJed Brown PetscFunctionReturn(0); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Two segments 41c4762a1bSJed Brown 42c4762a1bSJed Brown 2-------0-------3-------1-------4 43c4762a1bSJed Brown 44c4762a1bSJed Brown become 45c4762a1bSJed Brown 46c4762a1bSJed Brown 4---0---7---1---5---2---8---3---6 47c4762a1bSJed Brown 48c4762a1bSJed Brown */ 49*9371c9d4SSatish Balay PetscErrorCode CreateSimplex_1D(MPI_Comm comm, DM *dm) { 50c4762a1bSJed Brown PetscInt depth = 1; 51c4762a1bSJed Brown PetscMPIInt rank; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 55dd400576SPatrick Sanan if (rank == 0) { 56c4762a1bSJed Brown PetscInt numPoints[2] = {3, 2}; 57c4762a1bSJed Brown PetscInt coneSize[5] = {2, 2, 0, 0, 0}; 58c4762a1bSJed Brown PetscInt cones[4] = {2, 3, 3, 4}; 59c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0}; 60c4762a1bSJed Brown PetscScalar vertexCoords[3] = {-1.0, 0.0, 1.0}; 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 63c4762a1bSJed Brown } else { 64c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown PetscFunctionReturn(0); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Two triangles 72c4762a1bSJed Brown 4 73c4762a1bSJed Brown / | \ 74c4762a1bSJed Brown 8 | 10 75c4762a1bSJed Brown / | \ 76c4762a1bSJed Brown 2 0 7 1 5 77c4762a1bSJed Brown \ | / 78c4762a1bSJed Brown 6 | 9 79c4762a1bSJed Brown \ | / 80c4762a1bSJed Brown 3 81c4762a1bSJed Brown 82c4762a1bSJed Brown Becomes 83c4762a1bSJed Brown 10 84c4762a1bSJed Brown / | \ 85c4762a1bSJed Brown 21 | 26 86c4762a1bSJed Brown / | \ 87c4762a1bSJed Brown 14 2 20 4 16 88c4762a1bSJed Brown /|\ | /|\ 89c4762a1bSJed Brown 22 | 28 | 32 | 25 90c4762a1bSJed Brown / | \ | / | 6\ 91c4762a1bSJed Brown 8 29 3 13 7 31 11 92c4762a1bSJed Brown \0 | / | \ | / 93c4762a1bSJed Brown 17 | 27 | 30 | 24 94c4762a1bSJed Brown \|/ | \|/ 95c4762a1bSJed Brown 12 1 19 5 15 96c4762a1bSJed Brown \ | / 97c4762a1bSJed Brown 18 | 23 98c4762a1bSJed Brown \ | / 99c4762a1bSJed Brown 9 100c4762a1bSJed Brown */ 101*9371c9d4SSatish Balay PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *dm) { 102c4762a1bSJed Brown PetscInt depth = 2; 103c4762a1bSJed Brown PetscMPIInt rank; 104c4762a1bSJed Brown 105c4762a1bSJed Brown PetscFunctionBegin; 1069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 107dd400576SPatrick Sanan if (rank == 0) { 108c4762a1bSJed Brown PetscInt numPoints[3] = {4, 5, 2}; 109c4762a1bSJed Brown PetscInt coneSize[11] = {3, 3, 0, 0, 0, 0, 2, 2, 2, 2, 2}; 110c4762a1bSJed Brown PetscInt cones[16] = {6, 7, 8, 7, 9, 10, 2, 3, 3, 4, 4, 2, 3, 5, 5, 4}; 111b5a892a1SMatthew G. Knepley PetscInt coneOrientations[16] = {0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 112c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-0.5, 0.0, 0.0, -0.5, 0.0, 0.5, 0.5, 0.0}; 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 115c4762a1bSJed Brown } else { 116c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 117c4762a1bSJed Brown 1189566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown PetscFunctionReturn(0); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* Two triangles separated by a zero-volume cell with 4 vertices/2 edges 124c4762a1bSJed Brown 5--16--8 125c4762a1bSJed Brown / | | \ 126c4762a1bSJed Brown 11 | | 12 127c4762a1bSJed Brown / | | \ 128c4762a1bSJed Brown 3 0 10 2 14 1 6 129c4762a1bSJed Brown \ | | / 130c4762a1bSJed Brown 9 | | 13 131c4762a1bSJed Brown \ | | / 132c4762a1bSJed Brown 4--15--7 133c4762a1bSJed Brown */ 134*9371c9d4SSatish Balay PetscErrorCode CreateSimplexHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) { 135c4762a1bSJed Brown DM idm, hdm = NULL; 136c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 137c4762a1bSJed Brown PetscInt p; 138c4762a1bSJed Brown PetscMPIInt rank; 139c4762a1bSJed Brown 140c4762a1bSJed Brown PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 142dd400576SPatrick Sanan if (rank == 0) { 143c4762a1bSJed Brown switch (testNum) { 144*9371c9d4SSatish Balay case 0: { 145c4762a1bSJed Brown PetscInt numPoints[2] = {4, 2}; 146c4762a1bSJed Brown PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 147c4762a1bSJed Brown PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 148c4762a1bSJed Brown PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 149c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, 1.0, 0.5}; 150c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 1539566063dSJacob Faibussowitsch for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 154*9371c9d4SSatish Balay } break; 155*9371c9d4SSatish Balay case 1: { 156c4762a1bSJed Brown PetscInt numPoints[2] = {5, 4}; 157c4762a1bSJed Brown PetscInt coneSize[9] = {3, 3, 3, 3, 0, 0, 0, 0, 0}; 158c4762a1bSJed Brown PetscInt cones[12] = {4, 5, 6, 6, 7, 4, 6, 5, 8, 6, 8, 7}; 159c4762a1bSJed Brown PetscInt coneOrientations[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 160c4762a1bSJed Brown PetscScalar vertexCoords[10] = {-1.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0}; 161c4762a1bSJed Brown PetscInt faultPoints[3] = {5, 6, 7}; 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 1649566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 165*9371c9d4SSatish Balay } break; 16663a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 167c4762a1bSJed Brown } 1689566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 1709566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 1719566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 1729566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 1739566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault", &faultLabel)); 174caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm)); 1759566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&hybridLabel)); 1769566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 1779566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 178c4762a1bSJed Brown *dm = hdm; 179c4762a1bSJed Brown } else { 180c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 181c4762a1bSJed Brown 1829566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 1839566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 1859566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 1869566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 1879566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 188caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm)); 1899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 1909566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 191c4762a1bSJed Brown *dm = hdm; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* Two quadrilaterals 197c4762a1bSJed Brown 198c4762a1bSJed Brown 5----10-----4----14-----7 199c4762a1bSJed Brown | | | 200c4762a1bSJed Brown | | | 201c4762a1bSJed Brown | | | 202c4762a1bSJed Brown 11 0 9 1 13 203c4762a1bSJed Brown | | | 204c4762a1bSJed Brown | | | 205c4762a1bSJed Brown | | | 206c4762a1bSJed Brown 2-----8-----3----12-----6 207c4762a1bSJed Brown */ 208*9371c9d4SSatish Balay PetscErrorCode CreateTensorProduct_2D(MPI_Comm comm, PetscInt testNum, DM *dm) { 209c4762a1bSJed Brown PetscInt depth = 2; 210c4762a1bSJed Brown PetscMPIInt rank; 211c4762a1bSJed Brown 212c4762a1bSJed Brown PetscFunctionBegin; 2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 214dd400576SPatrick Sanan if (rank == 0) { 215c4762a1bSJed Brown PetscInt numPoints[3] = {6, 7, 2}; 216c4762a1bSJed Brown PetscInt coneSize[15] = {4, 4, 0, 0, 0, 0, 0, 0, 2, 2, 2, 2, 2, 2, 2}; 217c4762a1bSJed Brown PetscInt cones[22] = {8, 9, 10, 11, 12, 13, 14, 9, 2, 3, 3, 4, 4, 5, 5, 2, 3, 6, 6, 7, 7, 4}; 218b5a892a1SMatthew G. Knepley PetscInt coneOrientations[22] = {0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 219c4762a1bSJed Brown PetscScalar vertexCoords[12] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, -1.0, 0.5, 1.0, -0.5, 1.0, 0.5}; 220c4762a1bSJed Brown 2219566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 222c4762a1bSJed Brown } else { 223c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 224c4762a1bSJed Brown 2259566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown PetscFunctionReturn(0); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230*9371c9d4SSatish Balay PetscErrorCode CreateTensorProductHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) { 231c4762a1bSJed Brown DM idm, hdm = NULL; 232c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 233c4762a1bSJed Brown PetscInt p; 234c4762a1bSJed Brown PetscMPIInt rank; 235c4762a1bSJed Brown 236c4762a1bSJed Brown PetscFunctionBegin; 2379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 238dd400576SPatrick Sanan if (rank == 0) { 239c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 240c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 241*9371c9d4SSatish Balay PetscInt cones[8] = { 242*9371c9d4SSatish Balay 2, 3, 4, 5, 3, 6, 7, 4, 243*9371c9d4SSatish Balay }; 244c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 245c4762a1bSJed Brown PetscScalar vertexCoords[12] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, -1.0, 0.5, 1.0, -0.5, 1.0, 0.5}; 246c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 247c4762a1bSJed Brown 2489566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 2499566063dSJacob Faibussowitsch for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 2509566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 2529566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 2539566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 2549566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 2559566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault", &faultLabel)); 256caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm)); 2579566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&hybridLabel)); 258c4762a1bSJed Brown } else { 259c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 2629566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 2649566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 2659566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 2669566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 267caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm)); 268c4762a1bSJed Brown } 2699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 2709566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 271c4762a1bSJed Brown *dm = hdm; 272c4762a1bSJed Brown PetscFunctionReturn(0); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* Two tetrahedrons 276c4762a1bSJed Brown 277c4762a1bSJed Brown cell 5 5______ cell 278c4762a1bSJed Brown 0 / | \ |\ \ 1 279c4762a1bSJed Brown 17 | 18 | 18 13 21 280c4762a1bSJed Brown /8 19 10\ 19 \ \ 281c4762a1bSJed Brown 2-14-|----4 | 4--22--6 282c4762a1bSJed Brown \ 9 | 7 / |10 / / 283c4762a1bSJed Brown 16 | 15 | 15 12 20 284c4762a1bSJed Brown \ | / |/ / 285c4762a1bSJed Brown 3 3------ 286c4762a1bSJed Brown */ 287*9371c9d4SSatish Balay PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscInt testNum, DM *dm) { 288c4762a1bSJed Brown DM idm; 289c4762a1bSJed Brown PetscInt depth = 3; 290c4762a1bSJed Brown PetscMPIInt rank; 291c4762a1bSJed Brown 292c4762a1bSJed Brown PetscFunctionBegin; 2939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 294dd400576SPatrick Sanan if (rank == 0) { 295c4762a1bSJed Brown switch (testNum) { 296*9371c9d4SSatish Balay case 0: { 297c4762a1bSJed Brown PetscInt numPoints[4] = {5, 9, 7, 2}; 298c4762a1bSJed Brown PetscInt coneSize[23] = {4, 4, 0, 0, 0, 0, 0, 3, 3, 3, 3, 3, 3, 3, 2, 2, 2, 2, 2, 2, 2, 2, 2}; 299c4762a1bSJed Brown PetscInt cones[47] = {7, 8, 9, 10, 10, 11, 12, 13, 14, 15, 16, 17, 18, 14, 16, 19, 17, 15, 18, 19, 20, 21, 19, 15, 22, 20, 18, 21, 22, 2, 4, 4, 3, 3, 2, 2, 5, 5, 4, 3, 5, 3, 6, 6, 5, 4, 6}; 300b5a892a1SMatthew G. Knepley PetscInt coneOrientations[47] = {0, 0, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, -1, -1, 0, -1, -1, -1, -1, 0, 0, -1, -1, 0, -1, -1, -1, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 301c4762a1bSJed Brown PetscScalar vertexCoords[15] = {0.0, 0.0, -0.5, 0.0, -0.5, 0.0, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5}; 302c4762a1bSJed Brown 3039566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 304*9371c9d4SSatish Balay } break; 305*9371c9d4SSatish Balay case 1: { 306c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 307c4762a1bSJed Brown PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 308c4762a1bSJed Brown PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 309c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 310c4762a1bSJed Brown PetscScalar vertexCoords[15] = {-1.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0}; 311c4762a1bSJed Brown 312c4762a1bSJed Brown depth = 1; 3139566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3149566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 3159566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 3169566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 317c4762a1bSJed Brown *dm = idm; 318*9371c9d4SSatish Balay } break; 319*9371c9d4SSatish Balay case 2: { 320c4762a1bSJed Brown PetscInt numPoints[2] = {4, 1}; 321c4762a1bSJed Brown PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 322c4762a1bSJed Brown PetscInt cones[4] = {2, 3, 4, 1}; 323c4762a1bSJed Brown PetscInt coneOrientations[4] = {0, 0, 0, 0}; 324c4762a1bSJed Brown PetscScalar vertexCoords[12] = {0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 325c4762a1bSJed Brown 326c4762a1bSJed Brown depth = 1; 3279566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3289566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 3299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 3309566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 331c4762a1bSJed Brown *dm = idm; 332*9371c9d4SSatish Balay } break; 33363a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown } else { 336c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 337c4762a1bSJed Brown 3389566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 339c4762a1bSJed Brown switch (testNum) { 340c4762a1bSJed Brown case 1: 3419566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 3429566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 3439566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 344c4762a1bSJed Brown *dm = idm; 345c4762a1bSJed Brown break; 346c4762a1bSJed Brown } 347c4762a1bSJed Brown } 348c4762a1bSJed Brown PetscFunctionReturn(0); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* Two tetrahedrons separated by a zero-volume cell with 6 vertices 352c4762a1bSJed Brown 353c4762a1bSJed Brown cell 6 ___33___10______ cell 354c4762a1bSJed Brown 0 / | \ |\ \ 1 355c4762a1bSJed Brown 21 | 23 | 29 27 356c4762a1bSJed Brown /12 24 14\ 30 \ \ 357c4762a1bSJed Brown 3-20-|----5--32-|---9--26--7 358c4762a1bSJed Brown \ 13| 11/ |18 / / 359c4762a1bSJed Brown 19 | 22 | 28 25 360c4762a1bSJed Brown \ | / |/ / 361c4762a1bSJed Brown 4----31----8------ 362c4762a1bSJed Brown cell 2 363c4762a1bSJed Brown */ 364*9371c9d4SSatish Balay PetscErrorCode CreateSimplexHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) { 365c4762a1bSJed Brown DM idm, hdm = NULL; 366c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 367c4762a1bSJed Brown PetscInt p; 368c4762a1bSJed Brown PetscMPIInt rank; 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionBegin; 3719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 372dd400576SPatrick Sanan if (rank == 0) { 373c4762a1bSJed Brown switch (testNum) { 374*9371c9d4SSatish Balay case 0: { 375c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 376c4762a1bSJed Brown PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 377c4762a1bSJed Brown PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 378c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 379c4762a1bSJed Brown PetscScalar vertexCoords[15] = {-1.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0}; 380c4762a1bSJed Brown PetscInt faultPoints[3] = {3, 4, 5}; 381c4762a1bSJed Brown 3829566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3839566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 384*9371c9d4SSatish Balay } break; 385*9371c9d4SSatish Balay case 1: { 386c4762a1bSJed Brown /* Tets 0,3,5 and 1,2,4 */ 387c4762a1bSJed Brown PetscInt numPoints[2] = {9, 6}; 388c4762a1bSJed Brown PetscInt coneSize[15] = {4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 389*9371c9d4SSatish Balay PetscInt cones[24] = {7, 8, 9, 6, 11, 13, 9, 14, 10, 11, 13, 9, 9, 10, 11, 7, 9, 14, 13, 12, 7, 8, 11, 9}; 390*9371c9d4SSatish Balay PetscInt coneOrientations[24] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 391*9371c9d4SSatish Balay PetscScalar vertexCoords[27] = {-2.0, -1.0, 0.0, -2.0, 0.0, 0.0, -2.0, 0.0, 1.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 2.0, -1.0, 0.0, 2.0, 0.0, 0.0, 2.0, 0.0, 1.0}; 392c4762a1bSJed Brown PetscInt faultPoints[3] = {9, 10, 11}; 393c4762a1bSJed Brown 3949566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3959566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 396*9371c9d4SSatish Balay } break; 39763a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 398c4762a1bSJed Brown } 3999566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 4019566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 4029566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 4039566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 4049566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault", &faultLabel)); 405caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm)); 4069566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&hybridLabel)); 4079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 4089566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 409c4762a1bSJed Brown *dm = hdm; 410c4762a1bSJed Brown } else { 411c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 412c4762a1bSJed Brown 4139566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 4149566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 4159566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 4169566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 4179566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 4189566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 419caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm)); 4209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 4219566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 422c4762a1bSJed Brown *dm = hdm; 423c4762a1bSJed Brown } 424c4762a1bSJed Brown PetscFunctionReturn(0); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427*9371c9d4SSatish Balay PetscErrorCode CreateTensorProduct_3D(MPI_Comm comm, PetscInt testNum, DM *dm) { 428c4762a1bSJed Brown DM idm; 429c4762a1bSJed Brown PetscMPIInt rank; 430c4762a1bSJed Brown 431c4762a1bSJed Brown PetscFunctionBegin; 4329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 433dd400576SPatrick Sanan if (rank == 0) { 434c4762a1bSJed Brown switch (testNum) { 435*9371c9d4SSatish Balay case 0: { 436c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 437c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 438c4762a1bSJed Brown PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 439c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 440*9371c9d4SSatish Balay PetscScalar vertexCoords[36] = {-1.0, -0.5, -0.5, -1.0, 0.5, -0.5, 0.0, 0.5, -0.5, 0.0, -0.5, -0.5, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 441c4762a1bSJed Brown 4429566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 443*9371c9d4SSatish Balay } break; 444*9371c9d4SSatish Balay case 1: { 445c4762a1bSJed Brown PetscInt numPoints[2] = {8, 1}; 446c4762a1bSJed Brown PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 447c4762a1bSJed Brown PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8}; 448c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 449*9371c9d4SSatish Balay PetscScalar vertexCoords[24] = {-1.0, -1.0, -1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, -1.0, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 450c4762a1bSJed Brown 4519566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 452*9371c9d4SSatish Balay } break; 45363a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown } else { 456c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 457c4762a1bSJed Brown 4589566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 459c4762a1bSJed Brown } 4609566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 4619566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 4629566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 463c4762a1bSJed Brown *dm = idm; 464c4762a1bSJed Brown PetscFunctionReturn(0); 465c4762a1bSJed Brown } 466c4762a1bSJed Brown 467*9371c9d4SSatish Balay PetscErrorCode CreateTensorProductHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) { 468c4762a1bSJed Brown DM idm, hdm = NULL; 469c4762a1bSJed Brown DMLabel faultLabel; 470c4762a1bSJed Brown PetscInt p; 471c4762a1bSJed Brown PetscMPIInt rank; 472c4762a1bSJed Brown 473c4762a1bSJed Brown PetscFunctionBegin; 4749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 475dd400576SPatrick Sanan if (rank == 0) { 476c4762a1bSJed Brown switch (testNum) { 477*9371c9d4SSatish Balay case 0: { 478c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 479c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 480c4762a1bSJed Brown PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 481c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 482*9371c9d4SSatish Balay PetscScalar vertexCoords[36] = {-1.0, -0.5, -0.5, -1.0, 0.5, -0.5, 0.0, 0.5, -0.5, 0.0, -0.5, -0.5, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 483c4762a1bSJed Brown PetscInt faultPoints[4] = {2, 3, 5, 6}; 484c4762a1bSJed Brown 4859566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4869566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 487*9371c9d4SSatish Balay } break; 488*9371c9d4SSatish Balay case 1: { 489c4762a1bSJed Brown PetscInt numPoints[2] = {30, 7}; 490c4762a1bSJed Brown PetscInt coneSize[37] = {8, 8, 8, 8, 8, 8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 491*9371c9d4SSatish Balay PetscInt cones[56] = {8, 21, 20, 7, 13, 12, 23, 24, 14, 15, 10, 9, 13, 8, 21, 24, 15, 16, 11, 10, 24, 21, 22, 25, 30, 29, 28, 21, 35, 24, 33, 34, 24, 21, 30, 35, 25, 36, 31, 22, 27, 20, 21, 28, 32, 33, 24, 23, 15, 24, 13, 14, 19, 18, 17, 26}; 492c4762a1bSJed Brown PetscInt coneOrientations[56] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 493*9371c9d4SSatish Balay PetscScalar vertexCoords[90] = {-2.0, -2.0, -2.0, -2.0, -1.0, -2.0, -3.0, 0.0, -2.0, -2.0, 1.0, -2.0, -2.0, 2.0, -2.0, -2.0, -2.0, 0.0, -2.0, -1.0, 0.0, -3.0, 0.0, 0.0, -2.0, 1.0, 0.0, -2.0, 2.0, 0.0, 494*9371c9d4SSatish Balay -2.0, -1.0, 2.0, -3.0, 0.0, 2.0, -2.0, 1.0, 2.0, 0.0, -2.0, -2.0, 0.0, 0.0, -2.0, 0.0, 2.0, -2.0, 0.0, -2.0, 0.0, 0.0, 0.0, 0.0, 0.0, 2.0, 0.0, 0.0, 0.0, 2.0, 495*9371c9d4SSatish Balay 2.0, -2.0, -2.0, 2.0, -1.0, -2.0, 3.0, 0.0, -2.0, 2.0, 1.0, -2.0, 2.0, 2.0, -2.0, 2.0, -2.0, 0.0, 2.0, -1.0, 0.0, 3.0, 0.0, 0.0, 2.0, 1.0, 0.0, 2.0, 2.0, 0.0}; 496c4762a1bSJed Brown PetscInt faultPoints[6] = {20, 21, 22, 23, 24, 25}; 497c4762a1bSJed Brown 4989566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4999566063dSJacob Faibussowitsch for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 500*9371c9d4SSatish Balay } break; 50163a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 502c4762a1bSJed Brown } 5039566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 5049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 5059566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 5069566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 5079566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 5089566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault", &faultLabel)); 509caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, NULL, NULL, NULL, &hdm)); 5109566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 5119566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 512c4762a1bSJed Brown *dm = hdm; 513c4762a1bSJed Brown } else { 514c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 515c4762a1bSJed Brown 5169566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 5179566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 5189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_")); 5199566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 5209566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(idm)); 5219566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-dm_view")); 522caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm)); 5239566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 5249566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 525c4762a1bSJed Brown *dm = hdm; 526c4762a1bSJed Brown } 527c4762a1bSJed Brown PetscFunctionReturn(0); 528c4762a1bSJed Brown } 529c4762a1bSJed Brown 530*9371c9d4SSatish Balay PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 531c4762a1bSJed Brown PetscInt dim = user->dim; 532c4762a1bSJed Brown PetscBool cellHybrid = user->cellHybrid; 533c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 534c4762a1bSJed Brown PetscMPIInt rank, size; 535c4762a1bSJed Brown 536c4762a1bSJed Brown PetscFunctionBegin; 5379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 5399566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 5409566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 5419566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 542c4762a1bSJed Brown switch (dim) { 543c4762a1bSJed Brown case 1: 54463a3b9bcSJacob Faibussowitsch PetscCheck(!cellHybrid, comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim); 5459566063dSJacob Faibussowitsch PetscCall(CreateSimplex_1D(comm, dm)); 546c4762a1bSJed Brown break; 547c4762a1bSJed Brown case 2: 548c4762a1bSJed Brown if (cellSimplex) { 549c4762a1bSJed Brown if (cellHybrid) { 5509566063dSJacob Faibussowitsch PetscCall(CreateSimplexHybrid_2D(comm, user->testNum, dm)); 551c4762a1bSJed Brown } else { 5529566063dSJacob Faibussowitsch PetscCall(CreateSimplex_2D(comm, dm)); 553c4762a1bSJed Brown } 554c4762a1bSJed Brown } else { 555c4762a1bSJed Brown if (cellHybrid) { 5569566063dSJacob Faibussowitsch PetscCall(CreateTensorProductHybrid_2D(comm, user->testNum, dm)); 557c4762a1bSJed Brown } else { 5589566063dSJacob Faibussowitsch PetscCall(CreateTensorProduct_2D(comm, user->testNum, dm)); 559c4762a1bSJed Brown } 560c4762a1bSJed Brown } 561c4762a1bSJed Brown break; 562c4762a1bSJed Brown case 3: 563c4762a1bSJed Brown if (cellSimplex) { 564c4762a1bSJed Brown if (cellHybrid) { 5659566063dSJacob Faibussowitsch PetscCall(CreateSimplexHybrid_3D(comm, user->testNum, dm)); 566c4762a1bSJed Brown } else { 5679566063dSJacob Faibussowitsch PetscCall(CreateSimplex_3D(comm, user->testNum, dm)); 568c4762a1bSJed Brown } 569c4762a1bSJed Brown } else { 570c4762a1bSJed Brown if (cellHybrid) { 5719566063dSJacob Faibussowitsch PetscCall(CreateTensorProductHybrid_3D(comm, user->testNum, dm)); 572c4762a1bSJed Brown } else { 5739566063dSJacob Faibussowitsch PetscCall(CreateTensorProduct_3D(comm, user->testNum, dm)); 574c4762a1bSJed Brown } 575c4762a1bSJed Brown } 576c4762a1bSJed Brown break; 577*9371c9d4SSatish Balay default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim); 578c4762a1bSJed Brown } 579c4762a1bSJed Brown if (user->testPartition && size > 1) { 580c4762a1bSJed Brown PetscPartitioner part; 581c4762a1bSJed Brown PetscInt *sizes = NULL; 582c4762a1bSJed Brown PetscInt *points = NULL; 583c4762a1bSJed Brown 584dd400576SPatrick Sanan if (rank == 0) { 585c4762a1bSJed Brown if (dim == 2 && cellSimplex && !cellHybrid && size == 2) { 586c4762a1bSJed Brown switch (user->testNum) { 587c4762a1bSJed Brown case 0: { 588c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 1}; 589c4762a1bSJed Brown PetscInt triPoints_p2[2] = {0, 1}; 590c4762a1bSJed Brown 5919566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 5929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, triSizes_p2, 2)); 593*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, triPoints_p2, 2)); 594*9371c9d4SSatish Balay break; 595*9371c9d4SSatish Balay } 596*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum); 597c4762a1bSJed Brown } 598c4762a1bSJed Brown } else if (dim == 2 && cellSimplex && cellHybrid && size == 2) { 599c4762a1bSJed Brown switch (user->testNum) { 600c4762a1bSJed Brown case 0: { 601c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 2}; 602c4762a1bSJed Brown PetscInt triPoints_p2[3] = {0, 1, 2}; 603c4762a1bSJed Brown 6049566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 6059566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, triSizes_p2, 2)); 606*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, triPoints_p2, 3)); 607*9371c9d4SSatish Balay break; 608*9371c9d4SSatish Balay } 609*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular hybrid mesh on 2 procs", user->testNum); 610c4762a1bSJed Brown } 611c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && !cellHybrid && size == 2) { 612c4762a1bSJed Brown switch (user->testNum) { 613c4762a1bSJed Brown case 0: { 614c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 1}; 615c4762a1bSJed Brown PetscInt quadPoints_p2[2] = {0, 1}; 616c4762a1bSJed Brown 6179566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 6189566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2)); 619*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, quadPoints_p2, 2)); 620*9371c9d4SSatish Balay break; 621*9371c9d4SSatish Balay } 622*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral mesh on 2 procs", user->testNum); 623c4762a1bSJed Brown } 624c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && cellHybrid && size == 2) { 625c4762a1bSJed Brown switch (user->testNum) { 626c4762a1bSJed Brown case 0: { 627c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 2}; 628c4762a1bSJed Brown PetscInt quadPoints_p2[3] = {0, 1, 2}; 629c4762a1bSJed Brown 6309566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 6319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2)); 632*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, quadPoints_p2, 3)); 633*9371c9d4SSatish Balay break; 634*9371c9d4SSatish Balay } 635*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral hybrid mesh on 2 procs", user->testNum); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && !cellHybrid && size == 2) { 638c4762a1bSJed Brown switch (user->testNum) { 639c4762a1bSJed Brown case 0: { 640c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 1}; 641c4762a1bSJed Brown PetscInt tetPoints_p2[2] = {0, 1}; 642c4762a1bSJed Brown 6439566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 6449566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2)); 645*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, tetPoints_p2, 2)); 646*9371c9d4SSatish Balay break; 647*9371c9d4SSatish Balay } 648c4762a1bSJed Brown case 1: { 649c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 1}; 650c4762a1bSJed Brown PetscInt tetPoints_p2[2] = {0, 1}; 651c4762a1bSJed Brown 6529566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 6539566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2)); 654*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, tetPoints_p2, 2)); 655*9371c9d4SSatish Balay break; 656*9371c9d4SSatish Balay } 657*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for tetrahedral mesh on 2 procs", user->testNum); 658c4762a1bSJed Brown } 659c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && cellHybrid && size == 2) { 660c4762a1bSJed Brown switch (user->testNum) { 661c4762a1bSJed Brown case 0: { 662c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 2}; 663c4762a1bSJed Brown PetscInt tetPoints_p2[3] = {0, 1, 2}; 664c4762a1bSJed Brown 6659566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 6669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2)); 667*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, tetPoints_p2, 3)); 668*9371c9d4SSatish Balay break; 669*9371c9d4SSatish Balay } 670c4762a1bSJed Brown case 1: { 671c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {3, 4}; 672c4762a1bSJed Brown PetscInt tetPoints_p2[7] = {0, 3, 5, 1, 2, 4, 6}; 673c4762a1bSJed Brown 6749566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 7, &points)); 6759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2)); 676*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, tetPoints_p2, 7)); 677*9371c9d4SSatish Balay break; 678*9371c9d4SSatish Balay } 679*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for tetrahedral hybrid mesh on 2 procs", user->testNum); 680c4762a1bSJed Brown } 681c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && !cellHybrid && size == 2) { 682c4762a1bSJed Brown switch (user->testNum) { 683c4762a1bSJed Brown case 0: { 684c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 1}; 685c4762a1bSJed Brown PetscInt hexPoints_p2[2] = {0, 1}; 686c4762a1bSJed Brown 6879566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 6889566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2)); 689*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, hexPoints_p2, 2)); 690*9371c9d4SSatish Balay break; 691*9371c9d4SSatish Balay } 692*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral mesh on 2 procs", user->testNum); 693c4762a1bSJed Brown } 694c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && cellHybrid && size == 2) { 695c4762a1bSJed Brown switch (user->testNum) { 696c4762a1bSJed Brown case 0: { 697c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 1}; 698c4762a1bSJed Brown PetscInt hexPoints_p2[2] = {0, 1}; 699c4762a1bSJed Brown 7009566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 7019566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2)); 702*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, hexPoints_p2, 2)); 703*9371c9d4SSatish Balay break; 704*9371c9d4SSatish Balay } 705c4762a1bSJed Brown case 1: { 706c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {5, 4}; 707c4762a1bSJed Brown PetscInt hexPoints_p2[9] = {3, 4, 5, 7, 8, 0, 1, 2, 6}; 708c4762a1bSJed Brown 7099566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 9, &points)); 7109566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2)); 711*9371c9d4SSatish Balay PetscCall(PetscArraycpy(points, hexPoints_p2, 9)); 712*9371c9d4SSatish Balay break; 713*9371c9d4SSatish Balay } 714*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral hybrid mesh on 2 procs", user->testNum); 715c4762a1bSJed Brown } 716c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition"); 717c4762a1bSJed Brown } 7189566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 7199566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 7209566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points)); 7219566063dSJacob Faibussowitsch PetscCall(PetscFree2(sizes, points)); 722c4762a1bSJed Brown } else { 723c4762a1bSJed Brown PetscPartitioner part; 724c4762a1bSJed Brown 7259566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 7269566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 727c4762a1bSJed Brown } 728c4762a1bSJed Brown { 729c4762a1bSJed Brown DM pdm = NULL; 730c4762a1bSJed Brown 7319566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm)); 732c4762a1bSJed Brown if (pdm) { 7339566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view")); 7349566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 735c4762a1bSJed Brown *dm = pdm; 736c4762a1bSJed Brown } 737c4762a1bSJed Brown } 7389566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 7399566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre")); 7409566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 741c4762a1bSJed Brown if (user->uninterpolate || user->reinterpolate) { 742c4762a1bSJed Brown DM udm = NULL; 743c4762a1bSJed Brown 7449566063dSJacob Faibussowitsch PetscCall(DMPlexUninterpolate(*dm, &udm)); 7459566063dSJacob Faibussowitsch PetscCall(DMPlexCopyCoordinates(*dm, udm)); 7469566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 747c4762a1bSJed Brown *dm = udm; 748c4762a1bSJed Brown } 749c4762a1bSJed Brown if (user->reinterpolate) { 750c4762a1bSJed Brown DM idm = NULL; 751c4762a1bSJed Brown 7529566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 7539566063dSJacob Faibussowitsch PetscCall(DMPlexCopyCoordinates(*dm, idm)); 7549566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 755c4762a1bSJed Brown *dm = idm; 756c4762a1bSJed Brown } 7579566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 7589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Hybrid Mesh")); 7599566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 7609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "hyb_")); 7619566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 762c4762a1bSJed Brown PetscFunctionReturn(0); 763c4762a1bSJed Brown } 764c4762a1bSJed Brown 765*9371c9d4SSatish Balay int main(int argc, char **argv) { 766c4762a1bSJed Brown DM dm; 767c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 768c4762a1bSJed Brown 769327415f7SBarry Smith PetscFunctionBeginUser; 7709566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 7719566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 7729566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 7739566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 7749566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 775b122ec5aSJacob Faibussowitsch return 0; 776c4762a1bSJed Brown } 777c4762a1bSJed Brown 778c4762a1bSJed Brown /*TEST 779c4762a1bSJed Brown 780c4762a1bSJed Brown # 1D Simplex 29-31 781c4762a1bSJed Brown testset: 78254fcfd0cSMatthew G. Knepley args: -dim 1 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 783c4762a1bSJed Brown test: 784c4762a1bSJed Brown suffix: 29 785c4762a1bSJed Brown test: 786c4762a1bSJed Brown suffix: 30 787c4762a1bSJed Brown args: -dm_refine 1 788c4762a1bSJed Brown test: 789c4762a1bSJed Brown suffix: 31 790c4762a1bSJed Brown args: -dm_refine 5 791c4762a1bSJed Brown 792c4762a1bSJed Brown # 2D Simplex 0-3 793c4762a1bSJed Brown testset: 79454fcfd0cSMatthew G. Knepley args: -dim 2 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 795c4762a1bSJed Brown test: 796c4762a1bSJed Brown suffix: 0 797c4762a1bSJed Brown test: 798c4762a1bSJed Brown suffix: 1 7997f9d8d6cSVaclav Hapla args: -dm_refine 1 800c4762a1bSJed Brown test: 801c4762a1bSJed Brown suffix: 2 802c4762a1bSJed Brown nsize: 2 803c4762a1bSJed Brown test: 804c4762a1bSJed Brown suffix: 3 805c4762a1bSJed Brown nsize: 2 8067f9d8d6cSVaclav Hapla args: -dm_refine 1 807c4762a1bSJed Brown test: 808c4762a1bSJed Brown suffix: 32 809c4762a1bSJed Brown args: -dm_refine 1 -uninterpolate 810c4762a1bSJed Brown test: 811c4762a1bSJed Brown suffix: 33 812c4762a1bSJed Brown nsize: 2 813c4762a1bSJed Brown args: -dm_refine 1 -uninterpolate 814c4762a1bSJed Brown test: 815c4762a1bSJed Brown suffix: 34 816c4762a1bSJed Brown nsize: 2 817c4762a1bSJed Brown args: -dm_refine 3 -uninterpolate 818c4762a1bSJed Brown 819c4762a1bSJed Brown # 2D Hybrid Simplex 4-7 820c4762a1bSJed Brown testset: 82154fcfd0cSMatthew G. Knepley args: -dim 2 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 822c4762a1bSJed Brown test: 823c4762a1bSJed Brown suffix: 4 824c4762a1bSJed Brown test: 825c4762a1bSJed Brown suffix: 5 826c4762a1bSJed Brown args: -dm_refine 1 827c4762a1bSJed Brown test: 828c4762a1bSJed Brown suffix: 6 829c4762a1bSJed Brown nsize: 2 830c4762a1bSJed Brown test: 831c4762a1bSJed Brown suffix: 7 832c4762a1bSJed Brown nsize: 2 833c4762a1bSJed Brown args: -dm_refine 1 834c4762a1bSJed Brown test: 835c4762a1bSJed Brown suffix: 24 836c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 837c4762a1bSJed Brown 838c4762a1bSJed Brown # 2D Quad 12-13 839c4762a1bSJed Brown testset: 84054fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 841c4762a1bSJed Brown test: 842c4762a1bSJed Brown suffix: 12 843c4762a1bSJed Brown args: -dm_refine 1 844c4762a1bSJed Brown test: 845c4762a1bSJed Brown suffix: 13 846c4762a1bSJed Brown nsize: 2 847c4762a1bSJed Brown args: -dm_refine 1 848c4762a1bSJed Brown 849c4762a1bSJed Brown # 2D Hybrid Quad 27-28 850c4762a1bSJed Brown testset: 85154fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 852c4762a1bSJed Brown test: 853c4762a1bSJed Brown suffix: 27 854c4762a1bSJed Brown args: -dm_refine 1 855c4762a1bSJed Brown test: 856c4762a1bSJed Brown suffix: 28 857c4762a1bSJed Brown nsize: 2 858c4762a1bSJed Brown args: -dm_refine 1 859c4762a1bSJed Brown 860c4762a1bSJed Brown # 3D Simplex 8-11 861c4762a1bSJed Brown testset: 86254fcfd0cSMatthew G. Knepley args: -dim 3 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 863c4762a1bSJed Brown test: 864c4762a1bSJed Brown suffix: 8 865c4762a1bSJed Brown args: -dm_refine 1 866c4762a1bSJed Brown test: 867c4762a1bSJed Brown suffix: 9 868c4762a1bSJed Brown nsize: 2 869c4762a1bSJed Brown args: -dm_refine 1 870c4762a1bSJed Brown test: 871c4762a1bSJed Brown suffix: 10 872c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown suffix: 11 875c4762a1bSJed Brown nsize: 2 876c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 877c4762a1bSJed Brown test: 878c4762a1bSJed Brown suffix: 25 879c4762a1bSJed Brown args: -test_num 2 -dm_refine 2 880c4762a1bSJed Brown 881c4762a1bSJed Brown # 3D Hybrid Simplex 16-19 882c4762a1bSJed Brown testset: 88354fcfd0cSMatthew G. Knepley args: -dim 3 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 884c4762a1bSJed Brown test: 885c4762a1bSJed Brown suffix: 16 886c4762a1bSJed Brown args: -dm_refine 1 887c4762a1bSJed Brown test: 888c4762a1bSJed Brown suffix: 17 889c4762a1bSJed Brown nsize: 2 890c4762a1bSJed Brown args: -dm_refine 1 891c4762a1bSJed Brown test: 892c4762a1bSJed Brown suffix: 18 893c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 894c4762a1bSJed Brown test: 895c4762a1bSJed Brown suffix: 19 896c4762a1bSJed Brown nsize: 2 897c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 898c4762a1bSJed Brown 899c4762a1bSJed Brown # 3D Hex 14-15 900c4762a1bSJed Brown testset: 90154fcfd0cSMatthew G. Knepley args: -dim 3 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 902c4762a1bSJed Brown test: 903c4762a1bSJed Brown suffix: 14 904c4762a1bSJed Brown args: -dm_refine 1 905c4762a1bSJed Brown test: 906c4762a1bSJed Brown suffix: 15 907c4762a1bSJed Brown nsize: 2 908c4762a1bSJed Brown args: -dm_refine 1 909c4762a1bSJed Brown test: 910c4762a1bSJed Brown suffix: 26 911c4762a1bSJed Brown args: -test_num 1 -dm_refine 2 912c4762a1bSJed Brown 913c4762a1bSJed Brown # 3D Hybrid Hex 20-23 914c4762a1bSJed Brown testset: 91554fcfd0cSMatthew G. Knepley args: -dim 3 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 916c4762a1bSJed Brown test: 917c4762a1bSJed Brown suffix: 20 918c4762a1bSJed Brown args: -dm_refine 1 919c4762a1bSJed Brown test: 920c4762a1bSJed Brown suffix: 21 921c4762a1bSJed Brown nsize: 2 922c4762a1bSJed Brown args: -dm_refine 1 923c4762a1bSJed Brown test: 924c4762a1bSJed Brown suffix: 22 925c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 926c4762a1bSJed Brown test: 927c4762a1bSJed Brown suffix: 23 928c4762a1bSJed Brown nsize: 2 929c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 930c4762a1bSJed Brown 931c4762a1bSJed Brown # Hybrid interpolation 93254fcfd0cSMatthew G. Knepley # TODO Setup new tests (like -reinterpolate) that interpolate hybrid cells 933c4762a1bSJed Brown testset: 934c4762a1bSJed Brown nsize: 2 93554fcfd0cSMatthew G. Knepley args: -test_partition 0 -petscpartitioner_type simple -dm_view -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 936c4762a1bSJed Brown test: 937c4762a1bSJed Brown suffix: hybint_2d_0 938c4762a1bSJed Brown args: -dim 2 -dm_refine 2 939c4762a1bSJed Brown test: 940c4762a1bSJed Brown suffix: hybint_2d_1 941c4762a1bSJed Brown args: -dim 2 -dm_refine 2 -test_num 1 942c4762a1bSJed Brown test: 943c4762a1bSJed Brown suffix: hybint_3d_0 944c4762a1bSJed Brown args: -dim 3 -dm_refine 1 945c4762a1bSJed Brown test: 946c4762a1bSJed Brown suffix: hybint_3d_1 947c4762a1bSJed Brown args: -dim 3 -dm_refine 1 -test_num 1 948c4762a1bSJed Brown 949c4762a1bSJed Brown TEST*/ 950