1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3730356f6SLisandro Dalcin #include <petsc/private/hashmapi.h> 4331830f3SMatthew G. Knepley 57d282ae0SMichael Lange /*@C 67d282ae0SMichael Lange DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 77d282ae0SMichael Lange 87d282ae0SMichael Lange + comm - The MPI communicator 97d282ae0SMichael Lange . filename - Name of the Gmsh file 107d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh 117d282ae0SMichael Lange 127d282ae0SMichael Lange Output Parameter: 137d282ae0SMichael Lange . dm - The DM object representing the mesh 147d282ae0SMichael Lange 157d282ae0SMichael Lange Level: beginner 167d282ae0SMichael Lange 177d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 187d282ae0SMichael Lange @*/ 197d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 207d282ae0SMichael Lange { 2111c56cb0SLisandro Dalcin PetscViewer viewer; 22abc86ac4SMichael Lange PetscMPIInt rank; 2311c56cb0SLisandro Dalcin int fileType; 247d282ae0SMichael Lange PetscViewerType vtype; 257d282ae0SMichael Lange PetscErrorCode ierr; 267d282ae0SMichael Lange 277d282ae0SMichael Lange PetscFunctionBegin; 28abc86ac4SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2911c56cb0SLisandro Dalcin 307d282ae0SMichael Lange /* Determine Gmsh file type (ASCII or binary) from file header */ 3111c56cb0SLisandro Dalcin if (!rank) { 3211c56cb0SLisandro Dalcin PetscViewer vheader; 3311c56cb0SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 3411c56cb0SLisandro Dalcin PetscBool match; 3511c56cb0SLisandro Dalcin int snum; 3611c56cb0SLisandro Dalcin float version; 3711c56cb0SLisandro Dalcin 3811c56cb0SLisandro Dalcin ierr = PetscViewerCreate(PETSC_COMM_SELF, &vheader);CHKERRQ(ierr); 397d282ae0SMichael Lange ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr); 407d282ae0SMichael Lange ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr); 417d282ae0SMichael Lange ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr); 427d282ae0SMichael Lange /* Read only the first two lines of the Gmsh file */ 43060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 4411c56cb0SLisandro Dalcin ierr = PetscStrncmp(line, "$MeshFormat", sizeof(line), &match);CHKERRQ(ierr); 457d282ae0SMichael Lange if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 46060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 4711c56cb0SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 48f6ac7a6aSMichael Lange if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 49*698a79b9SLisandro Dalcin if (version < 2.2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f must be at least 2.2", (double)version); 50*698a79b9SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version); 51*698a79b9SLisandro Dalcin if (version > 4.1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f must be at most 4.1", (double)version); 5211c56cb0SLisandro Dalcin ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 53abc86ac4SMichael Lange } 5411c56cb0SLisandro Dalcin ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 5511c56cb0SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 5611c56cb0SLisandro Dalcin 577d282ae0SMichael Lange /* Create appropriate viewer and build plex */ 587d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 597d282ae0SMichael Lange ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 607d282ae0SMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 617d282ae0SMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 627d282ae0SMichael Lange ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 637d282ae0SMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 647d282ae0SMichael Lange PetscFunctionReturn(0); 657d282ae0SMichael Lange } 667d282ae0SMichael Lange 67de78e4feSLisandro Dalcin typedef struct { 68*698a79b9SLisandro Dalcin PetscViewer viewer; 69*698a79b9SLisandro Dalcin int fileFormat; 70*698a79b9SLisandro Dalcin int dataSize; 71*698a79b9SLisandro Dalcin PetscBool binary; 72*698a79b9SLisandro Dalcin PetscBool byteSwap; 73*698a79b9SLisandro Dalcin size_t wlen; 74*698a79b9SLisandro Dalcin void *wbuf; 75*698a79b9SLisandro Dalcin size_t slen; 76*698a79b9SLisandro Dalcin void *sbuf; 77*698a79b9SLisandro Dalcin } GmshFile; 78de78e4feSLisandro Dalcin 79*698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf) 80de78e4feSLisandro Dalcin { 81de78e4feSLisandro Dalcin size_t size = count * eltsize; 8211c56cb0SLisandro Dalcin PetscErrorCode ierr; 8311c56cb0SLisandro Dalcin 8411c56cb0SLisandro Dalcin PetscFunctionBegin; 85*698a79b9SLisandro Dalcin if (gmsh->wlen < size) { 86*698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 87*698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->wbuf);CHKERRQ(ierr); 88*698a79b9SLisandro Dalcin gmsh->wlen = size; 89de78e4feSLisandro Dalcin } 90*698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->wbuf : NULL; 91de78e4feSLisandro Dalcin PetscFunctionReturn(0); 92de78e4feSLisandro Dalcin } 93de78e4feSLisandro Dalcin 94*698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf) 95*698a79b9SLisandro Dalcin { 96*698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 97*698a79b9SLisandro Dalcin size_t size = count * dataSize; 98*698a79b9SLisandro Dalcin PetscErrorCode ierr; 99*698a79b9SLisandro Dalcin 100*698a79b9SLisandro Dalcin PetscFunctionBegin; 101*698a79b9SLisandro Dalcin if (gmsh->slen < size) { 102*698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 103*698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->sbuf);CHKERRQ(ierr); 104*698a79b9SLisandro Dalcin gmsh->slen = size; 105*698a79b9SLisandro Dalcin } 106*698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->sbuf : NULL; 107*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 108*698a79b9SLisandro Dalcin } 109*698a79b9SLisandro Dalcin 110*698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype) 111de78e4feSLisandro Dalcin { 112de78e4feSLisandro Dalcin PetscErrorCode ierr; 113de78e4feSLisandro Dalcin PetscFunctionBegin; 114*698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype);CHKERRQ(ierr); 115*698a79b9SLisandro Dalcin if (gmsh->byteSwap) {ierr = PetscByteSwap(buf, dtype, count);CHKERRQ(ierr);} 116*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 117*698a79b9SLisandro Dalcin } 118*698a79b9SLisandro Dalcin 119*698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count) 120*698a79b9SLisandro Dalcin { 121*698a79b9SLisandro Dalcin PetscErrorCode ierr; 122*698a79b9SLisandro Dalcin PetscFunctionBegin; 123*698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING);CHKERRQ(ierr); 124*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 125*698a79b9SLisandro Dalcin } 126*698a79b9SLisandro Dalcin 127*698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count) 128*698a79b9SLisandro Dalcin { 129*698a79b9SLisandro Dalcin PetscInt i; 130*698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 131*698a79b9SLisandro Dalcin PetscErrorCode ierr; 132*698a79b9SLisandro Dalcin 133*698a79b9SLisandro Dalcin PetscFunctionBegin; 134*698a79b9SLisandro Dalcin if (dataSize == sizeof(PetscInt)) { 135*698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_INT);CHKERRQ(ierr); 136*698a79b9SLisandro Dalcin } else if (dataSize == sizeof(int)) { 137*698a79b9SLisandro Dalcin int *ibuf = NULL; 138*698a79b9SLisandro Dalcin ierr = GmshBufferSizeGet(gmsh, count, &ibuf); 139*698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); 140*698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 141*698a79b9SLisandro Dalcin } else if (dataSize == sizeof(long)) { 142*698a79b9SLisandro Dalcin long *ibuf = NULL; 143*698a79b9SLisandro Dalcin ierr = GmshBufferSizeGet(gmsh, count, &ibuf); 144*698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_LONG);CHKERRQ(ierr); 145*698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 146*698a79b9SLisandro Dalcin } else if (dataSize == sizeof(PetscInt64)) { 147*698a79b9SLisandro Dalcin PetscInt64 *ibuf = NULL; 148*698a79b9SLisandro Dalcin ierr = GmshBufferSizeGet(gmsh, count, &ibuf); 149*698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_INT64);CHKERRQ(ierr); 150*698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 151*698a79b9SLisandro Dalcin } 152*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 153*698a79b9SLisandro Dalcin } 154*698a79b9SLisandro Dalcin 155*698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count) 156*698a79b9SLisandro Dalcin { 157*698a79b9SLisandro Dalcin PetscErrorCode ierr; 158*698a79b9SLisandro Dalcin PetscFunctionBegin; 159*698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_ENUM);CHKERRQ(ierr); 160*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 161*698a79b9SLisandro Dalcin } 162*698a79b9SLisandro Dalcin 163*698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count) 164*698a79b9SLisandro Dalcin { 165*698a79b9SLisandro Dalcin PetscErrorCode ierr; 166*698a79b9SLisandro Dalcin PetscFunctionBegin; 167*698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_DOUBLE);CHKERRQ(ierr); 168de78e4feSLisandro Dalcin PetscFunctionReturn(0); 169de78e4feSLisandro Dalcin } 170de78e4feSLisandro Dalcin 171de78e4feSLisandro Dalcin typedef struct { 172de78e4feSLisandro Dalcin PetscInt id; /* Entity number */ 173*698a79b9SLisandro Dalcin PetscInt dim; /* Entity dimension */ 174de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 175de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 176de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 177de78e4feSLisandro Dalcin } GmshEntity; 178de78e4feSLisandro Dalcin 179de78e4feSLisandro Dalcin typedef struct { 180730356f6SLisandro Dalcin GmshEntity *entity[4]; 181730356f6SLisandro Dalcin PetscHMapI entityMap[4]; 182730356f6SLisandro Dalcin } GmshEntities; 183730356f6SLisandro Dalcin 184*698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities) 185730356f6SLisandro Dalcin { 186730356f6SLisandro Dalcin PetscInt dim; 187730356f6SLisandro Dalcin PetscErrorCode ierr; 188730356f6SLisandro Dalcin 189730356f6SLisandro Dalcin PetscFunctionBegin; 190730356f6SLisandro Dalcin ierr = PetscNew(entities);CHKERRQ(ierr); 191730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 192730356f6SLisandro Dalcin ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr); 193730356f6SLisandro Dalcin ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 194730356f6SLisandro Dalcin } 195730356f6SLisandro Dalcin PetscFunctionReturn(0); 196730356f6SLisandro Dalcin } 197730356f6SLisandro Dalcin 198730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity) 199730356f6SLisandro Dalcin { 200730356f6SLisandro Dalcin PetscErrorCode ierr; 201730356f6SLisandro Dalcin PetscFunctionBegin; 202730356f6SLisandro Dalcin ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr); 203730356f6SLisandro Dalcin entities->entity[dim][index].dim = dim; 204730356f6SLisandro Dalcin entities->entity[dim][index].id = eid; 205730356f6SLisandro Dalcin if (entity) *entity = &entities->entity[dim][index]; 206730356f6SLisandro Dalcin PetscFunctionReturn(0); 207730356f6SLisandro Dalcin } 208730356f6SLisandro Dalcin 209730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity) 210730356f6SLisandro Dalcin { 211730356f6SLisandro Dalcin PetscInt index; 212730356f6SLisandro Dalcin PetscErrorCode ierr; 213730356f6SLisandro Dalcin 214730356f6SLisandro Dalcin PetscFunctionBegin; 215730356f6SLisandro Dalcin ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr); 216730356f6SLisandro Dalcin *entity = &entities->entity[dim][index]; 217730356f6SLisandro Dalcin PetscFunctionReturn(0); 218730356f6SLisandro Dalcin } 219730356f6SLisandro Dalcin 220730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) 221730356f6SLisandro Dalcin { 222730356f6SLisandro Dalcin PetscInt dim; 223730356f6SLisandro Dalcin PetscErrorCode ierr; 224730356f6SLisandro Dalcin 225730356f6SLisandro Dalcin PetscFunctionBegin; 226730356f6SLisandro Dalcin if (!*entities) PetscFunctionReturn(0); 227730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 228730356f6SLisandro Dalcin ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr); 229730356f6SLisandro Dalcin ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 230730356f6SLisandro Dalcin } 231730356f6SLisandro Dalcin ierr = PetscFree((*entities));CHKERRQ(ierr); 232730356f6SLisandro Dalcin PetscFunctionReturn(0); 233730356f6SLisandro Dalcin } 234730356f6SLisandro Dalcin 235730356f6SLisandro Dalcin typedef struct { 236*698a79b9SLisandro Dalcin PetscInt id; /* Entity number */ 237de78e4feSLisandro Dalcin PetscInt dim; /* Entity dimension */ 238de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 239de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 240*698a79b9SLisandro Dalcin PetscInt nodes[8]; /* Node array */ 241de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 242de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 243de78e4feSLisandro Dalcin } GmshElement; 244de78e4feSLisandro Dalcin 245*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v22(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 246de78e4feSLisandro Dalcin { 247*698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 248*698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 249de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 250*698a79b9SLisandro Dalcin int v, num, nid, snum; 251*698a79b9SLisandro Dalcin double *coordinates; 252de78e4feSLisandro Dalcin PetscErrorCode ierr; 253de78e4feSLisandro Dalcin 254de78e4feSLisandro Dalcin PetscFunctionBegin; 255de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 256*698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 257de78e4feSLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 258*698a79b9SLisandro Dalcin ierr = PetscMalloc1(num*3, &coordinates);CHKERRQ(ierr); 259*698a79b9SLisandro Dalcin *numVertices = num; 260*698a79b9SLisandro Dalcin *gmsh_nodes = coordinates; 261*698a79b9SLisandro Dalcin for (v = 0; v < num; ++v) { 262*698a79b9SLisandro Dalcin double *xyz = coordinates + v*3; 26311c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 26411c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 26511c56cb0SLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 26611c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 26711c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 26811c56cb0SLisandro Dalcin } 26911c56cb0SLisandro Dalcin PetscFunctionReturn(0); 27011c56cb0SLisandro Dalcin } 27111c56cb0SLisandro Dalcin 272de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 273de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 274de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 275de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 276*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v22(GmshFile* gmsh, PETSC_UNUSED int shift, PetscInt *numCells, GmshElement **gmsh_elems) 27711c56cb0SLisandro Dalcin { 278*698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 279*698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 280*698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 281de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 282df032b82SMatthew G. Knepley GmshElement *elements; 283*698a79b9SLisandro Dalcin int i, c, p, num, ibuf[1+4+512], snum; 28411c56cb0SLisandro Dalcin int cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 285df032b82SMatthew G. Knepley PetscErrorCode ierr; 286df032b82SMatthew G. Knepley 287df032b82SMatthew G. Knepley PetscFunctionBegin; 288de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 289*698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 290de78e4feSLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 291*698a79b9SLisandro Dalcin ierr = PetscMalloc1(num, &elements);CHKERRQ(ierr); 292*698a79b9SLisandro Dalcin *numCells = num; 293*698a79b9SLisandro Dalcin *gmsh_elems = elements; 294*698a79b9SLisandro Dalcin for (c = 0; c < num;) { 29511c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 29611c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 297df032b82SMatthew G. Knepley if (binary) { 298df032b82SMatthew G. Knepley cellType = ibuf[0]; 299df032b82SMatthew G. Knepley numElem = ibuf[1]; 300df032b82SMatthew G. Knepley numTags = ibuf[2]; 301df032b82SMatthew G. Knepley } else { 302df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 303df032b82SMatthew G. Knepley cellType = ibuf[1]; 304df032b82SMatthew G. Knepley numTags = ibuf[2]; 305df032b82SMatthew G. Knepley numElem = 1; 306df032b82SMatthew G. Knepley } 307bf6ba3a3SMatthew G. Knepley /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 308bf6ba3a3SMatthew G. Knepley numNodesIgnore = 0; 309df032b82SMatthew G. Knepley switch (cellType) { 310df032b82SMatthew G. Knepley case 1: /* 2-node line */ 311df032b82SMatthew G. Knepley dim = 1; 312df032b82SMatthew G. Knepley numNodes = 2; 313df032b82SMatthew G. Knepley break; 314df032b82SMatthew G. Knepley case 2: /* 3-node triangle */ 315df032b82SMatthew G. Knepley dim = 2; 316df032b82SMatthew G. Knepley numNodes = 3; 317df032b82SMatthew G. Knepley break; 318df032b82SMatthew G. Knepley case 3: /* 4-node quadrangle */ 319df032b82SMatthew G. Knepley dim = 2; 320df032b82SMatthew G. Knepley numNodes = 4; 321df032b82SMatthew G. Knepley break; 322df032b82SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 323df032b82SMatthew G. Knepley dim = 3; 324df032b82SMatthew G. Knepley numNodes = 4; 325df032b82SMatthew G. Knepley break; 326df032b82SMatthew G. Knepley case 5: /* 8-node hexahedron */ 327df032b82SMatthew G. Knepley dim = 3; 328df032b82SMatthew G. Knepley numNodes = 8; 329df032b82SMatthew G. Knepley break; 330dea2a3c7SStefano Zampini case 6: /* 6-node wedge */ 331dea2a3c7SStefano Zampini dim = 3; 332dea2a3c7SStefano Zampini numNodes = 6; 333dea2a3c7SStefano Zampini break; 334bf6ba3a3SMatthew G. Knepley case 8: /* 3-node 2nd order line */ 335bf6ba3a3SMatthew G. Knepley dim = 1; 336bf6ba3a3SMatthew G. Knepley numNodes = 2; 337bf6ba3a3SMatthew G. Knepley numNodesIgnore = 1; 338bf6ba3a3SMatthew G. Knepley break; 339bf6ba3a3SMatthew G. Knepley case 9: /* 6-node 2nd order triangle */ 340bf6ba3a3SMatthew G. Knepley dim = 2; 341bf6ba3a3SMatthew G. Knepley numNodes = 3; 342bf6ba3a3SMatthew G. Knepley numNodesIgnore = 3; 343bf6ba3a3SMatthew G. Knepley break; 344dea2a3c7SStefano Zampini case 13: /* 18-node 2nd wedge */ 345dea2a3c7SStefano Zampini dim = 3; 346dea2a3c7SStefano Zampini numNodes = 6; 347dea2a3c7SStefano Zampini numNodesIgnore = 12; 348dea2a3c7SStefano Zampini break; 349df032b82SMatthew G. Knepley case 15: /* 1-node vertex */ 350df032b82SMatthew G. Knepley dim = 0; 351df032b82SMatthew G. Knepley numNodes = 1; 352df032b82SMatthew G. Knepley break; 353bf6ba3a3SMatthew G. Knepley case 7: /* 5-node pyramid */ 354bf6ba3a3SMatthew G. Knepley case 10: /* 9-node 2nd order quadrangle */ 355bf6ba3a3SMatthew G. Knepley case 11: /* 10-node 2nd order tetrahedron */ 356bf6ba3a3SMatthew G. Knepley case 12: /* 27-node 2nd order hexhedron */ 357bf6ba3a3SMatthew G. Knepley case 14: /* 14-node 2nd order pyramid */ 358df032b82SMatthew G. Knepley default: 359df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 360df032b82SMatthew G. Knepley } 361df032b82SMatthew G. Knepley if (binary) { 36211c56cb0SLisandro Dalcin const int nint = 1 + numTags + numNodes + numNodesIgnore; 36311c56cb0SLisandro Dalcin /* Loop over element blocks */ 364df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 36511c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 36611c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 367df032b82SMatthew G. Knepley elements[c].dim = dim; 368df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 369df032b82SMatthew G. Knepley elements[c].numTags = numTags; 370df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 371dea2a3c7SStefano Zampini elements[c].cellType = cellType; 372df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 373df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 374df032b82SMatthew G. Knepley } 375df032b82SMatthew G. Knepley } else { 376*698a79b9SLisandro Dalcin const int nint = numTags + numNodes + numNodesIgnore; 377df032b82SMatthew G. Knepley elements[c].dim = dim; 378df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 379df032b82SMatthew G. Knepley elements[c].numTags = numTags; 380dea2a3c7SStefano Zampini elements[c].cellType = cellType; 381*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 382*698a79b9SLisandro Dalcin for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[p]; 383*698a79b9SLisandro Dalcin for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[numTags + p]; 384df032b82SMatthew G. Knepley c++; 385df032b82SMatthew G. Knepley } 386df032b82SMatthew G. Knepley } 387df032b82SMatthew G. Knepley PetscFunctionReturn(0); 388df032b82SMatthew G. Knepley } 3897d282ae0SMichael Lange 390de78e4feSLisandro Dalcin /* 391de78e4feSLisandro Dalcin $Entities 392de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 393de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 394de78e4feSLisandro Dalcin // points 395de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 396de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 397de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 398de78e4feSLisandro Dalcin ... 399de78e4feSLisandro Dalcin // curves 400de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 401de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 402de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 403de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 404de78e4feSLisandro Dalcin ... 405de78e4feSLisandro Dalcin // surfaces 406de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 407de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 408de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 409de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 410de78e4feSLisandro Dalcin ... 411de78e4feSLisandro Dalcin // volumes 412de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 413de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 414de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 415de78e4feSLisandro Dalcin numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 416de78e4feSLisandro Dalcin ... 417de78e4feSLisandro Dalcin $EndEntities 418de78e4feSLisandro Dalcin */ 419*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v40(GmshFile *gmsh, GmshEntities **entities) 420de78e4feSLisandro Dalcin { 421*698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 422*698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 423*698a79b9SLisandro Dalcin long index, num, lbuf[4]; 424730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 425*698a79b9SLisandro Dalcin PetscInt count[4], i; 426a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 427de78e4feSLisandro Dalcin PetscErrorCode ierr; 428de78e4feSLisandro Dalcin 429de78e4feSLisandro Dalcin PetscFunctionBegin; 430*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 431*698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);} 432*698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 433730356f6SLisandro Dalcin ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 434de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 435730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 436730356f6SLisandro Dalcin ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 437730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);} 438730356f6SLisandro Dalcin ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 439de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 440de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 441de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 442de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 443*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 444de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 445730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 446de78e4feSLisandro Dalcin entity->numTags = numTags = (int) PetscMin(num, 4); 447de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 448de78e4feSLisandro Dalcin if (dim == 0) continue; 449de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 450de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 451*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 452de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 453730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 454de78e4feSLisandro Dalcin } 455de78e4feSLisandro Dalcin } 456de78e4feSLisandro Dalcin PetscFunctionReturn(0); 457de78e4feSLisandro Dalcin } 458de78e4feSLisandro Dalcin 459de78e4feSLisandro Dalcin /* 460de78e4feSLisandro Dalcin $Nodes 461de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 462de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 463de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 464de78e4feSLisandro Dalcin ... 465de78e4feSLisandro Dalcin ... 466de78e4feSLisandro Dalcin $EndNodes 467de78e4feSLisandro Dalcin */ 468*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v40(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 469de78e4feSLisandro Dalcin { 470*698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 471*698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 472de78e4feSLisandro Dalcin long block, node, v, numEntityBlocks, numTotalNodes, numNodes; 473de78e4feSLisandro Dalcin int info[3], nid; 474de78e4feSLisandro Dalcin double *coordinates; 475de78e4feSLisandro Dalcin char *cbuf; 476de78e4feSLisandro Dalcin PetscErrorCode ierr; 477de78e4feSLisandro Dalcin 478de78e4feSLisandro Dalcin PetscFunctionBegin; 479de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 480de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 481de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 482de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 483de78e4feSLisandro Dalcin ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr); 484*698a79b9SLisandro Dalcin *numVertices = numTotalNodes; 485de78e4feSLisandro Dalcin *gmsh_nodes = coordinates; 486de78e4feSLisandro Dalcin for (v = 0, block = 0; block < numEntityBlocks; ++block) { 487de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 488de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 489de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 490*698a79b9SLisandro Dalcin if (gmsh->binary) { 491de78e4feSLisandro Dalcin int nbytes = sizeof(int) + 3*sizeof(double); 492*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 493de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 494de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 495de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 496de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 497de78e4feSLisandro Dalcin #if !defined(PETSC_WORDS_BIGENDIAN) 498de78e4feSLisandro Dalcin ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr); 499de78e4feSLisandro Dalcin ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr); 500de78e4feSLisandro Dalcin #endif 501de78e4feSLisandro Dalcin ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 502de78e4feSLisandro Dalcin ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 503de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 504de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 505de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 506de78e4feSLisandro Dalcin } 507de78e4feSLisandro Dalcin } else { 508de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 509de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 510de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 511de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 512de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 513de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 514de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 515de78e4feSLisandro Dalcin } 516de78e4feSLisandro Dalcin } 517de78e4feSLisandro Dalcin } 518de78e4feSLisandro Dalcin PetscFunctionReturn(0); 519de78e4feSLisandro Dalcin } 520de78e4feSLisandro Dalcin 521de78e4feSLisandro Dalcin /* 522de78e4feSLisandro Dalcin $Elements 523de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 524de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 525de78e4feSLisandro Dalcin tag(int) numVert[...](int) 526de78e4feSLisandro Dalcin ... 527de78e4feSLisandro Dalcin ... 528de78e4feSLisandro Dalcin $EndElements 529de78e4feSLisandro Dalcin */ 530*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v40(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems) 531de78e4feSLisandro Dalcin { 532*698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 533*698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 534de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 535de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 536de78e4feSLisandro Dalcin int eid, dim, numTags, *tags, cellType, numNodes; 537a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 538de78e4feSLisandro Dalcin GmshElement *elements; 539de78e4feSLisandro Dalcin PetscErrorCode ierr; 540de78e4feSLisandro Dalcin 541de78e4feSLisandro Dalcin PetscFunctionBegin; 542de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 543de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 544de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 545de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 546de78e4feSLisandro Dalcin ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr); 547*698a79b9SLisandro Dalcin *numCells = numTotalElements; 548de78e4feSLisandro Dalcin *gmsh_elems = elements; 549de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 550de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 551de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 552de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 553730356f6SLisandro Dalcin ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 554730356f6SLisandro Dalcin numTags = entity->numTags; 555730356f6SLisandro Dalcin tags = entity->tags; 556de78e4feSLisandro Dalcin switch (cellType) { 557de78e4feSLisandro Dalcin case 1: /* 2-node line */ 558de78e4feSLisandro Dalcin numNodes = 2; 559de78e4feSLisandro Dalcin break; 560de78e4feSLisandro Dalcin case 2: /* 3-node triangle */ 561*698a79b9SLisandro Dalcin numNodes = 3; 562*698a79b9SLisandro Dalcin break; 563de78e4feSLisandro Dalcin case 3: /* 4-node quadrangle */ 564de78e4feSLisandro Dalcin numNodes = 4; 565de78e4feSLisandro Dalcin break; 566de78e4feSLisandro Dalcin case 4: /* 4-node tetrahedron */ 567de78e4feSLisandro Dalcin numNodes = 4; 568de78e4feSLisandro Dalcin break; 569de78e4feSLisandro Dalcin case 5: /* 8-node hexahedron */ 570de78e4feSLisandro Dalcin numNodes = 8; 571de78e4feSLisandro Dalcin break; 572de78e4feSLisandro Dalcin case 6: /* 6-node wedge */ 573de78e4feSLisandro Dalcin numNodes = 6; 574de78e4feSLisandro Dalcin break; 575de78e4feSLisandro Dalcin case 15: /* 1-node vertex */ 576de78e4feSLisandro Dalcin numNodes = 1; 577de78e4feSLisandro Dalcin break; 578de78e4feSLisandro Dalcin default: 579de78e4feSLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 580de78e4feSLisandro Dalcin } 581de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 582de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 583*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 584de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 585de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 586de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 587de78e4feSLisandro Dalcin int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 588de78e4feSLisandro Dalcin GmshElement *element = elements + c; 589de78e4feSLisandro Dalcin element->dim = dim; 590de78e4feSLisandro Dalcin element->cellType = cellType; 591de78e4feSLisandro Dalcin element->numNodes = numNodes; 592de78e4feSLisandro Dalcin element->numTags = numTags; 593de78e4feSLisandro Dalcin element->id = id[0]; 594de78e4feSLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 595de78e4feSLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 596de78e4feSLisandro Dalcin } 597de78e4feSLisandro Dalcin } 598*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 599*698a79b9SLisandro Dalcin } 600*698a79b9SLisandro Dalcin 601*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v40(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt) 602*698a79b9SLisandro Dalcin { 603*698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 604*698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 605*698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 606*698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 607*698a79b9SLisandro Dalcin int numPeriodic, snum, i; 608*698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 609*698a79b9SLisandro Dalcin PetscErrorCode ierr; 610*698a79b9SLisandro Dalcin 611*698a79b9SLisandro Dalcin PetscFunctionBegin; 612*698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 613*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 614*698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 615*698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 616*698a79b9SLisandro Dalcin } else { 617*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 618*698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 619*698a79b9SLisandro Dalcin } 620*698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 621*698a79b9SLisandro Dalcin int ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode; 622*698a79b9SLisandro Dalcin long j, nNodes; 623*698a79b9SLisandro Dalcin double affine[16]; 624*698a79b9SLisandro Dalcin 625*698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 626*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 627*698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag); 628*698a79b9SLisandro Dalcin if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 629*698a79b9SLisandro Dalcin } else { 630*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 631*698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 632*698a79b9SLisandro Dalcin slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2]; 633*698a79b9SLisandro Dalcin } 634*698a79b9SLisandro Dalcin (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */ 635*698a79b9SLisandro Dalcin 636*698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 637*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 638*698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 639*698a79b9SLisandro Dalcin if (snum != 1) { /* discard tranformation and try again */ 640*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 641*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 642*698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 643*698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 644*698a79b9SLisandro Dalcin } 645*698a79b9SLisandro Dalcin } else { 646*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 647*698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 648*698a79b9SLisandro Dalcin if (nNodes == -1) { /* discard tranformation and try again */ 649*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 650*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 651*698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 652*698a79b9SLisandro Dalcin } 653*698a79b9SLisandro Dalcin } 654*698a79b9SLisandro Dalcin 655*698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 656*698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 657*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 658*698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 659*698a79b9SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 660*698a79b9SLisandro Dalcin } else { 661*698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 662*698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 663*698a79b9SLisandro Dalcin slaveNode = ibuf[0]; masterNode = ibuf[1]; 664*698a79b9SLisandro Dalcin } 665*698a79b9SLisandro Dalcin slaveMap[slaveNode - shift] = masterNode - shift; 666*698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, slaveNode - shift);CHKERRQ(ierr); 667*698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, masterNode - shift);CHKERRQ(ierr); 668*698a79b9SLisandro Dalcin } 669*698a79b9SLisandro Dalcin } 670*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 671*698a79b9SLisandro Dalcin } 672*698a79b9SLisandro Dalcin 673*698a79b9SLisandro Dalcin /* 674*698a79b9SLisandro Dalcin $Entities 675*698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 676*698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 677*698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 678*698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 679*698a79b9SLisandro Dalcin ... 680*698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 681*698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 682*698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 683*698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 684*698a79b9SLisandro Dalcin ... 685*698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 686*698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 687*698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 688*698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 689*698a79b9SLisandro Dalcin ... 690*698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 691*698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 692*698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 693*698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 694*698a79b9SLisandro Dalcin ... 695*698a79b9SLisandro Dalcin $EndEntities 696*698a79b9SLisandro Dalcin */ 697*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v41(GmshFile *gmsh, GmshEntities **entities) 698*698a79b9SLisandro Dalcin { 699*698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 700*698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 701*698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 702*698a79b9SLisandro Dalcin PetscErrorCode ierr; 703*698a79b9SLisandro Dalcin 704*698a79b9SLisandro Dalcin PetscFunctionBegin; 705*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr); 706*698a79b9SLisandro Dalcin ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 707*698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 708*698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 709*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr); 710*698a79b9SLisandro Dalcin ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 711*698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr); 712*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 713*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 714*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 715*698a79b9SLisandro Dalcin entity->numTags = PetscMin(numTags, 4); 716*698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 717*698a79b9SLisandro Dalcin if (dim == 0) continue; 718*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 719*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 720*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 721*698a79b9SLisandro Dalcin } 722*698a79b9SLisandro Dalcin } 723*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 724*698a79b9SLisandro Dalcin } 725*698a79b9SLisandro Dalcin 726*698a79b9SLisandro Dalcin /* 727*698a79b9SLisandro Dalcin $Nodes 728*698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 729*698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 730*698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 731*698a79b9SLisandro Dalcin nodeTag(size_t) 732*698a79b9SLisandro Dalcin ... 733*698a79b9SLisandro Dalcin x(double) y(double) z(double) 734*698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 735*698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 736*698a79b9SLisandro Dalcin ... 737*698a79b9SLisandro Dalcin ... 738*698a79b9SLisandro Dalcin $EndNodes 739*698a79b9SLisandro Dalcin */ 740*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v41(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 741*698a79b9SLisandro Dalcin { 742*698a79b9SLisandro Dalcin int info[3]; 743*698a79b9SLisandro Dalcin PetscInt sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, *nodeTag = NULL, block, node, i; 744*698a79b9SLisandro Dalcin double *coordinates; 745*698a79b9SLisandro Dalcin PetscErrorCode ierr; 746*698a79b9SLisandro Dalcin 747*698a79b9SLisandro Dalcin PetscFunctionBegin; 748*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 749*698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numNodes = sizes[1]; 750*698a79b9SLisandro Dalcin ierr = PetscMalloc1(numNodes*3, &coordinates);CHKERRQ(ierr); 751*698a79b9SLisandro Dalcin *numVertices = numNodes; 752*698a79b9SLisandro Dalcin *gmsh_nodes = coordinates; 753*698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 754*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 755*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr); 756*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodesBlock, sizeof(PetscInt), &nodeTag);CHKERRQ(ierr); 757*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTag, numNodesBlock);CHKERRQ(ierr); 758*698a79b9SLisandro Dalcin for (i = 0; i < numNodesBlock; ++i) if (nodeTag[i] != node+i+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nodeTag[i], node+i+shift); 759*698a79b9SLisandro Dalcin if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 760*698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, coordinates+node*3, numNodesBlock*3);CHKERRQ(ierr); 761*698a79b9SLisandro Dalcin } 762*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 763*698a79b9SLisandro Dalcin } 764*698a79b9SLisandro Dalcin 765*698a79b9SLisandro Dalcin /* 766*698a79b9SLisandro Dalcin $Elements 767*698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 768*698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 769*698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 770*698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 771*698a79b9SLisandro Dalcin ... 772*698a79b9SLisandro Dalcin ... 773*698a79b9SLisandro Dalcin $EndElements 774*698a79b9SLisandro Dalcin */ 775*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v41(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems) 776*698a79b9SLisandro Dalcin { 777*698a79b9SLisandro Dalcin int info[3], eid, dim, cellType, *tags; 778*698a79b9SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numNodes, numTags, block, elem, c, p; 779*698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 780*698a79b9SLisandro Dalcin GmshElement *elements; 781*698a79b9SLisandro Dalcin PetscErrorCode ierr; 782*698a79b9SLisandro Dalcin 783*698a79b9SLisandro Dalcin PetscFunctionBegin; 784*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 785*698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numElements = sizes[1]; 786*698a79b9SLisandro Dalcin ierr = PetscCalloc1(numElements, &elements);CHKERRQ(ierr); 787*698a79b9SLisandro Dalcin *numCells = numElements; 788*698a79b9SLisandro Dalcin *gmsh_elems = elements; 789*698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 790*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 791*698a79b9SLisandro Dalcin dim = info[0]; eid = info[1]; cellType = info[2]; 792*698a79b9SLisandro Dalcin ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 793*698a79b9SLisandro Dalcin numTags = entity->numTags; 794*698a79b9SLisandro Dalcin tags = entity->tags; 795*698a79b9SLisandro Dalcin switch (cellType) { 796*698a79b9SLisandro Dalcin case 1: /* 2-node line */ 797*698a79b9SLisandro Dalcin numNodes = 2; 798*698a79b9SLisandro Dalcin break; 799*698a79b9SLisandro Dalcin case 2: /* 3-node triangle */ 800*698a79b9SLisandro Dalcin numNodes = 3; 801*698a79b9SLisandro Dalcin break; 802*698a79b9SLisandro Dalcin case 3: /* 4-node quadrangle */ 803*698a79b9SLisandro Dalcin numNodes = 4; 804*698a79b9SLisandro Dalcin break; 805*698a79b9SLisandro Dalcin case 4: /* 4-node tetrahedron */ 806*698a79b9SLisandro Dalcin numNodes = 4; 807*698a79b9SLisandro Dalcin break; 808*698a79b9SLisandro Dalcin case 5: /* 8-node hexahedron */ 809*698a79b9SLisandro Dalcin numNodes = 8; 810*698a79b9SLisandro Dalcin break; 811*698a79b9SLisandro Dalcin case 6: /* 6-node wedge */ 812*698a79b9SLisandro Dalcin numNodes = 6; 813*698a79b9SLisandro Dalcin break; 814*698a79b9SLisandro Dalcin case 15: /* 1-node vertex */ 815*698a79b9SLisandro Dalcin numNodes = 1; 816*698a79b9SLisandro Dalcin break; 817*698a79b9SLisandro Dalcin default: 818*698a79b9SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 819*698a79b9SLisandro Dalcin } 820*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr); 821*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr); 822*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr); 823*698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 824*698a79b9SLisandro Dalcin GmshElement *element = elements + c; 825*698a79b9SLisandro Dalcin PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 826*698a79b9SLisandro Dalcin element->id = id[0]; 827*698a79b9SLisandro Dalcin element->dim = dim; 828*698a79b9SLisandro Dalcin element->cellType = cellType; 829*698a79b9SLisandro Dalcin element->numNodes = numNodes; 830*698a79b9SLisandro Dalcin element->numTags = numTags; 831*698a79b9SLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 832*698a79b9SLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 833*698a79b9SLisandro Dalcin } 834*698a79b9SLisandro Dalcin } 835*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 836*698a79b9SLisandro Dalcin } 837*698a79b9SLisandro Dalcin 838*698a79b9SLisandro Dalcin /* 839*698a79b9SLisandro Dalcin $Periodic 840*698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 841*698a79b9SLisandro Dalcin entityDim(int) entityTag(int) entityTagMaster(int) 842*698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 843*698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 844*698a79b9SLisandro Dalcin nodeTag(size_t) nodeTagMaster(size_t) 845*698a79b9SLisandro Dalcin ... 846*698a79b9SLisandro Dalcin ... 847*698a79b9SLisandro Dalcin $EndPeriodic 848*698a79b9SLisandro Dalcin */ 849*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v41(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt) 850*698a79b9SLisandro Dalcin { 851*698a79b9SLisandro Dalcin int info[3]; 852*698a79b9SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 853*698a79b9SLisandro Dalcin double dbuf[16]; 854*698a79b9SLisandro Dalcin PetscErrorCode ierr; 855*698a79b9SLisandro Dalcin 856*698a79b9SLisandro Dalcin PetscFunctionBegin; 857*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr); 858*698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 859*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 860*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr); 861*698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr); 862*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr); 863*698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr); 864*698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr); 865*698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 866*698a79b9SLisandro Dalcin PetscInt slaveNode = nodeTags[node*2+0] - shift; 867*698a79b9SLisandro Dalcin PetscInt masterNode = nodeTags[node*2+1] - shift; 868*698a79b9SLisandro Dalcin slaveMap[slaveNode] = masterNode; 869*698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, slaveNode);CHKERRQ(ierr); 870*698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, masterNode);CHKERRQ(ierr); 871*698a79b9SLisandro Dalcin } 872*698a79b9SLisandro Dalcin } 873*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 874*698a79b9SLisandro Dalcin } 875*698a79b9SLisandro Dalcin 876*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadMeshFormat(GmshFile *gmsh) 877*698a79b9SLisandro Dalcin { 878*698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 879*698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 880*698a79b9SLisandro Dalcin float version; 881*698a79b9SLisandro Dalcin PetscErrorCode ierr; 882*698a79b9SLisandro Dalcin 883*698a79b9SLisandro Dalcin PetscFunctionBegin; 884*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr); 885*698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 886*698a79b9SLisandro Dalcin if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 887*698a79b9SLisandro Dalcin if (version < 2.2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f must be at least 2.2", (double)version); 888*698a79b9SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version); 889*698a79b9SLisandro Dalcin if (version > 4.1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f must be at most 4.1", (double)version); 890*698a79b9SLisandro Dalcin if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is binary but Gmsh file is ASCII"); 891*698a79b9SLisandro Dalcin if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is ASCII but Gmsh file is binary"); 892*698a79b9SLisandro Dalcin fileFormat = (int)(version*10); /* XXX Should use (int)roundf(version*10) ? */ 893*698a79b9SLisandro Dalcin if (fileFormat <= 40 && dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 894*698a79b9SLisandro Dalcin if (fileFormat >= 41 && dataSize != sizeof(int) && dataSize != sizeof(PetscInt64)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 895*698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 896*698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 897*698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 898*698a79b9SLisandro Dalcin if (gmsh->binary) { 899*698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr); 900*698a79b9SLisandro Dalcin if (checkEndian != 1) { 901*698a79b9SLisandro Dalcin ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr); 902*698a79b9SLisandro Dalcin if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to detect endianness in Gmsh file header: %s", line); 903*698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 904*698a79b9SLisandro Dalcin } 905*698a79b9SLisandro Dalcin } 906*698a79b9SLisandro Dalcin PetscFunctionReturn(0); 907*698a79b9SLisandro Dalcin } 908*698a79b9SLisandro Dalcin 909*698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPhysicalNames(GmshFile *gmsh) 910*698a79b9SLisandro Dalcin { 911*698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 912*698a79b9SLisandro Dalcin int snum, numRegions, region; 913*698a79b9SLisandro Dalcin PetscErrorCode ierr; 914*698a79b9SLisandro Dalcin 915*698a79b9SLisandro Dalcin PetscFunctionBegin; 916*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 917*698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numRegions); 918*698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 919*698a79b9SLisandro Dalcin for (region = 0; region < numRegions; ++region) { 920*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr); 921*698a79b9SLisandro Dalcin } 922de78e4feSLisandro Dalcin PetscFunctionReturn(0); 923de78e4feSLisandro Dalcin } 924de78e4feSLisandro Dalcin 925331830f3SMatthew G. Knepley /*@ 9267d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 927331830f3SMatthew G. Knepley 928331830f3SMatthew G. Knepley Collective on comm 929331830f3SMatthew G. Knepley 930331830f3SMatthew G. Knepley Input Parameters: 931331830f3SMatthew G. Knepley + comm - The MPI communicator 932331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 933331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 934331830f3SMatthew G. Knepley 935331830f3SMatthew G. Knepley Output Parameter: 936331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 937331830f3SMatthew G. Knepley 938*698a79b9SLisandro Dalcin Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 939331830f3SMatthew G. Knepley 940331830f3SMatthew G. Knepley Level: beginner 941331830f3SMatthew G. Knepley 942331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 943331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 944331830f3SMatthew G. Knepley @*/ 945331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 946331830f3SMatthew G. Knepley { 94711c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 94811c56cb0SLisandro Dalcin double *coordsIn = NULL; 949730356f6SLisandro Dalcin GmshEntities *entities = NULL; 9503c67d5adSSatish Balay GmshElement *gmsh_elem = NULL; 951331830f3SMatthew G. Knepley PetscSection coordSection; 952331830f3SMatthew G. Knepley Vec coordinates; 9536fbe17bfSStefano Zampini PetscBT periodicV = NULL, periodicC = NULL; 95484572febSToby Isaac PetscScalar *coords; 955*698a79b9SLisandro Dalcin PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 956*698a79b9SLisandro Dalcin PetscInt numVertices = 0, numCells = 0, trueNumCells = 0; 957*698a79b9SLisandro Dalcin int i, shift = 1; 95811c56cb0SLisandro Dalcin PetscMPIInt rank; 959*698a79b9SLisandro Dalcin PetscBool binary, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 960dea2a3c7SStefano Zampini PetscBool enable_hybrid = PETSC_FALSE; 961331830f3SMatthew G. Knepley PetscErrorCode ierr; 962331830f3SMatthew G. Knepley 963331830f3SMatthew G. Knepley PetscFunctionBegin; 964331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 965dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr); 966dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 967dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 968dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 969dea2a3c7SStefano Zampini ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr); 97011c56cb0SLisandro Dalcin if (zerobase) shift = 0; 97111c56cb0SLisandro Dalcin 972331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 973331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 9743b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 97511c56cb0SLisandro Dalcin 97611c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 97711c56cb0SLisandro Dalcin 97811c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 9793b46f529SLisandro Dalcin if (binary) { 98011c56cb0SLisandro Dalcin parentviewer = viewer; 98111c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 98211c56cb0SLisandro Dalcin } 98311c56cb0SLisandro Dalcin 9843b46f529SLisandro Dalcin if (!rank) { 985*698a79b9SLisandro Dalcin GmshFile gmsh_, *gmsh = &gmsh_; 986*698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 987*698a79b9SLisandro Dalcin PetscBool match; 988331830f3SMatthew G. Knepley 989*698a79b9SLisandro Dalcin ierr = PetscMemzero(gmsh,sizeof(GmshFile));CHKERRQ(ierr); 990*698a79b9SLisandro Dalcin gmsh->viewer = viewer; 991*698a79b9SLisandro Dalcin gmsh->binary = binary; 992*698a79b9SLisandro Dalcin 993*698a79b9SLisandro Dalcin /* Read mesh format */ 994*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 99504d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 996331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 997*698a79b9SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadMeshFormat(gmsh);CHKERRQ(ierr); 998*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 99904d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1000331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 100111c56cb0SLisandro Dalcin 1002dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 1003*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1004dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1005dc0ede3bSMatthew G. Knepley if (match) { 1006*698a79b9SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadPhysicalNames(gmsh);CHKERRQ(ierr); 1007*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1008dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1009dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1010*698a79b9SLisandro Dalcin /* Initial read for entity section */ 1011*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1012dc0ede3bSMatthew G. Knepley } 101311c56cb0SLisandro Dalcin 1014de78e4feSLisandro Dalcin /* Read entities */ 1015*698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 1016de78e4feSLisandro Dalcin ierr = PetscStrncmp(line, "$Entities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1017de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1018*698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1019*698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadEntities_v41(gmsh, &entities);CHKERRQ(ierr); break; 1020*698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadEntities_v40(gmsh, &entities);CHKERRQ(ierr); break; 1021*698a79b9SLisandro Dalcin } 1022*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1023de78e4feSLisandro Dalcin ierr = PetscStrncmp(line, "$EndEntities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1024de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1025*698a79b9SLisandro Dalcin /* Initial read for nodes section */ 1026*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1027de78e4feSLisandro Dalcin } 1028de78e4feSLisandro Dalcin 1029*698a79b9SLisandro Dalcin /* Read nodes */ 103004d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1031331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1032*698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1033*698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadNodes_v41(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1034*698a79b9SLisandro Dalcin case 40: ierr = DMPlexCreateGmsh_ReadNodes_v40(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1035*698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadNodes_v22(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1036de78e4feSLisandro Dalcin } 1037*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 103804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1039331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 104011c56cb0SLisandro Dalcin 1041*698a79b9SLisandro Dalcin /* Read elements */ 1042*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);; 104304d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1044331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1045*698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1046*698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadElements_v41(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1047*698a79b9SLisandro Dalcin case 40: ierr = DMPlexCreateGmsh_ReadElements_v40(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1048*698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadElements_v22(gmsh, shift, &numCells, &gmsh_elem);CHKERRQ(ierr); 1049de78e4feSLisandro Dalcin } 1050*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1051*698a79b9SLisandro Dalcin ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1052de78e4feSLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1053de78e4feSLisandro Dalcin 1054*698a79b9SLisandro Dalcin /* OPTIONAL Read periodic section */ 1055*698a79b9SLisandro Dalcin if (periodic) { 1056*698a79b9SLisandro Dalcin PetscInt pVert, *periodicMapT, *aux; 1057de78e4feSLisandro Dalcin 1058*698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 1059*698a79b9SLisandro Dalcin ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 1060*698a79b9SLisandro Dalcin for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 1061*698a79b9SLisandro Dalcin 1062*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1063*698a79b9SLisandro Dalcin ierr = PetscStrncmp(line, "$Periodic", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1064*698a79b9SLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1065*698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1066*698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadPeriodic_v41(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break; 1067*698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadPeriodic_v40(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break; 1068*698a79b9SLisandro Dalcin } 1069*698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1070*698a79b9SLisandro Dalcin ierr = PetscStrncmp(line, "$EndPeriodic", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 1071*698a79b9SLisandro Dalcin if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 1072*698a79b9SLisandro Dalcin 1073*698a79b9SLisandro Dalcin /* we may have slaves of slaves */ 1074*698a79b9SLisandro Dalcin for (i = 0; i < numVertices; i++) { 1075*698a79b9SLisandro Dalcin while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 1076*698a79b9SLisandro Dalcin periodicMapT[i] = periodicMapT[periodicMapT[i]]; 1077*698a79b9SLisandro Dalcin } 1078*698a79b9SLisandro Dalcin } 1079*698a79b9SLisandro Dalcin /* periodicMap : from old to new numbering (periodic vertices excluded) 1080*698a79b9SLisandro Dalcin periodicMapI: from new to old numbering */ 1081*698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 1082*698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 1083*698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 1084*698a79b9SLisandro Dalcin for (i = 0, pVert = 0; i < numVertices; i++) { 1085*698a79b9SLisandro Dalcin if (periodicMapT[i] != i) { 1086*698a79b9SLisandro Dalcin pVert++; 1087*698a79b9SLisandro Dalcin } else { 1088*698a79b9SLisandro Dalcin aux[i] = i - pVert; 1089*698a79b9SLisandro Dalcin periodicMapI[i - pVert] = i; 1090*698a79b9SLisandro Dalcin } 1091*698a79b9SLisandro Dalcin } 1092*698a79b9SLisandro Dalcin for (i = 0 ; i < numVertices; i++) { 1093*698a79b9SLisandro Dalcin periodicMap[i] = aux[periodicMapT[i]]; 1094*698a79b9SLisandro Dalcin } 1095*698a79b9SLisandro Dalcin ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 1096*698a79b9SLisandro Dalcin ierr = PetscFree(aux);CHKERRQ(ierr); 1097*698a79b9SLisandro Dalcin /* remove periodic vertices */ 1098*698a79b9SLisandro Dalcin numVertices = numVertices - pVert; 1099*698a79b9SLisandro Dalcin } 1100*698a79b9SLisandro Dalcin 1101*698a79b9SLisandro Dalcin ierr = GmshEntitiesDestroy(&entities);CHKERRQ(ierr); 1102*698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 1103*698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 1104*698a79b9SLisandro Dalcin } 1105*698a79b9SLisandro Dalcin 1106*698a79b9SLisandro Dalcin if (parentviewer) { 1107*698a79b9SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 1108*698a79b9SLisandro Dalcin } 1109*698a79b9SLisandro Dalcin 1110*698a79b9SLisandro Dalcin if (!rank) { 1111*698a79b9SLisandro Dalcin PetscBool hybrid = PETSC_FALSE; 1112*698a79b9SLisandro Dalcin 1113a4bb7517SMichael Lange for (trueNumCells = 0, c = 0; c < numCells; ++c) { 1114dea2a3c7SStefano Zampini int on = -1; 1115a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 1116dea2a3c7SStefano 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++;} 1117dea2a3c7SStefano Zampini /* wedges always indicate an hybrid mesh in PLEX */ 1118dea2a3c7SStefano Zampini if (on == 6 || on == 13) hybrid = PETSC_TRUE; 1119db397164SMichael Lange } 1120dea2a3c7SStefano Zampini /* Renumber cells for hybrid grids */ 1121dea2a3c7SStefano Zampini if (hybrid && enable_hybrid) { 1122dea2a3c7SStefano Zampini PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 1123dea2a3c7SStefano Zampini PetscInt cell, tn, *tp; 1124dea2a3c7SStefano Zampini int n1 = 0,n2 = 0; 1125dea2a3c7SStefano Zampini 1126dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 1127dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 1128dea2a3c7SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1129dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) { 1130dea2a3c7SStefano Zampini if (!n1) n1 = gmsh_elem[c].cellType; 1131dea2a3c7SStefano Zampini else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 1132dea2a3c7SStefano Zampini 1133dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 1134dea2a3c7SStefano Zampini else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 1135dea2a3c7SStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 1136dea2a3c7SStefano Zampini cell++; 1137dea2a3c7SStefano Zampini } 1138dea2a3c7SStefano Zampini } 1139dea2a3c7SStefano Zampini 1140dea2a3c7SStefano Zampini switch (n1) { 1141dea2a3c7SStefano Zampini case 2: /* triangles */ 1142dea2a3c7SStefano Zampini case 9: 1143dea2a3c7SStefano Zampini switch (n2) { 1144dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1145dea2a3c7SStefano Zampini case 3: /* quads */ 1146dea2a3c7SStefano Zampini case 10: 1147dea2a3c7SStefano Zampini break; 1148dea2a3c7SStefano Zampini default: 1149dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1150dea2a3c7SStefano Zampini } 1151dea2a3c7SStefano Zampini break; 1152dea2a3c7SStefano Zampini case 3: 1153dea2a3c7SStefano Zampini case 10: 1154dea2a3c7SStefano Zampini switch (n2) { 1155dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1156dea2a3c7SStefano Zampini case 2: /* swap since we list simplices first */ 1157dea2a3c7SStefano Zampini case 9: 1158dea2a3c7SStefano Zampini tn = hc1; 1159dea2a3c7SStefano Zampini hc1 = hc2; 1160dea2a3c7SStefano Zampini hc2 = tn; 1161dea2a3c7SStefano Zampini 1162dea2a3c7SStefano Zampini tp = hybridCells1; 1163dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 1164dea2a3c7SStefano Zampini hybridCells2 = tp; 1165dea2a3c7SStefano Zampini break; 1166dea2a3c7SStefano Zampini default: 1167dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1168dea2a3c7SStefano Zampini } 1169dea2a3c7SStefano Zampini break; 1170dea2a3c7SStefano Zampini case 4: /* tetrahedra */ 1171dea2a3c7SStefano Zampini case 11: 1172dea2a3c7SStefano Zampini switch (n2) { 1173dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1174dea2a3c7SStefano Zampini case 6: /* wedges */ 1175dea2a3c7SStefano Zampini case 13: 1176dea2a3c7SStefano Zampini break; 1177dea2a3c7SStefano Zampini default: 1178dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1179dea2a3c7SStefano Zampini } 1180dea2a3c7SStefano Zampini break; 1181dea2a3c7SStefano Zampini case 6: 1182dea2a3c7SStefano Zampini case 13: 1183dea2a3c7SStefano Zampini switch (n2) { 1184dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1185dea2a3c7SStefano Zampini case 4: /* swap since we list simplices first */ 1186dea2a3c7SStefano Zampini case 11: 1187dea2a3c7SStefano Zampini tn = hc1; 1188dea2a3c7SStefano Zampini hc1 = hc2; 1189dea2a3c7SStefano Zampini hc2 = tn; 1190dea2a3c7SStefano Zampini 1191dea2a3c7SStefano Zampini tp = hybridCells1; 1192dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 1193dea2a3c7SStefano Zampini hybridCells2 = tp; 1194dea2a3c7SStefano Zampini break; 1195dea2a3c7SStefano Zampini default: 1196dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1197dea2a3c7SStefano Zampini } 1198dea2a3c7SStefano Zampini break; 1199dea2a3c7SStefano Zampini default: 1200dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1201dea2a3c7SStefano Zampini } 1202dea2a3c7SStefano Zampini cMax = hc1; 1203dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 1204dea2a3c7SStefano Zampini for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 1205dea2a3c7SStefano Zampini for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 1206dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 1207dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 1208dea2a3c7SStefano Zampini } 1209dea2a3c7SStefano Zampini 121011c56cb0SLisandro Dalcin } 121111c56cb0SLisandro Dalcin 1212a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 1213db397164SMichael Lange ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 1214a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 1215a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 1216dea2a3c7SStefano Zampini ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 1217a4bb7517SMichael Lange cell++; 1218db397164SMichael Lange } 1219331830f3SMatthew G. Knepley } 1220331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 1221a4bb7517SMichael Lange /* Add cell-vertex connections */ 1222a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 1223a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 122411c56cb0SLisandro Dalcin PetscInt pcone[8], corner; 1225a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1226dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 1227917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 1228db397164SMichael Lange } 122997ed6be6Ssarens if (dim == 3) { 123097ed6be6Ssarens /* Tetrahedra are inverted */ 1231dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 123297ed6be6Ssarens PetscInt tmp = pcone[0]; 123397ed6be6Ssarens pcone[0] = pcone[1]; 123497ed6be6Ssarens pcone[1] = tmp; 123597ed6be6Ssarens } 123697ed6be6Ssarens /* Hexahedra are inverted */ 1237dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 123897ed6be6Ssarens PetscInt tmp = pcone[1]; 123997ed6be6Ssarens pcone[1] = pcone[3]; 124097ed6be6Ssarens pcone[3] = tmp; 124197ed6be6Ssarens } 1242dea2a3c7SStefano Zampini /* Prisms are inverted */ 1243dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 1244dea2a3c7SStefano Zampini PetscInt tmp; 1245dea2a3c7SStefano Zampini 1246dea2a3c7SStefano Zampini tmp = pcone[1]; 1247dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1248dea2a3c7SStefano Zampini pcone[2] = tmp; 1249dea2a3c7SStefano Zampini tmp = pcone[4]; 1250dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1251dea2a3c7SStefano Zampini pcone[5] = tmp; 125297ed6be6Ssarens } 1253dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1254dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1255dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 1256dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1257dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1258dea2a3c7SStefano Zampini pcone[3] = tmp; 1259dea2a3c7SStefano Zampini } 1260dea2a3c7SStefano Zampini } 1261dea2a3c7SStefano Zampini ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 1262a4bb7517SMichael Lange cell++; 1263331830f3SMatthew G. Knepley } 1264a4bb7517SMichael Lange } 12656228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1266c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1267dea2a3c7SStefano Zampini ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1268331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1269331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1270331830f3SMatthew G. Knepley if (interpolate) { 12715fd9971aSMatthew G. Knepley DM idm; 1272331830f3SMatthew G. Knepley 1273331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1274331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 1275331830f3SMatthew G. Knepley *dm = idm; 1276331830f3SMatthew G. Knepley } 1277917266a4SLisandro Dalcin 1278f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 1279917266a4SLisandro Dalcin if (!rank && usemarker) { 1280d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 1281d3f73514SStefano Zampini 1282d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 1283d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1284d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 1285d3f73514SStefano Zampini PetscInt suppSize; 1286d3f73514SStefano Zampini 1287d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 1288d3f73514SStefano Zampini if (suppSize == 1) { 1289d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 1290d3f73514SStefano Zampini 1291d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1292d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 1293d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 1294d3f73514SStefano Zampini } 1295d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1296d3f73514SStefano Zampini } 1297d3f73514SStefano Zampini } 1298d3f73514SStefano Zampini } 129916ad7e67SMichael Lange 130016ad7e67SMichael Lange if (!rank) { 130111c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1302d1a54cefSMatthew G. Knepley 130316ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 130411c56cb0SLisandro Dalcin for (cell = 0, c = 0; c < numCells; ++c) { 130511c56cb0SLisandro Dalcin 130611c56cb0SLisandro Dalcin /* Create face sets */ 13075964b3dcSLisandro Dalcin if (interpolate && gmsh_elem[c].dim == dim-1) { 130816ad7e67SMichael Lange const PetscInt *join; 130911c56cb0SLisandro Dalcin PetscInt joinSize, pcone[8], corner; 131011c56cb0SLisandro Dalcin /* Find the relevant facet with vertex joins */ 1311a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1312dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 1313917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 131416ad7e67SMichael Lange } 131511c56cb0SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 1316f10f1c67SMatthew 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); 1317c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 1318a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 131916ad7e67SMichael Lange } 13206e1dd89aSLawrence Mitchell 13216e1dd89aSLawrence Mitchell /* Create cell sets */ 13226e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 13236e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 1324dea2a3c7SStefano Zampini ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 13256e1dd89aSLawrence Mitchell } 1326917266a4SLisandro Dalcin cell++; 13276e1dd89aSLawrence Mitchell } 13280c070f12SSander Arens 13290c070f12SSander Arens /* Create vertex sets */ 13300c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 13310c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 1332917266a4SLisandro Dalcin const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 133311c56cb0SLisandro Dalcin const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 1334d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 13350c070f12SSander Arens } 13360c070f12SSander Arens } 13370c070f12SSander Arens } 133816ad7e67SMichael Lange } 133916ad7e67SMichael Lange 134011c56cb0SLisandro Dalcin /* Create coordinates */ 1341dea2a3c7SStefano Zampini if (embedDim < 0) embedDim = dim; 1342dea2a3c7SStefano Zampini ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 1343979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1344331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 13451d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 1346f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 1347f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 1348f45c9589SStefano Zampini } else { 134975b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 1350f45c9589SStefano Zampini } 135175b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 13521d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 13531d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 1354331830f3SMatthew G. Knepley } 135511c56cb0SLisandro Dalcin if (periodicMap) { 1356437bdd5bSStefano Zampini ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 1357f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1358f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 1359437bdd5bSStefano Zampini PetscInt corner; 136011c56cb0SLisandro Dalcin PetscBool pc = PETSC_FALSE; 1361437bdd5bSStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1362917266a4SLisandro Dalcin pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 1363437bdd5bSStefano Zampini } 1364437bdd5bSStefano Zampini if (pc) { 136511c56cb0SLisandro Dalcin PetscInt dof = gmsh_elem[c].numNodes*embedDim; 1366dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1367dea2a3c7SStefano Zampini ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 1368dea2a3c7SStefano Zampini ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 1369dea2a3c7SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 13706fbe17bfSStefano Zampini } 1371f45c9589SStefano Zampini cell++; 1372f45c9589SStefano Zampini } 1373f45c9589SStefano Zampini } 1374f45c9589SStefano Zampini } 1375331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1376331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 13778b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1378331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1379331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 13801d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 1381331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 1382331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1383f45c9589SStefano Zampini if (periodicMap) { 1384f45c9589SStefano Zampini PetscInt off; 1385f45c9589SStefano Zampini 1386f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1387f45c9589SStefano Zampini PetscInt pcone[8], corner; 1388f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 1389dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1390dea2a3c7SStefano Zampini if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 1391f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1392917266a4SLisandro Dalcin pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 1393f45c9589SStefano Zampini } 1394f45c9589SStefano Zampini if (dim == 3) { 1395f45c9589SStefano Zampini /* Tetrahedra are inverted */ 1396dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 1397f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 1398f45c9589SStefano Zampini pcone[0] = pcone[1]; 1399f45c9589SStefano Zampini pcone[1] = tmp; 1400f45c9589SStefano Zampini } 1401f45c9589SStefano Zampini /* Hexahedra are inverted */ 1402dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 1403f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 1404f45c9589SStefano Zampini pcone[1] = pcone[3]; 1405f45c9589SStefano Zampini pcone[3] = tmp; 1406f45c9589SStefano Zampini } 1407dea2a3c7SStefano Zampini /* Prisms are inverted */ 1408dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 1409dea2a3c7SStefano Zampini PetscInt tmp; 1410dea2a3c7SStefano Zampini 1411dea2a3c7SStefano Zampini tmp = pcone[1]; 1412dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1413dea2a3c7SStefano Zampini pcone[2] = tmp; 1414dea2a3c7SStefano Zampini tmp = pcone[4]; 1415dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1416dea2a3c7SStefano Zampini pcone[5] = tmp; 1417f45c9589SStefano Zampini } 1418dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1419dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1420dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 1421dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1422dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1423dea2a3c7SStefano Zampini pcone[3] = tmp; 1424dea2a3c7SStefano Zampini } 1425dea2a3c7SStefano Zampini } 1426dea2a3c7SStefano Zampini ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 1427f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 142811c56cb0SLisandro Dalcin v = pcone[corner]; 1429dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 143011c56cb0SLisandro Dalcin coords[off++] = (PetscReal) coordsIn[v*3+d]; 1431f45c9589SStefano Zampini } 1432f45c9589SStefano Zampini } 14336fbe17bfSStefano Zampini } 1434f45c9589SStefano Zampini cell++; 1435f45c9589SStefano Zampini } 1436f45c9589SStefano Zampini } 1437f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 1438f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 1439dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 144011c56cb0SLisandro Dalcin coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 1441f45c9589SStefano Zampini } 1442f45c9589SStefano Zampini } 1443f45c9589SStefano Zampini } else { 1444331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 14451d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 144611c56cb0SLisandro Dalcin coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 1447331830f3SMatthew G. Knepley } 1448331830f3SMatthew G. Knepley } 1449331830f3SMatthew G. Knepley } 1450331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1451331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 145211c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 145311c56cb0SLisandro Dalcin 145411c56cb0SLisandro Dalcin ierr = PetscFree(coordsIn);CHKERRQ(ierr); 145511c56cb0SLisandro Dalcin ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 1456dea2a3c7SStefano Zampini ierr = PetscFree(hybridMap);CHKERRQ(ierr); 1457d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 1458fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 14596fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 14606fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 146111c56cb0SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 146211c56cb0SLisandro Dalcin 14633b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 1464331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1465331830f3SMatthew G. Knepley } 1466