1c4762a1bSJed Brown static char help[] = "Tests interpolation and output of hybrid meshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 7c4762a1bSJed Brown PetscBool interpolate; /* Interpolate the mesh */ 8c4762a1bSJed Brown PetscInt meshNum; /* Which mesh we should construct */ 9c4762a1bSJed Brown } AppCtx; 10c4762a1bSJed Brown 11*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 12c4762a1bSJed Brown PetscFunctionBeginUser; 13c4762a1bSJed Brown options->filename[0] = '\0'; 14c4762a1bSJed Brown options->interpolate = PETSC_FALSE; 15c4762a1bSJed Brown options->meshNum = 0; 16c4762a1bSJed Brown 17d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Hybrid Output Test Options", "DMPLEX"); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-filename", "The mesh file", "ex8.c", options->filename, options->filename, sizeof(options->filename), NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-interpolate", "Interpolate the mesh", "ex8.c", options->interpolate, &options->interpolate, NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-mesh_num", "The mesh we should construct", "ex8.c", options->meshNum, &options->meshNum, NULL, 0)); 21d0609cedSBarry Smith PetscOptionsEnd(); 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionReturn(0); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 26*9371c9d4SSatish Balay static PetscErrorCode CreateHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) { 27c4762a1bSJed Brown PetscInt dim; 28c4762a1bSJed Brown 29c4762a1bSJed Brown PetscFunctionBegin; 30c4762a1bSJed Brown dim = 3; 319566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 329566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Simple Hybrid Mesh")); 339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 349566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 35c4762a1bSJed Brown { 36c4762a1bSJed Brown /* Simple mesh with 2 tets and 1 wedge */ 37c4762a1bSJed Brown PetscInt numPoints[2] = {8, 3}; 38c4762a1bSJed Brown PetscInt coneSize[11] = {4, 4, 6, 0, 0, 0, 0, 0, 0, 0, 0}; 39c4762a1bSJed Brown PetscInt cones[14] = {4, 5, 6, 3, 7, 9, 8, 10, 4, 5, 6, 7, 8, 9}; 40c4762a1bSJed Brown PetscInt coneOrientations[14] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 41*9371c9d4SSatish Balay PetscScalar vertexCoords[48] = {-1.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, -1.0, 0.0, 1.0, 1.0, 1.0, 0.0, 0.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 2.0, 1.0, 0.0}; 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 44c4762a1bSJed Brown if (interpolate) { 45c4762a1bSJed Brown DM idm; 46c4762a1bSJed Brown 479566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 489566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 49c4762a1bSJed Brown *dm = idm; 50c4762a1bSJed Brown } 519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 52c4762a1bSJed Brown } 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56b5a892a1SMatthew G. Knepley /* 57b5a892a1SMatthew G. Knepley This is not a valid mesh. We need to either change to tensor quad prisms or regular triangular prisms. 58b5a892a1SMatthew G. Knepley 59b5a892a1SMatthew G. Knepley 10-------16--------20 60b5a892a1SMatthew G. Knepley /| | 61b5a892a1SMatthew G. Knepley / | | 62b5a892a1SMatthew G. Knepley / | | 63b5a892a1SMatthew G. Knepley 9---|---15 | 64b5a892a1SMatthew G. Knepley /| 7 | 13--------18 65b5a892a1SMatthew G. Knepley / | / | / ____/ 66b5a892a1SMatthew G. Knepley / | / | /____/ 67b5a892a1SMatthew G. Knepley 8 |/ 14---|//---19 68b5a892a1SMatthew G. Knepley | 6 | 12 69b5a892a1SMatthew G. Knepley | / | / \ 70b5a892a1SMatthew G. Knepley | / | / \__ 71b5a892a1SMatthew G. Knepley |/ |/ \ 72b5a892a1SMatthew G. Knepley 5--------11--------17 73b5a892a1SMatthew G. Knepley */ 74*9371c9d4SSatish Balay static PetscErrorCode CreateReverseHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) { 75c4762a1bSJed Brown PetscInt dim; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBegin; 78c4762a1bSJed Brown dim = 3; 799566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Reverse Hybrid Mesh")); 819566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 829566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 83c4762a1bSJed Brown { 84c4762a1bSJed Brown /* Simple mesh with 2 hexes and 3 wedges */ 85c4762a1bSJed Brown PetscInt numPoints[2] = {16, 5}; 86*9371c9d4SSatish Balay PetscInt coneSize[21] = {8, 8, 6, 6, 6, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 87*9371c9d4SSatish Balay PetscInt cones[34] = {5, 6, 12, 11, 8, 14, 15, 9, 6, 7, 13, 12, 9, 15, 16, 10, 11, 17, 12, 14, 19, 15, 12, 18, 13, 15, 20, 16, 12, 17, 18, 15, 19, 20}; 88*9371c9d4SSatish Balay PetscInt coneOrientations[34] = {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}; 89*9371c9d4SSatish Balay PetscScalar vertexCoords[48] = {-1.0, -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 1.0, 0.0, -1.0, -1.0, 1.0, -1.0, 0.0, 1.0, -1.0, 1.0, 1.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 90*9371c9d4SSatish Balay 0.0, 1.0, 0.0, 0.0, -1.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, 1.0, -1.0, 0.0, 1.0, 1.0, 0.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0}; 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 93c4762a1bSJed Brown if (interpolate) { 94c4762a1bSJed Brown DM idm; 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 979566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 98c4762a1bSJed Brown *dm = idm; 99c4762a1bSJed Brown } 1009566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown PetscFunctionReturn(0); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105*9371c9d4SSatish Balay static PetscErrorCode OrderHybridMesh(DM *dm) { 106c4762a1bSJed Brown DM pdm; 107c4762a1bSJed Brown IS perm; 108c4762a1bSJed Brown PetscInt *ind; 109c4762a1bSJed Brown PetscInt dim, pStart, pEnd, p, cStart, cEnd, c, Nhyb = 0, off[2]; 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscFunctionBegin; 1129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 11363a3b9bcSJacob Faibussowitsch PetscCheck(dim == 3, PetscObjectComm((PetscObject)*dm), PETSC_ERR_SUP, "No support for dimension %" PetscInt_FMT, dim); 1149566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(*dm, &pStart, &pEnd)); 1159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd - pStart, &ind)); 116c4762a1bSJed Brown for (p = 0; p < pEnd - pStart; ++p) ind[p] = p; 1179566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 118c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 119c4762a1bSJed Brown PetscInt coneSize; 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(*dm, c, &coneSize)); 122c4762a1bSJed Brown if (coneSize == 6) ++Nhyb; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown off[0] = 0; 125c4762a1bSJed Brown off[1] = cEnd - Nhyb; 126c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 127c4762a1bSJed Brown PetscInt coneSize; 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(*dm, c, &coneSize)); 130c4762a1bSJed Brown if (coneSize == 6) ind[c] = off[1]++; 131c4762a1bSJed Brown else ind[c] = off[0]++; 132c4762a1bSJed Brown } 1331dca8a05SBarry Smith PetscCheck(off[0] == cEnd - Nhyb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of normal cells %" PetscInt_FMT " should be %" PetscInt_FMT, off[0], cEnd - Nhyb); 13463a3b9bcSJacob Faibussowitsch PetscCheck(off[1] == cEnd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of hybrid cells %" PetscInt_FMT " should be %" PetscInt_FMT, off[1] - off[0], Nhyb); 1359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, ind, PETSC_OWN_POINTER, &perm)); 1369566063dSJacob Faibussowitsch PetscCall(DMPlexPermute(*dm, perm, &pdm)); 1379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1389566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 139c4762a1bSJed Brown *dm = pdm; 140c4762a1bSJed Brown PetscFunctionReturn(0); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 144c4762a1bSJed Brown const char *filename = user->filename; 145c4762a1bSJed Brown PetscBool interpolate = user->interpolate; 146c4762a1bSJed Brown PetscInt meshNum = user->meshNum; 147c4762a1bSJed Brown size_t len; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBegin; 1509566063dSJacob Faibussowitsch PetscCall(PetscStrlen(filename, &len)); 151c4762a1bSJed Brown if (len) { 1529566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex34_plex", PETSC_FALSE, dm)); 1539566063dSJacob Faibussowitsch PetscCall(OrderHybridMesh(dm)); 154c4762a1bSJed Brown if (interpolate) { 155c4762a1bSJed Brown DM idm; 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1589566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 159c4762a1bSJed Brown *dm = idm; 160c4762a1bSJed Brown } 1619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Input Mesh")); 1629566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 163c4762a1bSJed Brown } else { 164c4762a1bSJed Brown switch (meshNum) { 165*9371c9d4SSatish Balay case 0: PetscCall(CreateHybridMesh(comm, interpolate, dm)); break; 166*9371c9d4SSatish Balay case 1: PetscCall(CreateReverseHybridMesh(comm, interpolate, dm)); break; 16763a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown mesh number %" PetscInt_FMT, user->meshNum); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown } 170c4762a1bSJed Brown PetscFunctionReturn(0); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173*9371c9d4SSatish Balay int main(int argc, char **argv) { 174c4762a1bSJed Brown DM dm; 175c4762a1bSJed Brown AppCtx user; 176c4762a1bSJed Brown 177327415f7SBarry Smith PetscFunctionBeginUser; 1789566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1799566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1809566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1829566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 183b122ec5aSJacob Faibussowitsch return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown /*TEST 187c4762a1bSJed Brown 188c4762a1bSJed Brown test: 189c4762a1bSJed Brown suffix: 0 190c4762a1bSJed Brown args: -interpolate -dm_view ascii::ascii_info_detail 191c4762a1bSJed Brown 192b5a892a1SMatthew G. Knepley # Test needs to be reworked 193c4762a1bSJed Brown test: 194b5a892a1SMatthew G. Knepley requires: BROKEN 195c4762a1bSJed Brown suffix: 1 196c4762a1bSJed Brown args: -mesh_num 1 -interpolate -dm_view ascii::ascii_info_detail 197c4762a1bSJed Brown 198c4762a1bSJed Brown TEST*/ 199