1c4762a1bSJed Brown static char help[] = "Tests for mesh interpolation\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* List of test meshes 6c4762a1bSJed Brown 7c4762a1bSJed Brown Triangle 8c4762a1bSJed Brown -------- 9c4762a1bSJed Brown Test 0: 10c4762a1bSJed Brown Two triangles sharing a face 11c4762a1bSJed Brown 12c4762a1bSJed Brown 4 13c4762a1bSJed Brown / | \ 14c4762a1bSJed Brown / | \ 15c4762a1bSJed Brown / | \ 16c4762a1bSJed Brown 2 0 | 1 5 17c4762a1bSJed Brown \ | / 18c4762a1bSJed Brown \ | / 19c4762a1bSJed Brown \ | / 20c4762a1bSJed Brown 3 21c4762a1bSJed Brown 22c4762a1bSJed Brown should become 23c4762a1bSJed Brown 24c4762a1bSJed Brown 4 25c4762a1bSJed Brown / | \ 26c4762a1bSJed Brown 8 | 9 27c4762a1bSJed Brown / | \ 28c4762a1bSJed Brown 2 0 7 1 5 29c4762a1bSJed Brown \ | / 30c4762a1bSJed Brown 6 | 10 31c4762a1bSJed Brown \ | / 32c4762a1bSJed Brown 3 33c4762a1bSJed Brown 34c4762a1bSJed Brown Tetrahedron 35c4762a1bSJed Brown ----------- 36c4762a1bSJed Brown Test 0: 37c4762a1bSJed Brown Two tets sharing a face 38c4762a1bSJed Brown 39c4762a1bSJed Brown cell 5 _______ cell 40c4762a1bSJed Brown 0 / | \ \ 1 41c4762a1bSJed Brown / | \ \ 42c4762a1bSJed Brown / | \ \ 43c4762a1bSJed Brown 2----|----4-----6 44c4762a1bSJed Brown \ | / / 45c4762a1bSJed Brown \ | / / 46c4762a1bSJed Brown \ | / / 47c4762a1bSJed Brown 3------- 48c4762a1bSJed Brown 49c4762a1bSJed Brown should become 50c4762a1bSJed Brown 51c4762a1bSJed Brown cell 5 _______ cell 52c4762a1bSJed Brown 0 / | \ \ 1 53c4762a1bSJed Brown / | \ \ 54c4762a1bSJed Brown 17 | 18 13 22 55c4762a1bSJed Brown / 8 19 10 \ \ 56c4762a1bSJed Brown / | \ \ 57c4762a1bSJed Brown 2---14-|------4--20--6 58c4762a1bSJed Brown \ | / / 59c4762a1bSJed Brown \ 9 | 7 / / 60c4762a1bSJed Brown 16 | 15 11 21 61c4762a1bSJed Brown \ | / / 62c4762a1bSJed Brown \ | / / 63c4762a1bSJed Brown 3------- 64c4762a1bSJed Brown 65c4762a1bSJed Brown Quadrilateral 66c4762a1bSJed Brown ------------- 67c4762a1bSJed Brown Test 0: 68c4762a1bSJed Brown Two quads sharing a face 69c4762a1bSJed Brown 70c4762a1bSJed Brown 5-------4-------7 71c4762a1bSJed Brown | | | 72c4762a1bSJed Brown | 0 | 1 | 73c4762a1bSJed Brown | | | 74c4762a1bSJed Brown 2-------3-------6 75c4762a1bSJed Brown 76c4762a1bSJed Brown should become 77c4762a1bSJed Brown 78c4762a1bSJed Brown 5--10---4--14---7 79c4762a1bSJed Brown | | | 80c4762a1bSJed Brown 11 0 9 1 13 81c4762a1bSJed Brown | | | 82c4762a1bSJed Brown 2---8---3--12---6 83c4762a1bSJed Brown 84c4762a1bSJed Brown Test 1: 85c4762a1bSJed Brown A quad and a triangle sharing a face 86c4762a1bSJed Brown 87c4762a1bSJed Brown 5-------4 88c4762a1bSJed Brown | | \ 89c4762a1bSJed Brown | 0 | \ 90c4762a1bSJed Brown | | 1 \ 91c4762a1bSJed Brown 2-------3----6 92c4762a1bSJed Brown 93c4762a1bSJed Brown should become 94c4762a1bSJed Brown 95c4762a1bSJed Brown 5---9---4 96c4762a1bSJed Brown | | \ 97c4762a1bSJed Brown 10 0 8 12 98c4762a1bSJed Brown | | 1 \ 99c4762a1bSJed Brown 2---7---3-11-6 100c4762a1bSJed Brown 101c4762a1bSJed Brown Hexahedron 102c4762a1bSJed Brown ---------- 103c4762a1bSJed Brown Test 0: 104c4762a1bSJed Brown Two hexes sharing a face 105c4762a1bSJed Brown 106c4762a1bSJed Brown cell 9-------------8-------------13 cell 107c4762a1bSJed Brown 0 /| /| /| 1 108c4762a1bSJed Brown / | 15 / | 21 / | 109c4762a1bSJed Brown / | / | / | 110c4762a1bSJed Brown 6-------------7-------------12 | 111c4762a1bSJed Brown | | 18 | | 24 | | 112c4762a1bSJed Brown | | | | | | 113c4762a1bSJed Brown |19 | |17 | |23 | 114c4762a1bSJed Brown | | 16 | | 22 | | 115c4762a1bSJed Brown | 5---------|---4---------|---11 116c4762a1bSJed Brown | / | / | / 117c4762a1bSJed Brown | / 14 | / 20 | / 118c4762a1bSJed Brown |/ |/ |/ 119c4762a1bSJed Brown 2-------------3-------------10 120c4762a1bSJed Brown 121c4762a1bSJed Brown should become 122c4762a1bSJed Brown 123c4762a1bSJed Brown cell 9-----31------8-----42------13 cell 124c4762a1bSJed Brown 0 /| /| /| 1 125c4762a1bSJed Brown 32 | 15 30| 21 41| 126c4762a1bSJed Brown / | / | / | 127c4762a1bSJed Brown 6-----29------7-----40------12 | 128c4762a1bSJed Brown | | 17 | | 23 | | 129c4762a1bSJed Brown | 35 | 36 | 44 130c4762a1bSJed Brown |19 | |18 | |24 | 131c4762a1bSJed Brown 34 | 16 33 | 22 43 | 132c4762a1bSJed Brown | 5-----26--|---4-----37--|---11 133c4762a1bSJed Brown | / | / | / 134c4762a1bSJed Brown | 25 14 | 27 20 | 38 135c4762a1bSJed Brown |/ |/ |/ 136c4762a1bSJed Brown 2-----28------3-----39------10 137c4762a1bSJed Brown 138c4762a1bSJed Brown */ 139c4762a1bSJed Brown 140c4762a1bSJed Brown typedef struct { 141c4762a1bSJed Brown PetscInt testNum; /* Indicates the mesh to create */ 142c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 14330602db0SMatthew G. Knepley PetscBool simplex; /* Use simplices or hexes */ 144c4762a1bSJed Brown PetscBool useGenerator; /* Construct mesh with a mesh generator */ 145c4762a1bSJed Brown } AppCtx; 146c4762a1bSJed Brown 147d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 148d71ae5a4SJacob Faibussowitsch { 149c4762a1bSJed Brown PetscFunctionBegin; 150c4762a1bSJed Brown options->testNum = 0; 151c4762a1bSJed Brown options->dim = 2; 15230602db0SMatthew G. Knepley options->simplex = PETSC_TRUE; 153c4762a1bSJed Brown options->useGenerator = PETSC_FALSE; 154c4762a1bSJed Brown 155d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 1569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex7.c", options->testNum, &options->testNum, NULL, 0)); 1579566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex7.c", options->dim, &options->dim, NULL, 1, 3)); 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex7.c", options->simplex, &options->simplex, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex7.c", options->useGenerator, &options->useGenerator, NULL)); 160d0609cedSBarry Smith PetscOptionsEnd(); 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM dm) 165d71ae5a4SJacob Faibussowitsch { 166c4762a1bSJed Brown PetscInt depth = 1, testNum = 0, p; 167c4762a1bSJed Brown PetscMPIInt rank; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBegin; 1709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 171dd400576SPatrick Sanan if (rank == 0) { 172c4762a1bSJed Brown switch (testNum) { 1739371c9d4SSatish Balay case 0: { 174c4762a1bSJed Brown PetscInt numPoints[2] = {4, 2}; 175c4762a1bSJed Brown PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 176c4762a1bSJed Brown PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 177c4762a1bSJed Brown PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 178c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 179c4762a1bSJed Brown PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 180c4762a1bSJed Brown 1819566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 18248a46eb9SPierre Jolivet for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 1839371c9d4SSatish Balay } break; 184d71ae5a4SJacob Faibussowitsch default: 185d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } else { 188c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL)); 1919566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "marker")); 192c4762a1bSJed Brown } 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM dm) 197d71ae5a4SJacob Faibussowitsch { 198c4762a1bSJed Brown PetscInt depth = 1, testNum = 0, p; 199c4762a1bSJed Brown PetscMPIInt rank; 200c4762a1bSJed Brown 201c4762a1bSJed Brown PetscFunctionBegin; 2029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 203dd400576SPatrick Sanan if (rank == 0) { 204c4762a1bSJed Brown switch (testNum) { 2059371c9d4SSatish Balay case 0: { 206c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 207c4762a1bSJed Brown PetscInt coneSize[23] = {4, 4, 0, 0, 0, 0, 0}; 208c4762a1bSJed Brown PetscInt cones[8] = {2, 4, 3, 5, 3, 4, 6, 5}; 209c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 210c4762a1bSJed 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}; 211c4762a1bSJed Brown PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 21448a46eb9SPierre Jolivet for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 2159371c9d4SSatish Balay } break; 216d71ae5a4SJacob Faibussowitsch default: 217d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } else { 220c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 221c4762a1bSJed Brown 2229566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL)); 2239566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "marker")); 224c4762a1bSJed Brown } 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM dm) 229d71ae5a4SJacob Faibussowitsch { 230c4762a1bSJed Brown PetscInt depth = 1, p; 231c4762a1bSJed Brown PetscMPIInt rank; 232c4762a1bSJed Brown 233c4762a1bSJed Brown PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 235dd400576SPatrick Sanan if (rank == 0) { 236c4762a1bSJed Brown switch (testNum) { 2379371c9d4SSatish Balay case 0: { 238c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 239c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 240c4762a1bSJed Brown PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 241c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 242c4762a1bSJed Brown PetscScalar vertexCoords[12] = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0, 0.5, 1.0}; 243c4762a1bSJed Brown PetscInt markerPoints[12] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1}; 244c4762a1bSJed Brown 2459566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 24648a46eb9SPierre Jolivet for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 2479371c9d4SSatish Balay } break; 2489371c9d4SSatish Balay case 1: { 249c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 250c4762a1bSJed Brown PetscInt coneSize[7] = {4, 3, 0, 0, 0, 0, 0}; 251c4762a1bSJed Brown PetscInt cones[7] = {2, 3, 4, 5, 3, 6, 4}; 252c4762a1bSJed Brown PetscInt coneOrientations[7] = {0, 0, 0, 0, 0, 0, 0}; 253c4762a1bSJed Brown PetscScalar vertexCoords[10] = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0}; 254c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 255c4762a1bSJed Brown 2569566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 25748a46eb9SPierre Jolivet for (p = 0; p < 5; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 2589371c9d4SSatish Balay } break; 259d71ae5a4SJacob Faibussowitsch default: 260d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown } else { 263c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL)); 2669566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "marker")); 267c4762a1bSJed Brown } 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateHex_3D(MPI_Comm comm, DM dm) 272d71ae5a4SJacob Faibussowitsch { 273c4762a1bSJed Brown PetscInt depth = 1, testNum = 0, p; 274c4762a1bSJed Brown PetscMPIInt rank; 275c4762a1bSJed Brown 276c4762a1bSJed Brown PetscFunctionBegin; 2779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 278dd400576SPatrick Sanan if (rank == 0) { 279c4762a1bSJed Brown switch (testNum) { 2809371c9d4SSatish Balay case 0: { 281c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 282c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 283c4762a1bSJed Brown PetscInt cones[16] = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8}; 284c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 2859371c9d4SSatish Balay PetscScalar vertexCoords[36] = {-0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, -0.5, 1.0, 0.0, -0.5, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, -0.5, 1.0, 1.0, 0.5, 0.0, 0.0, 0.5, 1.0, 0.0, 0.5, 0.0, 1.0, 0.5, 1.0, 1.0}; 286c4762a1bSJed Brown PetscInt markerPoints[16] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1}; 287c4762a1bSJed Brown 2889566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 28948a46eb9SPierre Jolivet for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 2909371c9d4SSatish Balay } break; 291d71ae5a4SJacob Faibussowitsch default: 292d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 293c4762a1bSJed Brown } 294c4762a1bSJed Brown } else { 295c4762a1bSJed Brown PetscInt numPoints[2] = {0, 0}; 296c4762a1bSJed Brown 2979566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL)); 2989566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "marker")); 299c4762a1bSJed Brown } 3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, PetscInt testNum, AppCtx *user, DM *dm) 304d71ae5a4SJacob Faibussowitsch { 305c4762a1bSJed Brown PetscBool useGenerator = user->useGenerator; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBegin; 3089566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3099566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 310b5a892a1SMatthew G. Knepley if (!useGenerator) { 31130602db0SMatthew G. Knepley PetscInt dim = user->dim; 31230602db0SMatthew G. Knepley PetscBool simplex = user->simplex; 31330602db0SMatthew G. Knepley 3149566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 315c4762a1bSJed Brown switch (dim) { 316c4762a1bSJed Brown case 2: 31730602db0SMatthew G. Knepley if (simplex) { 3189566063dSJacob Faibussowitsch PetscCall(CreateSimplex_2D(comm, *dm)); 319c4762a1bSJed Brown } else { 3209566063dSJacob Faibussowitsch PetscCall(CreateQuad_2D(comm, testNum, *dm)); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown break; 323c4762a1bSJed Brown case 3: 32430602db0SMatthew G. Knepley if (simplex) { 3259566063dSJacob Faibussowitsch PetscCall(CreateSimplex_3D(comm, *dm)); 326c4762a1bSJed Brown } else { 3279566063dSJacob Faibussowitsch PetscCall(CreateHex_3D(comm, *dm)); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown break; 330d71ae5a4SJacob Faibussowitsch default: 331d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown } 3349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 3359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Mesh")); 3369566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 3373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 341d71ae5a4SJacob Faibussowitsch { 34230602db0SMatthew G. Knepley DM dm; 343c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 344c4762a1bSJed Brown 345327415f7SBarry Smith PetscFunctionBeginUser; 3469566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3479566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 3489566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, user.testNum, &user, &dm)); 3499566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3509566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 351b122ec5aSJacob Faibussowitsch return 0; 352c4762a1bSJed Brown } 353c4762a1bSJed Brown 354c4762a1bSJed Brown /*TEST 355b5a892a1SMatthew G. Knepley testset: 356b5a892a1SMatthew G. Knepley args: -dm_plex_interpolate_pre -dm_plex_check_all 357c4762a1bSJed Brown 358c4762a1bSJed Brown # Two cell test meshes 0-7 359c4762a1bSJed Brown test: 360c4762a1bSJed Brown suffix: 0 36167d101c4SVaclav Hapla args: -dim 2 -dm_view ascii::ascii_info_detail 362c4762a1bSJed Brown test: 363c4762a1bSJed Brown suffix: 1 364c4762a1bSJed Brown nsize: 2 365e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dim 2 -dm_view ascii::ascii_info_detail 366c4762a1bSJed Brown test: 367c4762a1bSJed Brown suffix: 2 36830602db0SMatthew G. Knepley args: -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail 369c4762a1bSJed Brown test: 370c4762a1bSJed Brown suffix: 3 371c4762a1bSJed Brown nsize: 2 372e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail 373c4762a1bSJed Brown test: 374c4762a1bSJed Brown suffix: 4 375c4762a1bSJed Brown args: -dim 3 -dm_view ascii::ascii_info_detail 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown suffix: 5 378c4762a1bSJed Brown nsize: 2 379e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dim 3 -dm_view ascii::ascii_info_detail 380c4762a1bSJed Brown test: 381c4762a1bSJed Brown suffix: 6 38230602db0SMatthew G. Knepley args: -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail 383c4762a1bSJed Brown test: 384c4762a1bSJed Brown suffix: 7 385c4762a1bSJed Brown nsize: 2 386e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail 387c4762a1bSJed Brown # 2D Hybrid Mesh 8 388c4762a1bSJed Brown test: 389c4762a1bSJed Brown suffix: 8 39030602db0SMatthew G. Knepley args: -dim 2 -simplex 0 -testnum 1 -dm_view ascii::ascii_info_detail 391b5a892a1SMatthew G. Knepley 392b5a892a1SMatthew G. Knepley testset: 393b5a892a1SMatthew G. Knepley args: -dm_plex_check_all 394b5a892a1SMatthew G. Knepley 395da9060c4SMatthew G. Knepley # Reference cells 396da9060c4SMatthew G. Knepley test: 397da9060c4SMatthew G. Knepley suffix: 12 398b5a892a1SMatthew G. Knepley args: -use_generator -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_plex_check_all 399*e0008caeSPierre Jolivet output_file: output/empty.out 400c4762a1bSJed Brown # TetGen meshes 9-10 401c4762a1bSJed Brown test: 402c4762a1bSJed Brown suffix: 9 403c4762a1bSJed Brown requires: triangle 40430602db0SMatthew G. Knepley args: -use_generator -dm_view ascii::ascii_info_detail -dm_coord_space 0 405c4762a1bSJed Brown test: 406c4762a1bSJed Brown suffix: 10 407c4762a1bSJed Brown requires: ctetgen 40830602db0SMatthew G. Knepley args: -use_generator -dm_plex_dim 3 -dm_view ascii::ascii_info_detail -dm_coord_space 0 409c4762a1bSJed Brown # Cubit meshes 11 410c4762a1bSJed Brown test: 411c4762a1bSJed Brown suffix: 11 412c4762a1bSJed Brown requires: exodusii 41346ac1a18SMatthew G. Knepley args: -use_generator -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -dm_view ascii::ascii_info_detail -dm_coord_space 0 414c4762a1bSJed Brown 415c4762a1bSJed Brown TEST*/ 416