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 65*d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match) 66*d6689ca9SLisandro Dalcin { 67*d6689ca9SLisandro Dalcin PetscErrorCode ierr; 68*d6689ca9SLisandro Dalcin PetscFunctionBegin; 69*d6689ca9SLisandro Dalcin ierr = PetscStrcmp(line, Section, match);CHKERRQ(ierr); 70*d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 71*d6689ca9SLisandro Dalcin } 72*d6689ca9SLisandro Dalcin 73*d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN]) 74*d6689ca9SLisandro Dalcin { 75*d6689ca9SLisandro Dalcin PetscBool match; 76*d6689ca9SLisandro Dalcin PetscErrorCode ierr; 77*d6689ca9SLisandro Dalcin 78*d6689ca9SLisandro Dalcin PetscFunctionBegin; 79*d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, Section, line, &match);CHKERRQ(ierr); 80*d6689ca9SLisandro Dalcin if (!match) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file, expecting %s",Section); 81*d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 82*d6689ca9SLisandro Dalcin } 83*d6689ca9SLisandro Dalcin 84*d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN]) 85*d6689ca9SLisandro Dalcin { 86*d6689ca9SLisandro Dalcin PetscBool match; 87*d6689ca9SLisandro Dalcin PetscErrorCode ierr; 88*d6689ca9SLisandro Dalcin 89*d6689ca9SLisandro Dalcin PetscFunctionBegin; 90*d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 91*d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 92*d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$Comments", line, &match);CHKERRQ(ierr); 93*d6689ca9SLisandro Dalcin if (!match) break; 94*d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 95*d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 96*d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$EndComments", line, &match);CHKERRQ(ierr); 97*d6689ca9SLisandro Dalcin if (match) break; 98*d6689ca9SLisandro Dalcin } 99*d6689ca9SLisandro Dalcin } 100*d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 101*d6689ca9SLisandro Dalcin } 102*d6689ca9SLisandro Dalcin 103*d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN]) 104*d6689ca9SLisandro Dalcin { 105*d6689ca9SLisandro Dalcin PetscErrorCode ierr; 106*d6689ca9SLisandro Dalcin PetscFunctionBegin; 107*d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 108*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, EndSection, line);CHKERRQ(ierr); 109*d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 110*d6689ca9SLisandro Dalcin } 111*d6689ca9SLisandro 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; 123698a79b9SLisandro Dalcin ierr = GmshBufferSizeGet(gmsh, count, &ibuf); 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; 128698a79b9SLisandro Dalcin ierr = GmshBufferSizeGet(gmsh, count, &ibuf); 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; 133698a79b9SLisandro Dalcin ierr = GmshBufferSizeGet(gmsh, count, &ibuf); 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; 329dea2a3c7SStefano Zampini case 13: /* 18-node 2nd wedge */ 330dea2a3c7SStefano Zampini dim = 3; 331dea2a3c7SStefano Zampini numNodes = 6; 332dea2a3c7SStefano Zampini numNodesIgnore = 12; 333dea2a3c7SStefano Zampini break; 334df032b82SMatthew G. Knepley case 15: /* 1-node vertex */ 335df032b82SMatthew G. Knepley dim = 0; 336df032b82SMatthew G. Knepley numNodes = 1; 337df032b82SMatthew G. Knepley break; 338bf6ba3a3SMatthew G. Knepley case 7: /* 5-node pyramid */ 339bf6ba3a3SMatthew G. Knepley case 10: /* 9-node 2nd order quadrangle */ 340bf6ba3a3SMatthew G. Knepley case 11: /* 10-node 2nd order tetrahedron */ 341bf6ba3a3SMatthew G. Knepley case 12: /* 27-node 2nd order hexhedron */ 342bf6ba3a3SMatthew G. Knepley case 14: /* 14-node 2nd order pyramid */ 343df032b82SMatthew G. Knepley default: 344df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 345df032b82SMatthew G. Knepley } 346df032b82SMatthew G. Knepley if (binary) { 34711c56cb0SLisandro Dalcin const int nint = 1 + numTags + numNodes + numNodesIgnore; 34811c56cb0SLisandro Dalcin /* Loop over element blocks */ 349df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 35011c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 35111c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 352df032b82SMatthew G. Knepley elements[c].dim = dim; 353df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 354df032b82SMatthew G. Knepley elements[c].numTags = numTags; 355df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 356dea2a3c7SStefano Zampini elements[c].cellType = cellType; 357df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 358df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 359df032b82SMatthew G. Knepley } 360df032b82SMatthew G. Knepley } else { 361698a79b9SLisandro Dalcin const int nint = numTags + numNodes + numNodesIgnore; 362df032b82SMatthew G. Knepley elements[c].dim = dim; 363df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 364df032b82SMatthew G. Knepley elements[c].numTags = numTags; 365dea2a3c7SStefano Zampini elements[c].cellType = cellType; 366698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 367698a79b9SLisandro Dalcin for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[p]; 368698a79b9SLisandro Dalcin for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[numTags + p]; 369df032b82SMatthew G. Knepley c++; 370df032b82SMatthew G. Knepley } 371df032b82SMatthew G. Knepley } 372df032b82SMatthew G. Knepley PetscFunctionReturn(0); 373df032b82SMatthew G. Knepley } 3747d282ae0SMichael Lange 375de78e4feSLisandro Dalcin /* 376de78e4feSLisandro Dalcin $Entities 377de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 378de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 379de78e4feSLisandro Dalcin // points 380de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 381de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 382de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 383de78e4feSLisandro Dalcin ... 384de78e4feSLisandro Dalcin // curves 385de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 386de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 387de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 388de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 389de78e4feSLisandro Dalcin ... 390de78e4feSLisandro Dalcin // surfaces 391de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 392de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 393de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 394de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 395de78e4feSLisandro Dalcin ... 396de78e4feSLisandro Dalcin // volumes 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 numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 401de78e4feSLisandro Dalcin ... 402de78e4feSLisandro Dalcin $EndEntities 403de78e4feSLisandro Dalcin */ 404698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v40(GmshFile *gmsh, GmshEntities **entities) 405de78e4feSLisandro Dalcin { 406698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 407698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 408698a79b9SLisandro Dalcin long index, num, lbuf[4]; 409730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 410698a79b9SLisandro Dalcin PetscInt count[4], i; 411a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 412de78e4feSLisandro Dalcin PetscErrorCode ierr; 413de78e4feSLisandro Dalcin 414de78e4feSLisandro Dalcin PetscFunctionBegin; 415698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 416698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);} 417698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 418730356f6SLisandro Dalcin ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 419de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 420730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 421730356f6SLisandro Dalcin ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 422730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);} 423730356f6SLisandro Dalcin ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 424de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 425de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 426de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 427de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 428698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 429de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 430730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 431de78e4feSLisandro Dalcin entity->numTags = numTags = (int) PetscMin(num, 4); 432de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 433de78e4feSLisandro Dalcin if (dim == 0) continue; 434de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 435de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 436698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 437de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 438730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 439de78e4feSLisandro Dalcin } 440de78e4feSLisandro Dalcin } 441de78e4feSLisandro Dalcin PetscFunctionReturn(0); 442de78e4feSLisandro Dalcin } 443de78e4feSLisandro Dalcin 444de78e4feSLisandro Dalcin /* 445de78e4feSLisandro Dalcin $Nodes 446de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 447de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 448de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 449de78e4feSLisandro Dalcin ... 450de78e4feSLisandro Dalcin ... 451de78e4feSLisandro Dalcin $EndNodes 452de78e4feSLisandro Dalcin */ 453698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v40(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 454de78e4feSLisandro Dalcin { 455698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 456698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 457de78e4feSLisandro Dalcin long block, node, v, numEntityBlocks, numTotalNodes, numNodes; 458de78e4feSLisandro Dalcin int info[3], nid; 459de78e4feSLisandro Dalcin double *coordinates; 460de78e4feSLisandro Dalcin char *cbuf; 461de78e4feSLisandro Dalcin PetscErrorCode ierr; 462de78e4feSLisandro Dalcin 463de78e4feSLisandro Dalcin PetscFunctionBegin; 464de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 465de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 466de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 467de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 468de78e4feSLisandro Dalcin ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr); 469698a79b9SLisandro Dalcin *numVertices = numTotalNodes; 470de78e4feSLisandro Dalcin *gmsh_nodes = coordinates; 471de78e4feSLisandro Dalcin for (v = 0, block = 0; block < numEntityBlocks; ++block) { 472de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 473de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 474de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 475698a79b9SLisandro Dalcin if (gmsh->binary) { 476de78e4feSLisandro Dalcin int nbytes = sizeof(int) + 3*sizeof(double); 477698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 478de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 479de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 480de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 481de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 482de78e4feSLisandro Dalcin #if !defined(PETSC_WORDS_BIGENDIAN) 483de78e4feSLisandro Dalcin ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr); 484de78e4feSLisandro Dalcin ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr); 485de78e4feSLisandro Dalcin #endif 486de78e4feSLisandro Dalcin ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 487de78e4feSLisandro Dalcin ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 488de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 489de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 490de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 491de78e4feSLisandro Dalcin } 492de78e4feSLisandro Dalcin } else { 493de78e4feSLisandro Dalcin for (node = 0; node < numNodes; ++node, ++v) { 494de78e4feSLisandro Dalcin double *xyz = coordinates + v*3; 495de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 496de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 497de78e4feSLisandro Dalcin if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 498de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 499de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 500de78e4feSLisandro Dalcin } 501de78e4feSLisandro Dalcin } 502de78e4feSLisandro Dalcin } 503de78e4feSLisandro Dalcin PetscFunctionReturn(0); 504de78e4feSLisandro Dalcin } 505de78e4feSLisandro Dalcin 506de78e4feSLisandro Dalcin /* 507de78e4feSLisandro Dalcin $Elements 508de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 509de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 510de78e4feSLisandro Dalcin tag(int) numVert[...](int) 511de78e4feSLisandro Dalcin ... 512de78e4feSLisandro Dalcin ... 513de78e4feSLisandro Dalcin $EndElements 514de78e4feSLisandro Dalcin */ 515698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v40(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems) 516de78e4feSLisandro Dalcin { 517698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 518698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 519de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 520de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 521de78e4feSLisandro Dalcin int eid, dim, numTags, *tags, cellType, numNodes; 522a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 523de78e4feSLisandro Dalcin GmshElement *elements; 524de78e4feSLisandro Dalcin PetscErrorCode ierr; 525de78e4feSLisandro Dalcin 526de78e4feSLisandro Dalcin PetscFunctionBegin; 527de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 528de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 529de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 530de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 531de78e4feSLisandro Dalcin ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr); 532698a79b9SLisandro Dalcin *numCells = numTotalElements; 533de78e4feSLisandro Dalcin *gmsh_elems = elements; 534de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 535de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 536de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 537de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 538730356f6SLisandro Dalcin ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 539730356f6SLisandro Dalcin numTags = entity->numTags; 540730356f6SLisandro Dalcin tags = entity->tags; 541de78e4feSLisandro Dalcin switch (cellType) { 542de78e4feSLisandro Dalcin case 1: /* 2-node line */ 543de78e4feSLisandro Dalcin numNodes = 2; 544de78e4feSLisandro Dalcin break; 545de78e4feSLisandro Dalcin case 2: /* 3-node triangle */ 546698a79b9SLisandro Dalcin numNodes = 3; 547698a79b9SLisandro Dalcin break; 548de78e4feSLisandro Dalcin case 3: /* 4-node quadrangle */ 549de78e4feSLisandro Dalcin numNodes = 4; 550de78e4feSLisandro Dalcin break; 551de78e4feSLisandro Dalcin case 4: /* 4-node tetrahedron */ 552de78e4feSLisandro Dalcin numNodes = 4; 553de78e4feSLisandro Dalcin break; 554de78e4feSLisandro Dalcin case 5: /* 8-node hexahedron */ 555de78e4feSLisandro Dalcin numNodes = 8; 556de78e4feSLisandro Dalcin break; 557de78e4feSLisandro Dalcin case 6: /* 6-node wedge */ 558de78e4feSLisandro Dalcin numNodes = 6; 559de78e4feSLisandro Dalcin break; 560de78e4feSLisandro Dalcin case 15: /* 1-node vertex */ 561de78e4feSLisandro Dalcin numNodes = 1; 562de78e4feSLisandro Dalcin break; 563de78e4feSLisandro Dalcin default: 564de78e4feSLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 565de78e4feSLisandro Dalcin } 566de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 567de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 568698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 569de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 570de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 571de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 572de78e4feSLisandro Dalcin int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 573de78e4feSLisandro Dalcin GmshElement *element = elements + c; 574de78e4feSLisandro Dalcin element->dim = dim; 575de78e4feSLisandro Dalcin element->cellType = cellType; 576de78e4feSLisandro Dalcin element->numNodes = numNodes; 577de78e4feSLisandro Dalcin element->numTags = numTags; 578de78e4feSLisandro Dalcin element->id = id[0]; 579de78e4feSLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 580de78e4feSLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 581de78e4feSLisandro Dalcin } 582de78e4feSLisandro Dalcin } 583698a79b9SLisandro Dalcin PetscFunctionReturn(0); 584698a79b9SLisandro Dalcin } 585698a79b9SLisandro Dalcin 586698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v40(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt) 587698a79b9SLisandro Dalcin { 588698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 589698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 590698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 591698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 592698a79b9SLisandro Dalcin int numPeriodic, snum, i; 593698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 594698a79b9SLisandro Dalcin PetscErrorCode ierr; 595698a79b9SLisandro Dalcin 596698a79b9SLisandro Dalcin PetscFunctionBegin; 597698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 598698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 599698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 600698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 601698a79b9SLisandro Dalcin } else { 602698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 603698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 604698a79b9SLisandro Dalcin } 605698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 606698a79b9SLisandro Dalcin int ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode; 607698a79b9SLisandro Dalcin long j, nNodes; 608698a79b9SLisandro Dalcin double affine[16]; 609698a79b9SLisandro Dalcin 610698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 611698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 612698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag); 613698a79b9SLisandro Dalcin if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 614698a79b9SLisandro Dalcin } else { 615698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 616698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 617698a79b9SLisandro Dalcin slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2]; 618698a79b9SLisandro Dalcin } 619698a79b9SLisandro Dalcin (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */ 620698a79b9SLisandro Dalcin 621698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 622698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 623698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 624698a79b9SLisandro Dalcin if (snum != 1) { /* discard tranformation and try again */ 625698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 626698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 627698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 628698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 629698a79b9SLisandro Dalcin } 630698a79b9SLisandro Dalcin } else { 631698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 632698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 633698a79b9SLisandro Dalcin if (nNodes == -1) { /* discard tranformation and try again */ 634698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 635698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 636698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 637698a79b9SLisandro Dalcin } 638698a79b9SLisandro Dalcin } 639698a79b9SLisandro Dalcin 640698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 641698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 642698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 643698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 644698a79b9SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 645698a79b9SLisandro Dalcin } else { 646698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 647698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 648698a79b9SLisandro Dalcin slaveNode = ibuf[0]; masterNode = ibuf[1]; 649698a79b9SLisandro Dalcin } 650698a79b9SLisandro Dalcin slaveMap[slaveNode - shift] = masterNode - shift; 651698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, slaveNode - shift);CHKERRQ(ierr); 652698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, masterNode - shift);CHKERRQ(ierr); 653698a79b9SLisandro Dalcin } 654698a79b9SLisandro Dalcin } 655698a79b9SLisandro Dalcin PetscFunctionReturn(0); 656698a79b9SLisandro Dalcin } 657698a79b9SLisandro Dalcin 658698a79b9SLisandro Dalcin /* 659698a79b9SLisandro Dalcin $Entities 660698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 661698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 662698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 663698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 664698a79b9SLisandro Dalcin ... 665698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 666698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 667698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 668698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 669698a79b9SLisandro Dalcin ... 670698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 671698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 672698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 673698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 674698a79b9SLisandro Dalcin ... 675698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 676698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 677698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 678698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 679698a79b9SLisandro Dalcin ... 680698a79b9SLisandro Dalcin $EndEntities 681698a79b9SLisandro Dalcin */ 682698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v41(GmshFile *gmsh, GmshEntities **entities) 683698a79b9SLisandro Dalcin { 684698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 685698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 686698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 687698a79b9SLisandro Dalcin PetscErrorCode ierr; 688698a79b9SLisandro Dalcin 689698a79b9SLisandro Dalcin PetscFunctionBegin; 690698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr); 691698a79b9SLisandro Dalcin ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 692698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 693698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 694698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr); 695698a79b9SLisandro Dalcin ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 696698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr); 697698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 698698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 699698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 700698a79b9SLisandro Dalcin entity->numTags = PetscMin(numTags, 4); 701698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 702698a79b9SLisandro Dalcin if (dim == 0) continue; 703698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 704698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 705698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 706698a79b9SLisandro Dalcin } 707698a79b9SLisandro Dalcin } 708698a79b9SLisandro Dalcin PetscFunctionReturn(0); 709698a79b9SLisandro Dalcin } 710698a79b9SLisandro Dalcin 711698a79b9SLisandro Dalcin /* 712698a79b9SLisandro Dalcin $Nodes 713698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 714698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 715698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 716698a79b9SLisandro Dalcin nodeTag(size_t) 717698a79b9SLisandro Dalcin ... 718698a79b9SLisandro Dalcin x(double) y(double) z(double) 719698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 720698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 721698a79b9SLisandro Dalcin ... 722698a79b9SLisandro Dalcin ... 723698a79b9SLisandro Dalcin $EndNodes 724698a79b9SLisandro Dalcin */ 725698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v41(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes) 726698a79b9SLisandro Dalcin { 727698a79b9SLisandro Dalcin int info[3]; 728698a79b9SLisandro Dalcin PetscInt sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, *nodeTag = NULL, block, node, i; 729698a79b9SLisandro Dalcin double *coordinates; 730698a79b9SLisandro Dalcin PetscErrorCode ierr; 731698a79b9SLisandro Dalcin 732698a79b9SLisandro Dalcin PetscFunctionBegin; 733698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 734698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numNodes = sizes[1]; 735698a79b9SLisandro Dalcin ierr = PetscMalloc1(numNodes*3, &coordinates);CHKERRQ(ierr); 736698a79b9SLisandro Dalcin *numVertices = numNodes; 737698a79b9SLisandro Dalcin *gmsh_nodes = coordinates; 738698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 739698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 740698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr); 741698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodesBlock, sizeof(PetscInt), &nodeTag);CHKERRQ(ierr); 742698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTag, numNodesBlock);CHKERRQ(ierr); 743698a79b9SLisandro 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); 744698a79b9SLisandro Dalcin if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 745698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, coordinates+node*3, numNodesBlock*3);CHKERRQ(ierr); 746698a79b9SLisandro Dalcin } 747698a79b9SLisandro Dalcin PetscFunctionReturn(0); 748698a79b9SLisandro Dalcin } 749698a79b9SLisandro Dalcin 750698a79b9SLisandro Dalcin /* 751698a79b9SLisandro Dalcin $Elements 752698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 753698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 754698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 755698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 756698a79b9SLisandro Dalcin ... 757698a79b9SLisandro Dalcin ... 758698a79b9SLisandro Dalcin $EndElements 759698a79b9SLisandro Dalcin */ 760698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v41(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems) 761698a79b9SLisandro Dalcin { 762698a79b9SLisandro Dalcin int info[3], eid, dim, cellType, *tags; 763698a79b9SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numNodes, numTags, block, elem, c, p; 764698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 765698a79b9SLisandro Dalcin GmshElement *elements; 766698a79b9SLisandro Dalcin PetscErrorCode ierr; 767698a79b9SLisandro Dalcin 768698a79b9SLisandro Dalcin PetscFunctionBegin; 769698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 770698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numElements = sizes[1]; 771698a79b9SLisandro Dalcin ierr = PetscCalloc1(numElements, &elements);CHKERRQ(ierr); 772698a79b9SLisandro Dalcin *numCells = numElements; 773698a79b9SLisandro Dalcin *gmsh_elems = elements; 774698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 775698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 776698a79b9SLisandro Dalcin dim = info[0]; eid = info[1]; cellType = info[2]; 777698a79b9SLisandro Dalcin ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 778698a79b9SLisandro Dalcin numTags = entity->numTags; 779698a79b9SLisandro Dalcin tags = entity->tags; 780698a79b9SLisandro Dalcin switch (cellType) { 781698a79b9SLisandro Dalcin case 1: /* 2-node line */ 782698a79b9SLisandro Dalcin numNodes = 2; 783698a79b9SLisandro Dalcin break; 784698a79b9SLisandro Dalcin case 2: /* 3-node triangle */ 785698a79b9SLisandro Dalcin numNodes = 3; 786698a79b9SLisandro Dalcin break; 787698a79b9SLisandro Dalcin case 3: /* 4-node quadrangle */ 788698a79b9SLisandro Dalcin numNodes = 4; 789698a79b9SLisandro Dalcin break; 790698a79b9SLisandro Dalcin case 4: /* 4-node tetrahedron */ 791698a79b9SLisandro Dalcin numNodes = 4; 792698a79b9SLisandro Dalcin break; 793698a79b9SLisandro Dalcin case 5: /* 8-node hexahedron */ 794698a79b9SLisandro Dalcin numNodes = 8; 795698a79b9SLisandro Dalcin break; 796698a79b9SLisandro Dalcin case 6: /* 6-node wedge */ 797698a79b9SLisandro Dalcin numNodes = 6; 798698a79b9SLisandro Dalcin break; 799698a79b9SLisandro Dalcin case 15: /* 1-node vertex */ 800698a79b9SLisandro Dalcin numNodes = 1; 801698a79b9SLisandro Dalcin break; 802698a79b9SLisandro Dalcin default: 803698a79b9SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 804698a79b9SLisandro Dalcin } 805698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr); 806698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr); 807698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr); 808698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 809698a79b9SLisandro Dalcin GmshElement *element = elements + c; 810698a79b9SLisandro Dalcin PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 811698a79b9SLisandro Dalcin element->id = id[0]; 812698a79b9SLisandro Dalcin element->dim = dim; 813698a79b9SLisandro Dalcin element->cellType = cellType; 814698a79b9SLisandro Dalcin element->numNodes = numNodes; 815698a79b9SLisandro Dalcin element->numTags = numTags; 816698a79b9SLisandro Dalcin for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 817698a79b9SLisandro Dalcin for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 818698a79b9SLisandro Dalcin } 819698a79b9SLisandro Dalcin } 820698a79b9SLisandro Dalcin PetscFunctionReturn(0); 821698a79b9SLisandro Dalcin } 822698a79b9SLisandro Dalcin 823698a79b9SLisandro Dalcin /* 824698a79b9SLisandro Dalcin $Periodic 825698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 826698a79b9SLisandro Dalcin entityDim(int) entityTag(int) entityTagMaster(int) 827698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 828698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 829698a79b9SLisandro Dalcin nodeTag(size_t) nodeTagMaster(size_t) 830698a79b9SLisandro Dalcin ... 831698a79b9SLisandro Dalcin ... 832698a79b9SLisandro Dalcin $EndPeriodic 833698a79b9SLisandro Dalcin */ 834698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v41(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt) 835698a79b9SLisandro Dalcin { 836698a79b9SLisandro Dalcin int info[3]; 837698a79b9SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 838698a79b9SLisandro Dalcin double dbuf[16]; 839698a79b9SLisandro Dalcin PetscErrorCode ierr; 840698a79b9SLisandro Dalcin 841698a79b9SLisandro Dalcin PetscFunctionBegin; 842698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr); 843698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 844698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 845698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr); 846698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr); 847698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr); 848698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr); 849698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr); 850698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 851698a79b9SLisandro Dalcin PetscInt slaveNode = nodeTags[node*2+0] - shift; 852698a79b9SLisandro Dalcin PetscInt masterNode = nodeTags[node*2+1] - shift; 853698a79b9SLisandro Dalcin slaveMap[slaveNode] = masterNode; 854698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, slaveNode);CHKERRQ(ierr); 855698a79b9SLisandro Dalcin ierr = PetscBTSet(bt, masterNode);CHKERRQ(ierr); 856698a79b9SLisandro Dalcin } 857698a79b9SLisandro Dalcin } 858698a79b9SLisandro Dalcin PetscFunctionReturn(0); 859698a79b9SLisandro Dalcin } 860698a79b9SLisandro Dalcin 861*d6689ca9SLisandro Dalcin /* 862*d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2 863*d6689ca9SLisandro Dalcin version(ASCII double; currently 4.1) 864*d6689ca9SLisandro Dalcin file-type(ASCII int; 0 for ASCII mode, 1 for binary mode) 865*d6689ca9SLisandro Dalcin data-size(ASCII int; sizeof(size_t)) 866*d6689ca9SLisandro Dalcin < int with value one; only in binary mode, to detect endianness > 867*d6689ca9SLisandro Dalcin $EndMeshFormat 868*d6689ca9SLisandro Dalcin */ 869698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadMeshFormat(GmshFile *gmsh) 870698a79b9SLisandro Dalcin { 871698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 872698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 873698a79b9SLisandro Dalcin float version; 874698a79b9SLisandro Dalcin PetscErrorCode ierr; 875698a79b9SLisandro Dalcin 876698a79b9SLisandro Dalcin PetscFunctionBegin; 877698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr); 878698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 879698a79b9SLisandro Dalcin if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 880698a79b9SLisandro 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); 881698a79b9SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version); 882698a79b9SLisandro 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); 883698a79b9SLisandro Dalcin if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is binary but Gmsh file is ASCII"); 884698a79b9SLisandro Dalcin if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is ASCII but Gmsh file is binary"); 885698a79b9SLisandro Dalcin fileFormat = (int)(version*10); /* XXX Should use (int)roundf(version*10) ? */ 886698a79b9SLisandro 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); 887698a79b9SLisandro 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); 888698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 889698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 890698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 891698a79b9SLisandro Dalcin if (gmsh->binary) { 892698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr); 893698a79b9SLisandro Dalcin if (checkEndian != 1) { 894698a79b9SLisandro Dalcin ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr); 895698a79b9SLisandro Dalcin if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to detect endianness in Gmsh file header: %s", line); 896698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 897698a79b9SLisandro Dalcin } 898698a79b9SLisandro Dalcin } 899698a79b9SLisandro Dalcin PetscFunctionReturn(0); 900698a79b9SLisandro Dalcin } 901698a79b9SLisandro Dalcin 9021f643af3SLisandro Dalcin /* 9031f643af3SLisandro Dalcin PhysicalNames 9041f643af3SLisandro Dalcin numPhysicalNames(ASCII int) 9051f643af3SLisandro Dalcin dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max) 9061f643af3SLisandro Dalcin ... 9071f643af3SLisandro Dalcin $EndPhysicalNames 9081f643af3SLisandro Dalcin */ 909698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPhysicalNames(GmshFile *gmsh) 910698a79b9SLisandro Dalcin { 9111f643af3SLisandro Dalcin char line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q; 9121f643af3SLisandro Dalcin int snum, numRegions, region, dim, tag; 913698a79b9SLisandro Dalcin PetscErrorCode ierr; 914698a79b9SLisandro Dalcin 915698a79b9SLisandro Dalcin PetscFunctionBegin; 916698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 917698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numRegions); 918698a79b9SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 919698a79b9SLisandro Dalcin for (region = 0; region < numRegions; ++region) { 9201f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 9211f643af3SLisandro Dalcin snum = sscanf(line, "%d %d", &dim, &tag); 9221f643af3SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 9231f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, -(PetscInt)sizeof(line));CHKERRQ(ierr); 9241f643af3SLisandro Dalcin ierr = PetscStrchr(line, '"', &p);CHKERRQ(ierr); 9251f643af3SLisandro Dalcin if (!p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 9261f643af3SLisandro Dalcin ierr = PetscStrrchr(line, '"', &q);CHKERRQ(ierr); 9271f643af3SLisandro Dalcin if (q == p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 9281f643af3SLisandro Dalcin ierr = PetscStrncpy(name, p+1, (size_t)(q-p-1));CHKERRQ(ierr); 929698a79b9SLisandro Dalcin } 930de78e4feSLisandro Dalcin PetscFunctionReturn(0); 931de78e4feSLisandro Dalcin } 932de78e4feSLisandro Dalcin 933*d6689ca9SLisandro Dalcin /*@C 934*d6689ca9SLisandro Dalcin DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 935*d6689ca9SLisandro Dalcin 936*d6689ca9SLisandro Dalcin + comm - The MPI communicator 937*d6689ca9SLisandro Dalcin . filename - Name of the Gmsh file 938*d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh 939*d6689ca9SLisandro Dalcin 940*d6689ca9SLisandro Dalcin Output Parameter: 941*d6689ca9SLisandro Dalcin . dm - The DM object representing the mesh 942*d6689ca9SLisandro Dalcin 943*d6689ca9SLisandro Dalcin Level: beginner 944*d6689ca9SLisandro Dalcin 945*d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 946*d6689ca9SLisandro Dalcin @*/ 947*d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 948*d6689ca9SLisandro Dalcin { 949*d6689ca9SLisandro Dalcin PetscViewer viewer; 950*d6689ca9SLisandro Dalcin PetscMPIInt rank; 951*d6689ca9SLisandro Dalcin int fileType; 952*d6689ca9SLisandro Dalcin PetscViewerType vtype; 953*d6689ca9SLisandro Dalcin PetscErrorCode ierr; 954*d6689ca9SLisandro Dalcin 955*d6689ca9SLisandro Dalcin PetscFunctionBegin; 956*d6689ca9SLisandro Dalcin ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 957*d6689ca9SLisandro Dalcin 958*d6689ca9SLisandro Dalcin /* Determine Gmsh file type (ASCII or binary) from file header */ 959*d6689ca9SLisandro Dalcin if (!rank) { 960*d6689ca9SLisandro Dalcin GmshFile gmsh_, *gmsh = &gmsh_; 961*d6689ca9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 962*d6689ca9SLisandro Dalcin int snum; 963*d6689ca9SLisandro Dalcin float version; 964*d6689ca9SLisandro Dalcin 965*d6689ca9SLisandro Dalcin ierr = PetscMemzero(gmsh,sizeof(GmshFile));CHKERRQ(ierr); 966*d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer);CHKERRQ(ierr); 967*d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 968*d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ);CHKERRQ(ierr); 969*d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(gmsh->viewer, filename);CHKERRQ(ierr); 970*d6689ca9SLisandro Dalcin /* Read only the first two lines of the Gmsh file */ 971*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 972*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 973*d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 974*d6689ca9SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 975*d6689ca9SLisandro Dalcin if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 976*d6689ca9SLisandro 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); 977*d6689ca9SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version); 978*d6689ca9SLisandro 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); 979*d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&gmsh->viewer);CHKERRQ(ierr); 980*d6689ca9SLisandro Dalcin } 981*d6689ca9SLisandro Dalcin ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 982*d6689ca9SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 983*d6689ca9SLisandro Dalcin 984*d6689ca9SLisandro Dalcin /* Create appropriate viewer and build plex */ 985*d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 986*d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 987*d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 988*d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 989*d6689ca9SLisandro Dalcin ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 990*d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 991*d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 992*d6689ca9SLisandro Dalcin } 993*d6689ca9SLisandro Dalcin 994331830f3SMatthew G. Knepley /*@ 9957d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 996331830f3SMatthew G. Knepley 997331830f3SMatthew G. Knepley Collective on comm 998331830f3SMatthew G. Knepley 999331830f3SMatthew G. Knepley Input Parameters: 1000331830f3SMatthew G. Knepley + comm - The MPI communicator 1001331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 1002331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 1003331830f3SMatthew G. Knepley 1004331830f3SMatthew G. Knepley Output Parameter: 1005331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 1006331830f3SMatthew G. Knepley 1007698a79b9SLisandro Dalcin Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 1008331830f3SMatthew G. Knepley 1009331830f3SMatthew G. Knepley Level: beginner 1010331830f3SMatthew G. Knepley 1011331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 1012331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 1013331830f3SMatthew G. Knepley @*/ 1014331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1015331830f3SMatthew G. Knepley { 101611c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 101711c56cb0SLisandro Dalcin double *coordsIn = NULL; 1018730356f6SLisandro Dalcin GmshEntities *entities = NULL; 10193c67d5adSSatish Balay GmshElement *gmsh_elem = NULL; 1020331830f3SMatthew G. Knepley PetscSection coordSection; 1021331830f3SMatthew G. Knepley Vec coordinates; 10226fbe17bfSStefano Zampini PetscBT periodicV = NULL, periodicC = NULL; 102384572febSToby Isaac PetscScalar *coords; 1024698a79b9SLisandro Dalcin PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 1025698a79b9SLisandro Dalcin PetscInt numVertices = 0, numCells = 0, trueNumCells = 0; 1026698a79b9SLisandro Dalcin int i, shift = 1; 102711c56cb0SLisandro Dalcin PetscMPIInt rank; 1028698a79b9SLisandro Dalcin PetscBool binary, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 1029dea2a3c7SStefano Zampini PetscBool enable_hybrid = PETSC_FALSE; 1030331830f3SMatthew G. Knepley PetscErrorCode ierr; 1031331830f3SMatthew G. Knepley 1032331830f3SMatthew G. Knepley PetscFunctionBegin; 1033331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1034dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr); 1035dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 1036dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 1037dea2a3c7SStefano Zampini ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 1038dea2a3c7SStefano Zampini ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr); 103911c56cb0SLisandro Dalcin if (zerobase) shift = 0; 104011c56cb0SLisandro Dalcin 1041331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1042331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 10433b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 104411c56cb0SLisandro Dalcin 104511c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 104611c56cb0SLisandro Dalcin 104711c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 10483b46f529SLisandro Dalcin if (binary) { 104911c56cb0SLisandro Dalcin parentviewer = viewer; 105011c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 105111c56cb0SLisandro Dalcin } 105211c56cb0SLisandro Dalcin 10533b46f529SLisandro Dalcin if (!rank) { 1054698a79b9SLisandro Dalcin GmshFile gmsh_, *gmsh = &gmsh_; 1055698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1056698a79b9SLisandro Dalcin PetscBool match; 1057331830f3SMatthew G. Knepley 1058698a79b9SLisandro Dalcin ierr = PetscMemzero(gmsh,sizeof(GmshFile));CHKERRQ(ierr); 1059698a79b9SLisandro Dalcin gmsh->viewer = viewer; 1060698a79b9SLisandro Dalcin gmsh->binary = binary; 1061698a79b9SLisandro Dalcin 1062698a79b9SLisandro Dalcin /* Read mesh format */ 1063*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1064*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 1065698a79b9SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadMeshFormat(gmsh);CHKERRQ(ierr); 1066*d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndMeshFormat", line);CHKERRQ(ierr); 106711c56cb0SLisandro Dalcin 1068dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 1069*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1070*d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh,"$PhysicalNames", line, &match);CHKERRQ(ierr); 1071dc0ede3bSMatthew G. Knepley if (match) { 1072698a79b9SLisandro Dalcin ierr = DMPlexCreateGmsh_ReadPhysicalNames(gmsh);CHKERRQ(ierr); 1073*d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPhysicalNames", line);CHKERRQ(ierr); 1074698a79b9SLisandro Dalcin /* Initial read for entity section */ 1075*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1076dc0ede3bSMatthew G. Knepley } 107711c56cb0SLisandro Dalcin 1078de78e4feSLisandro Dalcin /* Read entities */ 1079698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 1080*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Entities", line);CHKERRQ(ierr); 1081698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1082698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadEntities_v41(gmsh, &entities);CHKERRQ(ierr); break; 1083698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadEntities_v40(gmsh, &entities);CHKERRQ(ierr); break; 1084698a79b9SLisandro Dalcin } 1085*d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndEntities", line);CHKERRQ(ierr); 1086698a79b9SLisandro Dalcin /* Initial read for nodes section */ 1087*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1088de78e4feSLisandro Dalcin } 1089de78e4feSLisandro Dalcin 1090698a79b9SLisandro Dalcin /* Read nodes */ 1091*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Nodes", line);CHKERRQ(ierr); 1092698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1093698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadNodes_v41(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1094698a79b9SLisandro Dalcin case 40: ierr = DMPlexCreateGmsh_ReadNodes_v40(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1095698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadNodes_v22(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break; 1096de78e4feSLisandro Dalcin } 1097*d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndNodes", line);CHKERRQ(ierr); 109811c56cb0SLisandro Dalcin 1099698a79b9SLisandro Dalcin /* Read elements */ 1100*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);; 1101*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Elements", line);CHKERRQ(ierr); 1102698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1103698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadElements_v41(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1104698a79b9SLisandro Dalcin case 40: ierr = DMPlexCreateGmsh_ReadElements_v40(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1105*d6689ca9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadElements_v22(gmsh, shift, &numCells, &gmsh_elem);CHKERRQ(ierr); break; 1106de78e4feSLisandro Dalcin } 1107*d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndElements", line);CHKERRQ(ierr); 1108de78e4feSLisandro Dalcin 1109698a79b9SLisandro Dalcin /* OPTIONAL Read periodic section */ 1110698a79b9SLisandro Dalcin if (periodic) { 1111698a79b9SLisandro Dalcin PetscInt pVert, *periodicMapT, *aux; 1112de78e4feSLisandro Dalcin 1113698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 1114698a79b9SLisandro Dalcin ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 1115698a79b9SLisandro Dalcin for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 1116698a79b9SLisandro Dalcin 1117*d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1118*d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Periodic", line);CHKERRQ(ierr); 1119698a79b9SLisandro Dalcin switch (gmsh->fileFormat) { 1120698a79b9SLisandro Dalcin case 41: ierr = DMPlexCreateGmsh_ReadPeriodic_v41(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break; 1121698a79b9SLisandro Dalcin default: ierr = DMPlexCreateGmsh_ReadPeriodic_v40(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break; 1122698a79b9SLisandro Dalcin } 1123*d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPeriodic", line);CHKERRQ(ierr); 1124698a79b9SLisandro Dalcin 1125698a79b9SLisandro Dalcin /* we may have slaves of slaves */ 1126698a79b9SLisandro Dalcin for (i = 0; i < numVertices; i++) { 1127698a79b9SLisandro Dalcin while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 1128698a79b9SLisandro Dalcin periodicMapT[i] = periodicMapT[periodicMapT[i]]; 1129698a79b9SLisandro Dalcin } 1130698a79b9SLisandro Dalcin } 1131698a79b9SLisandro Dalcin /* periodicMap : from old to new numbering (periodic vertices excluded) 1132698a79b9SLisandro Dalcin periodicMapI: from new to old numbering */ 1133698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 1134698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 1135698a79b9SLisandro Dalcin ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 1136698a79b9SLisandro Dalcin for (i = 0, pVert = 0; i < numVertices; i++) { 1137698a79b9SLisandro Dalcin if (periodicMapT[i] != i) { 1138698a79b9SLisandro Dalcin pVert++; 1139698a79b9SLisandro Dalcin } else { 1140698a79b9SLisandro Dalcin aux[i] = i - pVert; 1141698a79b9SLisandro Dalcin periodicMapI[i - pVert] = i; 1142698a79b9SLisandro Dalcin } 1143698a79b9SLisandro Dalcin } 1144698a79b9SLisandro Dalcin for (i = 0 ; i < numVertices; i++) { 1145698a79b9SLisandro Dalcin periodicMap[i] = aux[periodicMapT[i]]; 1146698a79b9SLisandro Dalcin } 1147698a79b9SLisandro Dalcin ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 1148698a79b9SLisandro Dalcin ierr = PetscFree(aux);CHKERRQ(ierr); 1149698a79b9SLisandro Dalcin /* remove periodic vertices */ 1150698a79b9SLisandro Dalcin numVertices = numVertices - pVert; 1151698a79b9SLisandro Dalcin } 1152698a79b9SLisandro Dalcin 1153698a79b9SLisandro Dalcin ierr = GmshEntitiesDestroy(&entities);CHKERRQ(ierr); 1154698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 1155698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 1156698a79b9SLisandro Dalcin } 1157698a79b9SLisandro Dalcin 1158698a79b9SLisandro Dalcin if (parentviewer) { 1159698a79b9SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 1160698a79b9SLisandro Dalcin } 1161698a79b9SLisandro Dalcin 1162698a79b9SLisandro Dalcin if (!rank) { 1163698a79b9SLisandro Dalcin PetscBool hybrid = PETSC_FALSE; 1164698a79b9SLisandro Dalcin 1165a4bb7517SMichael Lange for (trueNumCells = 0, c = 0; c < numCells; ++c) { 1166dea2a3c7SStefano Zampini int on = -1; 1167a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 1168dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) {hybrid = (trueNumCells ? (on != gmsh_elem[c].cellType ? on = gmsh_elem[c].cellType,PETSC_TRUE : hybrid) : (on = gmsh_elem[c].cellType, PETSC_FALSE) ); trueNumCells++;} 1169dea2a3c7SStefano Zampini /* wedges always indicate an hybrid mesh in PLEX */ 1170dea2a3c7SStefano Zampini if (on == 6 || on == 13) hybrid = PETSC_TRUE; 1171db397164SMichael Lange } 1172dea2a3c7SStefano Zampini /* Renumber cells for hybrid grids */ 1173dea2a3c7SStefano Zampini if (hybrid && enable_hybrid) { 1174dea2a3c7SStefano Zampini PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 1175dea2a3c7SStefano Zampini PetscInt cell, tn, *tp; 1176dea2a3c7SStefano Zampini int n1 = 0,n2 = 0; 1177dea2a3c7SStefano Zampini 1178dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 1179dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 1180dea2a3c7SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1181dea2a3c7SStefano Zampini if (gmsh_elem[c].dim == dim) { 1182dea2a3c7SStefano Zampini if (!n1) n1 = gmsh_elem[c].cellType; 1183dea2a3c7SStefano Zampini else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 1184dea2a3c7SStefano Zampini 1185dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 1186dea2a3c7SStefano Zampini else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 1187dea2a3c7SStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 1188dea2a3c7SStefano Zampini cell++; 1189dea2a3c7SStefano Zampini } 1190dea2a3c7SStefano Zampini } 1191dea2a3c7SStefano Zampini 1192dea2a3c7SStefano Zampini switch (n1) { 1193dea2a3c7SStefano Zampini case 2: /* triangles */ 1194dea2a3c7SStefano Zampini case 9: 1195dea2a3c7SStefano Zampini switch (n2) { 1196dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1197dea2a3c7SStefano Zampini case 3: /* quads */ 1198dea2a3c7SStefano Zampini case 10: 1199dea2a3c7SStefano Zampini break; 1200dea2a3c7SStefano Zampini default: 1201dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1202dea2a3c7SStefano Zampini } 1203dea2a3c7SStefano Zampini break; 1204dea2a3c7SStefano Zampini case 3: 1205dea2a3c7SStefano Zampini case 10: 1206dea2a3c7SStefano Zampini switch (n2) { 1207dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1208dea2a3c7SStefano Zampini case 2: /* swap since we list simplices first */ 1209dea2a3c7SStefano Zampini case 9: 1210dea2a3c7SStefano Zampini tn = hc1; 1211dea2a3c7SStefano Zampini hc1 = hc2; 1212dea2a3c7SStefano Zampini hc2 = tn; 1213dea2a3c7SStefano Zampini 1214dea2a3c7SStefano Zampini tp = hybridCells1; 1215dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 1216dea2a3c7SStefano Zampini hybridCells2 = tp; 1217dea2a3c7SStefano Zampini break; 1218dea2a3c7SStefano Zampini default: 1219dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1220dea2a3c7SStefano Zampini } 1221dea2a3c7SStefano Zampini break; 1222dea2a3c7SStefano Zampini case 4: /* tetrahedra */ 1223dea2a3c7SStefano Zampini case 11: 1224dea2a3c7SStefano Zampini switch (n2) { 1225dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1226dea2a3c7SStefano Zampini case 6: /* wedges */ 1227dea2a3c7SStefano Zampini case 13: 1228dea2a3c7SStefano Zampini break; 1229dea2a3c7SStefano Zampini default: 1230dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1231dea2a3c7SStefano Zampini } 1232dea2a3c7SStefano Zampini break; 1233dea2a3c7SStefano Zampini case 6: 1234dea2a3c7SStefano Zampini case 13: 1235dea2a3c7SStefano Zampini switch (n2) { 1236dea2a3c7SStefano Zampini case 0: /* single type mesh */ 1237dea2a3c7SStefano Zampini case 4: /* swap since we list simplices first */ 1238dea2a3c7SStefano Zampini case 11: 1239dea2a3c7SStefano Zampini tn = hc1; 1240dea2a3c7SStefano Zampini hc1 = hc2; 1241dea2a3c7SStefano Zampini hc2 = tn; 1242dea2a3c7SStefano Zampini 1243dea2a3c7SStefano Zampini tp = hybridCells1; 1244dea2a3c7SStefano Zampini hybridCells1 = hybridCells2; 1245dea2a3c7SStefano Zampini hybridCells2 = tp; 1246dea2a3c7SStefano Zampini break; 1247dea2a3c7SStefano Zampini default: 1248dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1249dea2a3c7SStefano Zampini } 1250dea2a3c7SStefano Zampini break; 1251dea2a3c7SStefano Zampini default: 1252dea2a3c7SStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 1253dea2a3c7SStefano Zampini } 1254dea2a3c7SStefano Zampini cMax = hc1; 1255dea2a3c7SStefano Zampini ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 1256dea2a3c7SStefano Zampini for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 1257dea2a3c7SStefano Zampini for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 1258dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 1259dea2a3c7SStefano Zampini ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 1260dea2a3c7SStefano Zampini } 1261dea2a3c7SStefano Zampini 126211c56cb0SLisandro Dalcin } 126311c56cb0SLisandro Dalcin 1264a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 1265db397164SMichael Lange ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 1266a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 1267a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 1268dea2a3c7SStefano Zampini ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 1269a4bb7517SMichael Lange cell++; 1270db397164SMichael Lange } 1271331830f3SMatthew G. Knepley } 1272331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 1273a4bb7517SMichael Lange /* Add cell-vertex connections */ 1274a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 1275a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 127611c56cb0SLisandro Dalcin PetscInt pcone[8], corner; 1277a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1278dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 1279917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 1280db397164SMichael Lange } 128197ed6be6Ssarens if (dim == 3) { 128297ed6be6Ssarens /* Tetrahedra are inverted */ 1283dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 128497ed6be6Ssarens PetscInt tmp = pcone[0]; 128597ed6be6Ssarens pcone[0] = pcone[1]; 128697ed6be6Ssarens pcone[1] = tmp; 128797ed6be6Ssarens } 128897ed6be6Ssarens /* Hexahedra are inverted */ 1289dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 129097ed6be6Ssarens PetscInt tmp = pcone[1]; 129197ed6be6Ssarens pcone[1] = pcone[3]; 129297ed6be6Ssarens pcone[3] = tmp; 129397ed6be6Ssarens } 1294dea2a3c7SStefano Zampini /* Prisms are inverted */ 1295dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 1296dea2a3c7SStefano Zampini PetscInt tmp; 1297dea2a3c7SStefano Zampini 1298dea2a3c7SStefano Zampini tmp = pcone[1]; 1299dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1300dea2a3c7SStefano Zampini pcone[2] = tmp; 1301dea2a3c7SStefano Zampini tmp = pcone[4]; 1302dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1303dea2a3c7SStefano Zampini pcone[5] = tmp; 130497ed6be6Ssarens } 1305dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1306dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1307dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 1308dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1309dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1310dea2a3c7SStefano Zampini pcone[3] = tmp; 1311dea2a3c7SStefano Zampini } 1312dea2a3c7SStefano Zampini } 1313dea2a3c7SStefano Zampini ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 1314a4bb7517SMichael Lange cell++; 1315331830f3SMatthew G. Knepley } 1316a4bb7517SMichael Lange } 13176228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1318c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1319dea2a3c7SStefano Zampini ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1320331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1321331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1322331830f3SMatthew G. Knepley if (interpolate) { 13235fd9971aSMatthew G. Knepley DM idm; 1324331830f3SMatthew G. Knepley 1325331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1326331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 1327331830f3SMatthew G. Knepley *dm = idm; 1328331830f3SMatthew G. Knepley } 1329917266a4SLisandro Dalcin 1330f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 1331917266a4SLisandro Dalcin if (!rank && usemarker) { 1332d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 1333d3f73514SStefano Zampini 1334d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 1335d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1336d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 1337d3f73514SStefano Zampini PetscInt suppSize; 1338d3f73514SStefano Zampini 1339d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 1340d3f73514SStefano Zampini if (suppSize == 1) { 1341d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 1342d3f73514SStefano Zampini 1343d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1344d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 1345d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 1346d3f73514SStefano Zampini } 1347d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1348d3f73514SStefano Zampini } 1349d3f73514SStefano Zampini } 1350d3f73514SStefano Zampini } 135116ad7e67SMichael Lange 135216ad7e67SMichael Lange if (!rank) { 135311c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1354d1a54cefSMatthew G. Knepley 135516ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 135611c56cb0SLisandro Dalcin for (cell = 0, c = 0; c < numCells; ++c) { 135711c56cb0SLisandro Dalcin 135811c56cb0SLisandro Dalcin /* Create face sets */ 13595964b3dcSLisandro Dalcin if (interpolate && gmsh_elem[c].dim == dim-1) { 136016ad7e67SMichael Lange const PetscInt *join; 136111c56cb0SLisandro Dalcin PetscInt joinSize, pcone[8], corner; 136211c56cb0SLisandro Dalcin /* Find the relevant facet with vertex joins */ 1363a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1364dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 1365917266a4SLisandro Dalcin pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 136616ad7e67SMichael Lange } 136711c56cb0SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 1368f10f1c67SMatthew 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); 1369c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 1370a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 137116ad7e67SMichael Lange } 13726e1dd89aSLawrence Mitchell 13736e1dd89aSLawrence Mitchell /* Create cell sets */ 13746e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 13756e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 1376dea2a3c7SStefano Zampini ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 13776e1dd89aSLawrence Mitchell } 1378917266a4SLisandro Dalcin cell++; 13796e1dd89aSLawrence Mitchell } 13800c070f12SSander Arens 13810c070f12SSander Arens /* Create vertex sets */ 13820c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 13830c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 1384917266a4SLisandro Dalcin const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 138511c56cb0SLisandro Dalcin const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 1386d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 13870c070f12SSander Arens } 13880c070f12SSander Arens } 13890c070f12SSander Arens } 139016ad7e67SMichael Lange } 139116ad7e67SMichael Lange 139211c56cb0SLisandro Dalcin /* Create coordinates */ 1393dea2a3c7SStefano Zampini if (embedDim < 0) embedDim = dim; 1394dea2a3c7SStefano Zampini ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 1395979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1396331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 13971d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 1398f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 1399f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 1400f45c9589SStefano Zampini } else { 140175b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 1402f45c9589SStefano Zampini } 140375b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 14041d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 14051d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 1406331830f3SMatthew G. Knepley } 140711c56cb0SLisandro Dalcin if (periodicMap) { 1408437bdd5bSStefano Zampini ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 1409f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1410f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 1411437bdd5bSStefano Zampini PetscInt corner; 141211c56cb0SLisandro Dalcin PetscBool pc = PETSC_FALSE; 1413437bdd5bSStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1414917266a4SLisandro Dalcin pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 1415437bdd5bSStefano Zampini } 1416437bdd5bSStefano Zampini if (pc) { 141711c56cb0SLisandro Dalcin PetscInt dof = gmsh_elem[c].numNodes*embedDim; 1418dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1419dea2a3c7SStefano Zampini ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 1420dea2a3c7SStefano Zampini ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 1421dea2a3c7SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 14226fbe17bfSStefano Zampini } 1423f45c9589SStefano Zampini cell++; 1424f45c9589SStefano Zampini } 1425f45c9589SStefano Zampini } 1426f45c9589SStefano Zampini } 1427331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1428331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 14298b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1430331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1431331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 14321d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 1433331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 1434331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1435f45c9589SStefano Zampini if (periodicMap) { 1436f45c9589SStefano Zampini PetscInt off; 1437f45c9589SStefano Zampini 1438f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 1439f45c9589SStefano Zampini PetscInt pcone[8], corner; 1440f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 1441dea2a3c7SStefano Zampini PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1442dea2a3c7SStefano Zampini if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 1443f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1444917266a4SLisandro Dalcin pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 1445f45c9589SStefano Zampini } 1446f45c9589SStefano Zampini if (dim == 3) { 1447f45c9589SStefano Zampini /* Tetrahedra are inverted */ 1448dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 4) { 1449f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 1450f45c9589SStefano Zampini pcone[0] = pcone[1]; 1451f45c9589SStefano Zampini pcone[1] = tmp; 1452f45c9589SStefano Zampini } 1453f45c9589SStefano Zampini /* Hexahedra are inverted */ 1454dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 5) { 1455f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 1456f45c9589SStefano Zampini pcone[1] = pcone[3]; 1457f45c9589SStefano Zampini pcone[3] = tmp; 1458f45c9589SStefano Zampini } 1459dea2a3c7SStefano Zampini /* Prisms are inverted */ 1460dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 6) { 1461dea2a3c7SStefano Zampini PetscInt tmp; 1462dea2a3c7SStefano Zampini 1463dea2a3c7SStefano Zampini tmp = pcone[1]; 1464dea2a3c7SStefano Zampini pcone[1] = pcone[2]; 1465dea2a3c7SStefano Zampini pcone[2] = tmp; 1466dea2a3c7SStefano Zampini tmp = pcone[4]; 1467dea2a3c7SStefano Zampini pcone[4] = pcone[5]; 1468dea2a3c7SStefano Zampini pcone[5] = tmp; 1469f45c9589SStefano Zampini } 1470dea2a3c7SStefano Zampini } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1471dea2a3c7SStefano Zampini /* quads are input to PLEX as prisms */ 1472dea2a3c7SStefano Zampini if (gmsh_elem[c].cellType == 3) { 1473dea2a3c7SStefano Zampini PetscInt tmp = pcone[2]; 1474dea2a3c7SStefano Zampini pcone[2] = pcone[3]; 1475dea2a3c7SStefano Zampini pcone[3] = tmp; 1476dea2a3c7SStefano Zampini } 1477dea2a3c7SStefano Zampini } 1478dea2a3c7SStefano Zampini ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 1479f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 148011c56cb0SLisandro Dalcin v = pcone[corner]; 1481dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 148211c56cb0SLisandro Dalcin coords[off++] = (PetscReal) coordsIn[v*3+d]; 1483f45c9589SStefano Zampini } 1484f45c9589SStefano Zampini } 14856fbe17bfSStefano Zampini } 1486f45c9589SStefano Zampini cell++; 1487f45c9589SStefano Zampini } 1488f45c9589SStefano Zampini } 1489f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 1490f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 1491dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 149211c56cb0SLisandro Dalcin coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 1493f45c9589SStefano Zampini } 1494f45c9589SStefano Zampini } 1495f45c9589SStefano Zampini } else { 1496331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 14971d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 149811c56cb0SLisandro Dalcin coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 1499331830f3SMatthew G. Knepley } 1500331830f3SMatthew G. Knepley } 1501331830f3SMatthew G. Knepley } 1502331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1503331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 150411c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 150511c56cb0SLisandro Dalcin 150611c56cb0SLisandro Dalcin ierr = PetscFree(coordsIn);CHKERRQ(ierr); 150711c56cb0SLisandro Dalcin ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 1508dea2a3c7SStefano Zampini ierr = PetscFree(hybridMap);CHKERRQ(ierr); 1509d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 1510fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 15116fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 15126fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 151311c56cb0SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 151411c56cb0SLisandro Dalcin 15153b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 1516331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1517331830f3SMatthew G. Knepley } 1518