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 16c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscErrorCode ierr; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 21c4762a1bSJed Brown options->debug = 0; 22c4762a1bSJed Brown options->dim = 2; 23c4762a1bSJed Brown options->cellHybrid = PETSC_TRUE; 24c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 25c4762a1bSJed Brown options->testPartition = PETSC_TRUE; 26c4762a1bSJed Brown options->testNum = 0; 27c4762a1bSJed Brown options->uninterpolate = PETSC_FALSE; 28c4762a1bSJed Brown options->reinterpolate = PETSC_FALSE; 29c4762a1bSJed Brown 30c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-debug", "The debugging level", "ex4.c", options->debug, &options->debug, NULL,0)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL,1,3)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-cell_hybrid", "Use a hybrid mesh", "ex4.c", options->cellHybrid, &options->cellHybrid, NULL)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex4.c", options->cellSimplex, &options->cellSimplex, NULL)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex4.c", options->testPartition, &options->testPartition, NULL)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex4.c", options->testNum, &options->testNum, NULL,0)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex4.c", options->uninterpolate, &options->uninterpolate, NULL)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-reinterpolate", "Reinterpolate the mesh at the end", "ex4.c", options->reinterpolate, &options->reinterpolate, NULL)); 391e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 40c4762a1bSJed Brown PetscFunctionReturn(0); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Two segments 44c4762a1bSJed Brown 45c4762a1bSJed Brown 2-------0-------3-------1-------4 46c4762a1bSJed Brown 47c4762a1bSJed Brown become 48c4762a1bSJed Brown 49c4762a1bSJed Brown 4---0---7---1---5---2---8---3---6 50c4762a1bSJed Brown 51c4762a1bSJed Brown */ 52c4762a1bSJed Brown PetscErrorCode CreateSimplex_1D(MPI_Comm comm, DM *dm) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown PetscInt depth = 1; 55c4762a1bSJed Brown PetscMPIInt rank; 56c4762a1bSJed Brown 57c4762a1bSJed Brown PetscFunctionBegin; 585f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 59dd400576SPatrick Sanan if (rank == 0) { 60c4762a1bSJed Brown PetscInt numPoints[2] = {3, 2}; 61c4762a1bSJed Brown PetscInt coneSize[5] = {2, 2, 0, 0, 0}; 62c4762a1bSJed Brown PetscInt cones[4] = {2, 3, 3, 4}; 63c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0}; 64c4762a1bSJed Brown PetscScalar vertexCoords[3] = {-1.0, 0.0, 1.0}; 65c4762a1bSJed Brown 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 67c4762a1bSJed Brown } else { 68c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 69c4762a1bSJed Brown 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Two triangles 76c4762a1bSJed Brown 4 77c4762a1bSJed Brown / | \ 78c4762a1bSJed Brown 8 | 10 79c4762a1bSJed Brown / | \ 80c4762a1bSJed Brown 2 0 7 1 5 81c4762a1bSJed Brown \ | / 82c4762a1bSJed Brown 6 | 9 83c4762a1bSJed Brown \ | / 84c4762a1bSJed Brown 3 85c4762a1bSJed Brown 86c4762a1bSJed Brown Becomes 87c4762a1bSJed Brown 10 88c4762a1bSJed Brown / | \ 89c4762a1bSJed Brown 21 | 26 90c4762a1bSJed Brown / | \ 91c4762a1bSJed Brown 14 2 20 4 16 92c4762a1bSJed Brown /|\ | /|\ 93c4762a1bSJed Brown 22 | 28 | 32 | 25 94c4762a1bSJed Brown / | \ | / | 6\ 95c4762a1bSJed Brown 8 29 3 13 7 31 11 96c4762a1bSJed Brown \0 | / | \ | / 97c4762a1bSJed Brown 17 | 27 | 30 | 24 98c4762a1bSJed Brown \|/ | \|/ 99c4762a1bSJed Brown 12 1 19 5 15 100c4762a1bSJed Brown \ | / 101c4762a1bSJed Brown 18 | 23 102c4762a1bSJed Brown \ | / 103c4762a1bSJed Brown 9 104c4762a1bSJed Brown */ 105c4762a1bSJed Brown PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *dm) 106c4762a1bSJed Brown { 107c4762a1bSJed Brown PetscInt depth = 2; 108c4762a1bSJed Brown PetscMPIInt rank; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBegin; 1115f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 112dd400576SPatrick Sanan if (rank == 0) { 113c4762a1bSJed Brown PetscInt numPoints[3] = {4, 5, 2}; 114c4762a1bSJed Brown PetscInt coneSize[11] = {3, 3, 0, 0, 0, 0, 2, 2, 2, 2, 2}; 115c4762a1bSJed Brown PetscInt cones[16] = {6, 7, 8, 7, 9, 10, 2, 3, 3, 4, 4, 2, 3, 5, 5, 4}; 116b5a892a1SMatthew G. Knepley PetscInt coneOrientations[16] = {0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 117c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-0.5, 0.0, 0.0, -0.5, 0.0, 0.5, 0.5, 0.0}; 118c4762a1bSJed Brown 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 120c4762a1bSJed Brown } else { 121c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 122c4762a1bSJed Brown 1235f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown PetscFunctionReturn(0); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* Two triangles separated by a zero-volume cell with 4 vertices/2 edges 129c4762a1bSJed Brown 5--16--8 130c4762a1bSJed Brown / | | \ 131c4762a1bSJed Brown 11 | | 12 132c4762a1bSJed Brown / | | \ 133c4762a1bSJed Brown 3 0 10 2 14 1 6 134c4762a1bSJed Brown \ | | / 135c4762a1bSJed Brown 9 | | 13 136c4762a1bSJed Brown \ | | / 137c4762a1bSJed Brown 4--15--7 138c4762a1bSJed Brown */ 139c4762a1bSJed Brown PetscErrorCode CreateSimplexHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 140c4762a1bSJed Brown { 141c4762a1bSJed Brown DM idm, hdm = NULL; 142c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 143c4762a1bSJed Brown PetscInt p; 144c4762a1bSJed Brown PetscMPIInt rank; 145c4762a1bSJed Brown 146c4762a1bSJed Brown PetscFunctionBegin; 1475f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 148dd400576SPatrick Sanan if (rank == 0) { 149c4762a1bSJed Brown switch (testNum) { 150c4762a1bSJed Brown case 0: 151c4762a1bSJed Brown { 152c4762a1bSJed Brown PetscInt numPoints[2] = {4, 2}; 153c4762a1bSJed Brown PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 154c4762a1bSJed Brown PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 155c4762a1bSJed Brown PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 156c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, 1.0, 0.5}; 157c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 158c4762a1bSJed Brown 1595f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 1605f80ce2aSJacob Faibussowitsch for (p = 0; p < 2; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown break; 163c4762a1bSJed Brown case 1: 164c4762a1bSJed Brown { 165c4762a1bSJed Brown PetscInt numPoints[2] = {5, 4}; 166c4762a1bSJed Brown PetscInt coneSize[9] = {3, 3, 3, 3, 0, 0, 0, 0, 0}; 167c4762a1bSJed Brown PetscInt cones[12] = {4, 5, 6, 6, 7, 4, 6, 5, 8, 6, 8, 7}; 168c4762a1bSJed Brown PetscInt coneOrientations[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 169c4762a1bSJed Brown PetscScalar vertexCoords[10] = {-1.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0}; 170c4762a1bSJed Brown PetscInt faultPoints[3] = {5, 6, 7}; 171c4762a1bSJed Brown 1725f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 1735f80ce2aSJacob Faibussowitsch for (p = 0; p < 3; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown break; 17698921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 177c4762a1bSJed Brown } 1785f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, "fault", &faultLabel)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&hybridLabel)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 188c4762a1bSJed Brown *dm = hdm; 189c4762a1bSJed Brown } else { 190c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 191c4762a1bSJed Brown 1925f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 201c4762a1bSJed Brown *dm = hdm; 202c4762a1bSJed Brown } 203c4762a1bSJed Brown PetscFunctionReturn(0); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Two quadrilaterals 207c4762a1bSJed Brown 208c4762a1bSJed Brown 5----10-----4----14-----7 209c4762a1bSJed Brown | | | 210c4762a1bSJed Brown | | | 211c4762a1bSJed Brown | | | 212c4762a1bSJed Brown 11 0 9 1 13 213c4762a1bSJed Brown | | | 214c4762a1bSJed Brown | | | 215c4762a1bSJed Brown | | | 216c4762a1bSJed Brown 2-----8-----3----12-----6 217c4762a1bSJed Brown */ 218c4762a1bSJed Brown PetscErrorCode CreateTensorProduct_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 219c4762a1bSJed Brown { 220c4762a1bSJed Brown PetscInt depth = 2; 221c4762a1bSJed Brown PetscMPIInt rank; 222c4762a1bSJed Brown 223c4762a1bSJed Brown PetscFunctionBegin; 2245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 225dd400576SPatrick Sanan if (rank == 0) { 226c4762a1bSJed Brown PetscInt numPoints[3] = {6, 7, 2}; 227c4762a1bSJed Brown PetscInt coneSize[15] = {4, 4, 0, 0, 0, 0, 0, 0, 2, 2, 2, 2, 2, 2, 2}; 228c4762a1bSJed 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}; 229b5a892a1SMatthew 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}; 230c4762a1bSJed 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}; 231c4762a1bSJed Brown 2325f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 233c4762a1bSJed Brown } else { 234c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 235c4762a1bSJed Brown 2365f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown PetscFunctionReturn(0); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown PetscErrorCode CreateTensorProductHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 242c4762a1bSJed Brown { 243c4762a1bSJed Brown DM idm, hdm = NULL; 244c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 245c4762a1bSJed Brown PetscInt p; 246c4762a1bSJed Brown PetscMPIInt rank; 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionBegin; 2495f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 250dd400576SPatrick Sanan if (rank == 0) { 251c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 252c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 253c4762a1bSJed Brown PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4,}; 254c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 255c4762a1bSJed 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}; 256c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 257c4762a1bSJed Brown 2585f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 2595f80ce2aSJacob Faibussowitsch for (p = 0; p < 2; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, "fault", &faultLabel)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&hybridLabel)); 268c4762a1bSJed Brown } else { 269c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 270c4762a1bSJed Brown 2715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 2775f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm)); 278c4762a1bSJed Brown } 2795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 281c4762a1bSJed Brown *dm = hdm; 282c4762a1bSJed Brown PetscFunctionReturn(0); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* Two tetrahedrons 286c4762a1bSJed Brown 287c4762a1bSJed Brown cell 5 5______ cell 288c4762a1bSJed Brown 0 / | \ |\ \ 1 289c4762a1bSJed Brown 17 | 18 | 18 13 21 290c4762a1bSJed Brown /8 19 10\ 19 \ \ 291c4762a1bSJed Brown 2-14-|----4 | 4--22--6 292c4762a1bSJed Brown \ 9 | 7 / |10 / / 293c4762a1bSJed Brown 16 | 15 | 15 12 20 294c4762a1bSJed Brown \ | / |/ / 295c4762a1bSJed Brown 3 3------ 296c4762a1bSJed Brown */ 297c4762a1bSJed Brown PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 298c4762a1bSJed Brown { 299c4762a1bSJed Brown DM idm; 300c4762a1bSJed Brown PetscInt depth = 3; 301c4762a1bSJed Brown PetscMPIInt rank; 302c4762a1bSJed Brown 303c4762a1bSJed Brown PetscFunctionBegin; 3045f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 305dd400576SPatrick Sanan if (rank == 0) { 306c4762a1bSJed Brown switch (testNum) { 307c4762a1bSJed Brown case 0: 308c4762a1bSJed Brown { 309c4762a1bSJed Brown PetscInt numPoints[4] = {5, 9, 7, 2}; 310c4762a1bSJed 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}; 311c4762a1bSJed 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}; 312b5a892a1SMatthew 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}; 313c4762a1bSJed 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}; 314c4762a1bSJed Brown 3155f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown break; 318c4762a1bSJed Brown case 1: 319c4762a1bSJed Brown { 320c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 321c4762a1bSJed Brown PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 322c4762a1bSJed Brown PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 323c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 324c4762a1bSJed 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}; 325c4762a1bSJed Brown 326c4762a1bSJed Brown depth = 1; 3275f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-in_dm_view")); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 331c4762a1bSJed Brown *dm = idm; 332c4762a1bSJed Brown } 333c4762a1bSJed Brown break; 334c4762a1bSJed Brown case 2: 335c4762a1bSJed Brown { 336c4762a1bSJed Brown PetscInt numPoints[2] = {4, 1}; 337c4762a1bSJed Brown PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 338c4762a1bSJed Brown PetscInt cones[4] = {2, 3, 4, 1}; 339c4762a1bSJed Brown PetscInt coneOrientations[4] = {0, 0, 0, 0}; 340c4762a1bSJed 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}; 341c4762a1bSJed Brown 342c4762a1bSJed Brown depth = 1; 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-in_dm_view")); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 347c4762a1bSJed Brown *dm = idm; 348c4762a1bSJed Brown } 349c4762a1bSJed Brown break; 35098921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown } else { 353c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 354c4762a1bSJed Brown 3555f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 356c4762a1bSJed Brown switch (testNum) { 357c4762a1bSJed Brown case 1: 3585f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-in_dm_view")); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 361c4762a1bSJed Brown *dm = idm; 362c4762a1bSJed Brown break; 363c4762a1bSJed Brown } 364c4762a1bSJed Brown } 365c4762a1bSJed Brown PetscFunctionReturn(0); 366c4762a1bSJed Brown } 367c4762a1bSJed Brown 368c4762a1bSJed Brown /* Two tetrahedrons separated by a zero-volume cell with 6 vertices 369c4762a1bSJed Brown 370c4762a1bSJed Brown cell 6 ___33___10______ cell 371c4762a1bSJed Brown 0 / | \ |\ \ 1 372c4762a1bSJed Brown 21 | 23 | 29 27 373c4762a1bSJed Brown /12 24 14\ 30 \ \ 374c4762a1bSJed Brown 3-20-|----5--32-|---9--26--7 375c4762a1bSJed Brown \ 13| 11/ |18 / / 376c4762a1bSJed Brown 19 | 22 | 28 25 377c4762a1bSJed Brown \ | / |/ / 378c4762a1bSJed Brown 4----31----8------ 379c4762a1bSJed Brown cell 2 380c4762a1bSJed Brown */ 381c4762a1bSJed Brown PetscErrorCode CreateSimplexHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 382c4762a1bSJed Brown { 383c4762a1bSJed Brown DM idm, hdm = NULL; 384c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 385c4762a1bSJed Brown PetscInt p; 386c4762a1bSJed Brown PetscMPIInt rank; 387c4762a1bSJed Brown 388c4762a1bSJed Brown PetscFunctionBegin; 3895f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 390dd400576SPatrick Sanan if (rank == 0) { 391c4762a1bSJed Brown switch (testNum) { 392c4762a1bSJed Brown case 0: 393c4762a1bSJed Brown { 394c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 395c4762a1bSJed Brown PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 396c4762a1bSJed Brown PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 397c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 398c4762a1bSJed 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}; 399c4762a1bSJed Brown PetscInt faultPoints[3] = {3, 4, 5}; 400c4762a1bSJed Brown 4015f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4025f80ce2aSJacob Faibussowitsch for (p = 0; p < 3; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown break; 405c4762a1bSJed Brown case 1: 406c4762a1bSJed Brown { 407c4762a1bSJed Brown /* Tets 0,3,5 and 1,2,4 */ 408c4762a1bSJed Brown PetscInt numPoints[2] = {9, 6}; 409c4762a1bSJed Brown PetscInt coneSize[15] = {4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 41054fcfd0cSMatthew G. Knepley PetscInt cones[24] = { 7, 8, 9, 6, 11, 13, 9, 14, 10, 11, 13, 9, 41154fcfd0cSMatthew G. Knepley 9, 10, 11, 7, 9, 14, 13, 12, 7, 8, 11, 9}; 412c4762a1bSJed Brown PetscInt coneOrientations[24] = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 413c4762a1bSJed Brown 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 414c4762a1bSJed Brown PetscScalar vertexCoords[27] = {-2.0, -1.0, 0.0, -2.0, 0.0, 0.0, -2.0, 0.0, 1.0, 415c4762a1bSJed Brown 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 416c4762a1bSJed Brown 2.0, -1.0, 0.0, 2.0, 0.0, 0.0, 2.0, 0.0, 1.0}; 417c4762a1bSJed Brown PetscInt faultPoints[3] = {9, 10, 11}; 418c4762a1bSJed Brown 4195f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4205f80ce2aSJacob Faibussowitsch for (p = 0; p < 3; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 421c4762a1bSJed Brown } 422c4762a1bSJed Brown break; 42398921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 424c4762a1bSJed Brown } 4255f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 4265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 4275f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 4285f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, "fault", &faultLabel)); 4315f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm)); 4325f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&hybridLabel)); 4335f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 435c4762a1bSJed Brown *dm = hdm; 436c4762a1bSJed Brown } else { 437c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 438c4762a1bSJed Brown 4395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 4405f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 4415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 4425f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 4435f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 4445f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 4455f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm)); 4465f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 4475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 448c4762a1bSJed Brown *dm = hdm; 449c4762a1bSJed Brown } 450c4762a1bSJed Brown PetscFunctionReturn(0); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown PetscErrorCode CreateTensorProduct_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 454c4762a1bSJed Brown { 455c4762a1bSJed Brown DM idm; 456c4762a1bSJed Brown PetscMPIInt rank; 457c4762a1bSJed Brown 458c4762a1bSJed Brown PetscFunctionBegin; 4595f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 460dd400576SPatrick Sanan if (rank == 0) { 461c4762a1bSJed Brown switch (testNum) { 462c4762a1bSJed Brown case 0: 463c4762a1bSJed Brown { 464c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 465c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 466c4762a1bSJed Brown PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 467c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 468c4762a1bSJed Brown 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, 469c4762a1bSJed Brown -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 470c4762a1bSJed Brown 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 471c4762a1bSJed Brown 4725f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 473c4762a1bSJed Brown } 474c4762a1bSJed Brown break; 475c4762a1bSJed Brown case 1: 476c4762a1bSJed Brown { 477c4762a1bSJed Brown PetscInt numPoints[2] = {8, 1}; 478c4762a1bSJed Brown PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 479c4762a1bSJed Brown PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8}; 480c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 481c4762a1bSJed Brown 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, 482c4762a1bSJed Brown -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 483c4762a1bSJed Brown 4845f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 485c4762a1bSJed Brown } 486c4762a1bSJed Brown break; 48798921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown } else { 490c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 491c4762a1bSJed Brown 4925f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 493c4762a1bSJed Brown } 4945f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 4955f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-in_dm_view")); 4965f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 497c4762a1bSJed Brown *dm = idm; 498c4762a1bSJed Brown PetscFunctionReturn(0); 499c4762a1bSJed Brown } 500c4762a1bSJed Brown 501c4762a1bSJed Brown PetscErrorCode CreateTensorProductHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 502c4762a1bSJed Brown { 503c4762a1bSJed Brown DM idm, hdm = NULL; 504c4762a1bSJed Brown DMLabel faultLabel; 505c4762a1bSJed Brown PetscInt p; 506c4762a1bSJed Brown PetscMPIInt rank; 507c4762a1bSJed Brown 508c4762a1bSJed Brown PetscFunctionBegin; 5095f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 510dd400576SPatrick Sanan if (rank == 0) { 511c4762a1bSJed Brown switch (testNum) { 512c4762a1bSJed Brown case 0: 513c4762a1bSJed Brown { 514c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 515c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 516c4762a1bSJed Brown PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 517c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 518c4762a1bSJed Brown 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, 519c4762a1bSJed Brown -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 520c4762a1bSJed Brown 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 521c4762a1bSJed Brown PetscInt faultPoints[4] = {2, 3, 5, 6}; 522c4762a1bSJed Brown 5235f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 5245f80ce2aSJacob Faibussowitsch for (p = 0; p < 4; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 525c4762a1bSJed Brown } 526c4762a1bSJed Brown break; 527c4762a1bSJed Brown case 1: 528c4762a1bSJed Brown { 529c4762a1bSJed Brown PetscInt numPoints[2] = {30, 7}; 530c4762a1bSJed 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}; 531c4762a1bSJed Brown PetscInt cones[56] = { 8, 21, 20, 7, 13, 12, 23, 24, 532c4762a1bSJed Brown 14, 15, 10, 9, 13, 8, 21, 24, 533c4762a1bSJed Brown 15, 16, 11, 10, 24, 21, 22, 25, 534c4762a1bSJed Brown 30, 29, 28, 21, 35, 24, 33, 34, 535c4762a1bSJed Brown 24, 21, 30, 35, 25, 36, 31, 22, 536c4762a1bSJed Brown 27, 20, 21, 28, 32, 33, 24, 23, 537c4762a1bSJed Brown 15, 24, 13, 14, 19, 18, 17, 26}; 538c4762a1bSJed 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}; 539c4762a1bSJed Brown 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, 540c4762a1bSJed Brown -2.0, -1.0, 0.0, -3.0, 0.0, 0.0, -2.0, 1.0, 0.0, -2.0, 2.0, 0.0, -2.0, -1.0, 2.0, -3.0, 0.0, 2.0, 541c4762a1bSJed Brown -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, 542c4762a1bSJed Brown 0.0, 2.0, 0.0, 0.0, 0.0, 2.0, 2.0, -2.0, -2.0, 2.0, -1.0, -2.0, 3.0, 0.0, -2.0, 2.0, 1.0, -2.0, 543c4762a1bSJed Brown 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}; 544c4762a1bSJed Brown PetscInt faultPoints[6] = {20, 21, 22, 23, 24, 25}; 545c4762a1bSJed Brown 5465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 5475f80ce2aSJacob Faibussowitsch for (p = 0; p < 6; ++p) CHKERRQ(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 548c4762a1bSJed Brown } 549c4762a1bSJed Brown break; 55098921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 551c4762a1bSJed Brown } 5525f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 5545f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 5555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 5565f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, "fault", &faultLabel)); 5585f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, faultLabel, NULL, NULL, NULL, NULL, &hdm)); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 5605f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 561c4762a1bSJed Brown *dm = hdm; 562c4762a1bSJed Brown } else { 563c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 564c4762a1bSJed Brown 5655f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 5665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 5675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) idm, "in_")); 5685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(idm, PETSC_FALSE)); 5695f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(idm)); 5705f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(idm, NULL, "-dm_view")); 5715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm)); 5725f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&idm)); 5735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 574c4762a1bSJed Brown *dm = hdm; 575c4762a1bSJed Brown } 576c4762a1bSJed Brown PetscFunctionReturn(0); 577c4762a1bSJed Brown } 578c4762a1bSJed Brown 579c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 580c4762a1bSJed Brown { 581c4762a1bSJed Brown PetscInt dim = user->dim; 582c4762a1bSJed Brown PetscBool cellHybrid = user->cellHybrid; 583c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 584c4762a1bSJed Brown PetscMPIInt rank, size; 585c4762a1bSJed Brown 586c4762a1bSJed Brown PetscFunctionBegin; 5875f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 5885f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &size)); 5895f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 5905f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 5915f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*dm, dim)); 592c4762a1bSJed Brown switch (dim) { 593c4762a1bSJed Brown case 1: 594*28b400f6SJacob Faibussowitsch PetscCheck(!cellHybrid,comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %d", dim); 5955f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplex_1D(comm, dm)); 596c4762a1bSJed Brown break; 597c4762a1bSJed Brown case 2: 598c4762a1bSJed Brown if (cellSimplex) { 599c4762a1bSJed Brown if (cellHybrid) { 6005f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplexHybrid_2D(comm, user->testNum, dm)); 601c4762a1bSJed Brown } else { 6025f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplex_2D(comm, dm)); 603c4762a1bSJed Brown } 604c4762a1bSJed Brown } else { 605c4762a1bSJed Brown if (cellHybrid) { 6065f80ce2aSJacob Faibussowitsch CHKERRQ(CreateTensorProductHybrid_2D(comm, user->testNum, dm)); 607c4762a1bSJed Brown } else { 6085f80ce2aSJacob Faibussowitsch CHKERRQ(CreateTensorProduct_2D(comm, user->testNum, dm)); 609c4762a1bSJed Brown } 610c4762a1bSJed Brown } 611c4762a1bSJed Brown break; 612c4762a1bSJed Brown case 3: 613c4762a1bSJed Brown if (cellSimplex) { 614c4762a1bSJed Brown if (cellHybrid) { 6155f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplexHybrid_3D(comm, user->testNum, dm)); 616c4762a1bSJed Brown } else { 6175f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplex_3D(comm, user->testNum, dm)); 618c4762a1bSJed Brown } 619c4762a1bSJed Brown } else { 620c4762a1bSJed Brown if (cellHybrid) { 6215f80ce2aSJacob Faibussowitsch CHKERRQ(CreateTensorProductHybrid_3D(comm, user->testNum, dm)); 622c4762a1bSJed Brown } else { 6235f80ce2aSJacob Faibussowitsch CHKERRQ(CreateTensorProduct_3D(comm, user->testNum, dm)); 624c4762a1bSJed Brown } 625c4762a1bSJed Brown } 626c4762a1bSJed Brown break; 627c4762a1bSJed Brown default: 62898921bdaSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 629c4762a1bSJed Brown } 630c4762a1bSJed Brown if (user->testPartition && size > 1) { 631c4762a1bSJed Brown PetscPartitioner part; 632c4762a1bSJed Brown PetscInt *sizes = NULL; 633c4762a1bSJed Brown PetscInt *points = NULL; 634c4762a1bSJed Brown 635dd400576SPatrick Sanan if (rank == 0) { 636c4762a1bSJed Brown if (dim == 2 && cellSimplex && !cellHybrid && size == 2) { 637c4762a1bSJed Brown switch (user->testNum) { 638c4762a1bSJed Brown case 0: { 639c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 1}; 640c4762a1bSJed Brown PetscInt triPoints_p2[2] = {0, 1}; 641c4762a1bSJed Brown 6425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 2, &points)); 6435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, triSizes_p2, 2)); 6445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, triPoints_p2, 2));break;} 645c4762a1bSJed Brown default: 64698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum); 647c4762a1bSJed Brown } 648c4762a1bSJed Brown } else if (dim == 2 && cellSimplex && cellHybrid && size == 2) { 649c4762a1bSJed Brown switch (user->testNum) { 650c4762a1bSJed Brown case 0: { 651c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 2}; 652c4762a1bSJed Brown PetscInt triPoints_p2[3] = {0, 1, 2}; 653c4762a1bSJed Brown 6545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 3, &points)); 6555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, triSizes_p2, 2)); 6565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, triPoints_p2, 3));break;} 657c4762a1bSJed Brown default: 65898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular hybrid mesh on 2 procs", user->testNum); 659c4762a1bSJed Brown } 660c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && !cellHybrid && size == 2) { 661c4762a1bSJed Brown switch (user->testNum) { 662c4762a1bSJed Brown case 0: { 663c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 1}; 664c4762a1bSJed Brown PetscInt quadPoints_p2[2] = {0, 1}; 665c4762a1bSJed Brown 6665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 2, &points)); 6675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, quadSizes_p2, 2)); 6685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, quadPoints_p2, 2));break;} 669c4762a1bSJed Brown default: 67098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for quadrilateral mesh on 2 procs", user->testNum); 671c4762a1bSJed Brown } 672c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && cellHybrid && size == 2) { 673c4762a1bSJed Brown switch (user->testNum) { 674c4762a1bSJed Brown case 0: { 675c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 2}; 676c4762a1bSJed Brown PetscInt quadPoints_p2[3] = {0, 1, 2}; 677c4762a1bSJed Brown 6785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 3, &points)); 6795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, quadSizes_p2, 2)); 6805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, quadPoints_p2, 3));break;} 681c4762a1bSJed Brown default: 68298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for quadrilateral hybrid mesh on 2 procs", user->testNum); 683c4762a1bSJed Brown } 684c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && !cellHybrid && size == 2) { 685c4762a1bSJed Brown switch (user->testNum) { 686c4762a1bSJed Brown case 0: { 687c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 1}; 688c4762a1bSJed Brown PetscInt tetPoints_p2[2] = {0, 1}; 689c4762a1bSJed Brown 6905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 2, &points)); 6915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, tetSizes_p2, 2)); 6925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, tetPoints_p2, 2));break;} 693c4762a1bSJed Brown case 1: { 694c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 1}; 695c4762a1bSJed Brown PetscInt tetPoints_p2[2] = {0, 1}; 696c4762a1bSJed Brown 6975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 2, &points)); 6985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, tetSizes_p2, 2)); 6995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, tetPoints_p2, 2));break;} 700c4762a1bSJed Brown default: 70198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for tetrahedral mesh on 2 procs", user->testNum); 702c4762a1bSJed Brown } 703c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && cellHybrid && size == 2) { 704c4762a1bSJed Brown switch (user->testNum) { 705c4762a1bSJed Brown case 0: { 706c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 2}; 707c4762a1bSJed Brown PetscInt tetPoints_p2[3] = {0, 1, 2}; 708c4762a1bSJed Brown 7095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 3, &points)); 7105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, tetSizes_p2, 2)); 7115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, tetPoints_p2, 3));break;} 712c4762a1bSJed Brown case 1: { 713c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {3, 4}; 714c4762a1bSJed Brown PetscInt tetPoints_p2[7] = {0, 3, 5, 1, 2, 4, 6}; 715c4762a1bSJed Brown 7165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 7, &points)); 7175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, tetSizes_p2, 2)); 7185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, tetPoints_p2, 7));break;} 719c4762a1bSJed Brown default: 72098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for tetrahedral hybrid mesh on 2 procs", user->testNum); 721c4762a1bSJed Brown } 722c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && !cellHybrid && size == 2) { 723c4762a1bSJed Brown switch (user->testNum) { 724c4762a1bSJed Brown case 0: { 725c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 1}; 726c4762a1bSJed Brown PetscInt hexPoints_p2[2] = {0, 1}; 727c4762a1bSJed Brown 7285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 2, &points)); 7295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, hexSizes_p2, 2)); 7305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, hexPoints_p2, 2));break;} 731c4762a1bSJed Brown default: 73298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for hexahedral mesh on 2 procs", user->testNum); 733c4762a1bSJed Brown } 734c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && cellHybrid && size == 2) { 735c4762a1bSJed Brown switch (user->testNum) { 736c4762a1bSJed Brown case 0: { 737c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 1}; 738c4762a1bSJed Brown PetscInt hexPoints_p2[2] = {0, 1}; 739c4762a1bSJed Brown 7405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 2, &points)); 7415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, hexSizes_p2, 2)); 7425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, hexPoints_p2, 2));break;} 743c4762a1bSJed Brown case 1: { 744c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {5, 4}; 745c4762a1bSJed Brown PetscInt hexPoints_p2[9] = {3, 4, 5, 7, 8, 0, 1, 2, 6}; 746c4762a1bSJed Brown 7475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(2, &sizes, 9, &points)); 7485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(sizes, hexSizes_p2, 2)); 7495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(points, hexPoints_p2, 9));break;} 750c4762a1bSJed Brown default: 75198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for hexahedral hybrid mesh on 2 procs", user->testNum); 752c4762a1bSJed Brown } 753c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition"); 754c4762a1bSJed Brown } 7555f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(*dm, &part)); 7565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 7575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerShellSetPartition(part, size, sizes, points)); 7585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(sizes, points)); 759c4762a1bSJed Brown } else { 760c4762a1bSJed Brown PetscPartitioner part; 761c4762a1bSJed Brown 7625f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(*dm,&part)); 7635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 764c4762a1bSJed Brown } 765c4762a1bSJed Brown { 766c4762a1bSJed Brown DM pdm = NULL; 767c4762a1bSJed Brown 7685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(*dm, 0, NULL, &pdm)); 769c4762a1bSJed Brown if (pdm) { 7705f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(pdm, NULL, "-dm_view")); 7715f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 772c4762a1bSJed Brown *dm = pdm; 773c4762a1bSJed Brown } 774c4762a1bSJed Brown } 7755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 7765f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view_pre")); 7775f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 778c4762a1bSJed Brown if (user->uninterpolate || user->reinterpolate) { 779c4762a1bSJed Brown DM udm = NULL; 780c4762a1bSJed Brown 7815f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexUninterpolate(*dm, &udm)); 7825f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCopyCoordinates(*dm, udm)); 7835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 784c4762a1bSJed Brown *dm = udm; 785c4762a1bSJed Brown } 786c4762a1bSJed Brown if (user->reinterpolate) { 787c4762a1bSJed Brown DM idm = NULL; 788c4762a1bSJed Brown 7895f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 7905f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCopyCoordinates(*dm, idm)); 7915f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 792c4762a1bSJed Brown *dm = idm; 793c4762a1bSJed Brown } 7945f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 7955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Hybrid Mesh")); 7965f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 7975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "hyb_")); 7985f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 799c4762a1bSJed Brown PetscFunctionReturn(0); 800c4762a1bSJed Brown } 801c4762a1bSJed Brown 802c4762a1bSJed Brown int main(int argc, char **argv) 803c4762a1bSJed Brown { 804c4762a1bSJed Brown DM dm; 805c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 806c4762a1bSJed Brown PetscErrorCode ierr; 807c4762a1bSJed Brown 808c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 8095f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 8105f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 8115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 812c4762a1bSJed Brown ierr = PetscFinalize(); 813c4762a1bSJed Brown return ierr; 814c4762a1bSJed Brown } 815c4762a1bSJed Brown 816c4762a1bSJed Brown /*TEST 817c4762a1bSJed Brown 818c4762a1bSJed Brown # 1D Simplex 29-31 819c4762a1bSJed Brown testset: 82054fcfd0cSMatthew G. Knepley args: -dim 1 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 821c4762a1bSJed Brown test: 822c4762a1bSJed Brown suffix: 29 823c4762a1bSJed Brown test: 824c4762a1bSJed Brown suffix: 30 825c4762a1bSJed Brown args: -dm_refine 1 826c4762a1bSJed Brown test: 827c4762a1bSJed Brown suffix: 31 828c4762a1bSJed Brown args: -dm_refine 5 829c4762a1bSJed Brown 830c4762a1bSJed Brown # 2D Simplex 0-3 831c4762a1bSJed Brown testset: 83254fcfd0cSMatthew G. Knepley args: -dim 2 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 833c4762a1bSJed Brown test: 834c4762a1bSJed Brown suffix: 0 835c4762a1bSJed Brown args: -hyb_dm_plex_check_faces 836c4762a1bSJed Brown test: 837c4762a1bSJed Brown suffix: 1 838c4762a1bSJed Brown args: -dm_refine 1 -hyb_dm_plex_check_faces 839c4762a1bSJed Brown test: 840c4762a1bSJed Brown suffix: 2 841c4762a1bSJed Brown nsize: 2 842c4762a1bSJed Brown args: -hyb_dm_plex_check_faces 843c4762a1bSJed Brown test: 844c4762a1bSJed Brown suffix: 3 845c4762a1bSJed Brown nsize: 2 846c4762a1bSJed Brown args: -dm_refine 1 -hyb_dm_plex_check_faces 847c4762a1bSJed Brown test: 848c4762a1bSJed Brown suffix: 32 849c4762a1bSJed Brown args: -dm_refine 1 -uninterpolate 850c4762a1bSJed Brown test: 851c4762a1bSJed Brown suffix: 33 852c4762a1bSJed Brown nsize: 2 853c4762a1bSJed Brown args: -dm_refine 1 -uninterpolate 854c4762a1bSJed Brown test: 855c4762a1bSJed Brown suffix: 34 856c4762a1bSJed Brown nsize: 2 857c4762a1bSJed Brown args: -dm_refine 3 -uninterpolate 858c4762a1bSJed Brown 859c4762a1bSJed Brown # 2D Hybrid Simplex 4-7 860c4762a1bSJed Brown testset: 86154fcfd0cSMatthew G. Knepley args: -dim 2 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 862c4762a1bSJed Brown test: 863c4762a1bSJed Brown suffix: 4 864c4762a1bSJed Brown test: 865c4762a1bSJed Brown suffix: 5 866c4762a1bSJed Brown args: -dm_refine 1 867c4762a1bSJed Brown test: 868c4762a1bSJed Brown suffix: 6 869c4762a1bSJed Brown nsize: 2 870c4762a1bSJed Brown test: 871c4762a1bSJed Brown suffix: 7 872c4762a1bSJed Brown nsize: 2 873c4762a1bSJed Brown args: -dm_refine 1 874c4762a1bSJed Brown test: 875c4762a1bSJed Brown suffix: 24 876c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 877c4762a1bSJed Brown 878c4762a1bSJed Brown # 2D Quad 12-13 879c4762a1bSJed Brown testset: 88054fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 881c4762a1bSJed Brown test: 882c4762a1bSJed Brown suffix: 12 883c4762a1bSJed Brown args: -dm_refine 1 884c4762a1bSJed Brown test: 885c4762a1bSJed Brown suffix: 13 886c4762a1bSJed Brown nsize: 2 887c4762a1bSJed Brown args: -dm_refine 1 888c4762a1bSJed Brown 889c4762a1bSJed Brown # 2D Hybrid Quad 27-28 890c4762a1bSJed Brown testset: 89154fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 892c4762a1bSJed Brown test: 893c4762a1bSJed Brown suffix: 27 894c4762a1bSJed Brown args: -dm_refine 1 895c4762a1bSJed Brown test: 896c4762a1bSJed Brown suffix: 28 897c4762a1bSJed Brown nsize: 2 898c4762a1bSJed Brown args: -dm_refine 1 899c4762a1bSJed Brown 900c4762a1bSJed Brown # 3D Simplex 8-11 901c4762a1bSJed Brown testset: 90254fcfd0cSMatthew G. Knepley args: -dim 3 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 903c4762a1bSJed Brown test: 904c4762a1bSJed Brown suffix: 8 905c4762a1bSJed Brown args: -dm_refine 1 906c4762a1bSJed Brown test: 907c4762a1bSJed Brown suffix: 9 908c4762a1bSJed Brown nsize: 2 909c4762a1bSJed Brown args: -dm_refine 1 910c4762a1bSJed Brown test: 911c4762a1bSJed Brown suffix: 10 912c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 913c4762a1bSJed Brown test: 914c4762a1bSJed Brown suffix: 11 915c4762a1bSJed Brown nsize: 2 916c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 917c4762a1bSJed Brown test: 918c4762a1bSJed Brown suffix: 25 919c4762a1bSJed Brown args: -test_num 2 -dm_refine 2 920c4762a1bSJed Brown 921c4762a1bSJed Brown # 3D Hybrid Simplex 16-19 922c4762a1bSJed Brown testset: 92354fcfd0cSMatthew G. Knepley args: -dim 3 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 924c4762a1bSJed Brown test: 925c4762a1bSJed Brown suffix: 16 926c4762a1bSJed Brown args: -dm_refine 1 927c4762a1bSJed Brown test: 928c4762a1bSJed Brown suffix: 17 929c4762a1bSJed Brown nsize: 2 930c4762a1bSJed Brown args: -dm_refine 1 931c4762a1bSJed Brown test: 932c4762a1bSJed Brown suffix: 18 933c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 934c4762a1bSJed Brown test: 935c4762a1bSJed Brown suffix: 19 936c4762a1bSJed Brown nsize: 2 937c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 938c4762a1bSJed Brown 939c4762a1bSJed Brown # 3D Hex 14-15 940c4762a1bSJed Brown testset: 94154fcfd0cSMatthew G. Knepley args: -dim 3 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 942c4762a1bSJed Brown test: 943c4762a1bSJed Brown suffix: 14 944c4762a1bSJed Brown args: -dm_refine 1 945c4762a1bSJed Brown test: 946c4762a1bSJed Brown suffix: 15 947c4762a1bSJed Brown nsize: 2 948c4762a1bSJed Brown args: -dm_refine 1 949c4762a1bSJed Brown test: 950c4762a1bSJed Brown suffix: 26 951c4762a1bSJed Brown args: -test_num 1 -dm_refine 2 952c4762a1bSJed Brown 953c4762a1bSJed Brown # 3D Hybrid Hex 20-23 954c4762a1bSJed Brown testset: 95554fcfd0cSMatthew G. Knepley args: -dim 3 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 956c4762a1bSJed Brown test: 957c4762a1bSJed Brown suffix: 20 958c4762a1bSJed Brown args: -dm_refine 1 959c4762a1bSJed Brown test: 960c4762a1bSJed Brown suffix: 21 961c4762a1bSJed Brown nsize: 2 962c4762a1bSJed Brown args: -dm_refine 1 963c4762a1bSJed Brown test: 964c4762a1bSJed Brown suffix: 22 965c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 966c4762a1bSJed Brown test: 967c4762a1bSJed Brown suffix: 23 968c4762a1bSJed Brown nsize: 2 969c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 970c4762a1bSJed Brown 971c4762a1bSJed Brown # Hybrid interpolation 97254fcfd0cSMatthew G. Knepley # TODO Setup new tests (like -reinterpolate) that interpolate hybrid cells 973c4762a1bSJed Brown testset: 974c4762a1bSJed Brown nsize: 2 97554fcfd0cSMatthew 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 976c4762a1bSJed Brown test: 977c4762a1bSJed Brown suffix: hybint_2d_0 978c4762a1bSJed Brown args: -dim 2 -dm_refine 2 979c4762a1bSJed Brown test: 980c4762a1bSJed Brown suffix: hybint_2d_1 981c4762a1bSJed Brown args: -dim 2 -dm_refine 2 -test_num 1 982c4762a1bSJed Brown test: 983c4762a1bSJed Brown suffix: hybint_3d_0 984c4762a1bSJed Brown args: -dim 3 -dm_refine 1 985c4762a1bSJed Brown test: 986c4762a1bSJed Brown suffix: hybint_3d_1 987c4762a1bSJed Brown args: -dim 3 -dm_refine 1 -test_num 1 988c4762a1bSJed Brown 989c4762a1bSJed Brown TEST*/ 990