1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3331830f3SMatthew G. Knepley 47d282ae0SMichael Lange /*@C 57d282ae0SMichael Lange DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 67d282ae0SMichael Lange 77d282ae0SMichael Lange + comm - The MPI communicator 87d282ae0SMichael Lange . filename - Name of the Gmsh file 97d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh 107d282ae0SMichael Lange 117d282ae0SMichael Lange Output Parameter: 127d282ae0SMichael Lange . dm - The DM object representing the mesh 137d282ae0SMichael Lange 147d282ae0SMichael Lange Level: beginner 157d282ae0SMichael Lange 167d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 177d282ae0SMichael Lange @*/ 187d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 197d282ae0SMichael Lange { 2011c56cb0SLisandro Dalcin PetscViewer viewer; 21abc86ac4SMichael Lange PetscMPIInt rank; 2211c56cb0SLisandro Dalcin int fileType; 237d282ae0SMichael Lange PetscViewerType vtype; 247d282ae0SMichael Lange PetscErrorCode ierr; 257d282ae0SMichael Lange 267d282ae0SMichael Lange PetscFunctionBegin; 27abc86ac4SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2811c56cb0SLisandro Dalcin 297d282ae0SMichael Lange /* Determine Gmsh file type (ASCII or binary) from file header */ 3011c56cb0SLisandro Dalcin if (!rank) { 3111c56cb0SLisandro Dalcin PetscViewer vheader; 3211c56cb0SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 3311c56cb0SLisandro Dalcin PetscBool match; 3411c56cb0SLisandro Dalcin int snum; 3511c56cb0SLisandro Dalcin float version; 3611c56cb0SLisandro Dalcin 3711c56cb0SLisandro Dalcin ierr = PetscViewerCreate(PETSC_COMM_SELF, &vheader);CHKERRQ(ierr); 387d282ae0SMichael Lange ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr); 397d282ae0SMichael Lange ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr); 407d282ae0SMichael Lange ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr); 417d282ae0SMichael Lange /* Read only the first two lines of the Gmsh file */ 42060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 4311c56cb0SLisandro Dalcin ierr = PetscStrncmp(line, "$MeshFormat", sizeof(line), &match);CHKERRQ(ierr); 447d282ae0SMichael Lange if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 45060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 4611c56cb0SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 47f6ac7a6aSMichael Lange if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 48f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 4911c56cb0SLisandro Dalcin ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 50abc86ac4SMichael Lange } 5111c56cb0SLisandro Dalcin ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 5211c56cb0SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 5311c56cb0SLisandro Dalcin 547d282ae0SMichael Lange /* Create appropriate viewer and build plex */ 557d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 567d282ae0SMichael Lange ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 577d282ae0SMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 587d282ae0SMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 597d282ae0SMichael Lange ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 607d282ae0SMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 617d282ae0SMichael Lange PetscFunctionReturn(0); 627d282ae0SMichael Lange } 637d282ae0SMichael Lange 6411c56cb0SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int numVertices, double **coordinates) 65df032b82SMatthew G. Knepley { 6611c56cb0SLisandro Dalcin int v,nid; 6711c56cb0SLisandro Dalcin PetscErrorCode ierr; 6811c56cb0SLisandro Dalcin 6911c56cb0SLisandro Dalcin PetscFunctionBegin; 7011c56cb0SLisandro Dalcin ierr = PetscMalloc1(numVertices*3, coordinates);CHKERRQ(ierr); 7111c56cb0SLisandro Dalcin for (v = 0; v < numVertices; ++v) { 7211c56cb0SLisandro Dalcin double *xyz = *coordinates + v*3; 7311c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 7411c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 7511c56cb0SLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 7611c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 7711c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 7811c56cb0SLisandro Dalcin } 7911c56cb0SLisandro Dalcin PetscFunctionReturn(0); 8011c56cb0SLisandro Dalcin } 8111c56cb0SLisandro Dalcin 8211c56cb0SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, int numCells, GmshElement **gmsh_elems) 8311c56cb0SLisandro Dalcin { 84df032b82SMatthew G. Knepley GmshElement *elements; 8511c56cb0SLisandro Dalcin int i, c, p, ibuf[1+4+512]; 8611c56cb0SLisandro Dalcin int cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 87df032b82SMatthew G. Knepley PetscErrorCode ierr; 88df032b82SMatthew G. Knepley 89df032b82SMatthew G. Knepley PetscFunctionBegin; 90df032b82SMatthew G. Knepley ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); 91df032b82SMatthew G. Knepley for (c = 0; c < numCells;) { 9211c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 9311c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 94df032b82SMatthew G. Knepley if (binary) { 95df032b82SMatthew G. Knepley cellType = ibuf[0]; 96df032b82SMatthew G. Knepley numElem = ibuf[1]; 97df032b82SMatthew G. Knepley numTags = ibuf[2]; 98df032b82SMatthew G. Knepley } else { 99df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 100df032b82SMatthew G. Knepley cellType = ibuf[1]; 101df032b82SMatthew G. Knepley numTags = ibuf[2]; 102df032b82SMatthew G. Knepley numElem = 1; 103df032b82SMatthew G. Knepley } 104bf6ba3a3SMatthew G. Knepley /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 105bf6ba3a3SMatthew G. Knepley numNodesIgnore = 0; 106df032b82SMatthew G. Knepley switch (cellType) { 107df032b82SMatthew G. Knepley case 1: /* 2-node line */ 108df032b82SMatthew G. Knepley dim = 1; 109df032b82SMatthew G. Knepley numNodes = 2; 110df032b82SMatthew G. Knepley break; 111df032b82SMatthew G. Knepley case 2: /* 3-node triangle */ 112df032b82SMatthew G. Knepley dim = 2; 113df032b82SMatthew G. Knepley numNodes = 3; 114df032b82SMatthew G. Knepley break; 115df032b82SMatthew G. Knepley case 3: /* 4-node quadrangle */ 116df032b82SMatthew G. Knepley dim = 2; 117df032b82SMatthew G. Knepley numNodes = 4; 118df032b82SMatthew G. Knepley break; 119df032b82SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 120df032b82SMatthew G. Knepley dim = 3; 121df032b82SMatthew G. Knepley numNodes = 4; 122df032b82SMatthew G. Knepley break; 123df032b82SMatthew G. Knepley case 5: /* 8-node hexahedron */ 124df032b82SMatthew G. Knepley dim = 3; 125df032b82SMatthew G. Knepley numNodes = 8; 126df032b82SMatthew G. Knepley break; 127*dea2a3c7SStefano Zampini case 6: /* 6-node wedge */ 128*dea2a3c7SStefano Zampini dim = 3; 129*dea2a3c7SStefano Zampini numNodes = 6; 130*dea2a3c7SStefano Zampini break; 131bf6ba3a3SMatthew G. Knepley case 8: /* 3-node 2nd order line */ 132bf6ba3a3SMatthew G. Knepley dim = 1; 133bf6ba3a3SMatthew G. Knepley numNodes = 2; 134bf6ba3a3SMatthew G. Knepley numNodesIgnore = 1; 135bf6ba3a3SMatthew G. Knepley break; 136bf6ba3a3SMatthew G. Knepley case 9: /* 6-node 2nd order triangle */ 137bf6ba3a3SMatthew G. Knepley dim = 2; 138bf6ba3a3SMatthew G. Knepley numNodes = 3; 139bf6ba3a3SMatthew G. Knepley numNodesIgnore = 3; 140bf6ba3a3SMatthew G. Knepley break; 141*dea2a3c7SStefano Zampini case 13: /* 18-node 2nd wedge */ 142*dea2a3c7SStefano Zampini dim = 3; 143*dea2a3c7SStefano Zampini numNodes = 6; 144*dea2a3c7SStefano Zampini numNodesIgnore = 12; 145*dea2a3c7SStefano Zampini break; 146df032b82SMatthew G. Knepley case 15: /* 1-node vertex */ 147df032b82SMatthew G. Knepley dim = 0; 148df032b82SMatthew G. Knepley numNodes = 1; 149df032b82SMatthew G. Knepley break; 150bf6ba3a3SMatthew G. Knepley case 7: /* 5-node pyramid */ 151bf6ba3a3SMatthew G. Knepley case 10: /* 9-node 2nd order quadrangle */ 152bf6ba3a3SMatthew G. Knepley case 11: /* 10-node 2nd order tetrahedron */ 153bf6ba3a3SMatthew G. Knepley case 12: /* 27-node 2nd order hexhedron */ 154bf6ba3a3SMatthew G. Knepley case 14: /* 14-node 2nd order pyramid */ 155df032b82SMatthew G. Knepley default: 156df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 157df032b82SMatthew G. Knepley } 158df032b82SMatthew G. Knepley if (binary) { 15911c56cb0SLisandro Dalcin const int nint = 1 + numTags + numNodes + numNodesIgnore; 16011c56cb0SLisandro Dalcin /* Loop over element blocks */ 161df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 16211c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 16311c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 164df032b82SMatthew G. Knepley elements[c].dim = dim; 165df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 166df032b82SMatthew G. Knepley elements[c].numTags = numTags; 167df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 168*dea2a3c7SStefano Zampini elements[c].cellType = cellType; 169df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 170df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 171df032b82SMatthew G. Knepley } 172df032b82SMatthew G. Knepley } else { 173df032b82SMatthew G. Knepley elements[c].dim = dim; 174df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 175df032b82SMatthew G. Knepley elements[c].numTags = numTags; 176*dea2a3c7SStefano Zampini elements[c].cellType = cellType; 177df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 178df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 17911c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr); 180df032b82SMatthew G. Knepley c++; 181df032b82SMatthew G. Knepley } 182df032b82SMatthew G. Knepley } 183df032b82SMatthew G. Knepley *gmsh_elems = elements; 184df032b82SMatthew G. Knepley PetscFunctionReturn(0); 185df032b82SMatthew G. Knepley } 1867d282ae0SMichael Lange 187331830f3SMatthew G. Knepley /*@ 1887d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 189331830f3SMatthew G. Knepley 190331830f3SMatthew G. Knepley Collective on comm 191331830f3SMatthew G. Knepley 192331830f3SMatthew G. Knepley Input Parameters: 193331830f3SMatthew G. Knepley + comm - The MPI communicator 194331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 195331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 196331830f3SMatthew G. Knepley 197331830f3SMatthew G. Knepley Output Parameter: 198331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 199331830f3SMatthew G. Knepley 200331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 2013d51f2daSMichael Lange and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 202331830f3SMatthew G. Knepley 203331830f3SMatthew G. Knepley Level: beginner 204331830f3SMatthew G. Knepley 205331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 206331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 207331830f3SMatthew G. Knepley @*/ 208331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 209331830f3SMatthew G. Knepley { 21011c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 21111c56cb0SLisandro Dalcin double *coordsIn = NULL; 2123c67d5adSSatish Balay GmshElement *gmsh_elem = NULL; 213331830f3SMatthew G. Knepley PetscSection coordSection; 214331830f3SMatthew G. Knepley Vec coordinates; 2156fbe17bfSStefano Zampini PetscBT periodicV = NULL, periodicC = NULL; 21684572febSToby Isaac PetscScalar *coords; 217*dea2a3c7SStefano Zampini PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 2181d64f2ccSMatthew G. Knepley int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 21911c56cb0SLisandro Dalcin PetscMPIInt rank; 220331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 221*dea2a3c7SStefano Zampini PetscBool binary, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 222*dea2a3c7SStefano Zampini PetscBool enable_hybrid = PETSC_FALSE; 223331830f3SMatthew G. Knepley PetscErrorCode ierr; 224331830f3SMatthew G. Knepley 225331830f3SMatthew G. Knepley PetscFunctionBegin; 226331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 227*dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr); 228*dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 229*dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 230*dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 231*dea2a3c7SStefano Zampini ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr); 23211c56cb0SLisandro Dalcin if (zerobase) shift = 0; 23311c56cb0SLisandro Dalcin 234331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 235331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2363b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 23711c56cb0SLisandro Dalcin 23811c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 23911c56cb0SLisandro Dalcin 24011c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 2413b46f529SLisandro Dalcin if (binary) { 24211c56cb0SLisandro Dalcin parentviewer = viewer; 24311c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 24411c56cb0SLisandro Dalcin } 24511c56cb0SLisandro Dalcin 2463b46f529SLisandro Dalcin if (!rank) { 247*dea2a3c7SStefano Zampini PetscBool match, hybrid; 248331830f3SMatthew G. Knepley int fileType, dataSize; 249f6ac7a6aSMichael Lange float version; 250331830f3SMatthew G. Knepley 251331830f3SMatthew G. Knepley /* Read header */ 252060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 25304d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 254331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 255060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 256f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 257f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 258f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 259331830f3SMatthew G. Knepley if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 26004d1ad83SMichael Lange if (binary) { 261b9eae255SMichael Lange int checkInt; 262060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 26304d1ad83SMichael Lange if (checkInt != 1) { 264b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 26511c56cb0SLisandro Dalcin if (checkInt == 1) byteSwap = PETSC_TRUE; 26604d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 26704d1ad83SMichael Lange } 2680877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 269060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 27004d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 271331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 27211c56cb0SLisandro Dalcin 273dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 274060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 275dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 276dc0ede3bSMatthew G. Knepley if (match) { 277dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 278dc0ede3bSMatthew G. Knepley snum = sscanf(line, "%d", &numRegions); 279dc0ede3bSMatthew G. Knepley if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 28011c56cb0SLisandro Dalcin for (r = 0; r < numRegions; ++r) { 28111c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 28211c56cb0SLisandro Dalcin } 283dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 284dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 285dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 286dc0ede3bSMatthew G. Knepley /* Initial read for vertex section */ 287dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 288dc0ede3bSMatthew G. Knepley } 28911c56cb0SLisandro Dalcin 290dc0ede3bSMatthew G. Knepley /* Read vertices */ 29104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 292331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 293060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 2940877b519SMichael Lange snum = sscanf(line, "%d", &numVertices); 2950877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 29611c56cb0SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadNodes(viewer, binary, byteSwap, shift, numVertices, &coordsIn);CHKERRQ(ierr); 297060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 29804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 299331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 30011c56cb0SLisandro Dalcin 301331830f3SMatthew G. Knepley /* Read cells */ 302060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 30304d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 304331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 305060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 3060877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 3070877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 308a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 309a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 310a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 311a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 31211c56cb0SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadElements(viewer, binary, byteSwap, shift, numCells, &gmsh_elem);CHKERRQ(ierr); 313*dea2a3c7SStefano Zampini hybrid = PETSC_FALSE; 314a4bb7517SMichael Lange for (trueNumCells = 0, c = 0; c < numCells; ++c) { 315*dea2a3c7SStefano Zampini int on = -1; 316a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 317*dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) {hybrid = (trueNumCells ? (on != gmsh_elem[c].cellType ? on = gmsh_elem[c].cellType,PETSC_TRUE : hybrid) : (on = gmsh_elem[c].cellType, PETSC_FALSE) ); trueNumCells++;} 318*dea2a3c7SStefano Zampini /* wedges always indicate an hybrid mesh in PLEX */ 319*dea2a3c7SStefano Zampini if (on == 6 || on == 13) hybrid = PETSC_TRUE; 320db397164SMichael Lange } 321d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 3221b82c3ebSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 323331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 32411c56cb0SLisandro Dalcin 325*dea2a3c7SStefano Zampini /* Renumber cells for hybrid grids */ 326*dea2a3c7SStefano Zampini if (hybrid && enable_hybrid) { 327*dea2a3c7SStefano Zampini PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 328*dea2a3c7SStefano Zampini PetscInt cell, tn, *tp; 329*dea2a3c7SStefano Zampini int n1 = 0,n2 = 0; 330*dea2a3c7SStefano Zampini 331*dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 332*dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 333*dea2a3c7SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 334*dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) { 335*dea2a3c7SStefano Zampini if (!n1) n1 = gmsh_elem[c].cellType; 336*dea2a3c7SStefano Zampini else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 337*dea2a3c7SStefano Zampini 338*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 339*dea2a3c7SStefano Zampini else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 340*dea2a3c7SStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 341*dea2a3c7SStefano Zampini cell++; 342*dea2a3c7SStefano Zampini } 343*dea2a3c7SStefano Zampini } 344*dea2a3c7SStefano Zampini 345*dea2a3c7SStefano Zampini switch (n1) { 346*dea2a3c7SStefano Zampini case 2: /* triangles */ 347*dea2a3c7SStefano Zampini case 9: 348*dea2a3c7SStefano Zampini switch (n2) { 349*dea2a3c7SStefano Zampini case 0: /* single type mesh */ 350*dea2a3c7SStefano Zampini case 3: /* quads */ 351*dea2a3c7SStefano Zampini case 10: 352*dea2a3c7SStefano Zampini break; 353*dea2a3c7SStefano Zampini default: 354*dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 355*dea2a3c7SStefano Zampini } 356*dea2a3c7SStefano Zampini break; 357*dea2a3c7SStefano Zampini case 3: 358*dea2a3c7SStefano Zampini case 10: 359*dea2a3c7SStefano Zampini switch (n2) { 360*dea2a3c7SStefano Zampini case 0: /* single type mesh */ 361*dea2a3c7SStefano Zampini case 2: /* swap since we list simplices first */ 362*dea2a3c7SStefano Zampini case 9: 363*dea2a3c7SStefano Zampini tn = hc1; 364*dea2a3c7SStefano Zampini hc1 = hc2; 365*dea2a3c7SStefano Zampini hc2 = tn; 366*dea2a3c7SStefano Zampini 367*dea2a3c7SStefano Zampini tp = hybridCells1; 368*dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 369*dea2a3c7SStefano Zampini hybridCells2 = tp; 370*dea2a3c7SStefano Zampini break; 371*dea2a3c7SStefano Zampini default: 372*dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 373*dea2a3c7SStefano Zampini } 374*dea2a3c7SStefano Zampini break; 375*dea2a3c7SStefano Zampini case 4: /* tetrahedra */ 376*dea2a3c7SStefano Zampini case 11: 377*dea2a3c7SStefano Zampini switch (n2) { 378*dea2a3c7SStefano Zampini case 0: /* single type mesh */ 379*dea2a3c7SStefano Zampini case 6: /* wedges */ 380*dea2a3c7SStefano Zampini case 13: 381*dea2a3c7SStefano Zampini break; 382*dea2a3c7SStefano Zampini default: 383*dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 384*dea2a3c7SStefano Zampini } 385*dea2a3c7SStefano Zampini break; 386*dea2a3c7SStefano Zampini case 6: 387*dea2a3c7SStefano Zampini case 13: 388*dea2a3c7SStefano Zampini switch (n2) { 389*dea2a3c7SStefano Zampini case 0: /* single type mesh */ 390*dea2a3c7SStefano Zampini case 4: /* swap since we list simplices first */ 391*dea2a3c7SStefano Zampini case 11: 392*dea2a3c7SStefano Zampini tn = hc1; 393*dea2a3c7SStefano Zampini hc1 = hc2; 394*dea2a3c7SStefano Zampini hc2 = tn; 395*dea2a3c7SStefano Zampini 396*dea2a3c7SStefano Zampini tp = hybridCells1; 397*dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 398*dea2a3c7SStefano Zampini hybridCells2 = tp; 399*dea2a3c7SStefano Zampini break; 400*dea2a3c7SStefano Zampini default: 401*dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 402*dea2a3c7SStefano Zampini } 403*dea2a3c7SStefano Zampini break; 404*dea2a3c7SStefano Zampini default: 405*dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 406*dea2a3c7SStefano Zampini } 407*dea2a3c7SStefano Zampini cMax = hc1; 408*dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 409*dea2a3c7SStefano Zampini for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 410*dea2a3c7SStefano Zampini for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 411*dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 412*dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 413*dea2a3c7SStefano Zampini } 414*dea2a3c7SStefano Zampini 41511c56cb0SLisandro Dalcin /* OPTIONAL Read periodic section */ 416d08df55aSStefano Zampini if (periodic) { 417f45c9589SStefano Zampini PetscInt pVert, *periodicMapT, *aux; 418d08df55aSStefano Zampini int numPeriodic; 419d08df55aSStefano Zampini 420d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 421d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 422d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 423f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 4246fbe17bfSStefano Zampini ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 425f45c9589SStefano Zampini for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 426d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 427d08df55aSStefano Zampini snum = sscanf(line, "%d", &numPeriodic); 428d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 429d08df55aSStefano Zampini for (i = 0; i < numPeriodic; i++) { 430da83f57bSLisandro Dalcin int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes, slaveNode, masterNode; 431d08df55aSStefano Zampini 432d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 433d08df55aSStefano Zampini snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */ 434d08df55aSStefano Zampini if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 435da83f57bSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 436da83f57bSLisandro Dalcin snum = sscanf(line, "%d", &nNodes); 437da83f57bSLisandro Dalcin if (snum != 1) { /* discard tranformation and try again */ 43872ffbcc9SStefano Zampini ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 439d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 440d08df55aSStefano Zampini snum = sscanf(line, "%d", &nNodes); 441d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 442da83f57bSLisandro Dalcin } 443d08df55aSStefano Zampini for (j = 0; j < nNodes; j++) { 444d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 44511c56cb0SLisandro Dalcin snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 446d08df55aSStefano Zampini if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 447917266a4SLisandro Dalcin periodicMapT[slaveNode - shift] = masterNode - shift; 448917266a4SLisandro Dalcin ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr); 449917266a4SLisandro Dalcin ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr); 450d08df55aSStefano Zampini } 451d08df55aSStefano Zampini } 452d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 453d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 454d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 455d08df55aSStefano Zampini /* we may have slaves of slaves */ 456d08df55aSStefano Zampini for (i = 0; i < numVertices; i++) { 457f45c9589SStefano Zampini while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 458f45c9589SStefano Zampini periodicMapT[i] = periodicMapT[periodicMapT[i]]; 459d08df55aSStefano Zampini } 460d08df55aSStefano Zampini } 461f45c9589SStefano Zampini /* periodicMap : from old to new numbering (periodic vertices excluded) 462f45c9589SStefano Zampini periodicMapI: from new to old numbering */ 463f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 464f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 465f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 466f45c9589SStefano Zampini for (i = 0, pVert = 0; i < numVertices; i++) { 467f45c9589SStefano Zampini if (periodicMapT[i] != i) { 468f45c9589SStefano Zampini pVert++; 469f45c9589SStefano Zampini } else { 470f45c9589SStefano Zampini aux[i] = i - pVert; 471f45c9589SStefano Zampini periodicMapI[i - pVert] = i; 472f45c9589SStefano Zampini } 473f45c9589SStefano Zampini } 474f45c9589SStefano Zampini for (i = 0 ; i < numVertices; i++) { 475f45c9589SStefano Zampini periodicMap[i] = aux[periodicMapT[i]]; 476f45c9589SStefano Zampini } 477f45c9589SStefano Zampini ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 478f45c9589SStefano Zampini ierr = PetscFree(aux);CHKERRQ(ierr); 479f45c9589SStefano Zampini /* remove periodic vertices */ 480f45c9589SStefano Zampini numVertices = numVertices - pVert; 481d08df55aSStefano Zampini } 482331830f3SMatthew G. Knepley } 48311c56cb0SLisandro Dalcin 48411c56cb0SLisandro Dalcin if (parentviewer) { 48511c56cb0SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 48611c56cb0SLisandro Dalcin } 48711c56cb0SLisandro Dalcin 488a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 489331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 490a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 491a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 492*dea2a3c7SStefano Zampini ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 493a4bb7517SMichael Lange cell++; 494331830f3SMatthew G. Knepley } 495331830f3SMatthew G. Knepley } 496331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 497a4bb7517SMichael Lange /* Add cell-vertex connections */ 498a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 499a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 50011c56cb0SLisandro Dalcin PetscInt pcone[8], corner; 501a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 502dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 503917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 504331830f3SMatthew G. Knepley } 50597ed6be6Ssarens if (dim == 3) { 50697ed6be6Ssarens /* Tetrahedra are inverted */ 507*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 50897ed6be6Ssarens PetscInt tmp = pcone[0]; 50997ed6be6Ssarens pcone[0] = pcone[1]; 51097ed6be6Ssarens pcone[1] = tmp; 51197ed6be6Ssarens } 51297ed6be6Ssarens /* Hexahedra are inverted */ 513*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 51497ed6be6Ssarens PetscInt tmp = pcone[1]; 51597ed6be6Ssarens pcone[1] = pcone[3]; 51697ed6be6Ssarens pcone[3] = tmp; 51797ed6be6Ssarens } 518*dea2a3c7SStefano Zampini /* Prisms are inverted */ 519*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 520*dea2a3c7SStefano Zampini PetscInt tmp; 521*dea2a3c7SStefano Zampini 522*dea2a3c7SStefano Zampini tmp = pcone[1]; 523*dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 524*dea2a3c7SStefano Zampini pcone[2] = tmp; 525*dea2a3c7SStefano Zampini tmp = pcone[4]; 526*dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 527*dea2a3c7SStefano Zampini pcone[5] = tmp; 52897ed6be6Ssarens } 529*dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 530*dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 531*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 532*dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 533*dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 534*dea2a3c7SStefano Zampini pcone[3] = tmp; 535*dea2a3c7SStefano Zampini } 536*dea2a3c7SStefano Zampini } 537*dea2a3c7SStefano Zampini ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 538a4bb7517SMichael Lange cell++; 539331830f3SMatthew G. Knepley } 540a4bb7517SMichael Lange } 5416228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 542c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 543*dea2a3c7SStefano Zampini ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 544331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 545331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 546331830f3SMatthew G. Knepley if (interpolate) { 5475fd9971aSMatthew G. Knepley DM idm; 548331830f3SMatthew G. Knepley 549331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 550331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 551331830f3SMatthew G. Knepley *dm = idm; 552331830f3SMatthew G. Knepley } 553917266a4SLisandro Dalcin 554f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 555917266a4SLisandro Dalcin if (!rank && usemarker) { 556d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 557d3f73514SStefano Zampini 558d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 559d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 560d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 561d3f73514SStefano Zampini PetscInt suppSize; 562d3f73514SStefano Zampini 563d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 564d3f73514SStefano Zampini if (suppSize == 1) { 565d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 566d3f73514SStefano Zampini 567d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 568d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 569d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 570d3f73514SStefano Zampini } 571d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 572d3f73514SStefano Zampini } 573d3f73514SStefano Zampini } 574d3f73514SStefano Zampini } 57516ad7e67SMichael Lange 57616ad7e67SMichael Lange if (!rank) { 57711c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 578d1a54cefSMatthew G. Knepley 57916ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 58011c56cb0SLisandro Dalcin for (cell = 0, c = 0; c < numCells; ++c) { 58111c56cb0SLisandro Dalcin 58211c56cb0SLisandro Dalcin /* Create face sets */ 5835964b3dcSLisandro Dalcin if (interpolate && gmsh_elem[c].dim == dim-1) { 58416ad7e67SMichael Lange const PetscInt *join; 58511c56cb0SLisandro Dalcin PetscInt joinSize, pcone[8], corner; 58611c56cb0SLisandro Dalcin /* Find the relevant facet with vertex joins */ 587a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 588dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 589917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 59016ad7e67SMichael Lange } 59111c56cb0SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 592a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 593c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 594a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 59516ad7e67SMichael Lange } 5966e1dd89aSLawrence Mitchell 5976e1dd89aSLawrence Mitchell /* Create cell sets */ 5986e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 5996e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 600*dea2a3c7SStefano Zampini ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 6016e1dd89aSLawrence Mitchell } 602917266a4SLisandro Dalcin cell++; 6036e1dd89aSLawrence Mitchell } 6040c070f12SSander Arens 6050c070f12SSander Arens /* Create vertex sets */ 6060c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 6070c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 608917266a4SLisandro Dalcin const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 60911c56cb0SLisandro Dalcin const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 610d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 6110c070f12SSander Arens } 6120c070f12SSander Arens } 6130c070f12SSander Arens } 61416ad7e67SMichael Lange } 61516ad7e67SMichael Lange 61611c56cb0SLisandro Dalcin /* Create coordinates */ 617*dea2a3c7SStefano Zampini if (embedDim < 0) embedDim = dim; 618*dea2a3c7SStefano Zampini ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 619979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 620331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 6211d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 622f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 623f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 624f45c9589SStefano Zampini } else { 62575b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 626f45c9589SStefano Zampini } 62775b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 6281d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 6291d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 630331830f3SMatthew G. Knepley } 63111c56cb0SLisandro Dalcin if (periodicMap) { 632437bdd5bSStefano Zampini ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 633f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 634f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 635437bdd5bSStefano Zampini PetscInt corner; 63611c56cb0SLisandro Dalcin PetscBool pc = PETSC_FALSE; 637437bdd5bSStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 638917266a4SLisandro Dalcin pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 639437bdd5bSStefano Zampini } 640437bdd5bSStefano Zampini if (pc) { 64111c56cb0SLisandro Dalcin PetscInt dof = gmsh_elem[c].numNodes*embedDim; 642*dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 643*dea2a3c7SStefano Zampini ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 644*dea2a3c7SStefano Zampini ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 645*dea2a3c7SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 6466fbe17bfSStefano Zampini } 647f45c9589SStefano Zampini cell++; 648f45c9589SStefano Zampini } 649f45c9589SStefano Zampini } 650f45c9589SStefano Zampini } 651331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 652331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 6538b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 654331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 655331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 6561d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 657331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 658331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 659f45c9589SStefano Zampini if (periodicMap) { 660f45c9589SStefano Zampini PetscInt off; 661f45c9589SStefano Zampini 662f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 663f45c9589SStefano Zampini PetscInt pcone[8], corner; 664f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 665*dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 666*dea2a3c7SStefano Zampini if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 667f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 668917266a4SLisandro Dalcin pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 669f45c9589SStefano Zampini } 670f45c9589SStefano Zampini if (dim == 3) { 671f45c9589SStefano Zampini /* Tetrahedra are inverted */ 672*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 673f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 674f45c9589SStefano Zampini pcone[0] = pcone[1]; 675f45c9589SStefano Zampini pcone[1] = tmp; 676f45c9589SStefano Zampini } 677f45c9589SStefano Zampini /* Hexahedra are inverted */ 678*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 679f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 680f45c9589SStefano Zampini pcone[1] = pcone[3]; 681f45c9589SStefano Zampini pcone[3] = tmp; 682f45c9589SStefano Zampini } 683*dea2a3c7SStefano Zampini /* Prisms are inverted */ 684*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 685*dea2a3c7SStefano Zampini PetscInt tmp; 686*dea2a3c7SStefano Zampini 687*dea2a3c7SStefano Zampini tmp = pcone[1]; 688*dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 689*dea2a3c7SStefano Zampini pcone[2] = tmp; 690*dea2a3c7SStefano Zampini tmp = pcone[4]; 691*dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 692*dea2a3c7SStefano Zampini pcone[5] = tmp; 693f45c9589SStefano Zampini } 694*dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 695*dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 696*dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 697*dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 698*dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 699*dea2a3c7SStefano Zampini pcone[3] = tmp; 700*dea2a3c7SStefano Zampini } 701*dea2a3c7SStefano Zampini } 702*dea2a3c7SStefano Zampini ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 703f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 70411c56cb0SLisandro Dalcin v = pcone[corner]; 705dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 70611c56cb0SLisandro Dalcin coords[off++] = (PetscReal) coordsIn[v*3+d]; 707f45c9589SStefano Zampini } 708f45c9589SStefano Zampini } 7096fbe17bfSStefano Zampini } 710f45c9589SStefano Zampini cell++; 711f45c9589SStefano Zampini } 712f45c9589SStefano Zampini } 713f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 714f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 715dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 71611c56cb0SLisandro Dalcin coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 717f45c9589SStefano Zampini } 718f45c9589SStefano Zampini } 719f45c9589SStefano Zampini } else { 720331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 7211d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 72211c56cb0SLisandro Dalcin coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 723331830f3SMatthew G. Knepley } 724331830f3SMatthew G. Knepley } 725331830f3SMatthew G. Knepley } 726331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 727331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 72811c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 72911c56cb0SLisandro Dalcin 73011c56cb0SLisandro Dalcin ierr = PetscFree(coordsIn);CHKERRQ(ierr); 73111c56cb0SLisandro Dalcin ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 732*dea2a3c7SStefano Zampini ierr = PetscFree(hybridMap);CHKERRQ(ierr); 733d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 734fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 7356fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 7366fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 73711c56cb0SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 73811c56cb0SLisandro Dalcin 7393b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 740331830f3SMatthew G. Knepley PetscFunctionReturn(0); 741331830f3SMatthew G. Knepley } 742