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"); 49*de78e4feSLisandro Dalcin if ((int)version == 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version 3.0 not supported"); 50*de78e4feSLisandro Dalcin if (version > 4.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at most version 4.0"); 5111c56cb0SLisandro Dalcin ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 52abc86ac4SMichael Lange } 5311c56cb0SLisandro Dalcin ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 5411c56cb0SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 5511c56cb0SLisandro Dalcin 567d282ae0SMichael Lange /* Create appropriate viewer and build plex */ 577d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 587d282ae0SMichael Lange ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 597d282ae0SMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 607d282ae0SMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 617d282ae0SMichael Lange ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 627d282ae0SMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 637d282ae0SMichael Lange PetscFunctionReturn(0); 647d282ae0SMichael Lange } 657d282ae0SMichael Lange 66*de78e4feSLisandro Dalcin typedef struct { 67*de78e4feSLisandro Dalcin size_t size; 68*de78e4feSLisandro Dalcin void *buffer; 69*de78e4feSLisandro Dalcin } GmshWorkBuffer; 70*de78e4feSLisandro Dalcin 71*de78e4feSLisandro Dalcin static PetscErrorCode GmshWorkBufferInit(GmshWorkBuffer *ctx) 72df032b82SMatthew G. Knepley { 73*de78e4feSLisandro Dalcin PetscFunctionBegin; 74*de78e4feSLisandro Dalcin ctx->size = 0; 75*de78e4feSLisandro Dalcin ctx->buffer = NULL; 76*de78e4feSLisandro Dalcin PetscFunctionReturn(0); 77*de78e4feSLisandro Dalcin } 78*de78e4feSLisandro Dalcin 79*de78e4feSLisandro Dalcin static PetscErrorCode GmshWorkBufferGet(GmshWorkBuffer *ctx, size_t count, size_t eltsize, void *buf) 80*de78e4feSLisandro Dalcin { 81*de78e4feSLisandro Dalcin size_t size = count*eltsize; 8211c56cb0SLisandro Dalcin PetscErrorCode ierr; 8311c56cb0SLisandro Dalcin 8411c56cb0SLisandro Dalcin PetscFunctionBegin; 85*de78e4feSLisandro Dalcin if (ctx->size < size) { 86*de78e4feSLisandro Dalcin ierr = PetscFree(ctx->buffer);CHKERRQ(ierr); 87*de78e4feSLisandro Dalcin ierr = PetscMalloc(size, &ctx->buffer);CHKERRQ(ierr); 88*de78e4feSLisandro Dalcin ctx->size = size; 89*de78e4feSLisandro Dalcin } 90*de78e4feSLisandro Dalcin *(void**)buf = size ? ctx->buffer : NULL; 91*de78e4feSLisandro Dalcin PetscFunctionReturn(0); 92*de78e4feSLisandro Dalcin } 93*de78e4feSLisandro Dalcin 94*de78e4feSLisandro Dalcin static PetscErrorCode GmshWorkBufferFree(GmshWorkBuffer *ctx) 95*de78e4feSLisandro Dalcin { 96*de78e4feSLisandro Dalcin PetscErrorCode ierr; 97*de78e4feSLisandro Dalcin PetscFunctionBegin; 98*de78e4feSLisandro Dalcin ierr = PetscFree(ctx->buffer);CHKERRQ(ierr); 99*de78e4feSLisandro Dalcin PetscFunctionReturn(0); 100*de78e4feSLisandro Dalcin } 101*de78e4feSLisandro Dalcin 102*de78e4feSLisandro Dalcin typedef struct { 103*de78e4feSLisandro Dalcin PetscInt dim; /* Entity dimension */ 104*de78e4feSLisandro Dalcin PetscInt id; /* Entity number */ 105*de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 106*de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 107*de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 108*de78e4feSLisandro Dalcin } GmshEntity; 109*de78e4feSLisandro Dalcin 110*de78e4feSLisandro Dalcin typedef struct { 111*de78e4feSLisandro Dalcin PetscInt dim; /* Entity dimension */ 112*de78e4feSLisandro Dalcin PetscInt id; /* Element number */ 113*de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 114*de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 115*de78e4feSLisandro Dalcin int nodes[8]; /* Node array */ 116*de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 117*de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 118*de78e4feSLisandro Dalcin } GmshElement; 119*de78e4feSLisandro Dalcin 120*de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v2(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int *numVertices, double **coordinates) 121*de78e4feSLisandro Dalcin { 122*de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 123*de78e4feSLisandro Dalcin int v, nid, snum; 124*de78e4feSLisandro Dalcin PetscErrorCode ierr; 125*de78e4feSLisandro Dalcin 126*de78e4feSLisandro Dalcin PetscFunctionBegin; 127*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 128*de78e4feSLisandro Dalcin snum = sscanf(line, "%d", numVertices); 129*de78e4feSLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 130*de78e4feSLisandro Dalcin ierr = PetscMalloc1(*numVertices*3, coordinates);CHKERRQ(ierr); 131*de78e4feSLisandro Dalcin for (v = 0; v < *numVertices; ++v) { 13211c56cb0SLisandro Dalcin double *xyz = *coordinates + v*3; 13311c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 13411c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 13511c56cb0SLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 13611c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 13711c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 13811c56cb0SLisandro Dalcin } 13911c56cb0SLisandro Dalcin PetscFunctionReturn(0); 14011c56cb0SLisandro Dalcin } 14111c56cb0SLisandro Dalcin 142*de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 143*de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 144*de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 145*de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 146*de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v2(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, int *numCells, GmshElement **gmsh_elems) 14711c56cb0SLisandro Dalcin { 148*de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 149df032b82SMatthew G. Knepley GmshElement *elements; 150*de78e4feSLisandro Dalcin int i, c, p, ibuf[1+4+512], snum; 15111c56cb0SLisandro Dalcin int cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 152df032b82SMatthew G. Knepley PetscErrorCode ierr; 153df032b82SMatthew G. Knepley 154df032b82SMatthew G. Knepley PetscFunctionBegin; 155*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 156*de78e4feSLisandro Dalcin snum = sscanf(line, "%d", numCells); 157*de78e4feSLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 158*de78e4feSLisandro Dalcin ierr = PetscMalloc1(*numCells, &elements);CHKERRQ(ierr); 159*de78e4feSLisandro Dalcin for (c = 0; c < *numCells;) { 16011c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 16111c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 162df032b82SMatthew G. Knepley if (binary) { 163df032b82SMatthew G. Knepley cellType = ibuf[0]; 164df032b82SMatthew G. Knepley numElem = ibuf[1]; 165df032b82SMatthew G. Knepley numTags = ibuf[2]; 166df032b82SMatthew G. Knepley } else { 167df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 168df032b82SMatthew G. Knepley cellType = ibuf[1]; 169df032b82SMatthew G. Knepley numTags = ibuf[2]; 170df032b82SMatthew G. Knepley numElem = 1; 171df032b82SMatthew G. Knepley } 172bf6ba3a3SMatthew G. Knepley /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 173bf6ba3a3SMatthew G. Knepley numNodesIgnore = 0; 174df032b82SMatthew G. Knepley switch (cellType) { 175df032b82SMatthew G. Knepley case 1: /* 2-node line */ 176df032b82SMatthew G. Knepley dim = 1; 177df032b82SMatthew G. Knepley numNodes = 2; 178df032b82SMatthew G. Knepley break; 179df032b82SMatthew G. Knepley case 2: /* 3-node triangle */ 180df032b82SMatthew G. Knepley dim = 2; 181df032b82SMatthew G. Knepley numNodes = 3; 182df032b82SMatthew G. Knepley break; 183df032b82SMatthew G. Knepley case 3: /* 4-node quadrangle */ 184df032b82SMatthew G. Knepley dim = 2; 185df032b82SMatthew G. Knepley numNodes = 4; 186df032b82SMatthew G. Knepley break; 187df032b82SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 188df032b82SMatthew G. Knepley dim = 3; 189df032b82SMatthew G. Knepley numNodes = 4; 190df032b82SMatthew G. Knepley break; 191df032b82SMatthew G. Knepley case 5: /* 8-node hexahedron */ 192df032b82SMatthew G. Knepley dim = 3; 193df032b82SMatthew G. Knepley numNodes = 8; 194df032b82SMatthew G. Knepley break; 195dea2a3c7SStefano Zampini case 6: /* 6-node wedge */ 196dea2a3c7SStefano Zampini dim = 3; 197dea2a3c7SStefano Zampini numNodes = 6; 198dea2a3c7SStefano Zampini break; 199bf6ba3a3SMatthew G. Knepley case 8: /* 3-node 2nd order line */ 200bf6ba3a3SMatthew G. Knepley dim = 1; 201bf6ba3a3SMatthew G. Knepley numNodes = 2; 202bf6ba3a3SMatthew G. Knepley numNodesIgnore = 1; 203bf6ba3a3SMatthew G. Knepley break; 204bf6ba3a3SMatthew G. Knepley case 9: /* 6-node 2nd order triangle */ 205bf6ba3a3SMatthew G. Knepley dim = 2; 206bf6ba3a3SMatthew G. Knepley numNodes = 3; 207bf6ba3a3SMatthew G. Knepley numNodesIgnore = 3; 208bf6ba3a3SMatthew G. Knepley break; 209dea2a3c7SStefano Zampini case 13: /* 18-node 2nd wedge */ 210dea2a3c7SStefano Zampini dim = 3; 211dea2a3c7SStefano Zampini numNodes = 6; 212dea2a3c7SStefano Zampini numNodesIgnore = 12; 213dea2a3c7SStefano Zampini break; 214df032b82SMatthew G. Knepley case 15: /* 1-node vertex */ 215df032b82SMatthew G. Knepley dim = 0; 216df032b82SMatthew G. Knepley numNodes = 1; 217df032b82SMatthew G. Knepley break; 218bf6ba3a3SMatthew G. Knepley case 7: /* 5-node pyramid */ 219bf6ba3a3SMatthew G. Knepley case 10: /* 9-node 2nd order quadrangle */ 220bf6ba3a3SMatthew G. Knepley case 11: /* 10-node 2nd order tetrahedron */ 221bf6ba3a3SMatthew G. Knepley case 12: /* 27-node 2nd order hexhedron */ 222bf6ba3a3SMatthew G. Knepley case 14: /* 14-node 2nd order pyramid */ 223df032b82SMatthew G. Knepley default: 224df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 225df032b82SMatthew G. Knepley } 226df032b82SMatthew G. Knepley if (binary) { 22711c56cb0SLisandro Dalcin const int nint = 1 + numTags + numNodes + numNodesIgnore; 22811c56cb0SLisandro Dalcin /* Loop over element blocks */ 229df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 23011c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 23111c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 232df032b82SMatthew G. Knepley elements[c].dim = dim; 233df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 234df032b82SMatthew G. Knepley elements[c].numTags = numTags; 235df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 236dea2a3c7SStefano Zampini elements[c].cellType = cellType; 237df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 238df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 239df032b82SMatthew G. Knepley } 240df032b82SMatthew G. Knepley } else { 241df032b82SMatthew G. Knepley elements[c].dim = dim; 242df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 243df032b82SMatthew G. Knepley elements[c].numTags = numTags; 244dea2a3c7SStefano Zampini elements[c].cellType = cellType; 245df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 246df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 24711c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr); 248df032b82SMatthew G. Knepley c++; 249df032b82SMatthew G. Knepley } 250df032b82SMatthew G. Knepley } 251df032b82SMatthew G. Knepley *gmsh_elems = elements; 252df032b82SMatthew G. Knepley PetscFunctionReturn(0); 253df032b82SMatthew G. Knepley } 2547d282ae0SMichael Lange 255*de78e4feSLisandro Dalcin /* 256*de78e4feSLisandro Dalcin $Entities 257*de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 258*de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 259*de78e4feSLisandro Dalcin // points 260*de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 261*de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 262*de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 263*de78e4feSLisandro Dalcin ... 264*de78e4feSLisandro Dalcin // curves 265*de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 266*de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 267*de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 268*de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 269*de78e4feSLisandro Dalcin ... 270*de78e4feSLisandro Dalcin // surfaces 271*de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 272*de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 273*de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 274*de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 275*de78e4feSLisandro Dalcin ... 276*de78e4feSLisandro Dalcin // volumes 277*de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 278*de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 279*de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 280*de78e4feSLisandro Dalcin numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 281*de78e4feSLisandro Dalcin ... 282*de78e4feSLisandro Dalcin $EndEntities 283*de78e4feSLisandro Dalcin */ 284*de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PetscInt numEntities[4], GmshEntity *entities[4]) 285*de78e4feSLisandro Dalcin { 286*de78e4feSLisandro Dalcin long i, num, lbuf[4]; 287*de78e4feSLisandro Dalcin int t, dim, numTags, *ibuf; 288*de78e4feSLisandro Dalcin GmshWorkBuffer work; 289*de78e4feSLisandro Dalcin PetscErrorCode ierr; 290*de78e4feSLisandro Dalcin 291*de78e4feSLisandro Dalcin PetscFunctionBegin; 292*de78e4feSLisandro Dalcin ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr); 293*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 294*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);} 295*de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 296*de78e4feSLisandro Dalcin numEntities[dim] = (PetscInt) lbuf[dim]; 297*de78e4feSLisandro Dalcin ierr = PetscCalloc1(numEntities[dim], &entities[dim]);CHKERRQ(ierr); 298*de78e4feSLisandro Dalcin for (i = 0; i < numEntities[dim]; ++i) { 299*de78e4feSLisandro Dalcin GmshEntity *entity = &entities[dim][i]; 300*de78e4feSLisandro Dalcin entity->dim = dim; 301*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &entity->id, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 302*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&entity->id, PETSC_ENUM, 1);CHKERRQ(ierr);} 303*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 304*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 305*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 306*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 307*de78e4feSLisandro Dalcin ierr = GmshWorkBufferGet(&work, num, sizeof(int), &ibuf);CHKERRQ(ierr); 308*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 309*de78e4feSLisandro Dalcin entity->numTags = numTags = (int) PetscMin(num, 4); 310*de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 311*de78e4feSLisandro Dalcin if (dim == 0) continue; 312*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 313*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 314*de78e4feSLisandro Dalcin ierr = GmshWorkBufferGet(&work, num, sizeof(int), &ibuf);CHKERRQ(ierr); 315*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 316*de78e4feSLisandro Dalcin } 317*de78e4feSLisandro Dalcin } 318*de78e4feSLisandro Dalcin ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr); 319*de78e4feSLisandro Dalcin PetscFunctionReturn(0); 320*de78e4feSLisandro Dalcin } 321*de78e4feSLisandro Dalcin 322*de78e4feSLisandro Dalcin /* 323*de78e4feSLisandro Dalcin $Nodes 324*de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 325*de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 326*de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 327*de78e4feSLisandro Dalcin ... 328*de78e4feSLisandro Dalcin ... 329*de78e4feSLisandro Dalcin $EndNodes 330*de78e4feSLisandro Dalcin */ 331*de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int *numVertices, double **gmsh_nodes) 332*de78e4feSLisandro Dalcin { 333*de78e4feSLisandro Dalcin long block, node, v, numEntityBlocks, numTotalNodes, numNodes; 334*de78e4feSLisandro Dalcin int info[3], nid; 335*de78e4feSLisandro Dalcin double *coordinates; 336*de78e4feSLisandro Dalcin char *cbuf; 337*de78e4feSLisandro Dalcin GmshWorkBuffer work; 338*de78e4feSLisandro Dalcin PetscErrorCode ierr; 339*de78e4feSLisandro Dalcin 340*de78e4feSLisandro Dalcin PetscFunctionBegin; 341*de78e4feSLisandro Dalcin ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr); 342*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 343*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 344*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 345*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 346*de78e4feSLisandro Dalcin ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr); 347*de78e4feSLisandro Dalcin *numVertices = (int)numTotalNodes; 348*de78e4feSLisandro Dalcin *gmsh_nodes = coordinates; 349*de78e4feSLisandro Dalcin for (v = 0, block = 0; block < numEntityBlocks; ++block) { 350*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 351*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 352*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 353*de78e4feSLisandro Dalcin if (binary) { 354*de78e4feSLisandro Dalcin int nbytes = sizeof(int) + 3*sizeof(double); 355*de78e4feSLisandro Dalcin ierr = GmshWorkBufferGet(&work, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 356*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 357*de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 358*de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 359*de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 360*de78e4feSLisandro Dalcin #if !defined(PETSC_WORDS_BIGENDIAN) 361*de78e4feSLisandro Dalcin ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr); 362*de78e4feSLisandro Dalcin ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr); 363*de78e4feSLisandro Dalcin #endif 364*de78e4feSLisandro Dalcin ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 365*de78e4feSLisandro Dalcin ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 366*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 367*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 368*de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 369*de78e4feSLisandro Dalcin } 370*de78e4feSLisandro Dalcin } else { 371*de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 372*de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 373*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 374*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 375*de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 376*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 377*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 378*de78e4feSLisandro Dalcin } 379*de78e4feSLisandro Dalcin } 380*de78e4feSLisandro Dalcin } 381*de78e4feSLisandro Dalcin ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr); 382*de78e4feSLisandro Dalcin PetscFunctionReturn(0); 383*de78e4feSLisandro Dalcin } 384*de78e4feSLisandro Dalcin 385*de78e4feSLisandro Dalcin /* 386*de78e4feSLisandro Dalcin $Elements 387*de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 388*de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 389*de78e4feSLisandro Dalcin tag(int) numVert[...](int) 390*de78e4feSLisandro Dalcin ... 391*de78e4feSLisandro Dalcin ... 392*de78e4feSLisandro Dalcin $EndElements 393*de78e4feSLisandro Dalcin */ 394*de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, GmshEntity *entities[4], int *numCells, GmshElement **gmsh_elems) 395*de78e4feSLisandro Dalcin { 396*de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 397*de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 398*de78e4feSLisandro Dalcin int eid, dim, numTags, *tags, cellType, numNodes; 399*de78e4feSLisandro Dalcin GmshElement *elements; 400*de78e4feSLisandro Dalcin GmshWorkBuffer work; 401*de78e4feSLisandro Dalcin PetscErrorCode ierr; 402*de78e4feSLisandro Dalcin 403*de78e4feSLisandro Dalcin PetscFunctionBegin; 404*de78e4feSLisandro Dalcin ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr); 405*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 406*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 407*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 408*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 409*de78e4feSLisandro Dalcin ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr); 410*de78e4feSLisandro Dalcin *numCells = (int)numTotalElements; 411*de78e4feSLisandro Dalcin *gmsh_elems = elements; 412*de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 413*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 414*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 415*de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 416*de78e4feSLisandro Dalcin numTags = entities[dim][eid-1].numTags; 417*de78e4feSLisandro Dalcin tags = entities[dim][eid-1].tags; 418*de78e4feSLisandro Dalcin switch (cellType) { 419*de78e4feSLisandro Dalcin case 1: /* 2-node line */ 420*de78e4feSLisandro Dalcin numNodes = 2; 421*de78e4feSLisandro Dalcin break; 422*de78e4feSLisandro Dalcin case 2: /* 3-node triangle */ 423*de78e4feSLisandro Dalcin numNodes = 3; break; 424*de78e4feSLisandro Dalcin case 3: /* 4-node quadrangle */ 425*de78e4feSLisandro Dalcin numNodes = 4; 426*de78e4feSLisandro Dalcin break; 427*de78e4feSLisandro Dalcin case 4: /* 4-node tetrahedron */ 428*de78e4feSLisandro Dalcin numNodes = 4; 429*de78e4feSLisandro Dalcin break; 430*de78e4feSLisandro Dalcin case 5: /* 8-node hexahedron */ 431*de78e4feSLisandro Dalcin numNodes = 8; 432*de78e4feSLisandro Dalcin break; 433*de78e4feSLisandro Dalcin case 6: /* 6-node wedge */ 434*de78e4feSLisandro Dalcin numNodes = 6; 435*de78e4feSLisandro Dalcin break; 436*de78e4feSLisandro Dalcin case 15: /* 1-node vertex */ 437*de78e4feSLisandro Dalcin numNodes = 1; 438*de78e4feSLisandro Dalcin break; 439*de78e4feSLisandro Dalcin default: 440*de78e4feSLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 441*de78e4feSLisandro Dalcin } 442*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 443*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 444*de78e4feSLisandro Dalcin ierr = GmshWorkBufferGet(&work, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 445*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 446*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 447*de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 448*de78e4feSLisandro Dalcin int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 449*de78e4feSLisandro Dalcin GmshElement *element = elements + c; 450*de78e4feSLisandro Dalcin element->dim = dim; 451*de78e4feSLisandro Dalcin element->cellType = cellType; 452*de78e4feSLisandro Dalcin element->numNodes = numNodes; 453*de78e4feSLisandro Dalcin element->numTags = numTags; 454*de78e4feSLisandro Dalcin element->id = id[0]; 455*de78e4feSLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 456*de78e4feSLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 457*de78e4feSLisandro Dalcin } 458*de78e4feSLisandro Dalcin } 459*de78e4feSLisandro Dalcin ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr); 460*de78e4feSLisandro Dalcin PetscFunctionReturn(0); 461*de78e4feSLisandro Dalcin } 462*de78e4feSLisandro Dalcin 463331830f3SMatthew G. Knepley /*@ 4647d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 465331830f3SMatthew G. Knepley 466331830f3SMatthew G. Knepley Collective on comm 467331830f3SMatthew G. Knepley 468331830f3SMatthew G. Knepley Input Parameters: 469331830f3SMatthew G. Knepley + comm - The MPI communicator 470331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 471331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 472331830f3SMatthew G. Knepley 473331830f3SMatthew G. Knepley Output Parameter: 474331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 475331830f3SMatthew G. Knepley 476331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 4773d51f2daSMichael Lange and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 478331830f3SMatthew G. Knepley 479331830f3SMatthew G. Knepley Level: beginner 480331830f3SMatthew G. Knepley 481331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 482331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 483331830f3SMatthew G. Knepley @*/ 484331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 485331830f3SMatthew G. Knepley { 48611c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 48711c56cb0SLisandro Dalcin double *coordsIn = NULL; 488*de78e4feSLisandro Dalcin PetscInt numEntities[4] = {0, 0, 0, 0}; 489*de78e4feSLisandro Dalcin GmshEntity *entities[4] = {NULL, NULL, NULL, NULL}; 4903c67d5adSSatish Balay GmshElement *gmsh_elem = NULL; 491331830f3SMatthew G. Knepley PetscSection coordSection; 492331830f3SMatthew G. Knepley Vec coordinates; 4936fbe17bfSStefano Zampini PetscBT periodicV = NULL, periodicC = NULL; 49484572febSToby Isaac PetscScalar *coords; 495dea2a3c7SStefano Zampini PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 4961d64f2ccSMatthew G. Knepley int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 49711c56cb0SLisandro Dalcin PetscMPIInt rank; 498331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 499dea2a3c7SStefano Zampini PetscBool binary, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 500dea2a3c7SStefano Zampini PetscBool enable_hybrid = PETSC_FALSE; 501331830f3SMatthew G. Knepley PetscErrorCode ierr; 502331830f3SMatthew G. Knepley 503331830f3SMatthew G. Knepley PetscFunctionBegin; 504331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 505dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr); 506dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 507dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 508dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 509dea2a3c7SStefano Zampini ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr); 51011c56cb0SLisandro Dalcin if (zerobase) shift = 0; 51111c56cb0SLisandro Dalcin 512331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 513331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 5143b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 51511c56cb0SLisandro Dalcin 51611c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 51711c56cb0SLisandro Dalcin 51811c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 5193b46f529SLisandro Dalcin if (binary) { 52011c56cb0SLisandro Dalcin parentviewer = viewer; 52111c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 52211c56cb0SLisandro Dalcin } 52311c56cb0SLisandro Dalcin 5243b46f529SLisandro Dalcin if (!rank) { 525dea2a3c7SStefano Zampini PetscBool match, hybrid; 526*de78e4feSLisandro Dalcin int fileType, fileFormat, dataSize; 527f6ac7a6aSMichael Lange float version; 528331830f3SMatthew G. Knepley 529331830f3SMatthew G. Knepley /* Read header */ 530060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 53104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 532331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 533060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 534f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 535f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 536f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 537*de78e4feSLisandro Dalcin if ((int)version == 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version 3.0 not supported"); 538*de78e4feSLisandro Dalcin if (version > 4.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at most version 4.0"); 539331830f3SMatthew 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); 540*de78e4feSLisandro Dalcin fileFormat = (int)version; 54104d1ad83SMichael Lange if (binary) { 542b9eae255SMichael Lange int checkInt; 543060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 54404d1ad83SMichael Lange if (checkInt != 1) { 545b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 54611c56cb0SLisandro Dalcin if (checkInt == 1) byteSwap = PETSC_TRUE; 54704d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 54804d1ad83SMichael Lange } 5490877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 550060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 55104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 552331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 55311c56cb0SLisandro Dalcin 554dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 555060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 556dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 557dc0ede3bSMatthew G. Knepley if (match) { 558dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 559dc0ede3bSMatthew G. Knepley snum = sscanf(line, "%d", &numRegions); 560dc0ede3bSMatthew G. Knepley if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 56111c56cb0SLisandro Dalcin for (r = 0; r < numRegions; ++r) { 56211c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 56311c56cb0SLisandro Dalcin } 564dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 565dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 566dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 567dc0ede3bSMatthew G. Knepley /* Initial read for vertex section */ 568dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 569dc0ede3bSMatthew G. Knepley } 57011c56cb0SLisandro Dalcin 571*de78e4feSLisandro Dalcin /* Read entities */ 572*de78e4feSLisandro Dalcin if (fileFormat == 4) { 573*de78e4feSLisandro Dalcin ierr = PetscStrncmp(line, "$Entities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 574*de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 575*de78e4feSLisandro Dalcin ierr = DMPlexCreateGmsh_ReadEntities_v4(viewer, binary, byteSwap, numEntities, entities);CHKERRQ(ierr); 576*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 577*de78e4feSLisandro Dalcin ierr = PetscStrncmp(line, "$EndEntities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 578*de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 579*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 580*de78e4feSLisandro Dalcin } 581*de78e4feSLisandro Dalcin 582dc0ede3bSMatthew G. Knepley /* Read vertices */ 58304d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 584331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 585*de78e4feSLisandro Dalcin if (fileFormat == 2) { 586*de78e4feSLisandro Dalcin ierr = DMPlexCreateGmsh_ReadNodes_v2(viewer, binary, byteSwap, shift, &numVertices, &coordsIn);CHKERRQ(ierr); 587*de78e4feSLisandro Dalcin } else { 588*de78e4feSLisandro Dalcin ierr = DMPlexCreateGmsh_ReadNodes_v4(viewer, binary, byteSwap, shift, &numVertices, &coordsIn);CHKERRQ(ierr); 589*de78e4feSLisandro Dalcin } 590060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 59104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 592331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 59311c56cb0SLisandro Dalcin 594331830f3SMatthew G. Knepley /* Read cells */ 595060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 59604d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 597331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 598*de78e4feSLisandro Dalcin if (fileFormat == 2) { 599*de78e4feSLisandro Dalcin ierr = DMPlexCreateGmsh_ReadElements_v2(viewer, binary, byteSwap, shift, &numCells, &gmsh_elem);CHKERRQ(ierr); 600*de78e4feSLisandro Dalcin } else { 601*de78e4feSLisandro Dalcin ierr = DMPlexCreateGmsh_ReadElements_v4(viewer, binary, byteSwap, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); 602*de78e4feSLisandro Dalcin } 603*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 604*de78e4feSLisandro Dalcin ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 605*de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 606*de78e4feSLisandro Dalcin 607*de78e4feSLisandro Dalcin for (i = 0; i < 4; +i++) { 608*de78e4feSLisandro Dalcin ierr = PetscFree(entities[i]);CHKERRQ(ierr); 609*de78e4feSLisandro Dalcin } 610*de78e4feSLisandro Dalcin 611dea2a3c7SStefano Zampini hybrid = PETSC_FALSE; 612a4bb7517SMichael Lange for (trueNumCells = 0, c = 0; c < numCells; ++c) { 613dea2a3c7SStefano Zampini int on = -1; 614a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 615dea2a3c7SStefano 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++;} 616dea2a3c7SStefano Zampini /* wedges always indicate an hybrid mesh in PLEX */ 617dea2a3c7SStefano Zampini if (on == 6 || on == 13) hybrid = PETSC_TRUE; 618db397164SMichael Lange } 61911c56cb0SLisandro Dalcin 620dea2a3c7SStefano Zampini /* Renumber cells for hybrid grids */ 621dea2a3c7SStefano Zampini if (hybrid && enable_hybrid) { 622dea2a3c7SStefano Zampini PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 623dea2a3c7SStefano Zampini PetscInt cell, tn, *tp; 624dea2a3c7SStefano Zampini int n1 = 0,n2 = 0; 625dea2a3c7SStefano Zampini 626dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 627dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 628dea2a3c7SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 629dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) { 630dea2a3c7SStefano Zampini if (!n1) n1 = gmsh_elem[c].cellType; 631dea2a3c7SStefano Zampini else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 632dea2a3c7SStefano Zampini 633dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 634dea2a3c7SStefano Zampini else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 635dea2a3c7SStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 636dea2a3c7SStefano Zampini cell++; 637dea2a3c7SStefano Zampini } 638dea2a3c7SStefano Zampini } 639dea2a3c7SStefano Zampini 640dea2a3c7SStefano Zampini switch (n1) { 641dea2a3c7SStefano Zampini case 2: /* triangles */ 642dea2a3c7SStefano Zampini case 9: 643dea2a3c7SStefano Zampini switch (n2) { 644dea2a3c7SStefano Zampini case 0: /* single type mesh */ 645dea2a3c7SStefano Zampini case 3: /* quads */ 646dea2a3c7SStefano Zampini case 10: 647dea2a3c7SStefano Zampini break; 648dea2a3c7SStefano Zampini default: 649dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 650dea2a3c7SStefano Zampini } 651dea2a3c7SStefano Zampini break; 652dea2a3c7SStefano Zampini case 3: 653dea2a3c7SStefano Zampini case 10: 654dea2a3c7SStefano Zampini switch (n2) { 655dea2a3c7SStefano Zampini case 0: /* single type mesh */ 656dea2a3c7SStefano Zampini case 2: /* swap since we list simplices first */ 657dea2a3c7SStefano Zampini case 9: 658dea2a3c7SStefano Zampini tn = hc1; 659dea2a3c7SStefano Zampini hc1 = hc2; 660dea2a3c7SStefano Zampini hc2 = tn; 661dea2a3c7SStefano Zampini 662dea2a3c7SStefano Zampini tp = hybridCells1; 663dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 664dea2a3c7SStefano Zampini hybridCells2 = tp; 665dea2a3c7SStefano Zampini break; 666dea2a3c7SStefano Zampini default: 667dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 668dea2a3c7SStefano Zampini } 669dea2a3c7SStefano Zampini break; 670dea2a3c7SStefano Zampini case 4: /* tetrahedra */ 671dea2a3c7SStefano Zampini case 11: 672dea2a3c7SStefano Zampini switch (n2) { 673dea2a3c7SStefano Zampini case 0: /* single type mesh */ 674dea2a3c7SStefano Zampini case 6: /* wedges */ 675dea2a3c7SStefano Zampini case 13: 676dea2a3c7SStefano Zampini break; 677dea2a3c7SStefano Zampini default: 678dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 679dea2a3c7SStefano Zampini } 680dea2a3c7SStefano Zampini break; 681dea2a3c7SStefano Zampini case 6: 682dea2a3c7SStefano Zampini case 13: 683dea2a3c7SStefano Zampini switch (n2) { 684dea2a3c7SStefano Zampini case 0: /* single type mesh */ 685dea2a3c7SStefano Zampini case 4: /* swap since we list simplices first */ 686dea2a3c7SStefano Zampini case 11: 687dea2a3c7SStefano Zampini tn = hc1; 688dea2a3c7SStefano Zampini hc1 = hc2; 689dea2a3c7SStefano Zampini hc2 = tn; 690dea2a3c7SStefano Zampini 691dea2a3c7SStefano Zampini tp = hybridCells1; 692dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 693dea2a3c7SStefano Zampini hybridCells2 = tp; 694dea2a3c7SStefano Zampini break; 695dea2a3c7SStefano Zampini default: 696dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 697dea2a3c7SStefano Zampini } 698dea2a3c7SStefano Zampini break; 699dea2a3c7SStefano Zampini default: 700dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 701dea2a3c7SStefano Zampini } 702dea2a3c7SStefano Zampini cMax = hc1; 703dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 704dea2a3c7SStefano Zampini for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 705dea2a3c7SStefano Zampini for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 706dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 707dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 708dea2a3c7SStefano Zampini } 709dea2a3c7SStefano Zampini 71011c56cb0SLisandro Dalcin /* OPTIONAL Read periodic section */ 711d08df55aSStefano Zampini if (periodic) { 712f45c9589SStefano Zampini PetscInt pVert, *periodicMapT, *aux; 713d08df55aSStefano Zampini int numPeriodic; 714d08df55aSStefano Zampini 715f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 7166fbe17bfSStefano Zampini ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 717f45c9589SStefano Zampini for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 718*de78e4feSLisandro Dalcin 719*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 720*de78e4feSLisandro Dalcin ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 721*de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 722*de78e4feSLisandro Dalcin if (fileFormat == 2 || !binary) { 723d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 724d08df55aSStefano Zampini snum = sscanf(line, "%d", &numPeriodic); 725d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 726*de78e4feSLisandro Dalcin } else { 727*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 728*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 729*de78e4feSLisandro Dalcin } 730d08df55aSStefano Zampini for (i = 0; i < numPeriodic; i++) { 731*de78e4feSLisandro Dalcin int ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode; 732*de78e4feSLisandro Dalcin long j, nNodes; 733*de78e4feSLisandro Dalcin double affine[16]; 734d08df55aSStefano Zampini 735*de78e4feSLisandro Dalcin if (fileFormat == 2 || !binary) { 736d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 737*de78e4feSLisandro Dalcin snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag); 738d08df55aSStefano Zampini if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 739*de78e4feSLisandro Dalcin } else { 740*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 741*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 742*de78e4feSLisandro Dalcin slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2]; 743*de78e4feSLisandro Dalcin } 744*de78e4feSLisandro Dalcin (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */ 745*de78e4feSLisandro Dalcin 746*de78e4feSLisandro Dalcin if (fileFormat == 2 || !binary) { 747da83f57bSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 748*de78e4feSLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 749da83f57bSLisandro Dalcin if (snum != 1) { /* discard tranformation and try again */ 75072ffbcc9SStefano Zampini ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 751d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 752*de78e4feSLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 753d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 754da83f57bSLisandro Dalcin } 755*de78e4feSLisandro Dalcin } else { 756*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 757*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 758*de78e4feSLisandro Dalcin if (nNodes == -1) { /* discard tranformation and try again */ 759*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 760*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 761*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 762*de78e4feSLisandro Dalcin } 763*de78e4feSLisandro Dalcin } 764*de78e4feSLisandro Dalcin 765d08df55aSStefano Zampini for (j = 0; j < nNodes; j++) { 766*de78e4feSLisandro Dalcin if (fileFormat == 2 || !binary) { 767d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 76811c56cb0SLisandro Dalcin snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 769d08df55aSStefano Zampini if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 770*de78e4feSLisandro Dalcin } else { 771*de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 772*de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 773*de78e4feSLisandro Dalcin slaveNode = ibuf[0]; masterNode = ibuf[1]; 774*de78e4feSLisandro Dalcin } 775917266a4SLisandro Dalcin periodicMapT[slaveNode - shift] = masterNode - shift; 776917266a4SLisandro Dalcin ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr); 777917266a4SLisandro Dalcin ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr); 778d08df55aSStefano Zampini } 779d08df55aSStefano Zampini } 780d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 781d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 782d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 783*de78e4feSLisandro Dalcin 784d08df55aSStefano Zampini /* we may have slaves of slaves */ 785d08df55aSStefano Zampini for (i = 0; i < numVertices; i++) { 786f45c9589SStefano Zampini while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 787f45c9589SStefano Zampini periodicMapT[i] = periodicMapT[periodicMapT[i]]; 788d08df55aSStefano Zampini } 789d08df55aSStefano Zampini } 790f45c9589SStefano Zampini /* periodicMap : from old to new numbering (periodic vertices excluded) 791f45c9589SStefano Zampini periodicMapI: from new to old numbering */ 792f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 793f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 794f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 795f45c9589SStefano Zampini for (i = 0, pVert = 0; i < numVertices; i++) { 796f45c9589SStefano Zampini if (periodicMapT[i] != i) { 797f45c9589SStefano Zampini pVert++; 798f45c9589SStefano Zampini } else { 799f45c9589SStefano Zampini aux[i] = i - pVert; 800f45c9589SStefano Zampini periodicMapI[i - pVert] = i; 801f45c9589SStefano Zampini } 802f45c9589SStefano Zampini } 803f45c9589SStefano Zampini for (i = 0 ; i < numVertices; i++) { 804f45c9589SStefano Zampini periodicMap[i] = aux[periodicMapT[i]]; 805f45c9589SStefano Zampini } 806f45c9589SStefano Zampini ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 807f45c9589SStefano Zampini ierr = PetscFree(aux);CHKERRQ(ierr); 808f45c9589SStefano Zampini /* remove periodic vertices */ 809f45c9589SStefano Zampini numVertices = numVertices - pVert; 810d08df55aSStefano Zampini } 811db397164SMichael Lange } 81211c56cb0SLisandro Dalcin 81311c56cb0SLisandro Dalcin if (parentviewer) { 81411c56cb0SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 81511c56cb0SLisandro Dalcin } 81611c56cb0SLisandro Dalcin 817a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 818db397164SMichael Lange ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 819a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 820a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 821dea2a3c7SStefano Zampini ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 822a4bb7517SMichael Lange cell++; 823db397164SMichael Lange } 824331830f3SMatthew G. Knepley } 825331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 826a4bb7517SMichael Lange /* Add cell-vertex connections */ 827a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 828a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 82911c56cb0SLisandro Dalcin PetscInt pcone[8], corner; 830a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 831dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 832917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 833db397164SMichael Lange } 83497ed6be6Ssarens if (dim == 3) { 83597ed6be6Ssarens /* Tetrahedra are inverted */ 836dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 83797ed6be6Ssarens PetscInt tmp = pcone[0]; 83897ed6be6Ssarens pcone[0] = pcone[1]; 83997ed6be6Ssarens pcone[1] = tmp; 84097ed6be6Ssarens } 84197ed6be6Ssarens /* Hexahedra are inverted */ 842dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 84397ed6be6Ssarens PetscInt tmp = pcone[1]; 84497ed6be6Ssarens pcone[1] = pcone[3]; 84597ed6be6Ssarens pcone[3] = tmp; 84697ed6be6Ssarens } 847dea2a3c7SStefano Zampini /* Prisms are inverted */ 848dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 849dea2a3c7SStefano Zampini PetscInt tmp; 850dea2a3c7SStefano Zampini 851dea2a3c7SStefano Zampini tmp = pcone[1]; 852dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 853dea2a3c7SStefano Zampini pcone[2] = tmp; 854dea2a3c7SStefano Zampini tmp = pcone[4]; 855dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 856dea2a3c7SStefano Zampini pcone[5] = tmp; 85797ed6be6Ssarens } 858dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 859dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 860dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 861dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 862dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 863dea2a3c7SStefano Zampini pcone[3] = tmp; 864dea2a3c7SStefano Zampini } 865dea2a3c7SStefano Zampini } 866dea2a3c7SStefano Zampini ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 867a4bb7517SMichael Lange cell++; 868331830f3SMatthew G. Knepley } 869a4bb7517SMichael Lange } 8706228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 871c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 872dea2a3c7SStefano Zampini ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 873331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 874331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 875331830f3SMatthew G. Knepley if (interpolate) { 8765fd9971aSMatthew G. Knepley DM idm; 877331830f3SMatthew G. Knepley 878331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 879331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 880331830f3SMatthew G. Knepley *dm = idm; 881331830f3SMatthew G. Knepley } 882917266a4SLisandro Dalcin 883f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 884917266a4SLisandro Dalcin if (!rank && usemarker) { 885d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 886d3f73514SStefano Zampini 887d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 888d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 889d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 890d3f73514SStefano Zampini PetscInt suppSize; 891d3f73514SStefano Zampini 892d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 893d3f73514SStefano Zampini if (suppSize == 1) { 894d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 895d3f73514SStefano Zampini 896d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 897d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 898d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 899d3f73514SStefano Zampini } 900d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 901d3f73514SStefano Zampini } 902d3f73514SStefano Zampini } 903d3f73514SStefano Zampini } 90416ad7e67SMichael Lange 90516ad7e67SMichael Lange if (!rank) { 90611c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 907d1a54cefSMatthew G. Knepley 90816ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 90911c56cb0SLisandro Dalcin for (cell = 0, c = 0; c < numCells; ++c) { 91011c56cb0SLisandro Dalcin 91111c56cb0SLisandro Dalcin /* Create face sets */ 9125964b3dcSLisandro Dalcin if (interpolate && gmsh_elem[c].dim == dim-1) { 91316ad7e67SMichael Lange const PetscInt *join; 91411c56cb0SLisandro Dalcin PetscInt joinSize, pcone[8], corner; 91511c56cb0SLisandro Dalcin /* Find the relevant facet with vertex joins */ 916a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 917dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 918917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 91916ad7e67SMichael Lange } 92011c56cb0SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 921f10f1c67SMatthew G. Knepley if (joinSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for GMsh element %d (Plex cell %D)", gmsh_elem[c].id, c); 922c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 923a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 92416ad7e67SMichael Lange } 9256e1dd89aSLawrence Mitchell 9266e1dd89aSLawrence Mitchell /* Create cell sets */ 9276e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 9286e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 929dea2a3c7SStefano Zampini ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 9306e1dd89aSLawrence Mitchell } 931917266a4SLisandro Dalcin cell++; 9326e1dd89aSLawrence Mitchell } 9330c070f12SSander Arens 9340c070f12SSander Arens /* Create vertex sets */ 9350c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 9360c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 937917266a4SLisandro Dalcin const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 93811c56cb0SLisandro Dalcin const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 939d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 9400c070f12SSander Arens } 9410c070f12SSander Arens } 9420c070f12SSander Arens } 94316ad7e67SMichael Lange } 94416ad7e67SMichael Lange 94511c56cb0SLisandro Dalcin /* Create coordinates */ 946dea2a3c7SStefano Zampini if (embedDim < 0) embedDim = dim; 947dea2a3c7SStefano Zampini ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 948979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 949331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 9501d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 951f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 952f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 953f45c9589SStefano Zampini } else { 95475b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 955f45c9589SStefano Zampini } 95675b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 9571d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 9581d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 959331830f3SMatthew G. Knepley } 96011c56cb0SLisandro Dalcin if (periodicMap) { 961437bdd5bSStefano Zampini ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 962f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 963f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 964437bdd5bSStefano Zampini PetscInt corner; 96511c56cb0SLisandro Dalcin PetscBool pc = PETSC_FALSE; 966437bdd5bSStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 967917266a4SLisandro Dalcin pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 968437bdd5bSStefano Zampini } 969437bdd5bSStefano Zampini if (pc) { 97011c56cb0SLisandro Dalcin PetscInt dof = gmsh_elem[c].numNodes*embedDim; 971dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 972dea2a3c7SStefano Zampini ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 973dea2a3c7SStefano Zampini ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 974dea2a3c7SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 9756fbe17bfSStefano Zampini } 976f45c9589SStefano Zampini cell++; 977f45c9589SStefano Zampini } 978f45c9589SStefano Zampini } 979f45c9589SStefano Zampini } 980331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 981331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 9828b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 983331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 984331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 9851d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 986331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 987331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 988f45c9589SStefano Zampini if (periodicMap) { 989f45c9589SStefano Zampini PetscInt off; 990f45c9589SStefano Zampini 991f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 992f45c9589SStefano Zampini PetscInt pcone[8], corner; 993f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 994dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 995dea2a3c7SStefano Zampini if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 996f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 997917266a4SLisandro Dalcin pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 998f45c9589SStefano Zampini } 999f45c9589SStefano Zampini if (dim == 3) { 1000f45c9589SStefano Zampini /* Tetrahedra are inverted */ 1001dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 1002f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 1003f45c9589SStefano Zampini pcone[0] = pcone[1]; 1004f45c9589SStefano Zampini pcone[1] = tmp; 1005f45c9589SStefano Zampini } 1006f45c9589SStefano Zampini /* Hexahedra are inverted */ 1007dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 1008f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 1009f45c9589SStefano Zampini pcone[1] = pcone[3]; 1010f45c9589SStefano Zampini pcone[3] = tmp; 1011f45c9589SStefano Zampini } 1012dea2a3c7SStefano Zampini /* Prisms are inverted */ 1013dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 1014dea2a3c7SStefano Zampini PetscInt tmp; 1015dea2a3c7SStefano Zampini 1016dea2a3c7SStefano Zampini tmp = pcone[1]; 1017dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1018dea2a3c7SStefano Zampini pcone[2] = tmp; 1019dea2a3c7SStefano Zampini tmp = pcone[4]; 1020dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1021dea2a3c7SStefano Zampini pcone[5] = tmp; 1022f45c9589SStefano Zampini } 1023dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1024dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1025dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 1026dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1027dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1028dea2a3c7SStefano Zampini pcone[3] = tmp; 1029dea2a3c7SStefano Zampini } 1030dea2a3c7SStefano Zampini } 1031dea2a3c7SStefano Zampini ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 1032f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 103311c56cb0SLisandro Dalcin v = pcone[corner]; 1034dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 103511c56cb0SLisandro Dalcin coords[off++] = (PetscReal) coordsIn[v*3+d]; 1036f45c9589SStefano Zampini } 1037f45c9589SStefano Zampini } 10386fbe17bfSStefano Zampini } 1039f45c9589SStefano Zampini cell++; 1040f45c9589SStefano Zampini } 1041f45c9589SStefano Zampini } 1042f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 1043f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 1044dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 104511c56cb0SLisandro Dalcin coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 1046f45c9589SStefano Zampini } 1047f45c9589SStefano Zampini } 1048f45c9589SStefano Zampini } else { 1049331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 10501d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 105111c56cb0SLisandro Dalcin coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 1052331830f3SMatthew G. Knepley } 1053331830f3SMatthew G. Knepley } 1054331830f3SMatthew G. Knepley } 1055331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1056331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 105711c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 105811c56cb0SLisandro Dalcin 105911c56cb0SLisandro Dalcin ierr = PetscFree(coordsIn);CHKERRQ(ierr); 106011c56cb0SLisandro Dalcin ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 1061dea2a3c7SStefano Zampini ierr = PetscFree(hybridMap);CHKERRQ(ierr); 1062d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 1063fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 10646fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 10656fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 106611c56cb0SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 106711c56cb0SLisandro Dalcin 10683b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 1069331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1070331830f3SMatthew G. Knepley } 1071