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 11d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12d71ae5a4SJacob Faibussowitsch { 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 24*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) 28d71ae5a4SJacob Faibussowitsch { 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}; 439371c9d4SSatish 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}; 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 46c4762a1bSJed Brown if (interpolate) { 47c4762a1bSJed Brown DM idm; 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 509566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 51c4762a1bSJed Brown *dm = idm; 52c4762a1bSJed Brown } 539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 54c4762a1bSJed Brown } 55*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58b5a892a1SMatthew G. Knepley /* 59b5a892a1SMatthew G. Knepley This is not a valid mesh. We need to either change to tensor quad prisms or regular triangular prisms. 60b5a892a1SMatthew G. Knepley 61b5a892a1SMatthew G. Knepley 10-------16--------20 62b5a892a1SMatthew G. Knepley /| | 63b5a892a1SMatthew G. Knepley / | | 64b5a892a1SMatthew G. Knepley / | | 65b5a892a1SMatthew G. Knepley 9---|---15 | 66b5a892a1SMatthew G. Knepley /| 7 | 13--------18 67b5a892a1SMatthew G. Knepley / | / | / ____/ 68b5a892a1SMatthew G. Knepley / | / | /____/ 69b5a892a1SMatthew G. Knepley 8 |/ 14---|//---19 70b5a892a1SMatthew G. Knepley | 6 | 12 71b5a892a1SMatthew G. Knepley | / | / \ 72b5a892a1SMatthew G. Knepley | / | / \__ 73b5a892a1SMatthew G. Knepley |/ |/ \ 74b5a892a1SMatthew G. Knepley 5--------11--------17 75b5a892a1SMatthew G. Knepley */ 76d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateReverseHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm) 77d71ae5a4SJacob Faibussowitsch { 78c4762a1bSJed Brown PetscInt dim; 79c4762a1bSJed Brown 80c4762a1bSJed Brown PetscFunctionBegin; 81c4762a1bSJed Brown dim = 3; 829566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Reverse Hybrid Mesh")); 849566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 859566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 86c4762a1bSJed Brown { 87c4762a1bSJed Brown /* Simple mesh with 2 hexes and 3 wedges */ 88c4762a1bSJed Brown PetscInt numPoints[2] = {16, 5}; 899371c9d4SSatish 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}; 909371c9d4SSatish 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}; 919371c9d4SSatish 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}; 929371c9d4SSatish 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, 939371c9d4SSatish 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}; 94c4762a1bSJed Brown 959566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 96c4762a1bSJed Brown if (interpolate) { 97c4762a1bSJed Brown DM idm; 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1009566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 101c4762a1bSJed Brown *dm = idm; 102c4762a1bSJed Brown } 1039566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 104c4762a1bSJed Brown } 105*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108d71ae5a4SJacob Faibussowitsch static PetscErrorCode OrderHybridMesh(DM *dm) 109d71ae5a4SJacob Faibussowitsch { 110c4762a1bSJed Brown DM pdm; 111c4762a1bSJed Brown IS perm; 112c4762a1bSJed Brown PetscInt *ind; 113c4762a1bSJed Brown PetscInt dim, pStart, pEnd, p, cStart, cEnd, c, Nhyb = 0, off[2]; 114c4762a1bSJed Brown 115c4762a1bSJed Brown PetscFunctionBegin; 1169566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 11763a3b9bcSJacob Faibussowitsch PetscCheck(dim == 3, PetscObjectComm((PetscObject)*dm), PETSC_ERR_SUP, "No support for dimension %" PetscInt_FMT, dim); 1189566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(*dm, &pStart, &pEnd)); 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd - pStart, &ind)); 120c4762a1bSJed Brown for (p = 0; p < pEnd - pStart; ++p) ind[p] = p; 1219566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 122c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 123c4762a1bSJed Brown PetscInt coneSize; 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(*dm, c, &coneSize)); 126c4762a1bSJed Brown if (coneSize == 6) ++Nhyb; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown off[0] = 0; 129c4762a1bSJed Brown off[1] = cEnd - Nhyb; 130c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 131c4762a1bSJed Brown PetscInt coneSize; 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(*dm, c, &coneSize)); 134c4762a1bSJed Brown if (coneSize == 6) ind[c] = off[1]++; 135c4762a1bSJed Brown else ind[c] = off[0]++; 136c4762a1bSJed Brown } 1371dca8a05SBarry 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); 13863a3b9bcSJacob 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); 1399566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, ind, PETSC_OWN_POINTER, &perm)); 1409566063dSJacob Faibussowitsch PetscCall(DMPlexPermute(*dm, perm, &pdm)); 1419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1429566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 143c4762a1bSJed Brown *dm = pdm; 144*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 148d71ae5a4SJacob Faibussowitsch { 149c4762a1bSJed Brown const char *filename = user->filename; 150c4762a1bSJed Brown PetscBool interpolate = user->interpolate; 151c4762a1bSJed Brown PetscInt meshNum = user->meshNum; 152c4762a1bSJed Brown size_t len; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBegin; 1559566063dSJacob Faibussowitsch PetscCall(PetscStrlen(filename, &len)); 156c4762a1bSJed Brown if (len) { 1579566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex34_plex", PETSC_FALSE, dm)); 1589566063dSJacob Faibussowitsch PetscCall(OrderHybridMesh(dm)); 159c4762a1bSJed Brown if (interpolate) { 160c4762a1bSJed Brown DM idm; 161c4762a1bSJed Brown 1629566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 1639566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 164c4762a1bSJed Brown *dm = idm; 165c4762a1bSJed Brown } 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Input Mesh")); 1679566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 168c4762a1bSJed Brown } else { 169c4762a1bSJed Brown switch (meshNum) { 170d71ae5a4SJacob Faibussowitsch case 0: 171d71ae5a4SJacob Faibussowitsch PetscCall(CreateHybridMesh(comm, interpolate, dm)); 172d71ae5a4SJacob Faibussowitsch break; 173d71ae5a4SJacob Faibussowitsch case 1: 174d71ae5a4SJacob Faibussowitsch PetscCall(CreateReverseHybridMesh(comm, interpolate, dm)); 175d71ae5a4SJacob Faibussowitsch break; 176d71ae5a4SJacob Faibussowitsch default: 177d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown mesh number %" PetscInt_FMT, user->meshNum); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown } 180*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 184d71ae5a4SJacob Faibussowitsch { 185c4762a1bSJed Brown DM dm; 186c4762a1bSJed Brown AppCtx user; 187c4762a1bSJed Brown 188327415f7SBarry Smith PetscFunctionBeginUser; 1899566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1909566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1919566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1939566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 194b122ec5aSJacob Faibussowitsch return 0; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /*TEST 198c4762a1bSJed Brown 199c4762a1bSJed Brown test: 200c4762a1bSJed Brown suffix: 0 201c4762a1bSJed Brown args: -interpolate -dm_view ascii::ascii_info_detail 202c4762a1bSJed Brown 203b5a892a1SMatthew G. Knepley # Test needs to be reworked 204c4762a1bSJed Brown test: 205b5a892a1SMatthew G. Knepley requires: BROKEN 206c4762a1bSJed Brown suffix: 1 207c4762a1bSJed Brown args: -mesh_num 1 -interpolate -dm_view ascii::ascii_info_detail 208c4762a1bSJed Brown 209c4762a1bSJed Brown TEST*/ 210