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 11c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown PetscFunctionBeginUser; 14c4762a1bSJed Brown options->filename[0] = '\0'; 15c4762a1bSJed Brown options->interpolate = PETSC_FALSE; 16c4762a1bSJed Brown options->meshNum = 0; 17c4762a1bSJed Brown 18d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Hybrid Output Test Options", "DMPLEX"); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-filename", "The mesh file", "ex8.c", options->filename, options->filename, sizeof(options->filename), NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-interpolate", "Interpolate the mesh", "ex8.c", options->interpolate, &options->interpolate, NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-mesh_num", "The mesh we should construct", "ex8.c", options->meshNum, &options->meshNum, NULL,0)); 22d0609cedSBarry Smith PetscOptionsEnd(); 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscFunctionReturn(0); 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27c4762a1bSJed Brown static PetscErrorCode CreateHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown PetscInt dim; 30c4762a1bSJed Brown 31c4762a1bSJed Brown PetscFunctionBegin; 32c4762a1bSJed Brown dim = 3; 339566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Simple Hybrid Mesh")); 359566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 369566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 37c4762a1bSJed Brown { 38c4762a1bSJed Brown /* Simple mesh with 2 tets and 1 wedge */ 39c4762a1bSJed Brown PetscInt numPoints[2] = {8, 3}; 40c4762a1bSJed Brown PetscInt coneSize[11] = {4, 4, 6, 0, 0, 0, 0, 0, 0, 0, 0}; 41c4762a1bSJed Brown PetscInt cones[14] = {4, 5, 6, 3, 7, 9, 8, 10, 4, 5, 6, 7, 8, 9}; 42c4762a1bSJed Brown PetscInt coneOrientations[14] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 43c4762a1bSJed Brown PetscScalar vertexCoords[48] = {-1.0, 1.0, 0.0, 44c4762a1bSJed Brown 0.0, 0.0, 0.0, 0.0, 1.0, -1.0, 0.0, 1.0, 1.0, 45c4762a1bSJed Brown 1.0, 0.0, 0.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 46c4762a1bSJed Brown 2.0, 1.0, 0.0}; 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 49c4762a1bSJed Brown if (interpolate) { 50c4762a1bSJed Brown DM idm; 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 539566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 54c4762a1bSJed Brown *dm = idm; 55c4762a1bSJed Brown } 569566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61b5a892a1SMatthew G. Knepley /* 62b5a892a1SMatthew G. Knepley This is not a valid mesh. We need to either change to tensor quad prisms or regular triangular prisms. 63b5a892a1SMatthew G. Knepley 64b5a892a1SMatthew G. Knepley 10-------16--------20 65b5a892a1SMatthew G. Knepley /| | 66b5a892a1SMatthew G. Knepley / | | 67b5a892a1SMatthew G. Knepley / | | 68b5a892a1SMatthew G. Knepley 9---|---15 | 69b5a892a1SMatthew G. Knepley /| 7 | 13--------18 70b5a892a1SMatthew G. Knepley / | / | / ____/ 71b5a892a1SMatthew G. Knepley / | / | /____/ 72b5a892a1SMatthew G. Knepley 8 |/ 14---|//---19 73b5a892a1SMatthew G. Knepley | 6 | 12 74b5a892a1SMatthew G. Knepley | / | / \ 75b5a892a1SMatthew G. Knepley | / | / \__ 76b5a892a1SMatthew G. Knepley |/ |/ \ 77b5a892a1SMatthew G. Knepley 5--------11--------17 78b5a892a1SMatthew G. Knepley */ 79c4762a1bSJed Brown static PetscErrorCode CreateReverseHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt dim; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBegin; 84c4762a1bSJed Brown dim = 3; 859566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Reverse Hybrid Mesh")); 879566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 889566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 89c4762a1bSJed Brown { 90c4762a1bSJed Brown /* Simple mesh with 2 hexes and 3 wedges */ 91c4762a1bSJed Brown PetscInt numPoints[2] = {16, 5}; 92c4762a1bSJed Brown PetscInt coneSize[21] = {8, 8, 6, 6, 6, 93c4762a1bSJed Brown 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 94c4762a1bSJed Brown PetscInt cones[34] = { 5, 6, 12, 11, 8, 14, 15, 9, 95c4762a1bSJed Brown 6, 7, 13, 12, 9, 15, 16, 10, 96c4762a1bSJed Brown 11, 17, 12, 14, 19, 15, 97c4762a1bSJed Brown 12, 18, 13, 15, 20, 16, 98c4762a1bSJed Brown 12, 17, 18, 15, 19, 20}; 99c4762a1bSJed Brown PetscInt coneOrientations[34] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 100c4762a1bSJed Brown 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 101c4762a1bSJed Brown PetscScalar vertexCoords[48] = {-1.0, -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 1.0, 0.0, 102c4762a1bSJed Brown -1.0, -1.0, 1.0, -1.0, 0.0, 1.0, -1.0, 1.0, 1.0, 103c4762a1bSJed Brown 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 104c4762a1bSJed Brown 0.0, -1.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, 105c4762a1bSJed Brown 1.0, -1.0, 0.0, 1.0, 1.0, 0.0, 106c4762a1bSJed Brown 1.0, -1.0, 1.0, 1.0, 1.0, 1.0}; 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 109c4762a1bSJed Brown if (interpolate) { 110c4762a1bSJed Brown DM idm; 111c4762a1bSJed Brown 1129566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1139566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 114c4762a1bSJed Brown *dm = idm; 115c4762a1bSJed Brown } 1169566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown PetscFunctionReturn(0); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown static PetscErrorCode OrderHybridMesh(DM *dm) 122c4762a1bSJed Brown { 123c4762a1bSJed Brown DM pdm; 124c4762a1bSJed Brown IS perm; 125c4762a1bSJed Brown PetscInt *ind; 126c4762a1bSJed Brown PetscInt dim, pStart, pEnd, p, cStart, cEnd, c, Nhyb = 0, off[2]; 127c4762a1bSJed Brown 128c4762a1bSJed Brown PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 13063a3b9bcSJacob Faibussowitsch PetscCheck(dim == 3,PetscObjectComm((PetscObject) *dm), PETSC_ERR_SUP, "No support for dimension %" PetscInt_FMT, dim); 1319566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(*dm, &pStart, &pEnd)); 1329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd-pStart, &ind)); 133c4762a1bSJed Brown for (p = 0; p < pEnd-pStart; ++p) ind[p] = p; 1349566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 135c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 136c4762a1bSJed Brown PetscInt coneSize; 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(*dm, c, &coneSize)); 139c4762a1bSJed Brown if (coneSize == 6) ++Nhyb; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown off[0] = 0; 142c4762a1bSJed Brown off[1] = cEnd - Nhyb; 143c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 144c4762a1bSJed Brown PetscInt coneSize; 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(*dm, c, &coneSize)); 147c4762a1bSJed Brown if (coneSize == 6) ind[c] = off[1]++; 148c4762a1bSJed Brown else ind[c] = off[0]++; 149c4762a1bSJed Brown } 1501dca8a05SBarry 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); 15163a3b9bcSJacob 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); 1529566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd-pStart, ind, PETSC_OWN_POINTER, &perm)); 1539566063dSJacob Faibussowitsch PetscCall(DMPlexPermute(*dm, perm, &pdm)); 1549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1559566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 156c4762a1bSJed Brown *dm = pdm; 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown const char *filename = user->filename; 163c4762a1bSJed Brown PetscBool interpolate = user->interpolate; 164c4762a1bSJed Brown PetscInt meshNum = user->meshNum; 165c4762a1bSJed Brown size_t len; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(PetscStrlen(filename, &len)); 169c4762a1bSJed Brown if (len) { 1709566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex34_plex", PETSC_FALSE, dm)); 1719566063dSJacob Faibussowitsch PetscCall(OrderHybridMesh(dm)); 172c4762a1bSJed Brown if (interpolate) { 173c4762a1bSJed Brown DM idm; 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1769566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 177c4762a1bSJed Brown *dm = idm; 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Input Mesh")); 1809566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 181c4762a1bSJed Brown } else { 182c4762a1bSJed Brown switch (meshNum) { 183c4762a1bSJed Brown case 0: 1849566063dSJacob Faibussowitsch PetscCall(CreateHybridMesh(comm, interpolate, dm));break; 185c4762a1bSJed Brown case 1: 1869566063dSJacob Faibussowitsch PetscCall(CreateReverseHybridMesh(comm, interpolate, dm));break; 18763a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown mesh number %" PetscInt_FMT, user->meshNum); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown } 190c4762a1bSJed Brown PetscFunctionReturn(0); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown int main(int argc, char **argv) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown DM dm; 196c4762a1bSJed Brown AppCtx user; 197c4762a1bSJed Brown 198*327415f7SBarry Smith PetscFunctionBeginUser; 1999566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 2009566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2019566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2029566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 204b122ec5aSJacob Faibussowitsch return 0; 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /*TEST 208c4762a1bSJed Brown 209c4762a1bSJed Brown test: 210c4762a1bSJed Brown suffix: 0 211c4762a1bSJed Brown args: -interpolate -dm_view ascii::ascii_info_detail 212c4762a1bSJed Brown 213b5a892a1SMatthew G. Knepley # Test needs to be reworked 214c4762a1bSJed Brown test: 215b5a892a1SMatthew G. Knepley requires: BROKEN 216c4762a1bSJed Brown suffix: 1 217c4762a1bSJed Brown args: -mesh_num 1 -interpolate -dm_view ascii::ascii_info_detail 218c4762a1bSJed Brown 219c4762a1bSJed Brown TEST*/ 220