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 PetscErrorCode ierr; 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscFunctionBeginUser; 16c4762a1bSJed Brown options->filename[0] = '\0'; 17c4762a1bSJed Brown options->interpolate = PETSC_FALSE; 18c4762a1bSJed Brown options->meshNum = 0; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Hybrid Output Test Options", "DMPLEX");CHKERRQ(ierr); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-filename", "The mesh file", "ex8.c", options->filename, options->filename, sizeof(options->filename), NULL)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-interpolate", "Interpolate the mesh", "ex8.c", options->interpolate, &options->interpolate, NULL)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-mesh_num", "The mesh we should construct", "ex8.c", options->meshNum, &options->meshNum, NULL,0)); 241e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 25c4762a1bSJed Brown 26c4762a1bSJed Brown PetscFunctionReturn(0); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown 29c4762a1bSJed Brown static PetscErrorCode CreateHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) 30c4762a1bSJed Brown { 31c4762a1bSJed Brown PetscInt dim; 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscFunctionBegin; 34c4762a1bSJed Brown dim = 3; 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Simple Hybrid Mesh")); 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*dm, dim)); 39c4762a1bSJed Brown { 40c4762a1bSJed Brown /* Simple mesh with 2 tets and 1 wedge */ 41c4762a1bSJed Brown PetscInt numPoints[2] = {8, 3}; 42c4762a1bSJed Brown PetscInt coneSize[11] = {4, 4, 6, 0, 0, 0, 0, 0, 0, 0, 0}; 43c4762a1bSJed Brown PetscInt cones[14] = {4, 5, 6, 3, 7, 9, 8, 10, 4, 5, 6, 7, 8, 9}; 44c4762a1bSJed Brown PetscInt coneOrientations[14] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 45c4762a1bSJed Brown PetscScalar vertexCoords[48] = {-1.0, 1.0, 0.0, 46c4762a1bSJed Brown 0.0, 0.0, 0.0, 0.0, 1.0, -1.0, 0.0, 1.0, 1.0, 47c4762a1bSJed Brown 1.0, 0.0, 0.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 48c4762a1bSJed Brown 2.0, 1.0, 0.0}; 49c4762a1bSJed Brown 505f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 51c4762a1bSJed Brown if (interpolate) { 52c4762a1bSJed Brown DM idm; 53c4762a1bSJed Brown 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 56c4762a1bSJed Brown *dm = idm; 57c4762a1bSJed Brown } 585f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown PetscFunctionReturn(0); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63b5a892a1SMatthew G. Knepley /* 64b5a892a1SMatthew G. Knepley This is not a valid mesh. We need to either change to tensor quad prisms or regular triangular prisms. 65b5a892a1SMatthew G. Knepley 66b5a892a1SMatthew G. Knepley 10-------16--------20 67b5a892a1SMatthew G. Knepley /| | 68b5a892a1SMatthew G. Knepley / | | 69b5a892a1SMatthew G. Knepley / | | 70b5a892a1SMatthew G. Knepley 9---|---15 | 71b5a892a1SMatthew G. Knepley /| 7 | 13--------18 72b5a892a1SMatthew G. Knepley / | / | / ____/ 73b5a892a1SMatthew G. Knepley / | / | /____/ 74b5a892a1SMatthew G. Knepley 8 |/ 14---|//---19 75b5a892a1SMatthew G. Knepley | 6 | 12 76b5a892a1SMatthew G. Knepley | / | / \ 77b5a892a1SMatthew G. Knepley | / | / \__ 78b5a892a1SMatthew G. Knepley |/ |/ \ 79b5a892a1SMatthew G. Knepley 5--------11--------17 80b5a892a1SMatthew G. Knepley */ 81c4762a1bSJed Brown static PetscErrorCode CreateReverseHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) 82c4762a1bSJed Brown { 83c4762a1bSJed Brown PetscInt dim; 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscFunctionBegin; 86c4762a1bSJed Brown dim = 3; 875f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Reverse Hybrid Mesh")); 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*dm, dim)); 91c4762a1bSJed Brown { 92c4762a1bSJed Brown /* Simple mesh with 2 hexes and 3 wedges */ 93c4762a1bSJed Brown PetscInt numPoints[2] = {16, 5}; 94c4762a1bSJed Brown PetscInt coneSize[21] = {8, 8, 6, 6, 6, 95c4762a1bSJed Brown 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 96c4762a1bSJed Brown PetscInt cones[34] = { 5, 6, 12, 11, 8, 14, 15, 9, 97c4762a1bSJed Brown 6, 7, 13, 12, 9, 15, 16, 10, 98c4762a1bSJed Brown 11, 17, 12, 14, 19, 15, 99c4762a1bSJed Brown 12, 18, 13, 15, 20, 16, 100c4762a1bSJed Brown 12, 17, 18, 15, 19, 20}; 101c4762a1bSJed Brown PetscInt coneOrientations[34] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 102c4762a1bSJed Brown 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 103c4762a1bSJed Brown PetscScalar vertexCoords[48] = {-1.0, -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 1.0, 0.0, 104c4762a1bSJed Brown -1.0, -1.0, 1.0, -1.0, 0.0, 1.0, -1.0, 1.0, 1.0, 105c4762a1bSJed Brown 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 106c4762a1bSJed Brown 0.0, -1.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, 107c4762a1bSJed Brown 1.0, -1.0, 0.0, 1.0, 1.0, 0.0, 108c4762a1bSJed Brown 1.0, -1.0, 1.0, 1.0, 1.0, 1.0}; 109c4762a1bSJed Brown 1105f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 111c4762a1bSJed Brown if (interpolate) { 112c4762a1bSJed Brown DM idm; 113c4762a1bSJed Brown 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 116c4762a1bSJed Brown *dm = idm; 117c4762a1bSJed Brown } 1185f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown PetscFunctionReturn(0); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown static PetscErrorCode OrderHybridMesh(DM *dm) 124c4762a1bSJed Brown { 125c4762a1bSJed Brown DM pdm; 126c4762a1bSJed Brown IS perm; 127c4762a1bSJed Brown PetscInt *ind; 128c4762a1bSJed Brown PetscInt dim, pStart, pEnd, p, cStart, cEnd, c, Nhyb = 0, off[2]; 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscFunctionBegin; 1315f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*dm, &dim)); 1322c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim != 3,PetscObjectComm((PetscObject) *dm), PETSC_ERR_SUP, "No support for dimension %D", dim); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetChart(*dm, &pStart, &pEnd)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(pEnd-pStart, &ind)); 135c4762a1bSJed Brown for (p = 0; p < pEnd-pStart; ++p) ind[p] = p; 1365f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 137c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 138c4762a1bSJed Brown PetscInt coneSize; 139c4762a1bSJed Brown 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeSize(*dm, c, &coneSize)); 141c4762a1bSJed Brown if (coneSize == 6) ++Nhyb; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown off[0] = 0; 144c4762a1bSJed Brown off[1] = cEnd - Nhyb; 145c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 146c4762a1bSJed Brown PetscInt coneSize; 147c4762a1bSJed Brown 1485f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeSize(*dm, c, &coneSize)); 149c4762a1bSJed Brown if (coneSize == 6) ind[c] = off[1]++; 150c4762a1bSJed Brown else ind[c] = off[0]++; 151c4762a1bSJed Brown } 1522c71b3e2SJacob Faibussowitsch PetscCheckFalse(off[0] != cEnd - Nhyb,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of normal cells %D should be %D", off[0], cEnd - Nhyb); 1532c71b3e2SJacob Faibussowitsch PetscCheckFalse(off[1] != cEnd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of hybrid cells %D should be %D", off[1] - off[0], Nhyb); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF, pEnd-pStart, ind, PETSC_OWN_POINTER, &perm)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPermute(*dm, perm, &pdm)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 158c4762a1bSJed Brown *dm = pdm; 159c4762a1bSJed Brown PetscFunctionReturn(0); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 163c4762a1bSJed Brown { 164c4762a1bSJed Brown const char *filename = user->filename; 165c4762a1bSJed Brown PetscBool interpolate = user->interpolate; 166c4762a1bSJed Brown PetscInt meshNum = user->meshNum; 167c4762a1bSJed Brown size_t len; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBegin; 1705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlen(filename, &len)); 171c4762a1bSJed Brown if (len) { 1725f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromFile(comm, filename, "ex34_plex", PETSC_FALSE, dm)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(OrderHybridMesh(dm)); 174c4762a1bSJed Brown if (interpolate) { 175c4762a1bSJed Brown DM idm; 176c4762a1bSJed Brown 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 179c4762a1bSJed Brown *dm = idm; 180c4762a1bSJed Brown } 1815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Input Mesh")); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 183c4762a1bSJed Brown } else { 184c4762a1bSJed Brown switch (meshNum) { 185c4762a1bSJed Brown case 0: 1865f80ce2aSJacob Faibussowitsch CHKERRQ(CreateHybridMesh(comm, interpolate, dm));break; 187c4762a1bSJed Brown case 1: 1885f80ce2aSJacob Faibussowitsch CHKERRQ(CreateReverseHybridMesh(comm, interpolate, dm));break; 18998921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown mesh number %D", user->meshNum); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown } 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown int main(int argc, char **argv) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown DM dm; 198c4762a1bSJed Brown AppCtx user; 199c4762a1bSJed Brown 200*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 204*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 205*b122ec5aSJacob Faibussowitsch return 0; 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown /*TEST 209c4762a1bSJed Brown 210c4762a1bSJed Brown test: 211c4762a1bSJed Brown suffix: 0 212c4762a1bSJed Brown args: -interpolate -dm_view ascii::ascii_info_detail 213c4762a1bSJed Brown 214b5a892a1SMatthew G. Knepley # Test needs to be reworked 215c4762a1bSJed Brown test: 216b5a892a1SMatthew G. Knepley requires: BROKEN 217c4762a1bSJed Brown suffix: 1 218c4762a1bSJed Brown args: -mesh_num 1 -interpolate -dm_view ascii::ascii_info_detail 219c4762a1bSJed Brown 220c4762a1bSJed Brown TEST*/ 221