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 5de78e4feSLisandro Dalcin typedef struct { 6698a79b9SLisandro Dalcin PetscViewer viewer; 7698a79b9SLisandro Dalcin int fileFormat; 8698a79b9SLisandro Dalcin int dataSize; 9698a79b9SLisandro Dalcin PetscBool binary; 10698a79b9SLisandro Dalcin PetscBool byteSwap; 11698a79b9SLisandro Dalcin size_t wlen; 12698a79b9SLisandro Dalcin void *wbuf; 13698a79b9SLisandro Dalcin size_t slen; 14698a79b9SLisandro Dalcin void *sbuf; 15698a79b9SLisandro Dalcin } GmshFile; 16de78e4feSLisandro Dalcin 17698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf) 18de78e4feSLisandro Dalcin { 19de78e4feSLisandro Dalcin size_t size = count * eltsize; 2011c56cb0SLisandro Dalcin PetscErrorCode ierr; 2111c56cb0SLisandro Dalcin 2211c56cb0SLisandro Dalcin PetscFunctionBegin; 23698a79b9SLisandro Dalcin if (gmsh->wlen < size) { 24698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 25698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->wbuf);CHKERRQ(ierr); 26698a79b9SLisandro Dalcin gmsh->wlen = size; 27de78e4feSLisandro Dalcin } 28698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->wbuf : NULL; 29de78e4feSLisandro Dalcin PetscFunctionReturn(0); 30de78e4feSLisandro Dalcin } 31de78e4feSLisandro Dalcin 32698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf) 33698a79b9SLisandro Dalcin { 34698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 35698a79b9SLisandro Dalcin size_t size = count * dataSize; 36698a79b9SLisandro Dalcin PetscErrorCode ierr; 37698a79b9SLisandro Dalcin 38698a79b9SLisandro Dalcin PetscFunctionBegin; 39698a79b9SLisandro Dalcin if (gmsh->slen < size) { 40698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 41698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->sbuf);CHKERRQ(ierr); 42698a79b9SLisandro Dalcin gmsh->slen = size; 43698a79b9SLisandro Dalcin } 44698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->sbuf : NULL; 45698a79b9SLisandro Dalcin PetscFunctionReturn(0); 46698a79b9SLisandro Dalcin } 47698a79b9SLisandro Dalcin 48698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype) 49de78e4feSLisandro Dalcin { 50de78e4feSLisandro Dalcin PetscErrorCode ierr; 51de78e4feSLisandro Dalcin PetscFunctionBegin; 52698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype);CHKERRQ(ierr); 53698a79b9SLisandro Dalcin if (gmsh->byteSwap) {ierr = PetscByteSwap(buf, dtype, count);CHKERRQ(ierr);} 54698a79b9SLisandro Dalcin PetscFunctionReturn(0); 55698a79b9SLisandro Dalcin } 56698a79b9SLisandro Dalcin 57698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count) 58698a79b9SLisandro Dalcin { 59698a79b9SLisandro Dalcin PetscErrorCode ierr; 60698a79b9SLisandro Dalcin PetscFunctionBegin; 61698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING);CHKERRQ(ierr); 62698a79b9SLisandro Dalcin PetscFunctionReturn(0); 63698a79b9SLisandro Dalcin } 64698a79b9SLisandro Dalcin 65d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match) 66d6689ca9SLisandro Dalcin { 67d6689ca9SLisandro Dalcin PetscErrorCode ierr; 68d6689ca9SLisandro Dalcin PetscFunctionBegin; 69d6689ca9SLisandro Dalcin ierr = PetscStrcmp(line, Section, match);CHKERRQ(ierr); 70d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 71d6689ca9SLisandro Dalcin } 72d6689ca9SLisandro Dalcin 73d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN]) 74d6689ca9SLisandro Dalcin { 75d6689ca9SLisandro Dalcin PetscBool match; 76d6689ca9SLisandro Dalcin PetscErrorCode ierr; 77d6689ca9SLisandro Dalcin 78d6689ca9SLisandro Dalcin PetscFunctionBegin; 79d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, Section, line, &match);CHKERRQ(ierr); 80d6689ca9SLisandro Dalcin if (!match) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file, expecting %s",Section); 81d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 82d6689ca9SLisandro Dalcin } 83d6689ca9SLisandro Dalcin 84d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN]) 85d6689ca9SLisandro Dalcin { 86d6689ca9SLisandro Dalcin PetscBool match; 87d6689ca9SLisandro Dalcin PetscErrorCode ierr; 88d6689ca9SLisandro Dalcin 89d6689ca9SLisandro Dalcin PetscFunctionBegin; 90d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 91d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 92d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$Comments", line, &match);CHKERRQ(ierr); 93d6689ca9SLisandro Dalcin if (!match) break; 94d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 95d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 96d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$EndComments", line, &match);CHKERRQ(ierr); 97d6689ca9SLisandro Dalcin if (match) break; 98d6689ca9SLisandro Dalcin } 99d6689ca9SLisandro Dalcin } 100d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 101d6689ca9SLisandro Dalcin } 102d6689ca9SLisandro Dalcin 103d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN]) 104d6689ca9SLisandro Dalcin { 105d6689ca9SLisandro Dalcin PetscErrorCode ierr; 106d6689ca9SLisandro Dalcin PetscFunctionBegin; 107d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 108d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, EndSection, line);CHKERRQ(ierr); 109d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 110d6689ca9SLisandro Dalcin } 111d6689ca9SLisandro Dalcin 112698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count) 113698a79b9SLisandro Dalcin { 114698a79b9SLisandro Dalcin PetscInt i; 115698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 116698a79b9SLisandro Dalcin PetscErrorCode ierr; 117698a79b9SLisandro Dalcin 118698a79b9SLisandro Dalcin PetscFunctionBegin; 119698a79b9SLisandro Dalcin if (dataSize == sizeof(PetscInt)) { 120698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_INT);CHKERRQ(ierr); 121698a79b9SLisandro Dalcin } else if (dataSize == sizeof(int)) { 122698a79b9SLisandro Dalcin int *ibuf = NULL; 12389d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 124698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); 125698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 126698a79b9SLisandro Dalcin } else if (dataSize == sizeof(long)) { 127698a79b9SLisandro Dalcin long *ibuf = NULL; 12889d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 129698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_LONG);CHKERRQ(ierr); 130698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 131698a79b9SLisandro Dalcin } else if (dataSize == sizeof(PetscInt64)) { 132698a79b9SLisandro Dalcin PetscInt64 *ibuf = NULL; 13389d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 134698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_INT64);CHKERRQ(ierr); 135698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 136698a79b9SLisandro Dalcin } 137698a79b9SLisandro Dalcin PetscFunctionReturn(0); 138698a79b9SLisandro Dalcin } 139698a79b9SLisandro Dalcin 140698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count) 141698a79b9SLisandro Dalcin { 142698a79b9SLisandro Dalcin PetscErrorCode ierr; 143698a79b9SLisandro Dalcin PetscFunctionBegin; 144698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_ENUM);CHKERRQ(ierr); 145698a79b9SLisandro Dalcin PetscFunctionReturn(0); 146698a79b9SLisandro Dalcin } 147698a79b9SLisandro Dalcin 148698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count) 149698a79b9SLisandro Dalcin { 150698a79b9SLisandro Dalcin PetscErrorCode ierr; 151698a79b9SLisandro Dalcin PetscFunctionBegin; 152698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_DOUBLE);CHKERRQ(ierr); 153de78e4feSLisandro Dalcin PetscFunctionReturn(0); 154de78e4feSLisandro Dalcin } 155de78e4feSLisandro Dalcin 156de78e4feSLisandro Dalcin typedef struct { 157de78e4feSLisandro Dalcin PetscInt id; /* Entity number */ 158698a79b9SLisandro Dalcin PetscInt dim; /* Entity dimension */ 159de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 160de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 161de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 162de78e4feSLisandro Dalcin } GmshEntity; 163de78e4feSLisandro Dalcin 164de78e4feSLisandro Dalcin typedef struct { 165730356f6SLisandro Dalcin GmshEntity *entity[4]; 166730356f6SLisandro Dalcin PetscHMapI entityMap[4]; 167730356f6SLisandro Dalcin } GmshEntities; 168730356f6SLisandro Dalcin 169698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities) 170730356f6SLisandro Dalcin { 171730356f6SLisandro Dalcin PetscInt dim; 172730356f6SLisandro Dalcin PetscErrorCode ierr; 173730356f6SLisandro Dalcin 174730356f6SLisandro Dalcin PetscFunctionBegin; 175730356f6SLisandro Dalcin ierr = PetscNew(entities);CHKERRQ(ierr); 176730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 177730356f6SLisandro Dalcin ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr); 178730356f6SLisandro Dalcin ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 179730356f6SLisandro Dalcin } 180730356f6SLisandro Dalcin PetscFunctionReturn(0); 181730356f6SLisandro Dalcin } 182730356f6SLisandro Dalcin 183730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity) 184730356f6SLisandro Dalcin { 185730356f6SLisandro Dalcin PetscErrorCode ierr; 186730356f6SLisandro Dalcin PetscFunctionBegin; 187730356f6SLisandro Dalcin ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr); 188730356f6SLisandro Dalcin entities->entity[dim][index].dim = dim; 189730356f6SLisandro Dalcin entities->entity[dim][index].id = eid; 190730356f6SLisandro Dalcin if (entity) *entity = &entities->entity[dim][index]; 191730356f6SLisandro Dalcin PetscFunctionReturn(0); 192730356f6SLisandro Dalcin } 193730356f6SLisandro Dalcin 194730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity) 195730356f6SLisandro Dalcin { 196730356f6SLisandro Dalcin PetscInt index; 197730356f6SLisandro Dalcin PetscErrorCode ierr; 198730356f6SLisandro Dalcin 199730356f6SLisandro Dalcin PetscFunctionBegin; 200730356f6SLisandro Dalcin ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr); 201730356f6SLisandro Dalcin *entity = &entities->entity[dim][index]; 202730356f6SLisandro Dalcin PetscFunctionReturn(0); 203730356f6SLisandro Dalcin } 204730356f6SLisandro Dalcin 205730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) 206730356f6SLisandro Dalcin { 207730356f6SLisandro Dalcin PetscInt dim; 208730356f6SLisandro Dalcin PetscErrorCode ierr; 209730356f6SLisandro Dalcin 210730356f6SLisandro Dalcin PetscFunctionBegin; 211730356f6SLisandro Dalcin if (!*entities) PetscFunctionReturn(0); 212730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 213730356f6SLisandro Dalcin ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr); 214730356f6SLisandro Dalcin ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 215730356f6SLisandro Dalcin } 216730356f6SLisandro Dalcin ierr = PetscFree((*entities));CHKERRQ(ierr); 217730356f6SLisandro Dalcin PetscFunctionReturn(0); 218730356f6SLisandro Dalcin } 219730356f6SLisandro Dalcin 220730356f6SLisandro Dalcin typedef struct { 221698a79b9SLisandro Dalcin PetscInt id; /* Entity number */ 222de78e4feSLisandro Dalcin PetscInt dim; /* Entity dimension */ 223de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 224de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 225698a79b9SLisandro Dalcin PetscInt nodes[8]; /* Node array */ 226de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 227de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 228de78e4feSLisandro Dalcin } GmshElement; 229de78e4feSLisandro Dalcin 230698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v22(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 231de78e4feSLisandro Dalcin { 232698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 233698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 234de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 235698a79b9SLisandro Dalcin int v, num, nid, snum; 236698a79b9SLisandro Dalcin double *coordinates; 237de78e4feSLisandro Dalcin PetscErrorCode ierr; 238de78e4feSLisandro Dalcin 239de78e4feSLisandro Dalcin PetscFunctionBegin; 240de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 241698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 242de78e4feSLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 243698a79b9SLisandro Dalcin ierr = PetscMalloc1(num*3, &coordinates);CHKERRQ(ierr); 244698a79b9SLisandro Dalcin *numVertices = num; 245698a79b9SLisandro Dalcin *gmsh_nodes = coordinates; 246698a79b9SLisandro Dalcin for (v = 0; v < num; ++v) { 247698a79b9SLisandro Dalcin double *xyz = coordinates + v*3; 24811c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 24911c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 25011c56cb0SLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 25111c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 25211c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 25311c56cb0SLisandro Dalcin } 25411c56cb0SLisandro Dalcin PetscFunctionReturn(0); 25511c56cb0SLisandro Dalcin } 25611c56cb0SLisandro Dalcin 257de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 258de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 259de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 260de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 261698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v22(GmshFile* gmsh, PETSC_UNUSED int shift, PetscInt *numCells, GmshElement **gmsh_elems) 26211c56cb0SLisandro Dalcin { 263698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 264698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 265698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 266de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 267df032b82SMatthew G. Knepley GmshElement *elements; 268698a79b9SLisandro Dalcin int i, c, p, num, ibuf[1+4+512], snum; 26911c56cb0SLisandro Dalcin int cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 270df032b82SMatthew G. Knepley PetscErrorCode ierr; 271df032b82SMatthew G. Knepley 272df032b82SMatthew G. Knepley PetscFunctionBegin; 273de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 274698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 275de78e4feSLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 276698a79b9SLisandro Dalcin ierr = PetscMalloc1(num, &elements);CHKERRQ(ierr); 277698a79b9SLisandro Dalcin *numCells = num; 278698a79b9SLisandro Dalcin *gmsh_elems = elements; 279698a79b9SLisandro Dalcin for (c = 0; c < num;) { 28011c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 28111c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 282df032b82SMatthew G. Knepley if (binary) { 283df032b82SMatthew G. Knepley cellType = ibuf[0]; 284df032b82SMatthew G. Knepley numElem = ibuf[1]; 285df032b82SMatthew G. Knepley numTags = ibuf[2]; 286df032b82SMatthew G. Knepley } else { 287df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 288df032b82SMatthew G. Knepley cellType = ibuf[1]; 289df032b82SMatthew G. Knepley numTags = ibuf[2]; 290df032b82SMatthew G. Knepley numElem = 1; 291df032b82SMatthew G. Knepley } 292bf6ba3a3SMatthew G. Knepley /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 293bf6ba3a3SMatthew G. Knepley numNodesIgnore = 0; 294df032b82SMatthew G. Knepley switch (cellType) { 295df032b82SMatthew G. Knepley case 1: /* 2-node line */ 296df032b82SMatthew G. Knepley dim = 1; 297df032b82SMatthew G. Knepley numNodes = 2; 298df032b82SMatthew G. Knepley break; 299df032b82SMatthew G. Knepley case 2: /* 3-node triangle */ 300df032b82SMatthew G. Knepley dim = 2; 301df032b82SMatthew G. Knepley numNodes = 3; 302df032b82SMatthew G. Knepley break; 303df032b82SMatthew G. Knepley case 3: /* 4-node quadrangle */ 304df032b82SMatthew G. Knepley dim = 2; 305df032b82SMatthew G. Knepley numNodes = 4; 306df032b82SMatthew G. Knepley break; 307df032b82SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 308df032b82SMatthew G. Knepley dim = 3; 309df032b82SMatthew G. Knepley numNodes = 4; 310df032b82SMatthew G. Knepley break; 311df032b82SMatthew G. Knepley case 5: /* 8-node hexahedron */ 312df032b82SMatthew G. Knepley dim = 3; 313df032b82SMatthew G. Knepley numNodes = 8; 314df032b82SMatthew G. Knepley break; 315dea2a3c7SStefano Zampini case 6: /* 6-node wedge */ 316dea2a3c7SStefano Zampini dim = 3; 317dea2a3c7SStefano Zampini numNodes = 6; 318dea2a3c7SStefano Zampini break; 319bf6ba3a3SMatthew G. Knepley case 8: /* 3-node 2nd order line */ 320bf6ba3a3SMatthew G. Knepley dim = 1; 321bf6ba3a3SMatthew G. Knepley numNodes = 2; 322bf6ba3a3SMatthew G. Knepley numNodesIgnore = 1; 323bf6ba3a3SMatthew G. Knepley break; 324bf6ba3a3SMatthew G. Knepley case 9: /* 6-node 2nd order triangle */ 325bf6ba3a3SMatthew G. Knepley dim = 2; 326bf6ba3a3SMatthew G. Knepley numNodes = 3; 327bf6ba3a3SMatthew G. Knepley numNodesIgnore = 3; 328bf6ba3a3SMatthew G. Knepley break; 329a1e86c74SStefano Zampini case 10: /* 9-node 2nd order quadrangle */ 330a1e86c74SStefano Zampini dim = 2; 331a1e86c74SStefano Zampini numNodes = 4; 332a1e86c74SStefano Zampini numNodesIgnore = 5; 333a1e86c74SStefano Zampini break; 334a1e86c74SStefano Zampini case 11: /* 10-node 2nd order tetrahedron */ 335a1e86c74SStefano Zampini dim = 3; 336a1e86c74SStefano Zampini numNodes = 4; 337a1e86c74SStefano Zampini numNodesIgnore = 6; 338a1e86c74SStefano Zampini break; 339a1e86c74SStefano Zampini case 12: /* 27-node 2nd order hexhedron */ 340a1e86c74SStefano Zampini dim = 3; 341a1e86c74SStefano Zampini numNodes = 8; 342a1e86c74SStefano Zampini numNodesIgnore = 19; 343a1e86c74SStefano Zampini 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 14: /* 14-node 2nd order pyramid */ 355df032b82SMatthew G. Knepley default: 356df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 357df032b82SMatthew G. Knepley } 358df032b82SMatthew G. Knepley if (binary) { 35911c56cb0SLisandro Dalcin const int nint = 1 + numTags + numNodes + numNodesIgnore; 36011c56cb0SLisandro Dalcin /* Loop over element blocks */ 361df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 36211c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 36311c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 364df032b82SMatthew G. Knepley elements[c].dim = dim; 365df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 366df032b82SMatthew G. Knepley elements[c].numTags = numTags; 367df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 368dea2a3c7SStefano Zampini elements[c].cellType = cellType; 369df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 370df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 371df032b82SMatthew G. Knepley } 372df032b82SMatthew G. Knepley } else { 373698a79b9SLisandro Dalcin const int nint = numTags + numNodes + numNodesIgnore; 374df032b82SMatthew G. Knepley elements[c].dim = dim; 375df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 376df032b82SMatthew G. Knepley elements[c].numTags = numTags; 377dea2a3c7SStefano Zampini elements[c].cellType = cellType; 378698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 379698a79b9SLisandro Dalcin for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[p]; 380698a79b9SLisandro Dalcin for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[numTags + p]; 381df032b82SMatthew G. Knepley c++; 382df032b82SMatthew G. Knepley } 383df032b82SMatthew G. Knepley } 384df032b82SMatthew G. Knepley PetscFunctionReturn(0); 385df032b82SMatthew G. Knepley } 3867d282ae0SMichael Lange 387de78e4feSLisandro Dalcin /* 388de78e4feSLisandro Dalcin $Entities 389de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 390de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 391de78e4feSLisandro Dalcin // points 392de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 393de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 394de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 395de78e4feSLisandro Dalcin ... 396de78e4feSLisandro Dalcin // curves 397de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 398de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 399de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 400de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 401de78e4feSLisandro Dalcin ... 402de78e4feSLisandro Dalcin // surfaces 403de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 404de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 405de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 406de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 407de78e4feSLisandro Dalcin ... 408de78e4feSLisandro Dalcin // volumes 409de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 410de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 411de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 412de78e4feSLisandro Dalcin numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 413de78e4feSLisandro Dalcin ... 414de78e4feSLisandro Dalcin $EndEntities 415de78e4feSLisandro Dalcin */ 416698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v40(GmshFile *gmsh, GmshEntities **entities) 417de78e4feSLisandro Dalcin { 418698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 419698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 420698a79b9SLisandro Dalcin long index, num, lbuf[4]; 421730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 422698a79b9SLisandro Dalcin PetscInt count[4], i; 423a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 424de78e4feSLisandro Dalcin PetscErrorCode ierr; 425de78e4feSLisandro Dalcin 426de78e4feSLisandro Dalcin PetscFunctionBegin; 427698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 428698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);} 429698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 430730356f6SLisandro Dalcin ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 431de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 432730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 433730356f6SLisandro Dalcin ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 434730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);} 435730356f6SLisandro Dalcin ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 436de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 437de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 438de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 439de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 440698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 441de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 442730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 443de78e4feSLisandro Dalcin entity->numTags = numTags = (int) PetscMin(num, 4); 444de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 445de78e4feSLisandro Dalcin if (dim == 0) continue; 446de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 447de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 448698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 449de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 450730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 451de78e4feSLisandro Dalcin } 452de78e4feSLisandro Dalcin } 453de78e4feSLisandro Dalcin PetscFunctionReturn(0); 454de78e4feSLisandro Dalcin } 455de78e4feSLisandro Dalcin 456de78e4feSLisandro Dalcin /* 457de78e4feSLisandro Dalcin $Nodes 458de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 459de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 460de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 461de78e4feSLisandro Dalcin ... 462de78e4feSLisandro Dalcin ... 463de78e4feSLisandro Dalcin $EndNodes 464de78e4feSLisandro Dalcin */ 465698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v40(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 466de78e4feSLisandro Dalcin { 467698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 468698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 469de78e4feSLisandro Dalcin long block, node, v, numEntityBlocks, numTotalNodes, numNodes; 470de78e4feSLisandro Dalcin int info[3], nid; 471de78e4feSLisandro Dalcin double *coordinates; 472de78e4feSLisandro Dalcin char *cbuf; 473de78e4feSLisandro Dalcin PetscErrorCode ierr; 474de78e4feSLisandro Dalcin 475de78e4feSLisandro Dalcin PetscFunctionBegin; 476de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 477de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 478de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 479de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 480de78e4feSLisandro Dalcin ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr); 481698a79b9SLisandro Dalcin *numVertices = numTotalNodes; 482de78e4feSLisandro Dalcin *gmsh_nodes = coordinates; 483de78e4feSLisandro Dalcin for (v = 0, block = 0; block < numEntityBlocks; ++block) { 484de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 485de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 486de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 487698a79b9SLisandro Dalcin if (gmsh->binary) { 488de78e4feSLisandro Dalcin int nbytes = sizeof(int) + 3*sizeof(double); 489698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 490de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 491de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 492de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 493de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 49430815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr);} 49530815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 496de78e4feSLisandro Dalcin ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 497de78e4feSLisandro Dalcin ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 498de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 499de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 500de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 501de78e4feSLisandro Dalcin } 502de78e4feSLisandro Dalcin } else { 503de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 504de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 505de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 506de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 507de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 508de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 509de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 510de78e4feSLisandro Dalcin } 511de78e4feSLisandro Dalcin } 512de78e4feSLisandro Dalcin } 513de78e4feSLisandro Dalcin PetscFunctionReturn(0); 514de78e4feSLisandro Dalcin } 515de78e4feSLisandro Dalcin 516de78e4feSLisandro Dalcin /* 517de78e4feSLisandro Dalcin $Elements 518de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 519de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 520de78e4feSLisandro Dalcin tag(int) numVert[...](int) 521de78e4feSLisandro Dalcin ... 522de78e4feSLisandro Dalcin ... 523de78e4feSLisandro Dalcin $EndElements 524de78e4feSLisandro Dalcin */ 525698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v40(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems) 526de78e4feSLisandro Dalcin { 527698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 528698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 529de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 530de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 531de78e4feSLisandro Dalcin int eid, dim, numTags, *tags, cellType, numNodes; 532a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 533de78e4feSLisandro Dalcin GmshElement *elements; 534de78e4feSLisandro Dalcin PetscErrorCode ierr; 535de78e4feSLisandro Dalcin 536de78e4feSLisandro Dalcin PetscFunctionBegin; 537de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 538de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 539de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 540de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 541de78e4feSLisandro Dalcin ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr); 542698a79b9SLisandro Dalcin *numCells = numTotalElements; 543de78e4feSLisandro Dalcin *gmsh_elems = elements; 544de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 545de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 546de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 547de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 548730356f6SLisandro Dalcin ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 549730356f6SLisandro Dalcin numTags = entity->numTags; 550730356f6SLisandro Dalcin tags = entity->tags; 551de78e4feSLisandro Dalcin switch (cellType) { 552de78e4feSLisandro Dalcin case 1: /* 2-node line */ 553de78e4feSLisandro Dalcin numNodes = 2; 554de78e4feSLisandro Dalcin break; 555de78e4feSLisandro Dalcin case 2: /* 3-node triangle */ 556698a79b9SLisandro Dalcin numNodes = 3; 557698a79b9SLisandro Dalcin break; 558de78e4feSLisandro Dalcin case 3: /* 4-node quadrangle */ 559de78e4feSLisandro Dalcin numNodes = 4; 560de78e4feSLisandro Dalcin break; 561de78e4feSLisandro Dalcin case 4: /* 4-node tetrahedron */ 562de78e4feSLisandro Dalcin numNodes = 4; 563de78e4feSLisandro Dalcin break; 564de78e4feSLisandro Dalcin case 5: /* 8-node hexahedron */ 565de78e4feSLisandro Dalcin numNodes = 8; 566de78e4feSLisandro Dalcin break; 567de78e4feSLisandro Dalcin case 6: /* 6-node wedge */ 568de78e4feSLisandro Dalcin numNodes = 6; 569de78e4feSLisandro Dalcin break; 570de78e4feSLisandro Dalcin case 15: /* 1-node vertex */ 571de78e4feSLisandro Dalcin numNodes = 1; 572de78e4feSLisandro Dalcin break; 573de78e4feSLisandro Dalcin default: 574de78e4feSLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 575de78e4feSLisandro Dalcin } 576de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 577de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 578698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 579de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 580de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 581de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 582de78e4feSLisandro Dalcin int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 583de78e4feSLisandro Dalcin GmshElement *element = elements + c; 584de78e4feSLisandro Dalcin element->dim = dim; 585de78e4feSLisandro Dalcin element->cellType = cellType; 586de78e4feSLisandro Dalcin element->numNodes = numNodes; 587de78e4feSLisandro Dalcin element->numTags = numTags; 588de78e4feSLisandro Dalcin element->id = id[0]; 589de78e4feSLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 590de78e4feSLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 591de78e4feSLisandro Dalcin } 592de78e4feSLisandro Dalcin } 593698a79b9SLisandro Dalcin PetscFunctionReturn(0); 594698a79b9SLisandro Dalcin } 595698a79b9SLisandro Dalcin 596698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v40(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt) 597698a79b9SLisandro Dalcin { 598698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 599698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 600698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 601698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 602698a79b9SLisandro Dalcin int numPeriodic, snum, i; 603698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 604698a79b9SLisandro Dalcin PetscErrorCode ierr; 605698a79b9SLisandro Dalcin 606698a79b9SLisandro Dalcin PetscFunctionBegin; 607698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 608698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 609698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 610698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 611698a79b9SLisandro Dalcin } else { 612698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 613698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 614698a79b9SLisandro Dalcin } 615698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 616698a79b9SLisandro Dalcin int ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode; 617698a79b9SLisandro Dalcin long j, nNodes; 618698a79b9SLisandro Dalcin double affine[16]; 619698a79b9SLisandro Dalcin 620698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 621698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 622698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag); 623698a79b9SLisandro Dalcin if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 624698a79b9SLisandro Dalcin } else { 625698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 626698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 627698a79b9SLisandro Dalcin slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2]; 628698a79b9SLisandro Dalcin } 629698a79b9SLisandro Dalcin (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */ 630698a79b9SLisandro Dalcin 631698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 632698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 633698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 634*4c500f23SPierre Jolivet if (snum != 1) { /* discard transformation and try again */ 635698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 636698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 637698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 638698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 639698a79b9SLisandro Dalcin } 640698a79b9SLisandro Dalcin } else { 641698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 642698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 643*4c500f23SPierre Jolivet if (nNodes == -1) { /* discard transformation and try again */ 644698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 645698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 646698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 647698a79b9SLisandro Dalcin } 648698a79b9SLisandro Dalcin } 649698a79b9SLisandro Dalcin 650698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 651698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 652698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 653698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 654698a79b9SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 655698a79b9SLisandro Dalcin } else { 656698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 657698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 658698a79b9SLisandro Dalcin slaveNode = ibuf[0]; masterNode = ibuf[1]; 659698a79b9SLisandro Dalcin } 660698a79b9SLisandro Dalcin slaveMap[slaveNode - shift] = masterNode - shift; 661698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, slaveNode - shift);CHKERRQ(ierr); 662698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, masterNode - shift);CHKERRQ(ierr); 663698a79b9SLisandro Dalcin } 664698a79b9SLisandro Dalcin } 665698a79b9SLisandro Dalcin PetscFunctionReturn(0); 666698a79b9SLisandro Dalcin } 667698a79b9SLisandro Dalcin 668698a79b9SLisandro Dalcin /* 669698a79b9SLisandro Dalcin $Entities 670698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 671698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 672698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 673698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 674698a79b9SLisandro Dalcin ... 675698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 676698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 677698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 678698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 679698a79b9SLisandro Dalcin ... 680698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 681698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 682698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 683698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 684698a79b9SLisandro Dalcin ... 685698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 686698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 687698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 688698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 689698a79b9SLisandro Dalcin ... 690698a79b9SLisandro Dalcin $EndEntities 691698a79b9SLisandro Dalcin */ 692698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v41(GmshFile *gmsh, GmshEntities **entities) 693698a79b9SLisandro Dalcin { 694698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 695698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 696698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 697698a79b9SLisandro Dalcin PetscErrorCode ierr; 698698a79b9SLisandro Dalcin 699698a79b9SLisandro Dalcin PetscFunctionBegin; 700698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr); 701698a79b9SLisandro Dalcin ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 702698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 703698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 704698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr); 705698a79b9SLisandro Dalcin ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 706698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr); 707698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 708698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 709698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 710698a79b9SLisandro Dalcin entity->numTags = PetscMin(numTags, 4); 711698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 712698a79b9SLisandro Dalcin if (dim == 0) continue; 713698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 714698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 715698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 716698a79b9SLisandro Dalcin } 717698a79b9SLisandro Dalcin } 718698a79b9SLisandro Dalcin PetscFunctionReturn(0); 719698a79b9SLisandro Dalcin } 720698a79b9SLisandro Dalcin 721698a79b9SLisandro Dalcin /* 722698a79b9SLisandro Dalcin $Nodes 723698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 724698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 725698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 726698a79b9SLisandro Dalcin nodeTag(size_t) 727698a79b9SLisandro Dalcin ... 728698a79b9SLisandro Dalcin x(double) y(double) z(double) 729698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 730698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 731698a79b9SLisandro Dalcin ... 732698a79b9SLisandro Dalcin ... 733698a79b9SLisandro Dalcin $EndNodes 734698a79b9SLisandro Dalcin */ 735698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v41(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 736698a79b9SLisandro Dalcin { 737698a79b9SLisandro Dalcin int info[3]; 738698a79b9SLisandro Dalcin PetscInt sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, *nodeTag = NULL, block, node, i; 739698a79b9SLisandro Dalcin double *coordinates; 740698a79b9SLisandro Dalcin PetscErrorCode ierr; 741698a79b9SLisandro Dalcin 742698a79b9SLisandro Dalcin PetscFunctionBegin; 743698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 744698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numNodes = sizes[1]; 745698a79b9SLisandro Dalcin ierr = PetscMalloc1(numNodes*3, &coordinates);CHKERRQ(ierr); 746698a79b9SLisandro Dalcin *numVertices = numNodes; 747698a79b9SLisandro Dalcin *gmsh_nodes = coordinates; 748698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 749698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 750698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr); 751698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodesBlock, sizeof(PetscInt), &nodeTag);CHKERRQ(ierr); 752698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTag, numNodesBlock);CHKERRQ(ierr); 753698a79b9SLisandro 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); 754698a79b9SLisandro Dalcin if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 755698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, coordinates+node*3, numNodesBlock*3);CHKERRQ(ierr); 756698a79b9SLisandro Dalcin } 757698a79b9SLisandro Dalcin PetscFunctionReturn(0); 758698a79b9SLisandro Dalcin } 759698a79b9SLisandro Dalcin 760698a79b9SLisandro Dalcin /* 761698a79b9SLisandro Dalcin $Elements 762698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 763698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 764698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 765698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 766698a79b9SLisandro Dalcin ... 767698a79b9SLisandro Dalcin ... 768698a79b9SLisandro Dalcin $EndElements 769698a79b9SLisandro Dalcin */ 770698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v41(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems) 771698a79b9SLisandro Dalcin { 772698a79b9SLisandro Dalcin int info[3], eid, dim, cellType, *tags; 773698a79b9SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numNodes, numTags, block, elem, c, p; 774698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 775698a79b9SLisandro Dalcin GmshElement *elements; 776698a79b9SLisandro Dalcin PetscErrorCode ierr; 777698a79b9SLisandro Dalcin 778698a79b9SLisandro Dalcin PetscFunctionBegin; 779698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 780698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numElements = sizes[1]; 781698a79b9SLisandro Dalcin ierr = PetscCalloc1(numElements, &elements);CHKERRQ(ierr); 782698a79b9SLisandro Dalcin *numCells = numElements; 783698a79b9SLisandro Dalcin *gmsh_elems = elements; 784698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 785698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 786698a79b9SLisandro Dalcin dim = info[0]; eid = info[1]; cellType = info[2]; 787698a79b9SLisandro Dalcin ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 788698a79b9SLisandro Dalcin numTags = entity->numTags; 789698a79b9SLisandro Dalcin tags = entity->tags; 790698a79b9SLisandro Dalcin switch (cellType) { 791698a79b9SLisandro Dalcin case 1: /* 2-node line */ 792698a79b9SLisandro Dalcin numNodes = 2; 793698a79b9SLisandro Dalcin break; 794698a79b9SLisandro Dalcin case 2: /* 3-node triangle */ 795698a79b9SLisandro Dalcin numNodes = 3; 796698a79b9SLisandro Dalcin break; 797698a79b9SLisandro Dalcin case 3: /* 4-node quadrangle */ 798698a79b9SLisandro Dalcin numNodes = 4; 799698a79b9SLisandro Dalcin break; 800698a79b9SLisandro Dalcin case 4: /* 4-node tetrahedron */ 801698a79b9SLisandro Dalcin numNodes = 4; 802698a79b9SLisandro Dalcin break; 803698a79b9SLisandro Dalcin case 5: /* 8-node hexahedron */ 804698a79b9SLisandro Dalcin numNodes = 8; 805698a79b9SLisandro Dalcin break; 806698a79b9SLisandro Dalcin case 6: /* 6-node wedge */ 807698a79b9SLisandro Dalcin numNodes = 6; 808698a79b9SLisandro Dalcin break; 809698a79b9SLisandro Dalcin case 15: /* 1-node vertex */ 810698a79b9SLisandro Dalcin numNodes = 1; 811698a79b9SLisandro Dalcin break; 812698a79b9SLisandro Dalcin default: 813698a79b9SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 814698a79b9SLisandro Dalcin } 815698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr); 816698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr); 817698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr); 818698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 819698a79b9SLisandro Dalcin GmshElement *element = elements + c; 820698a79b9SLisandro Dalcin PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 821698a79b9SLisandro Dalcin element->id = id[0]; 822698a79b9SLisandro Dalcin element->dim = dim; 823698a79b9SLisandro Dalcin element->cellType = cellType; 824698a79b9SLisandro Dalcin element->numNodes = numNodes; 825698a79b9SLisandro Dalcin element->numTags = numTags; 826698a79b9SLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 827698a79b9SLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 828698a79b9SLisandro Dalcin } 829698a79b9SLisandro Dalcin } 830698a79b9SLisandro Dalcin PetscFunctionReturn(0); 831698a79b9SLisandro Dalcin } 832698a79b9SLisandro Dalcin 833698a79b9SLisandro Dalcin /* 834698a79b9SLisandro Dalcin $Periodic 835698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 836698a79b9SLisandro Dalcin entityDim(int) entityTag(int) entityTagMaster(int) 837698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 838698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 839698a79b9SLisandro Dalcin nodeTag(size_t) nodeTagMaster(size_t) 840698a79b9SLisandro Dalcin ... 841698a79b9SLisandro Dalcin ... 842698a79b9SLisandro Dalcin $EndPeriodic 843698a79b9SLisandro Dalcin */ 844698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v41(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt) 845698a79b9SLisandro Dalcin { 846698a79b9SLisandro Dalcin int info[3]; 847698a79b9SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 848698a79b9SLisandro Dalcin double dbuf[16]; 849698a79b9SLisandro Dalcin PetscErrorCode ierr; 850698a79b9SLisandro Dalcin 851698a79b9SLisandro Dalcin PetscFunctionBegin; 852698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr); 853698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 854698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 855698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr); 856698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr); 857698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr); 858698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr); 859698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr); 860698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 861698a79b9SLisandro Dalcin PetscInt slaveNode = nodeTags[node*2+0] - shift; 862698a79b9SLisandro Dalcin PetscInt masterNode = nodeTags[node*2+1] - shift; 863698a79b9SLisandro Dalcin slaveMap[slaveNode] = masterNode; 864698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, slaveNode);CHKERRQ(ierr); 865698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, masterNode);CHKERRQ(ierr); 866698a79b9SLisandro Dalcin } 867698a79b9SLisandro Dalcin } 868698a79b9SLisandro Dalcin PetscFunctionReturn(0); 869698a79b9SLisandro Dalcin } 870698a79b9SLisandro Dalcin 871d6689ca9SLisandro Dalcin /* 872d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2 873d6689ca9SLisandro Dalcin version(ASCII double; currently 4.1) 874d6689ca9SLisandro Dalcin file-type(ASCII int; 0 for ASCII mode, 1 for binary mode) 875d6689ca9SLisandro Dalcin data-size(ASCII int; sizeof(size_t)) 876d6689ca9SLisandro Dalcin < int with value one; only in binary mode, to detect endianness > 877d6689ca9SLisandro Dalcin $EndMeshFormat 878d6689ca9SLisandro Dalcin */ 879698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadMeshFormat(GmshFile *gmsh) 880698a79b9SLisandro Dalcin { 881698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 882698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 883698a79b9SLisandro Dalcin float version; 884698a79b9SLisandro Dalcin PetscErrorCode ierr; 885698a79b9SLisandro Dalcin 886698a79b9SLisandro Dalcin PetscFunctionBegin; 887698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr); 888698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 889698a79b9SLisandro Dalcin if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 890698a79b9SLisandro 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); 891698a79b9SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version); 892698a79b9SLisandro 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); 893698a79b9SLisandro Dalcin if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is binary but Gmsh file is ASCII"); 894698a79b9SLisandro Dalcin if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is ASCII but Gmsh file is binary"); 895698a79b9SLisandro Dalcin fileFormat = (int)(version*10); /* XXX Should use (int)roundf(version*10) ? */ 896698a79b9SLisandro 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); 897698a79b9SLisandro 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); 898698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 899698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 900698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 901698a79b9SLisandro Dalcin if (gmsh->binary) { 902698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr); 903698a79b9SLisandro Dalcin if (checkEndian != 1) { 904698a79b9SLisandro Dalcin ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr); 905698a79b9SLisandro Dalcin if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to detect endianness in Gmsh file header: %s", line); 906698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 907698a79b9SLisandro Dalcin } 908698a79b9SLisandro Dalcin } 909698a79b9SLisandro Dalcin PetscFunctionReturn(0); 910698a79b9SLisandro Dalcin } 911698a79b9SLisandro Dalcin 9121f643af3SLisandro Dalcin /* 9131f643af3SLisandro Dalcin PhysicalNames 9141f643af3SLisandro Dalcin numPhysicalNames(ASCII int) 9151f643af3SLisandro Dalcin dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max) 9161f643af3SLisandro Dalcin ... 9171f643af3SLisandro Dalcin $EndPhysicalNames 9181f643af3SLisandro Dalcin */ 919698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPhysicalNames(GmshFile *gmsh) 920698a79b9SLisandro Dalcin { 9211f643af3SLisandro Dalcin char line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q; 9221f643af3SLisandro Dalcin int snum, numRegions, region, dim, tag; 923698a79b9SLisandro Dalcin PetscErrorCode ierr; 924698a79b9SLisandro Dalcin 925698a79b9SLisandro Dalcin PetscFunctionBegin; 926698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 927698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numRegions); 928698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 929698a79b9SLisandro Dalcin for (region = 0; region < numRegions; ++region) { 9301f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 9311f643af3SLisandro Dalcin snum = sscanf(line, "%d %d", &dim, &tag); 9321f643af3SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 9331f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, -(PetscInt)sizeof(line));CHKERRQ(ierr); 9341f643af3SLisandro Dalcin ierr = PetscStrchr(line, '"', &p);CHKERRQ(ierr); 9351f643af3SLisandro Dalcin if (!p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 9361f643af3SLisandro Dalcin ierr = PetscStrrchr(line, '"', &q);CHKERRQ(ierr); 9371f643af3SLisandro Dalcin if (q == p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 9381f643af3SLisandro Dalcin ierr = PetscStrncpy(name, p+1, (size_t)(q-p-1));CHKERRQ(ierr); 939698a79b9SLisandro Dalcin } 940de78e4feSLisandro Dalcin PetscFunctionReturn(0); 941de78e4feSLisandro Dalcin } 942de78e4feSLisandro Dalcin 943d6689ca9SLisandro Dalcin /*@C 944d6689ca9SLisandro Dalcin DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 945d6689ca9SLisandro Dalcin 946d6689ca9SLisandro Dalcin + comm - The MPI communicator 947d6689ca9SLisandro Dalcin . filename - Name of the Gmsh file 948d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh 949d6689ca9SLisandro Dalcin 950d6689ca9SLisandro Dalcin Output Parameter: 951d6689ca9SLisandro Dalcin . dm - The DM object representing the mesh 952d6689ca9SLisandro Dalcin 953d6689ca9SLisandro Dalcin Level: beginner 954d6689ca9SLisandro Dalcin 955d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 956d6689ca9SLisandro Dalcin @*/ 957d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 958d6689ca9SLisandro Dalcin { 959d6689ca9SLisandro Dalcin PetscViewer viewer; 960d6689ca9SLisandro Dalcin PetscMPIInt rank; 961d6689ca9SLisandro Dalcin int fileType; 962d6689ca9SLisandro Dalcin PetscViewerType vtype; 963d6689ca9SLisandro Dalcin PetscErrorCode ierr; 964d6689ca9SLisandro Dalcin 965d6689ca9SLisandro Dalcin PetscFunctionBegin; 966d6689ca9SLisandro Dalcin ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 967d6689ca9SLisandro Dalcin 968d6689ca9SLisandro Dalcin /* Determine Gmsh file type (ASCII or binary) from file header */ 969d6689ca9SLisandro Dalcin if (!rank) { 970d6689ca9SLisandro Dalcin GmshFile gmsh_, *gmsh = &gmsh_; 971d6689ca9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 972d6689ca9SLisandro Dalcin int snum; 973d6689ca9SLisandro Dalcin float version; 974d6689ca9SLisandro Dalcin 975580bdb30SBarry Smith ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr); 976d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer);CHKERRQ(ierr); 977d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 978d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ);CHKERRQ(ierr); 979d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(gmsh->viewer, filename);CHKERRQ(ierr); 980d6689ca9SLisandro Dalcin /* Read only the first two lines of the Gmsh file */ 981d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 982d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 983d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 984d6689ca9SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 985d6689ca9SLisandro Dalcin if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 986d6689ca9SLisandro 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); 987d6689ca9SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version); 988d6689ca9SLisandro 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); 989d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&gmsh->viewer);CHKERRQ(ierr); 990d6689ca9SLisandro Dalcin } 991d6689ca9SLisandro Dalcin ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 992d6689ca9SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 993d6689ca9SLisandro Dalcin 994d6689ca9SLisandro Dalcin /* Create appropriate viewer and build plex */ 995d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 996d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 997d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 998d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 999d6689ca9SLisandro Dalcin ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 1000d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1001d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 1002d6689ca9SLisandro Dalcin } 1003d6689ca9SLisandro Dalcin 1004331830f3SMatthew G. Knepley /*@ 10057d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 1006331830f3SMatthew G. Knepley 1007d083f849SBarry Smith Collective 1008331830f3SMatthew G. Knepley 1009331830f3SMatthew G. Knepley Input Parameters: 1010331830f3SMatthew G. Knepley + comm - The MPI communicator 1011331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 1012331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 1013331830f3SMatthew G. Knepley 1014331830f3SMatthew G. Knepley Output Parameter: 1015331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 1016331830f3SMatthew G. Knepley 1017698a79b9SLisandro Dalcin Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 1018331830f3SMatthew G. Knepley 1019331830f3SMatthew G. Knepley Level: beginner 1020331830f3SMatthew G. Knepley 1021331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 1022331830f3SMatthew G. Knepley @*/ 1023331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1024331830f3SMatthew G. Knepley { 102511c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 102611c56cb0SLisandro Dalcin double *coordsIn = NULL; 1027730356f6SLisandro Dalcin GmshEntities *entities = NULL; 10283c67d5adSSatish Balay GmshElement *gmsh_elem = NULL; 1029331830f3SMatthew G. Knepley PetscSection coordSection; 1030331830f3SMatthew G. Knepley Vec coordinates; 10316fbe17bfSStefano Zampini PetscBT periodicV = NULL, periodicC = NULL; 103284572febSToby Isaac PetscScalar *coords; 1033698a79b9SLisandro Dalcin PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 1034698a79b9SLisandro Dalcin PetscInt numVertices = 0, numCells = 0, trueNumCells = 0; 1035698a79b9SLisandro Dalcin int i, shift = 1; 103611c56cb0SLisandro Dalcin PetscMPIInt rank; 1037ef12879bSLisandro Dalcin PetscBool binary, zerobase = PETSC_FALSE, usemarker = PETSC_FALSE; 10380db3fc9eSLisandro Dalcin PetscBool enable_hybrid = interpolate, periodic = PETSC_TRUE; 1039331830f3SMatthew G. Knepley PetscErrorCode ierr; 1040331830f3SMatthew G. Knepley 1041331830f3SMatthew G. Knepley PetscFunctionBegin; 1042331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1043ef12879bSLisandro Dalcin ierr = PetscObjectOptionsBegin((PetscObject)viewer);CHKERRQ(ierr); 1044ef12879bSLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Gmsh options");CHKERRQ(ierr); 1045ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", enable_hybrid, &enable_hybrid, NULL);CHKERRQ(ierr); 1046ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL);CHKERRQ(ierr); 1047ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL);CHKERRQ(ierr); 1048ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_zero_base", "Read Gmsh file with zero base indices", "DMPlexCreateGmsh", zerobase, &zerobase, NULL);CHKERRQ(ierr); 10495a856986SBarry Smith ierr = PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", embedDim, &embedDim, NULL,PETSC_DECIDE);CHKERRQ(ierr); 1050ef12879bSLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 1051ef12879bSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 105211c56cb0SLisandro Dalcin if (zerobase) shift = 0; 105311c56cb0SLisandro Dalcin 1054331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1055331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 10563b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 105711c56cb0SLisandro Dalcin 105811c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 105911c56cb0SLisandro Dalcin 106011c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 10613b46f529SLisandro Dalcin if (binary) { 106211c56cb0SLisandro Dalcin parentviewer = viewer; 106311c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 106411c56cb0SLisandro Dalcin } 106511c56cb0SLisandro Dalcin 10663b46f529SLisandro Dalcin if (!rank) { 1067698a79b9SLisandro Dalcin GmshFile gmsh_, *gmsh = &gmsh_; 1068698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1069698a79b9SLisandro Dalcin PetscBool match; 1070331830f3SMatthew G. Knepley 1071580bdb30SBarry Smith ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr); 1072698a79b9SLisandro Dalcin gmsh->viewer = viewer; 1073698a79b9SLisandro Dalcin gmsh->binary = binary; 1074698a79b9SLisandro Dalcin 1075698a79b9SLisandro Dalcin /* Read mesh format */ 1076d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1077d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 1078698a79b9SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadMeshFormat(gmsh);CHKERRQ(ierr); 1079d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndMeshFormat", line);CHKERRQ(ierr); 108011c56cb0SLisandro Dalcin 1081dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 1082d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1083d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh,"$PhysicalNames", line, &match);CHKERRQ(ierr); 1084dc0ede3bSMatthew G. Knepley if (match) { 1085698a79b9SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadPhysicalNames(gmsh);CHKERRQ(ierr); 1086d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPhysicalNames", line);CHKERRQ(ierr); 1087698a79b9SLisandro Dalcin /* Initial read for entity section */ 1088d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1089dc0ede3bSMatthew G. Knepley } 109011c56cb0SLisandro Dalcin 1091de78e4feSLisandro Dalcin /* Read entities */ 1092698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 1093d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Entities", line);CHKERRQ(ierr); 1094698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1095698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadEntities_v41(gmsh, &entities);CHKERRQ(ierr); break; 1096698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadEntities_v40(gmsh, &entities);CHKERRQ(ierr); break; 1097698a79b9SLisandro Dalcin } 1098d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndEntities", line);CHKERRQ(ierr); 1099698a79b9SLisandro Dalcin /* Initial read for nodes section */ 1100d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1101de78e4feSLisandro Dalcin } 1102de78e4feSLisandro Dalcin 1103698a79b9SLisandro Dalcin /* Read nodes */ 1104d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Nodes", line);CHKERRQ(ierr); 1105698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1106698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadNodes_v41(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1107698a79b9SLisandro Dalcin case 40: ierr = DMPlexCreateGmsh_ReadNodes_v40(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1108698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadNodes_v22(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1109de78e4feSLisandro Dalcin } 1110d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndNodes", line);CHKERRQ(ierr); 111111c56cb0SLisandro Dalcin 1112698a79b9SLisandro Dalcin /* Read elements */ 1113d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);; 1114d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Elements", line);CHKERRQ(ierr); 1115698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1116698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadElements_v41(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1117698a79b9SLisandro Dalcin case 40: ierr = DMPlexCreateGmsh_ReadElements_v40(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1118d6689ca9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadElements_v22(gmsh, shift, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1119de78e4feSLisandro Dalcin } 1120d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndElements", line);CHKERRQ(ierr); 1121de78e4feSLisandro Dalcin 1122698a79b9SLisandro Dalcin /* OPTIONAL Read periodic section */ 1123698a79b9SLisandro Dalcin if (periodic) { 1124ef12879bSLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1125ef12879bSLisandro Dalcin ierr = GmshMatch(gmsh, "$Periodic", line, &periodic);CHKERRQ(ierr); 1126ef12879bSLisandro Dalcin } 1127ef12879bSLisandro Dalcin if (periodic) { 1128698a79b9SLisandro Dalcin PetscInt pVert, *periodicMapT, *aux; 1129de78e4feSLisandro Dalcin 1130698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 1131698a79b9SLisandro Dalcin ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 1132698a79b9SLisandro Dalcin for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 1133698a79b9SLisandro Dalcin 1134d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Periodic", line);CHKERRQ(ierr); 1135698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1136698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadPeriodic_v41(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break; 1137698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadPeriodic_v40(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break; 1138698a79b9SLisandro Dalcin } 1139d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPeriodic", line);CHKERRQ(ierr); 1140698a79b9SLisandro Dalcin 1141698a79b9SLisandro Dalcin /* we may have slaves of slaves */ 1142698a79b9SLisandro Dalcin for (i = 0; i < numVertices; i++) { 1143698a79b9SLisandro Dalcin while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 1144698a79b9SLisandro Dalcin periodicMapT[i] = periodicMapT[periodicMapT[i]]; 1145698a79b9SLisandro Dalcin } 1146698a79b9SLisandro Dalcin } 1147698a79b9SLisandro Dalcin /* periodicMap : from old to new numbering (periodic vertices excluded) 1148698a79b9SLisandro Dalcin periodicMapI: from new to old numbering */ 1149698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 1150698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 1151698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 1152698a79b9SLisandro Dalcin for (i = 0, pVert = 0; i < numVertices; i++) { 1153698a79b9SLisandro Dalcin if (periodicMapT[i] != i) { 1154698a79b9SLisandro Dalcin pVert++; 1155698a79b9SLisandro Dalcin } else { 1156698a79b9SLisandro Dalcin aux[i] = i - pVert; 1157698a79b9SLisandro Dalcin periodicMapI[i - pVert] = i; 1158698a79b9SLisandro Dalcin } 1159698a79b9SLisandro Dalcin } 1160698a79b9SLisandro Dalcin for (i = 0 ; i < numVertices; i++) { 1161698a79b9SLisandro Dalcin periodicMap[i] = aux[periodicMapT[i]]; 1162698a79b9SLisandro Dalcin } 1163698a79b9SLisandro Dalcin ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 1164698a79b9SLisandro Dalcin ierr = PetscFree(aux);CHKERRQ(ierr); 1165698a79b9SLisandro Dalcin /* remove periodic vertices */ 1166698a79b9SLisandro Dalcin numVertices = numVertices - pVert; 1167698a79b9SLisandro Dalcin } 1168698a79b9SLisandro Dalcin 1169698a79b9SLisandro Dalcin ierr = GmshEntitiesDestroy(&entities);CHKERRQ(ierr); 1170698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 1171698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 1172698a79b9SLisandro Dalcin } 1173698a79b9SLisandro Dalcin 1174698a79b9SLisandro Dalcin if (parentviewer) { 1175698a79b9SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 1176698a79b9SLisandro Dalcin } 1177698a79b9SLisandro Dalcin 1178698a79b9SLisandro Dalcin if (!rank) { 1179698a79b9SLisandro Dalcin PetscBool hybrid = PETSC_FALSE; 11800db3fc9eSLisandro Dalcin PetscInt cellType = -1; 1181698a79b9SLisandro Dalcin 1182a4bb7517SMichael Lange for (trueNumCells = 0, c = 0; c < numCells; ++c) { 11830db3fc9eSLisandro Dalcin if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0; hybrid = PETSC_FALSE; cellType = -1;} 11840db3fc9eSLisandro Dalcin if (gmsh_elem[c].dim < dim) continue; 11850db3fc9eSLisandro Dalcin if (cellType == -1) cellType = gmsh_elem[c].cellType; 11860db3fc9eSLisandro Dalcin /* different cell type indicate an hybrid mesh in PLEX */ 11870db3fc9eSLisandro Dalcin if (cellType != gmsh_elem[c].cellType) hybrid = PETSC_TRUE; 1188dea2a3c7SStefano Zampini /* wedges always indicate an hybrid mesh in PLEX */ 11890db3fc9eSLisandro Dalcin if (cellType == 6 || cellType == 13) hybrid = PETSC_TRUE; 11900db3fc9eSLisandro Dalcin trueNumCells++; 1191db397164SMichael Lange } 1192dea2a3c7SStefano Zampini /* Renumber cells for hybrid grids */ 1193dea2a3c7SStefano Zampini if (hybrid && enable_hybrid) { 1194dea2a3c7SStefano Zampini PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 1195dea2a3c7SStefano Zampini PetscInt cell, tn, *tp; 1196dea2a3c7SStefano Zampini int n1 = 0,n2 = 0; 1197dea2a3c7SStefano Zampini 1198dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 1199dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 1200dea2a3c7SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1201dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) { 1202dea2a3c7SStefano Zampini if (!n1) n1 = gmsh_elem[c].cellType; 1203dea2a3c7SStefano Zampini else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 1204dea2a3c7SStefano Zampini 1205dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 1206dea2a3c7SStefano Zampini else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 1207dea2a3c7SStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 1208dea2a3c7SStefano Zampini cell++; 1209dea2a3c7SStefano Zampini } 1210dea2a3c7SStefano Zampini } 1211dea2a3c7SStefano Zampini 1212dea2a3c7SStefano Zampini switch (n1) { 1213dea2a3c7SStefano Zampini case 2: /* triangles */ 1214dea2a3c7SStefano Zampini case 9: 1215dea2a3c7SStefano Zampini switch (n2) { 1216dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1217dea2a3c7SStefano Zampini case 3: /* quads */ 1218dea2a3c7SStefano Zampini case 10: 1219dea2a3c7SStefano Zampini break; 1220dea2a3c7SStefano Zampini default: 1221dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1222dea2a3c7SStefano Zampini } 1223dea2a3c7SStefano Zampini break; 1224a5b208b6SMatthew G. Knepley case 3: /* quadrilateral */ 1225dea2a3c7SStefano Zampini case 10: 1226dea2a3c7SStefano Zampini switch (n2) { 1227dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1228dea2a3c7SStefano Zampini case 2: /* swap since we list simplices first */ 1229dea2a3c7SStefano Zampini case 9: 1230dea2a3c7SStefano Zampini tn = hc1; 1231dea2a3c7SStefano Zampini hc1 = hc2; 1232dea2a3c7SStefano Zampini hc2 = tn; 1233dea2a3c7SStefano Zampini 1234dea2a3c7SStefano Zampini tp = hybridCells1; 1235dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 1236dea2a3c7SStefano Zampini hybridCells2 = tp; 1237dea2a3c7SStefano Zampini break; 1238dea2a3c7SStefano Zampini default: 1239dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1240dea2a3c7SStefano Zampini } 1241dea2a3c7SStefano Zampini break; 1242dea2a3c7SStefano Zampini case 4: /* tetrahedra */ 1243dea2a3c7SStefano Zampini case 11: 1244dea2a3c7SStefano Zampini switch (n2) { 1245dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1246dea2a3c7SStefano Zampini case 6: /* wedges */ 1247dea2a3c7SStefano Zampini case 13: 1248dea2a3c7SStefano Zampini break; 1249dea2a3c7SStefano Zampini default: 1250dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1251dea2a3c7SStefano Zampini } 1252dea2a3c7SStefano Zampini break; 1253a5b208b6SMatthew G. Knepley case 5: /* hexahedra */ 1254a5b208b6SMatthew G. Knepley case 12: 1255a5b208b6SMatthew G. Knepley switch (n2) { 1256a5b208b6SMatthew G. Knepley case 0: /* single type mesh */ 1257a5b208b6SMatthew G. Knepley case 6: /* wedges */ 1258a5b208b6SMatthew G. Knepley case 13: 1259a5b208b6SMatthew G. Knepley break; 1260a5b208b6SMatthew G. Knepley default: 1261a5b208b6SMatthew G. Knepley SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1262a5b208b6SMatthew G. Knepley } 1263a5b208b6SMatthew G. Knepley break; 1264a5b208b6SMatthew G. Knepley case 6: /* wedge */ 1265dea2a3c7SStefano Zampini case 13: 1266dea2a3c7SStefano Zampini switch (n2) { 1267dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1268a5b208b6SMatthew G. Knepley case 4: /* tetrahedra: swap since we list simplices first */ 1269dea2a3c7SStefano Zampini case 11: 1270a5b208b6SMatthew G. Knepley case 5: /* hexahedra */ 1271a5b208b6SMatthew G. Knepley case 12: 1272dea2a3c7SStefano Zampini tn = hc1; 1273dea2a3c7SStefano Zampini hc1 = hc2; 1274dea2a3c7SStefano Zampini hc2 = tn; 1275dea2a3c7SStefano Zampini 1276dea2a3c7SStefano Zampini tp = hybridCells1; 1277dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 1278dea2a3c7SStefano Zampini hybridCells2 = tp; 1279dea2a3c7SStefano Zampini break; 1280dea2a3c7SStefano Zampini default: 1281dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1282dea2a3c7SStefano Zampini } 1283dea2a3c7SStefano Zampini break; 1284dea2a3c7SStefano Zampini default: 1285dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1286dea2a3c7SStefano Zampini } 1287dea2a3c7SStefano Zampini cMax = hc1; 1288dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 1289dea2a3c7SStefano Zampini for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 1290dea2a3c7SStefano Zampini for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 1291dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 1292dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 1293dea2a3c7SStefano Zampini } 1294dea2a3c7SStefano Zampini 129511c56cb0SLisandro Dalcin } 129611c56cb0SLisandro Dalcin 1297a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 1298db397164SMichael Lange ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 1299a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 1300a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 1301dea2a3c7SStefano Zampini ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 1302a4bb7517SMichael Lange cell++; 1303db397164SMichael Lange } 1304331830f3SMatthew G. Knepley } 1305331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 1306a4bb7517SMichael Lange /* Add cell-vertex connections */ 1307a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 1308a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 130911c56cb0SLisandro Dalcin PetscInt pcone[8], corner; 1310a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1311dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 1312917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 1313db397164SMichael Lange } 131497ed6be6Ssarens if (dim == 3) { 131597ed6be6Ssarens /* Tetrahedra are inverted */ 1316a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 4 || gmsh_elem[c].cellType == 11) { 131797ed6be6Ssarens PetscInt tmp = pcone[0]; 131897ed6be6Ssarens pcone[0] = pcone[1]; 131997ed6be6Ssarens pcone[1] = tmp; 132097ed6be6Ssarens } 132197ed6be6Ssarens /* Hexahedra are inverted */ 1322a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 5 || gmsh_elem[c].cellType == 12) { 132397ed6be6Ssarens PetscInt tmp = pcone[1]; 132497ed6be6Ssarens pcone[1] = pcone[3]; 132597ed6be6Ssarens pcone[3] = tmp; 132697ed6be6Ssarens } 1327dea2a3c7SStefano Zampini /* Prisms are inverted */ 1328a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 6 || gmsh_elem[c].cellType == 13) { 1329dea2a3c7SStefano Zampini PetscInt tmp; 1330dea2a3c7SStefano Zampini 1331dea2a3c7SStefano Zampini tmp = pcone[1]; 1332dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1333dea2a3c7SStefano Zampini pcone[2] = tmp; 1334dea2a3c7SStefano Zampini tmp = pcone[4]; 1335dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1336dea2a3c7SStefano Zampini pcone[5] = tmp; 133797ed6be6Ssarens } 1338dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1339dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1340a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 3 || gmsh_elem[c].cellType == 10) { 1341dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1342dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1343dea2a3c7SStefano Zampini pcone[3] = tmp; 1344dea2a3c7SStefano Zampini } 1345dea2a3c7SStefano Zampini } 1346dea2a3c7SStefano Zampini ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 1347a4bb7517SMichael Lange cell++; 1348331830f3SMatthew G. Knepley } 1349a4bb7517SMichael Lange } 13506228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1351c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1352dea2a3c7SStefano Zampini ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1353331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1354331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1355331830f3SMatthew G. Knepley if (interpolate) { 13565fd9971aSMatthew G. Knepley DM idm; 1357331830f3SMatthew G. Knepley 1358331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1359331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 1360331830f3SMatthew G. Knepley *dm = idm; 1361331830f3SMatthew G. Knepley } 1362917266a4SLisandro Dalcin 1363f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 1364917266a4SLisandro Dalcin if (!rank && usemarker) { 1365d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 1366d3f73514SStefano Zampini 1367d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 1368d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1369d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 1370d3f73514SStefano Zampini PetscInt suppSize; 1371d3f73514SStefano Zampini 1372d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 1373d3f73514SStefano Zampini if (suppSize == 1) { 1374d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 1375d3f73514SStefano Zampini 1376d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1377d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 1378d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 1379d3f73514SStefano Zampini } 1380d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1381d3f73514SStefano Zampini } 1382d3f73514SStefano Zampini } 1383d3f73514SStefano Zampini } 138416ad7e67SMichael Lange 138516ad7e67SMichael Lange if (!rank) { 138611c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1387d1a54cefSMatthew G. Knepley 138816ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 138911c56cb0SLisandro Dalcin for (cell = 0, c = 0; c < numCells; ++c) { 139011c56cb0SLisandro Dalcin 139111c56cb0SLisandro Dalcin /* Create face sets */ 13925964b3dcSLisandro Dalcin if (interpolate && gmsh_elem[c].dim == dim-1) { 139316ad7e67SMichael Lange const PetscInt *join; 139411c56cb0SLisandro Dalcin PetscInt joinSize, pcone[8], corner; 139511c56cb0SLisandro Dalcin /* Find the relevant facet with vertex joins */ 1396a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1397dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 1398917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 139916ad7e67SMichael Lange } 140011c56cb0SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 1401f10f1c67SMatthew 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); 1402c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 1403a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 140416ad7e67SMichael Lange } 14056e1dd89aSLawrence Mitchell 14066e1dd89aSLawrence Mitchell /* Create cell sets */ 14076e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 14086e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 1409dea2a3c7SStefano Zampini ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 14106e1dd89aSLawrence Mitchell } 1411917266a4SLisandro Dalcin cell++; 14126e1dd89aSLawrence Mitchell } 14130c070f12SSander Arens 14140c070f12SSander Arens /* Create vertex sets */ 14150c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 14160c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 1417917266a4SLisandro Dalcin const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 141811c56cb0SLisandro Dalcin const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 1419d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 14200c070f12SSander Arens } 14210c070f12SSander Arens } 14220c070f12SSander Arens } 142316ad7e67SMichael Lange } 142416ad7e67SMichael Lange 142511c56cb0SLisandro Dalcin /* Create coordinates */ 1426dea2a3c7SStefano Zampini if (embedDim < 0) embedDim = dim; 1427dea2a3c7SStefano Zampini ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 1428979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1429331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 14301d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 1431f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 1432f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 1433f45c9589SStefano Zampini } else { 143475b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 1435f45c9589SStefano Zampini } 143675b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 14371d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 14381d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 1439331830f3SMatthew G. Knepley } 144011c56cb0SLisandro Dalcin if (periodicMap) { 1441437bdd5bSStefano Zampini ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 1442f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1443f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 1444437bdd5bSStefano Zampini PetscInt corner; 144511c56cb0SLisandro Dalcin PetscBool pc = PETSC_FALSE; 1446437bdd5bSStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1447917266a4SLisandro Dalcin pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 1448437bdd5bSStefano Zampini } 1449437bdd5bSStefano Zampini if (pc) { 145011c56cb0SLisandro Dalcin PetscInt dof = gmsh_elem[c].numNodes*embedDim; 1451dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1452dea2a3c7SStefano Zampini ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 1453dea2a3c7SStefano Zampini ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 1454dea2a3c7SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 14556fbe17bfSStefano Zampini } 1456f45c9589SStefano Zampini cell++; 1457f45c9589SStefano Zampini } 1458f45c9589SStefano Zampini } 1459f45c9589SStefano Zampini } 1460331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1461331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 14628b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1463331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1464331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 14651d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 1466331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 1467331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1468f45c9589SStefano Zampini if (periodicMap) { 1469f45c9589SStefano Zampini PetscInt off; 1470f45c9589SStefano Zampini 1471f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1472f45c9589SStefano Zampini PetscInt pcone[8], corner; 1473f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 1474dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1475dea2a3c7SStefano Zampini if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 1476f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1477917266a4SLisandro Dalcin pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 1478f45c9589SStefano Zampini } 1479f45c9589SStefano Zampini if (dim == 3) { 1480f45c9589SStefano Zampini /* Tetrahedra are inverted */ 1481a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 4 || gmsh_elem[c].cellType == 11) { 1482f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 1483f45c9589SStefano Zampini pcone[0] = pcone[1]; 1484f45c9589SStefano Zampini pcone[1] = tmp; 1485f45c9589SStefano Zampini } 1486f45c9589SStefano Zampini /* Hexahedra are inverted */ 1487a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 5 || gmsh_elem[c].cellType == 12) { 1488f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 1489f45c9589SStefano Zampini pcone[1] = pcone[3]; 1490f45c9589SStefano Zampini pcone[3] = tmp; 1491f45c9589SStefano Zampini } 1492dea2a3c7SStefano Zampini /* Prisms are inverted */ 1493a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 6 || gmsh_elem[c].cellType == 13) { 1494dea2a3c7SStefano Zampini PetscInt tmp; 1495dea2a3c7SStefano Zampini 1496dea2a3c7SStefano Zampini tmp = pcone[1]; 1497dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1498dea2a3c7SStefano Zampini pcone[2] = tmp; 1499dea2a3c7SStefano Zampini tmp = pcone[4]; 1500dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1501dea2a3c7SStefano Zampini pcone[5] = tmp; 1502f45c9589SStefano Zampini } 1503dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1504dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1505a1e86c74SStefano Zampini if (gmsh_elem[c].cellType == 3 || gmsh_elem[c].cellType == 10) { 1506dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1507dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1508dea2a3c7SStefano Zampini pcone[3] = tmp; 1509dea2a3c7SStefano Zampini } 1510dea2a3c7SStefano Zampini } 1511dea2a3c7SStefano Zampini ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 1512f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 151311c56cb0SLisandro Dalcin v = pcone[corner]; 1514dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 151511c56cb0SLisandro Dalcin coords[off++] = (PetscReal) coordsIn[v*3+d]; 1516f45c9589SStefano Zampini } 1517f45c9589SStefano Zampini } 15186fbe17bfSStefano Zampini } 1519f45c9589SStefano Zampini cell++; 1520f45c9589SStefano Zampini } 1521f45c9589SStefano Zampini } 1522f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 1523f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 1524dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 152511c56cb0SLisandro Dalcin coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 1526f45c9589SStefano Zampini } 1527f45c9589SStefano Zampini } 1528f45c9589SStefano Zampini } else { 1529331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 15301d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 153111c56cb0SLisandro Dalcin coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 1532331830f3SMatthew G. Knepley } 1533331830f3SMatthew G. Knepley } 1534331830f3SMatthew G. Knepley } 1535331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1536331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1537ef12879bSLisandro Dalcin 1538ef12879bSLisandro Dalcin periodic = periodicMap ? PETSC_TRUE : PETSC_FALSE; 1539ef12879bSLisandro Dalcin ierr = MPI_Bcast(&periodic, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 154011c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 154111c56cb0SLisandro Dalcin 154211c56cb0SLisandro Dalcin ierr = PetscFree(coordsIn);CHKERRQ(ierr); 154311c56cb0SLisandro Dalcin ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 1544dea2a3c7SStefano Zampini ierr = PetscFree(hybridMap);CHKERRQ(ierr); 1545d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 1546fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 15476fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 15486fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 154911c56cb0SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 155011c56cb0SLisandro Dalcin 15513b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 1552331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1553331830f3SMatthew G. Knepley } 1554