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); 31c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex4.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = PetscOptionsBool("-cell_hybrid", "Use a hybrid mesh", "ex4.c", options->cellHybrid, &options->cellHybrid, NULL);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex4.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex4.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex4.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex4.c", options->uninterpolate, &options->uninterpolate, NULL);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscOptionsBool("-reinterpolate", "Reinterpolate the mesh at the end", "ex4.c", options->reinterpolate, &options->reinterpolate, NULL);CHKERRQ(ierr); 39*1e1ea65dSPierre 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 PetscErrorCode ierr; 57c4762a1bSJed Brown 58c4762a1bSJed Brown PetscFunctionBegin; 59ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 60c4762a1bSJed Brown if (!rank) { 61c4762a1bSJed Brown PetscInt numPoints[2] = {3, 2}; 62c4762a1bSJed Brown PetscInt coneSize[5] = {2, 2, 0, 0, 0}; 63c4762a1bSJed Brown PetscInt cones[4] = {2, 3, 3, 4}; 64c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0}; 65c4762a1bSJed Brown PetscScalar vertexCoords[3] = {-1.0, 0.0, 1.0}; 66c4762a1bSJed Brown 67c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 68c4762a1bSJed Brown } else { 69c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 70c4762a1bSJed Brown 71c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown PetscFunctionReturn(0); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Two triangles 77c4762a1bSJed Brown 4 78c4762a1bSJed Brown / | \ 79c4762a1bSJed Brown 8 | 10 80c4762a1bSJed Brown / | \ 81c4762a1bSJed Brown 2 0 7 1 5 82c4762a1bSJed Brown \ | / 83c4762a1bSJed Brown 6 | 9 84c4762a1bSJed Brown \ | / 85c4762a1bSJed Brown 3 86c4762a1bSJed Brown 87c4762a1bSJed Brown Becomes 88c4762a1bSJed Brown 10 89c4762a1bSJed Brown / | \ 90c4762a1bSJed Brown 21 | 26 91c4762a1bSJed Brown / | \ 92c4762a1bSJed Brown 14 2 20 4 16 93c4762a1bSJed Brown /|\ | /|\ 94c4762a1bSJed Brown 22 | 28 | 32 | 25 95c4762a1bSJed Brown / | \ | / | 6\ 96c4762a1bSJed Brown 8 29 3 13 7 31 11 97c4762a1bSJed Brown \0 | / | \ | / 98c4762a1bSJed Brown 17 | 27 | 30 | 24 99c4762a1bSJed Brown \|/ | \|/ 100c4762a1bSJed Brown 12 1 19 5 15 101c4762a1bSJed Brown \ | / 102c4762a1bSJed Brown 18 | 23 103c4762a1bSJed Brown \ | / 104c4762a1bSJed Brown 9 105c4762a1bSJed Brown */ 106c4762a1bSJed Brown PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *dm) 107c4762a1bSJed Brown { 108c4762a1bSJed Brown PetscInt depth = 2; 109c4762a1bSJed Brown PetscMPIInt rank; 110c4762a1bSJed Brown PetscErrorCode ierr; 111c4762a1bSJed Brown 112c4762a1bSJed Brown PetscFunctionBegin; 113ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 114c4762a1bSJed Brown if (!rank) { 115c4762a1bSJed Brown PetscInt numPoints[3] = {4, 5, 2}; 116c4762a1bSJed Brown PetscInt coneSize[11] = {3, 3, 0, 0, 0, 0, 2, 2, 2, 2, 2}; 117c4762a1bSJed Brown PetscInt cones[16] = {6, 7, 8, 7, 9, 10, 2, 3, 3, 4, 4, 2, 3, 5, 5, 4}; 118c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 119c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-0.5, 0.0, 0.0, -0.5, 0.0, 0.5, 0.5, 0.0}; 120c4762a1bSJed Brown 121c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 122c4762a1bSJed Brown } else { 123c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 124c4762a1bSJed Brown 125c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown PetscFunctionReturn(0); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Two triangles separated by a zero-volume cell with 4 vertices/2 edges 131c4762a1bSJed Brown 5--16--8 132c4762a1bSJed Brown / | | \ 133c4762a1bSJed Brown 11 | | 12 134c4762a1bSJed Brown / | | \ 135c4762a1bSJed Brown 3 0 10 2 14 1 6 136c4762a1bSJed Brown \ | | / 137c4762a1bSJed Brown 9 | | 13 138c4762a1bSJed Brown \ | | / 139c4762a1bSJed Brown 4--15--7 140c4762a1bSJed Brown */ 141c4762a1bSJed Brown PetscErrorCode CreateSimplexHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 142c4762a1bSJed Brown { 143c4762a1bSJed Brown DM idm, hdm = NULL; 144c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 145c4762a1bSJed Brown PetscInt p; 146c4762a1bSJed Brown PetscMPIInt rank; 147c4762a1bSJed Brown PetscErrorCode ierr; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBegin; 150ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 151c4762a1bSJed Brown if (!rank) { 152c4762a1bSJed Brown switch (testNum) { 153c4762a1bSJed Brown case 0: 154c4762a1bSJed Brown { 155c4762a1bSJed Brown PetscInt numPoints[2] = {4, 2}; 156c4762a1bSJed Brown PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 157c4762a1bSJed Brown PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 158c4762a1bSJed Brown PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 159c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, 1.0, 0.5}; 160c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 161c4762a1bSJed Brown 162c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 163c4762a1bSJed Brown for (p = 0; p < 2; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 164c4762a1bSJed Brown } 165c4762a1bSJed Brown break; 166c4762a1bSJed Brown case 1: 167c4762a1bSJed Brown { 168c4762a1bSJed Brown PetscInt numPoints[2] = {5, 4}; 169c4762a1bSJed Brown PetscInt coneSize[9] = {3, 3, 3, 3, 0, 0, 0, 0, 0}; 170c4762a1bSJed Brown PetscInt cones[12] = {4, 5, 6, 6, 7, 4, 6, 5, 8, 6, 8, 7}; 171c4762a1bSJed Brown PetscInt coneOrientations[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 172c4762a1bSJed Brown PetscScalar vertexCoords[10] = {-1.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0}; 173c4762a1bSJed Brown PetscInt faultPoints[3] = {5, 6, 7}; 174c4762a1bSJed Brown 175c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 176c4762a1bSJed Brown for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 177c4762a1bSJed Brown } 178c4762a1bSJed Brown break; 179c4762a1bSJed Brown default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = DMLabelDestroy(&hybridLabel);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 190c4762a1bSJed Brown *dm = hdm; 191c4762a1bSJed Brown } else { 192c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 193c4762a1bSJed Brown 194c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 19654fcfd0cSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 19754fcfd0cSMatthew G. Knepley ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 19854fcfd0cSMatthew G. Knepley ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 202c4762a1bSJed Brown *dm = hdm; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown PetscFunctionReturn(0); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Two quadrilaterals 208c4762a1bSJed Brown 209c4762a1bSJed Brown 5----10-----4----14-----7 210c4762a1bSJed Brown | | | 211c4762a1bSJed Brown | | | 212c4762a1bSJed Brown | | | 213c4762a1bSJed Brown 11 0 9 1 13 214c4762a1bSJed Brown | | | 215c4762a1bSJed Brown | | | 216c4762a1bSJed Brown | | | 217c4762a1bSJed Brown 2-----8-----3----12-----6 218c4762a1bSJed Brown */ 219c4762a1bSJed Brown PetscErrorCode CreateTensorProduct_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown PetscInt depth = 2; 222c4762a1bSJed Brown PetscMPIInt rank; 223c4762a1bSJed Brown PetscErrorCode ierr; 224c4762a1bSJed Brown 225c4762a1bSJed Brown PetscFunctionBegin; 226ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 227c4762a1bSJed Brown if (!rank) { 228c4762a1bSJed Brown PetscInt numPoints[3] = {6, 7, 2}; 229c4762a1bSJed Brown PetscInt coneSize[15] = {4, 4, 0, 0, 0, 0, 0, 0, 2, 2, 2, 2, 2, 2, 2}; 230c4762a1bSJed 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}; 231c4762a1bSJed Brown PetscInt coneOrientations[22] = {0, 0, 0, 0, 0, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 232c4762a1bSJed 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}; 233c4762a1bSJed Brown 234c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 235c4762a1bSJed Brown } else { 236c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 237c4762a1bSJed Brown 238c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown PetscFunctionReturn(0); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 243c4762a1bSJed Brown PetscErrorCode CreateTensorProductHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 244c4762a1bSJed Brown { 245c4762a1bSJed Brown DM idm, hdm = NULL; 246c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 247c4762a1bSJed Brown PetscInt p; 248c4762a1bSJed Brown PetscMPIInt rank; 249c4762a1bSJed Brown PetscErrorCode ierr; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBegin; 252ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 253c4762a1bSJed Brown if (!rank) { 254c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 255c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 256c4762a1bSJed Brown PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4,}; 257c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 258c4762a1bSJed 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}; 259c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 260c4762a1bSJed Brown 261c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 262c4762a1bSJed Brown for (p = 0; p < 2; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 263c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 264c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 265c4762a1bSJed Brown ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr); 268c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm);CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = DMLabelDestroy(&hybridLabel);CHKERRQ(ierr); 270c4762a1bSJed Brown } else { 271c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 272c4762a1bSJed Brown 273c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 27554fcfd0cSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 27654fcfd0cSMatthew G. Knepley ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 27754fcfd0cSMatthew G. Knepley ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 281c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 282c4762a1bSJed Brown *dm = hdm; 283c4762a1bSJed Brown PetscFunctionReturn(0); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* Two tetrahedrons 287c4762a1bSJed Brown 288c4762a1bSJed Brown cell 5 5______ cell 289c4762a1bSJed Brown 0 / | \ |\ \ 1 290c4762a1bSJed Brown 17 | 18 | 18 13 21 291c4762a1bSJed Brown /8 19 10\ 19 \ \ 292c4762a1bSJed Brown 2-14-|----4 | 4--22--6 293c4762a1bSJed Brown \ 9 | 7 / |10 / / 294c4762a1bSJed Brown 16 | 15 | 15 12 20 295c4762a1bSJed Brown \ | / |/ / 296c4762a1bSJed Brown 3 3------ 297c4762a1bSJed Brown */ 298c4762a1bSJed Brown PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 299c4762a1bSJed Brown { 300c4762a1bSJed Brown DM idm; 301c4762a1bSJed Brown PetscInt depth = 3; 302c4762a1bSJed Brown PetscMPIInt rank; 303c4762a1bSJed Brown PetscErrorCode ierr; 304c4762a1bSJed Brown 305c4762a1bSJed Brown PetscFunctionBegin; 306ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 307c4762a1bSJed Brown if (!rank) { 308c4762a1bSJed Brown switch (testNum) { 309c4762a1bSJed Brown case 0: 310c4762a1bSJed Brown { 311c4762a1bSJed Brown PetscInt numPoints[4] = {5, 9, 7, 2}; 312c4762a1bSJed 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}; 313c4762a1bSJed 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}; 314c4762a1bSJed Brown PetscInt coneOrientations[47] = { 0, 0, 0, 0, -3, 0, 0, 0, 0, 0, 0, 0, 0, -2, -2, 0, -2, -2, -2, -2, 0, 0, -2, -2, 0, -2, -2, -2, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 315c4762a1bSJed 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}; 316c4762a1bSJed Brown 317c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown break; 320c4762a1bSJed Brown case 1: 321c4762a1bSJed Brown { 322c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 323c4762a1bSJed Brown PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 324c4762a1bSJed Brown PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 325c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 326c4762a1bSJed 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}; 327c4762a1bSJed Brown 328c4762a1bSJed Brown depth = 1; 329c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 330c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 331c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr); 332c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 333c4762a1bSJed Brown *dm = idm; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown break; 336c4762a1bSJed Brown case 2: 337c4762a1bSJed Brown { 338c4762a1bSJed Brown PetscInt numPoints[2] = {4, 1}; 339c4762a1bSJed Brown PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 340c4762a1bSJed Brown PetscInt cones[4] = {2, 3, 4, 1}; 341c4762a1bSJed Brown PetscInt coneOrientations[4] = {0, 0, 0, 0}; 342c4762a1bSJed 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}; 343c4762a1bSJed Brown 344c4762a1bSJed Brown depth = 1; 345c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 347c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr); 348c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 349c4762a1bSJed Brown *dm = idm; 350c4762a1bSJed Brown } 351c4762a1bSJed Brown break; 352c4762a1bSJed Brown default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown } else { 355c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 356c4762a1bSJed Brown 357c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 358c4762a1bSJed Brown switch (testNum) { 359c4762a1bSJed Brown case 1: 360c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 361c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr); 362c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 363c4762a1bSJed Brown *dm = idm; 364c4762a1bSJed Brown break; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown } 367c4762a1bSJed Brown PetscFunctionReturn(0); 368c4762a1bSJed Brown } 369c4762a1bSJed Brown 370c4762a1bSJed Brown /* Two tetrahedrons separated by a zero-volume cell with 6 vertices 371c4762a1bSJed Brown 372c4762a1bSJed Brown cell 6 ___33___10______ cell 373c4762a1bSJed Brown 0 / | \ |\ \ 1 374c4762a1bSJed Brown 21 | 23 | 29 27 375c4762a1bSJed Brown /12 24 14\ 30 \ \ 376c4762a1bSJed Brown 3-20-|----5--32-|---9--26--7 377c4762a1bSJed Brown \ 13| 11/ |18 / / 378c4762a1bSJed Brown 19 | 22 | 28 25 379c4762a1bSJed Brown \ | / |/ / 380c4762a1bSJed Brown 4----31----8------ 381c4762a1bSJed Brown cell 2 382c4762a1bSJed Brown */ 383c4762a1bSJed Brown PetscErrorCode CreateSimplexHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 384c4762a1bSJed Brown { 385c4762a1bSJed Brown DM idm, hdm = NULL; 386c4762a1bSJed Brown DMLabel faultLabel, hybridLabel; 387c4762a1bSJed Brown PetscInt p; 388c4762a1bSJed Brown PetscMPIInt rank; 389c4762a1bSJed Brown PetscErrorCode ierr; 390c4762a1bSJed Brown 391c4762a1bSJed Brown PetscFunctionBegin; 392ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 393c4762a1bSJed Brown if (!rank) { 394c4762a1bSJed Brown switch (testNum) { 395c4762a1bSJed Brown case 0: 396c4762a1bSJed Brown { 397c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 398c4762a1bSJed Brown PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 399c4762a1bSJed Brown PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 400c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 401c4762a1bSJed 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}; 402c4762a1bSJed Brown PetscInt faultPoints[3] = {3, 4, 5}; 403c4762a1bSJed Brown 404c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 405c4762a1bSJed Brown for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 406c4762a1bSJed Brown } 407c4762a1bSJed Brown break; 408c4762a1bSJed Brown case 1: 409c4762a1bSJed Brown { 410c4762a1bSJed Brown /* Tets 0,3,5 and 1,2,4 */ 411c4762a1bSJed Brown PetscInt numPoints[2] = {9, 6}; 412c4762a1bSJed Brown PetscInt coneSize[15] = {4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 41354fcfd0cSMatthew G. Knepley PetscInt cones[24] = { 7, 8, 9, 6, 11, 13, 9, 14, 10, 11, 13, 9, 41454fcfd0cSMatthew G. Knepley 9, 10, 11, 7, 9, 14, 13, 12, 7, 8, 11, 9}; 415c4762a1bSJed Brown PetscInt coneOrientations[24] = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 416c4762a1bSJed Brown 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 417c4762a1bSJed Brown PetscScalar vertexCoords[27] = {-2.0, -1.0, 0.0, -2.0, 0.0, 0.0, -2.0, 0.0, 1.0, 418c4762a1bSJed Brown 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 419c4762a1bSJed Brown 2.0, -1.0, 0.0, 2.0, 0.0, 0.0, 2.0, 0.0, 1.0}; 420c4762a1bSJed Brown PetscInt faultPoints[3] = {9, 10, 11}; 421c4762a1bSJed Brown 422c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 423c4762a1bSJed Brown for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 424c4762a1bSJed Brown } 425c4762a1bSJed Brown break; 426c4762a1bSJed Brown default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 427c4762a1bSJed Brown } 428c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 429c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 430c4762a1bSJed Brown ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 431c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 432c4762a1bSJed Brown ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr); 433c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm);CHKERRQ(ierr); 434c4762a1bSJed Brown ierr = DMLabelDestroy(&hybridLabel);CHKERRQ(ierr); 435c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 436c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 437c4762a1bSJed Brown *dm = hdm; 438c4762a1bSJed Brown } else { 439c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 440c4762a1bSJed Brown 441c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 442c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 44354fcfd0cSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 44454fcfd0cSMatthew G. Knepley ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 44554fcfd0cSMatthew G. Knepley ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr); 447c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 448c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 449c4762a1bSJed Brown *dm = hdm; 450c4762a1bSJed Brown } 451c4762a1bSJed Brown PetscFunctionReturn(0); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454c4762a1bSJed Brown PetscErrorCode CreateTensorProduct_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 455c4762a1bSJed Brown { 456c4762a1bSJed Brown DM idm; 457c4762a1bSJed Brown PetscMPIInt rank; 458c4762a1bSJed Brown PetscErrorCode ierr; 459c4762a1bSJed Brown 460c4762a1bSJed Brown PetscFunctionBegin; 461ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 462c4762a1bSJed Brown if (!rank) { 463c4762a1bSJed Brown switch (testNum) { 464c4762a1bSJed Brown case 0: 465c4762a1bSJed Brown { 466c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 467c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 468c4762a1bSJed Brown PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 469c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 470c4762a1bSJed 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, 471c4762a1bSJed 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, 472c4762a1bSJed 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}; 473c4762a1bSJed Brown 474c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 475c4762a1bSJed Brown } 476c4762a1bSJed Brown break; 477c4762a1bSJed Brown case 1: 478c4762a1bSJed Brown { 479c4762a1bSJed Brown PetscInt numPoints[2] = {8, 1}; 480c4762a1bSJed Brown PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 481c4762a1bSJed Brown PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8}; 482c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 483c4762a1bSJed 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, 484c4762a1bSJed 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}; 485c4762a1bSJed Brown 486c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 487c4762a1bSJed Brown } 488c4762a1bSJed Brown break; 489c4762a1bSJed Brown default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown } else { 492c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 493c4762a1bSJed Brown 494c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 495c4762a1bSJed Brown } 496c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr); 498c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 499c4762a1bSJed Brown *dm = idm; 500c4762a1bSJed Brown PetscFunctionReturn(0); 501c4762a1bSJed Brown } 502c4762a1bSJed Brown 503c4762a1bSJed Brown PetscErrorCode CreateTensorProductHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 504c4762a1bSJed Brown { 505c4762a1bSJed Brown DM idm, hdm = NULL; 506c4762a1bSJed Brown DMLabel faultLabel; 507c4762a1bSJed Brown PetscInt p; 508c4762a1bSJed Brown PetscMPIInt rank; 509c4762a1bSJed Brown PetscErrorCode ierr; 510c4762a1bSJed Brown 511c4762a1bSJed Brown PetscFunctionBegin; 512ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 513c4762a1bSJed Brown if (!rank) { 514c4762a1bSJed Brown switch (testNum) { 515c4762a1bSJed Brown case 0: 516c4762a1bSJed Brown { 517c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 518c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 519c4762a1bSJed Brown PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 520c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 521c4762a1bSJed 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, 522c4762a1bSJed 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, 523c4762a1bSJed 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}; 524c4762a1bSJed Brown PetscInt faultPoints[4] = {2, 3, 5, 6}; 525c4762a1bSJed Brown 526c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 527c4762a1bSJed Brown for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 528c4762a1bSJed Brown } 529c4762a1bSJed Brown break; 530c4762a1bSJed Brown case 1: 531c4762a1bSJed Brown { 532c4762a1bSJed Brown PetscInt numPoints[2] = {30, 7}; 533c4762a1bSJed 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}; 534c4762a1bSJed Brown PetscInt cones[56] = { 8, 21, 20, 7, 13, 12, 23, 24, 535c4762a1bSJed Brown 14, 15, 10, 9, 13, 8, 21, 24, 536c4762a1bSJed Brown 15, 16, 11, 10, 24, 21, 22, 25, 537c4762a1bSJed Brown 30, 29, 28, 21, 35, 24, 33, 34, 538c4762a1bSJed Brown 24, 21, 30, 35, 25, 36, 31, 22, 539c4762a1bSJed Brown 27, 20, 21, 28, 32, 33, 24, 23, 540c4762a1bSJed Brown 15, 24, 13, 14, 19, 18, 17, 26}; 541c4762a1bSJed 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}; 542c4762a1bSJed 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, 543c4762a1bSJed 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, 544c4762a1bSJed 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, 545c4762a1bSJed 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, 546c4762a1bSJed 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}; 547c4762a1bSJed Brown PetscInt faultPoints[6] = {20, 21, 22, 23, 24, 25}; 548c4762a1bSJed Brown 549c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 550c4762a1bSJed Brown for (p = 0; p < 6; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);} 551c4762a1bSJed Brown } 552c4762a1bSJed Brown break; 553c4762a1bSJed Brown default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum); 554c4762a1bSJed Brown } 555c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 556c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 557c4762a1bSJed Brown ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 558c4762a1bSJed Brown ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 559c4762a1bSJed Brown ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr); 560c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr); 561c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 562c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 563c4762a1bSJed Brown *dm = hdm; 564c4762a1bSJed Brown } else { 565c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 566c4762a1bSJed Brown 567c4762a1bSJed Brown ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 568c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 56954fcfd0cSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr); 57054fcfd0cSMatthew G. Knepley ierr = DMSetFromOptions(idm);CHKERRQ(ierr); 57154fcfd0cSMatthew G. Knepley ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr); 572c4762a1bSJed Brown ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr); 573c4762a1bSJed Brown ierr = DMDestroy(&idm);CHKERRQ(ierr); 574c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 575c4762a1bSJed Brown *dm = hdm; 576c4762a1bSJed Brown } 577c4762a1bSJed Brown PetscFunctionReturn(0); 578c4762a1bSJed Brown } 579c4762a1bSJed Brown 580c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 581c4762a1bSJed Brown { 582c4762a1bSJed Brown PetscInt dim = user->dim; 583c4762a1bSJed Brown PetscBool cellHybrid = user->cellHybrid; 584c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 585c4762a1bSJed Brown PetscMPIInt rank, size; 586c4762a1bSJed Brown PetscErrorCode ierr; 587c4762a1bSJed Brown 588c4762a1bSJed Brown PetscFunctionBegin; 589ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 590ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 591c4762a1bSJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 592c4762a1bSJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 593c4762a1bSJed Brown ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 594c4762a1bSJed Brown switch (dim) { 595c4762a1bSJed Brown case 1: 596c4762a1bSJed Brown if (cellHybrid) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %d", dim); 597c4762a1bSJed Brown ierr = CreateSimplex_1D(comm, dm);CHKERRQ(ierr); 598c4762a1bSJed Brown break; 599c4762a1bSJed Brown case 2: 600c4762a1bSJed Brown if (cellSimplex) { 601c4762a1bSJed Brown if (cellHybrid) { 602c4762a1bSJed Brown ierr = CreateSimplexHybrid_2D(comm, user->testNum, dm);CHKERRQ(ierr); 603c4762a1bSJed Brown } else { 604c4762a1bSJed Brown ierr = CreateSimplex_2D(comm, dm);CHKERRQ(ierr); 605c4762a1bSJed Brown } 606c4762a1bSJed Brown } else { 607c4762a1bSJed Brown if (cellHybrid) { 608c4762a1bSJed Brown ierr = CreateTensorProductHybrid_2D(comm, user->testNum, dm);CHKERRQ(ierr); 609c4762a1bSJed Brown } else { 610c4762a1bSJed Brown ierr = CreateTensorProduct_2D(comm, user->testNum, dm);CHKERRQ(ierr); 611c4762a1bSJed Brown } 612c4762a1bSJed Brown } 613c4762a1bSJed Brown break; 614c4762a1bSJed Brown case 3: 615c4762a1bSJed Brown if (cellSimplex) { 616c4762a1bSJed Brown if (cellHybrid) { 617c4762a1bSJed Brown ierr = CreateSimplexHybrid_3D(comm, user->testNum, dm);CHKERRQ(ierr); 618c4762a1bSJed Brown } else { 619c4762a1bSJed Brown ierr = CreateSimplex_3D(comm, user->testNum, dm);CHKERRQ(ierr); 620c4762a1bSJed Brown } 621c4762a1bSJed Brown } else { 622c4762a1bSJed Brown if (cellHybrid) { 623c4762a1bSJed Brown ierr = CreateTensorProductHybrid_3D(comm, user->testNum, dm);CHKERRQ(ierr); 624c4762a1bSJed Brown } else { 625c4762a1bSJed Brown ierr = CreateTensorProduct_3D(comm, user->testNum, dm);CHKERRQ(ierr); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown } 628c4762a1bSJed Brown break; 629c4762a1bSJed Brown default: 630c4762a1bSJed Brown SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 631c4762a1bSJed Brown } 632c4762a1bSJed Brown if (user->testPartition && size > 1) { 633c4762a1bSJed Brown PetscPartitioner part; 634c4762a1bSJed Brown PetscInt *sizes = NULL; 635c4762a1bSJed Brown PetscInt *points = NULL; 636c4762a1bSJed Brown 637c4762a1bSJed Brown if (!rank) { 638c4762a1bSJed Brown if (dim == 2 && cellSimplex && !cellHybrid && size == 2) { 639c4762a1bSJed Brown switch (user->testNum) { 640c4762a1bSJed Brown case 0: { 641c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 1}; 642c4762a1bSJed Brown PetscInt triPoints_p2[2] = {0, 1}; 643c4762a1bSJed Brown 644c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr); 645c4762a1bSJed Brown ierr = PetscArraycpy(sizes, triSizes_p2, 2);CHKERRQ(ierr); 646c4762a1bSJed Brown ierr = PetscArraycpy(points, triPoints_p2, 2);CHKERRQ(ierr);break;} 647c4762a1bSJed Brown default: 648c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum); 649c4762a1bSJed Brown } 650c4762a1bSJed Brown } else if (dim == 2 && cellSimplex && cellHybrid && size == 2) { 651c4762a1bSJed Brown switch (user->testNum) { 652c4762a1bSJed Brown case 0: { 653c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 2}; 654c4762a1bSJed Brown PetscInt triPoints_p2[3] = {0, 1, 2}; 655c4762a1bSJed Brown 656c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 3, &points);CHKERRQ(ierr); 657c4762a1bSJed Brown ierr = PetscArraycpy(sizes, triSizes_p2, 2);CHKERRQ(ierr); 658c4762a1bSJed Brown ierr = PetscArraycpy(points, triPoints_p2, 3);CHKERRQ(ierr);break;} 659c4762a1bSJed Brown default: 660c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular hybrid mesh on 2 procs", user->testNum); 661c4762a1bSJed Brown } 662c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && !cellHybrid && size == 2) { 663c4762a1bSJed Brown switch (user->testNum) { 664c4762a1bSJed Brown case 0: { 665c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 1}; 666c4762a1bSJed Brown PetscInt quadPoints_p2[2] = {0, 1}; 667c4762a1bSJed Brown 668c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr); 669c4762a1bSJed Brown ierr = PetscArraycpy(sizes, quadSizes_p2, 2);CHKERRQ(ierr); 670c4762a1bSJed Brown ierr = PetscArraycpy(points, quadPoints_p2, 2);CHKERRQ(ierr);break;} 671c4762a1bSJed Brown default: 672c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for quadrilateral mesh on 2 procs", user->testNum); 673c4762a1bSJed Brown } 674c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && cellHybrid && size == 2) { 675c4762a1bSJed Brown switch (user->testNum) { 676c4762a1bSJed Brown case 0: { 677c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 2}; 678c4762a1bSJed Brown PetscInt quadPoints_p2[3] = {0, 1, 2}; 679c4762a1bSJed Brown 680c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 3, &points);CHKERRQ(ierr); 681c4762a1bSJed Brown ierr = PetscArraycpy(sizes, quadSizes_p2, 2);CHKERRQ(ierr); 682c4762a1bSJed Brown ierr = PetscArraycpy(points, quadPoints_p2, 3);CHKERRQ(ierr);break;} 683c4762a1bSJed Brown default: 684c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for quadrilateral hybrid mesh on 2 procs", user->testNum); 685c4762a1bSJed Brown } 686c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && !cellHybrid && size == 2) { 687c4762a1bSJed Brown switch (user->testNum) { 688c4762a1bSJed Brown case 0: { 689c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 1}; 690c4762a1bSJed Brown PetscInt tetPoints_p2[2] = {0, 1}; 691c4762a1bSJed Brown 692c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr); 693c4762a1bSJed Brown ierr = PetscArraycpy(sizes, tetSizes_p2, 2);CHKERRQ(ierr); 694c4762a1bSJed Brown ierr = PetscArraycpy(points, tetPoints_p2, 2);CHKERRQ(ierr);break;} 695c4762a1bSJed Brown case 1: { 696c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 1}; 697c4762a1bSJed Brown PetscInt tetPoints_p2[2] = {0, 1}; 698c4762a1bSJed Brown 699c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr); 700c4762a1bSJed Brown ierr = PetscArraycpy(sizes, tetSizes_p2, 2);CHKERRQ(ierr); 701c4762a1bSJed Brown ierr = PetscArraycpy(points, tetPoints_p2, 2);CHKERRQ(ierr);break;} 702c4762a1bSJed Brown default: 703c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for tetrahedral mesh on 2 procs", user->testNum); 704c4762a1bSJed Brown } 705c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && cellHybrid && size == 2) { 706c4762a1bSJed Brown switch (user->testNum) { 707c4762a1bSJed Brown case 0: { 708c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 2}; 709c4762a1bSJed Brown PetscInt tetPoints_p2[3] = {0, 1, 2}; 710c4762a1bSJed Brown 711c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 3, &points);CHKERRQ(ierr); 712c4762a1bSJed Brown ierr = PetscArraycpy(sizes, tetSizes_p2, 2);CHKERRQ(ierr); 713c4762a1bSJed Brown ierr = PetscArraycpy(points, tetPoints_p2, 3);CHKERRQ(ierr);break;} 714c4762a1bSJed Brown case 1: { 715c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {3, 4}; 716c4762a1bSJed Brown PetscInt tetPoints_p2[7] = {0, 3, 5, 1, 2, 4, 6}; 717c4762a1bSJed Brown 718c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 7, &points);CHKERRQ(ierr); 719c4762a1bSJed Brown ierr = PetscArraycpy(sizes, tetSizes_p2, 2);CHKERRQ(ierr); 720c4762a1bSJed Brown ierr = PetscArraycpy(points, tetPoints_p2, 7);CHKERRQ(ierr);break;} 721c4762a1bSJed Brown default: 722c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for tetrahedral hybrid mesh on 2 procs", user->testNum); 723c4762a1bSJed Brown } 724c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && !cellHybrid && size == 2) { 725c4762a1bSJed Brown switch (user->testNum) { 726c4762a1bSJed Brown case 0: { 727c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 1}; 728c4762a1bSJed Brown PetscInt hexPoints_p2[2] = {0, 1}; 729c4762a1bSJed Brown 730c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr); 731c4762a1bSJed Brown ierr = PetscArraycpy(sizes, hexSizes_p2, 2);CHKERRQ(ierr); 732c4762a1bSJed Brown ierr = PetscArraycpy(points, hexPoints_p2, 2);CHKERRQ(ierr);break;} 733c4762a1bSJed Brown default: 734c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for hexahedral mesh on 2 procs", user->testNum); 735c4762a1bSJed Brown } 736c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && cellHybrid && size == 2) { 737c4762a1bSJed Brown switch (user->testNum) { 738c4762a1bSJed Brown case 0: { 739c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 1}; 740c4762a1bSJed Brown PetscInt hexPoints_p2[2] = {0, 1}; 741c4762a1bSJed Brown 742c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr); 743c4762a1bSJed Brown ierr = PetscArraycpy(sizes, hexSizes_p2, 2);CHKERRQ(ierr); 744c4762a1bSJed Brown ierr = PetscArraycpy(points, hexPoints_p2, 2);CHKERRQ(ierr);break;} 745c4762a1bSJed Brown case 1: { 746c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {5, 4}; 747c4762a1bSJed Brown PetscInt hexPoints_p2[9] = {3, 4, 5, 7, 8, 0, 1, 2, 6}; 748c4762a1bSJed Brown 749c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 9, &points);CHKERRQ(ierr); 750c4762a1bSJed Brown ierr = PetscArraycpy(sizes, hexSizes_p2, 2);CHKERRQ(ierr); 751c4762a1bSJed Brown ierr = PetscArraycpy(points, hexPoints_p2, 9);CHKERRQ(ierr);break;} 752c4762a1bSJed Brown default: 753c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for hexahedral hybrid mesh on 2 procs", user->testNum); 754c4762a1bSJed Brown } 755c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition"); 756c4762a1bSJed Brown } 757c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 758c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 759c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr); 760c4762a1bSJed Brown ierr = PetscFree2(sizes, points);CHKERRQ(ierr); 761c4762a1bSJed Brown } else { 762c4762a1bSJed Brown PetscPartitioner part; 763c4762a1bSJed Brown 764c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm,&part);CHKERRQ(ierr); 765c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 766c4762a1bSJed Brown } 767c4762a1bSJed Brown { 768c4762a1bSJed Brown DM pdm = NULL; 769c4762a1bSJed Brown 770c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 771c4762a1bSJed Brown if (pdm) { 772c4762a1bSJed Brown ierr = DMViewFromOptions(pdm, NULL, "-dm_view");CHKERRQ(ierr); 773c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 774c4762a1bSJed Brown *dm = pdm; 775c4762a1bSJed Brown } 776c4762a1bSJed Brown } 777c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 778c4762a1bSJed Brown if (user->uninterpolate || user->reinterpolate) { 779c4762a1bSJed Brown DM udm = NULL; 780c4762a1bSJed Brown 781c4762a1bSJed Brown ierr = DMPlexUninterpolate(*dm, &udm);CHKERRQ(ierr); 782c4762a1bSJed Brown ierr = DMPlexCopyCoordinates(*dm, udm);CHKERRQ(ierr); 783c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 784c4762a1bSJed Brown *dm = udm; 785c4762a1bSJed Brown } 786c4762a1bSJed Brown if (user->reinterpolate) { 787c4762a1bSJed Brown DM idm = NULL; 788c4762a1bSJed Brown 789c4762a1bSJed Brown ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 790c4762a1bSJed Brown ierr = DMPlexCopyCoordinates(*dm, idm);CHKERRQ(ierr); 791c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 792c4762a1bSJed Brown *dm = idm; 793c4762a1bSJed Brown } 794c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Hybrid Mesh");CHKERRQ(ierr); 795c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 796c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "hyb_");CHKERRQ(ierr); 797c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 798c4762a1bSJed Brown PetscFunctionReturn(0); 799c4762a1bSJed Brown } 800c4762a1bSJed Brown 801c4762a1bSJed Brown int main(int argc, char **argv) 802c4762a1bSJed Brown { 803c4762a1bSJed Brown DM dm; 804c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 805c4762a1bSJed Brown PetscErrorCode ierr; 806c4762a1bSJed Brown 807c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 808c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 809c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 810c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 811c4762a1bSJed Brown ierr = PetscFinalize(); 812c4762a1bSJed Brown return ierr; 813c4762a1bSJed Brown } 814c4762a1bSJed Brown 815c4762a1bSJed Brown /*TEST 816c4762a1bSJed Brown 817c4762a1bSJed Brown # 1D Simplex 29-31 818c4762a1bSJed Brown testset: 81954fcfd0cSMatthew G. Knepley args: -dim 1 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 820c4762a1bSJed Brown test: 821c4762a1bSJed Brown suffix: 29 822c4762a1bSJed Brown test: 823c4762a1bSJed Brown suffix: 30 824c4762a1bSJed Brown args: -dm_refine 1 825c4762a1bSJed Brown test: 826c4762a1bSJed Brown suffix: 31 827c4762a1bSJed Brown args: -dm_refine 5 828c4762a1bSJed Brown 829c4762a1bSJed Brown # 2D Simplex 0-3 830c4762a1bSJed Brown testset: 83154fcfd0cSMatthew G. Knepley args: -dim 2 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 832c4762a1bSJed Brown test: 833c4762a1bSJed Brown suffix: 0 834c4762a1bSJed Brown args: -hyb_dm_plex_check_faces 835c4762a1bSJed Brown test: 836c4762a1bSJed Brown suffix: 1 837c4762a1bSJed Brown args: -dm_refine 1 -hyb_dm_plex_check_faces 838c4762a1bSJed Brown test: 839c4762a1bSJed Brown suffix: 2 840c4762a1bSJed Brown nsize: 2 841c4762a1bSJed Brown args: -hyb_dm_plex_check_faces 842c4762a1bSJed Brown test: 843c4762a1bSJed Brown suffix: 3 844c4762a1bSJed Brown nsize: 2 845c4762a1bSJed Brown args: -dm_refine 1 -hyb_dm_plex_check_faces 846c4762a1bSJed Brown test: 847c4762a1bSJed Brown suffix: 32 848c4762a1bSJed Brown args: -dm_refine 1 -uninterpolate 849c4762a1bSJed Brown test: 850c4762a1bSJed Brown suffix: 33 851c4762a1bSJed Brown nsize: 2 852c4762a1bSJed Brown args: -dm_refine 1 -uninterpolate 853c4762a1bSJed Brown test: 854c4762a1bSJed Brown suffix: 34 855c4762a1bSJed Brown nsize: 2 856c4762a1bSJed Brown args: -dm_refine 3 -uninterpolate 857c4762a1bSJed Brown 858c4762a1bSJed Brown # 2D Hybrid Simplex 4-7 859c4762a1bSJed Brown testset: 86054fcfd0cSMatthew G. Knepley args: -dim 2 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 861c4762a1bSJed Brown test: 862c4762a1bSJed Brown suffix: 4 863c4762a1bSJed Brown test: 864c4762a1bSJed Brown suffix: 5 865c4762a1bSJed Brown args: -dm_refine 1 866c4762a1bSJed Brown test: 867c4762a1bSJed Brown suffix: 6 868c4762a1bSJed Brown nsize: 2 869c4762a1bSJed Brown test: 870c4762a1bSJed Brown suffix: 7 871c4762a1bSJed Brown nsize: 2 872c4762a1bSJed Brown args: -dm_refine 1 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown suffix: 24 875c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 876c4762a1bSJed Brown 877c4762a1bSJed Brown # 2D Quad 12-13 878c4762a1bSJed Brown testset: 87954fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 880c4762a1bSJed Brown test: 881c4762a1bSJed Brown suffix: 12 882c4762a1bSJed Brown args: -dm_refine 1 883c4762a1bSJed Brown test: 884c4762a1bSJed Brown suffix: 13 885c4762a1bSJed Brown nsize: 2 886c4762a1bSJed Brown args: -dm_refine 1 887c4762a1bSJed Brown 888c4762a1bSJed Brown # 2D Hybrid Quad 27-28 889c4762a1bSJed Brown testset: 89054fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 891c4762a1bSJed Brown test: 892c4762a1bSJed Brown suffix: 27 893c4762a1bSJed Brown args: -dm_refine 1 894c4762a1bSJed Brown test: 895c4762a1bSJed Brown suffix: 28 896c4762a1bSJed Brown nsize: 2 897c4762a1bSJed Brown args: -dm_refine 1 898c4762a1bSJed Brown 899c4762a1bSJed Brown # 3D Simplex 8-11 900c4762a1bSJed Brown testset: 90154fcfd0cSMatthew G. Knepley args: -dim 3 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 902c4762a1bSJed Brown test: 903c4762a1bSJed Brown suffix: 8 904c4762a1bSJed Brown args: -dm_refine 1 905c4762a1bSJed Brown test: 906c4762a1bSJed Brown suffix: 9 907c4762a1bSJed Brown nsize: 2 908c4762a1bSJed Brown args: -dm_refine 1 909c4762a1bSJed Brown test: 910c4762a1bSJed Brown suffix: 10 911c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 912c4762a1bSJed Brown test: 913c4762a1bSJed Brown suffix: 11 914c4762a1bSJed Brown nsize: 2 915c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 916c4762a1bSJed Brown test: 917c4762a1bSJed Brown suffix: 25 918c4762a1bSJed Brown args: -test_num 2 -dm_refine 2 919c4762a1bSJed Brown 920c4762a1bSJed Brown # 3D Hybrid Simplex 16-19 921c4762a1bSJed Brown testset: 92254fcfd0cSMatthew G. Knepley args: -dim 3 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 923c4762a1bSJed Brown test: 924c4762a1bSJed Brown suffix: 16 925c4762a1bSJed Brown args: -dm_refine 1 926c4762a1bSJed Brown test: 927c4762a1bSJed Brown suffix: 17 928c4762a1bSJed Brown nsize: 2 929c4762a1bSJed Brown args: -dm_refine 1 930c4762a1bSJed Brown test: 931c4762a1bSJed Brown suffix: 18 932c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 933c4762a1bSJed Brown test: 934c4762a1bSJed Brown suffix: 19 935c4762a1bSJed Brown nsize: 2 936c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 937c4762a1bSJed Brown 938c4762a1bSJed Brown # 3D Hex 14-15 939c4762a1bSJed Brown testset: 94054fcfd0cSMatthew G. Knepley args: -dim 3 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all 941c4762a1bSJed Brown test: 942c4762a1bSJed Brown suffix: 14 943c4762a1bSJed Brown args: -dm_refine 1 944c4762a1bSJed Brown test: 945c4762a1bSJed Brown suffix: 15 946c4762a1bSJed Brown nsize: 2 947c4762a1bSJed Brown args: -dm_refine 1 948c4762a1bSJed Brown test: 949c4762a1bSJed Brown suffix: 26 950c4762a1bSJed Brown args: -test_num 1 -dm_refine 2 951c4762a1bSJed Brown 952c4762a1bSJed Brown # 3D Hybrid Hex 20-23 953c4762a1bSJed Brown testset: 95454fcfd0cSMatthew G. Knepley args: -dim 3 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all 955c4762a1bSJed Brown test: 956c4762a1bSJed Brown suffix: 20 957c4762a1bSJed Brown args: -dm_refine 1 958c4762a1bSJed Brown test: 959c4762a1bSJed Brown suffix: 21 960c4762a1bSJed Brown nsize: 2 961c4762a1bSJed Brown args: -dm_refine 1 962c4762a1bSJed Brown test: 963c4762a1bSJed Brown suffix: 22 964c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 965c4762a1bSJed Brown test: 966c4762a1bSJed Brown suffix: 23 967c4762a1bSJed Brown nsize: 2 968c4762a1bSJed Brown args: -test_num 1 -dm_refine 1 969c4762a1bSJed Brown 970c4762a1bSJed Brown # Hybrid interpolation 97154fcfd0cSMatthew G. Knepley # TODO Setup new tests (like -reinterpolate) that interpolate hybrid cells 972c4762a1bSJed Brown testset: 973c4762a1bSJed Brown nsize: 2 97454fcfd0cSMatthew 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 975c4762a1bSJed Brown test: 976c4762a1bSJed Brown suffix: hybint_2d_0 977c4762a1bSJed Brown args: -dim 2 -dm_refine 2 978c4762a1bSJed Brown test: 979c4762a1bSJed Brown suffix: hybint_2d_1 980c4762a1bSJed Brown args: -dim 2 -dm_refine 2 -test_num 1 981c4762a1bSJed Brown test: 982c4762a1bSJed Brown suffix: hybint_3d_0 983c4762a1bSJed Brown args: -dim 3 -dm_refine 1 984c4762a1bSJed Brown test: 985c4762a1bSJed Brown suffix: hybint_3d_1 986c4762a1bSJed Brown args: -dim 3 -dm_refine 1 -test_num 1 987c4762a1bSJed Brown 988c4762a1bSJed Brown TEST*/ 989