xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 33a088b6677de5fc5871673ed8de1d0bcd999ed1)
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;
15*33a088b6SMatthew G. Knepley   PetscInt     nodeTagMin, nodeTagMax;
16*33a088b6SMatthew G. Knepley   PetscInt    *nodeMap;
17698a79b9SLisandro Dalcin } GmshFile;
18de78e4feSLisandro Dalcin 
19698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf)
20de78e4feSLisandro Dalcin {
21de78e4feSLisandro Dalcin   size_t         size = count * eltsize;
2211c56cb0SLisandro Dalcin   PetscErrorCode ierr;
2311c56cb0SLisandro Dalcin 
2411c56cb0SLisandro Dalcin   PetscFunctionBegin;
25698a79b9SLisandro Dalcin   if (gmsh->wlen < size) {
26698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr);
27698a79b9SLisandro Dalcin     ierr = PetscMalloc(size, &gmsh->wbuf);CHKERRQ(ierr);
28698a79b9SLisandro Dalcin     gmsh->wlen = size;
29de78e4feSLisandro Dalcin   }
30698a79b9SLisandro Dalcin   *(void**)buf = size ? gmsh->wbuf : NULL;
31de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
32de78e4feSLisandro Dalcin }
33de78e4feSLisandro Dalcin 
34698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf)
35698a79b9SLisandro Dalcin {
36698a79b9SLisandro Dalcin   size_t         dataSize = (size_t)gmsh->dataSize;
37698a79b9SLisandro Dalcin   size_t         size = count * dataSize;
38698a79b9SLisandro Dalcin   PetscErrorCode ierr;
39698a79b9SLisandro Dalcin 
40698a79b9SLisandro Dalcin   PetscFunctionBegin;
41698a79b9SLisandro Dalcin   if (gmsh->slen < size) {
42698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr);
43698a79b9SLisandro Dalcin     ierr = PetscMalloc(size, &gmsh->sbuf);CHKERRQ(ierr);
44698a79b9SLisandro Dalcin     gmsh->slen = size;
45698a79b9SLisandro Dalcin   }
46698a79b9SLisandro Dalcin   *(void**)buf = size ? gmsh->sbuf : NULL;
47698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
48698a79b9SLisandro Dalcin }
49698a79b9SLisandro Dalcin 
50698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype)
51de78e4feSLisandro Dalcin {
52de78e4feSLisandro Dalcin   PetscErrorCode ierr;
53de78e4feSLisandro Dalcin   PetscFunctionBegin;
54698a79b9SLisandro Dalcin   ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype);CHKERRQ(ierr);
55698a79b9SLisandro Dalcin   if (gmsh->byteSwap) {ierr = PetscByteSwap(buf, dtype, count);CHKERRQ(ierr);}
56698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
57698a79b9SLisandro Dalcin }
58698a79b9SLisandro Dalcin 
59698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count)
60698a79b9SLisandro Dalcin {
61698a79b9SLisandro Dalcin   PetscErrorCode ierr;
62698a79b9SLisandro Dalcin   PetscFunctionBegin;
63698a79b9SLisandro Dalcin   ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING);CHKERRQ(ierr);
64698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
65698a79b9SLisandro Dalcin }
66698a79b9SLisandro Dalcin 
67d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match)
68d6689ca9SLisandro Dalcin {
69d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
70d6689ca9SLisandro Dalcin   PetscFunctionBegin;
71d6689ca9SLisandro Dalcin   ierr = PetscStrcmp(line, Section, match);CHKERRQ(ierr);
72d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
73d6689ca9SLisandro Dalcin }
74d6689ca9SLisandro Dalcin 
75d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN])
76d6689ca9SLisandro Dalcin {
77d6689ca9SLisandro Dalcin   PetscBool      match;
78d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
79d6689ca9SLisandro Dalcin 
80d6689ca9SLisandro Dalcin   PetscFunctionBegin;
81d6689ca9SLisandro Dalcin   ierr = GmshMatch(gmsh, Section, line, &match);CHKERRQ(ierr);
82d6689ca9SLisandro Dalcin   if (!match) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file, expecting %s",Section);
83d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
84d6689ca9SLisandro Dalcin }
85d6689ca9SLisandro Dalcin 
86d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN])
87d6689ca9SLisandro Dalcin {
88d6689ca9SLisandro Dalcin   PetscBool      match;
89d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
90d6689ca9SLisandro Dalcin 
91d6689ca9SLisandro Dalcin   PetscFunctionBegin;
92d6689ca9SLisandro Dalcin   while (PETSC_TRUE) {
93d6689ca9SLisandro Dalcin     ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
94d6689ca9SLisandro Dalcin     ierr = GmshMatch(gmsh, "$Comments", line, &match);CHKERRQ(ierr);
95d6689ca9SLisandro Dalcin     if (!match) break;
96d6689ca9SLisandro Dalcin     while (PETSC_TRUE) {
97d6689ca9SLisandro Dalcin       ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
98d6689ca9SLisandro Dalcin       ierr = GmshMatch(gmsh, "$EndComments", line, &match);CHKERRQ(ierr);
99d6689ca9SLisandro Dalcin       if (match) break;
100d6689ca9SLisandro Dalcin     }
101d6689ca9SLisandro Dalcin   }
102d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
103d6689ca9SLisandro Dalcin }
104d6689ca9SLisandro Dalcin 
105d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN])
106d6689ca9SLisandro Dalcin {
107d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
108d6689ca9SLisandro Dalcin   PetscFunctionBegin;
109d6689ca9SLisandro Dalcin   ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
110d6689ca9SLisandro Dalcin   ierr = GmshExpect(gmsh, EndSection, line);CHKERRQ(ierr);
111d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
112d6689ca9SLisandro Dalcin }
113d6689ca9SLisandro Dalcin 
114698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count)
115698a79b9SLisandro Dalcin {
116698a79b9SLisandro Dalcin   PetscInt       i;
117698a79b9SLisandro Dalcin   size_t         dataSize = (size_t)gmsh->dataSize;
118698a79b9SLisandro Dalcin   PetscErrorCode ierr;
119698a79b9SLisandro Dalcin 
120698a79b9SLisandro Dalcin   PetscFunctionBegin;
121698a79b9SLisandro Dalcin   if (dataSize == sizeof(PetscInt)) {
122698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, buf, count, PETSC_INT);CHKERRQ(ierr);
123698a79b9SLisandro Dalcin   } else  if (dataSize == sizeof(int)) {
124698a79b9SLisandro Dalcin     int *ibuf = NULL;
12589d5b078SSatish Balay     ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr);
126698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, ibuf, count, PETSC_ENUM);CHKERRQ(ierr);
127698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
128698a79b9SLisandro Dalcin   } else  if (dataSize == sizeof(long)) {
129698a79b9SLisandro Dalcin     long *ibuf = NULL;
13089d5b078SSatish Balay     ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr);
131698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, ibuf, count, PETSC_LONG);CHKERRQ(ierr);
132698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
133698a79b9SLisandro Dalcin   } else if (dataSize == sizeof(PetscInt64)) {
134698a79b9SLisandro Dalcin     PetscInt64 *ibuf = NULL;
13589d5b078SSatish Balay     ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr);
136698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, ibuf, count, PETSC_INT64);CHKERRQ(ierr);
137698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
138698a79b9SLisandro Dalcin   }
139698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
140698a79b9SLisandro Dalcin }
141698a79b9SLisandro Dalcin 
142698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count)
143698a79b9SLisandro Dalcin {
144698a79b9SLisandro Dalcin   PetscErrorCode ierr;
145698a79b9SLisandro Dalcin   PetscFunctionBegin;
146698a79b9SLisandro Dalcin   ierr = GmshRead(gmsh, buf, count, PETSC_ENUM);CHKERRQ(ierr);
147698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
148698a79b9SLisandro Dalcin }
149698a79b9SLisandro Dalcin 
150698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count)
151698a79b9SLisandro Dalcin {
152698a79b9SLisandro Dalcin   PetscErrorCode ierr;
153698a79b9SLisandro Dalcin   PetscFunctionBegin;
154698a79b9SLisandro Dalcin   ierr = GmshRead(gmsh, buf, count, PETSC_DOUBLE);CHKERRQ(ierr);
155de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
156de78e4feSLisandro Dalcin }
157de78e4feSLisandro Dalcin 
158de78e4feSLisandro Dalcin typedef struct {
159de78e4feSLisandro Dalcin   PetscInt id;       /* Entity number */
160698a79b9SLisandro Dalcin   PetscInt dim;      /* Entity dimension */
161de78e4feSLisandro Dalcin   double   bbox[6];  /* Bounding box */
162de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
163de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
164de78e4feSLisandro Dalcin } GmshEntity;
165de78e4feSLisandro Dalcin 
166de78e4feSLisandro Dalcin typedef struct {
167730356f6SLisandro Dalcin   GmshEntity *entity[4];
168730356f6SLisandro Dalcin   PetscHMapI  entityMap[4];
169730356f6SLisandro Dalcin } GmshEntities;
170730356f6SLisandro Dalcin 
171698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities)
172730356f6SLisandro Dalcin {
173730356f6SLisandro Dalcin   PetscInt       dim;
174730356f6SLisandro Dalcin   PetscErrorCode ierr;
175730356f6SLisandro Dalcin 
176730356f6SLisandro Dalcin   PetscFunctionBegin;
177730356f6SLisandro Dalcin   ierr = PetscNew(entities);CHKERRQ(ierr);
178730356f6SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
179730356f6SLisandro Dalcin     ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr);
180730356f6SLisandro Dalcin     ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr);
181730356f6SLisandro Dalcin   }
182730356f6SLisandro Dalcin   PetscFunctionReturn(0);
183730356f6SLisandro Dalcin }
184730356f6SLisandro Dalcin 
185730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity)
186730356f6SLisandro Dalcin {
187730356f6SLisandro Dalcin   PetscErrorCode ierr;
188730356f6SLisandro Dalcin   PetscFunctionBegin;
189730356f6SLisandro Dalcin   ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr);
190730356f6SLisandro Dalcin   entities->entity[dim][index].dim = dim;
191730356f6SLisandro Dalcin   entities->entity[dim][index].id  = eid;
192730356f6SLisandro Dalcin   if (entity) *entity = &entities->entity[dim][index];
193730356f6SLisandro Dalcin   PetscFunctionReturn(0);
194730356f6SLisandro Dalcin }
195730356f6SLisandro Dalcin 
196730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity)
197730356f6SLisandro Dalcin {
198730356f6SLisandro Dalcin   PetscInt       index;
199730356f6SLisandro Dalcin   PetscErrorCode ierr;
200730356f6SLisandro Dalcin 
201730356f6SLisandro Dalcin   PetscFunctionBegin;
202730356f6SLisandro Dalcin   ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr);
203730356f6SLisandro Dalcin   *entity = &entities->entity[dim][index];
204730356f6SLisandro Dalcin   PetscFunctionReturn(0);
205730356f6SLisandro Dalcin }
206730356f6SLisandro Dalcin 
207730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities)
208730356f6SLisandro Dalcin {
209730356f6SLisandro Dalcin   PetscInt       dim;
210730356f6SLisandro Dalcin   PetscErrorCode ierr;
211730356f6SLisandro Dalcin 
212730356f6SLisandro Dalcin   PetscFunctionBegin;
213730356f6SLisandro Dalcin   if (!*entities) PetscFunctionReturn(0);
214730356f6SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
215730356f6SLisandro Dalcin     ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr);
216730356f6SLisandro Dalcin     ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr);
217730356f6SLisandro Dalcin   }
218730356f6SLisandro Dalcin   ierr = PetscFree((*entities));CHKERRQ(ierr);
219730356f6SLisandro Dalcin   PetscFunctionReturn(0);
220730356f6SLisandro Dalcin }
221730356f6SLisandro Dalcin 
222730356f6SLisandro Dalcin typedef struct {
223698a79b9SLisandro Dalcin   PetscInt id;       /* Entity number */
224de78e4feSLisandro Dalcin   PetscInt dim;      /* Entity dimension */
225de78e4feSLisandro Dalcin   PetscInt cellType; /* Cell type */
226de78e4feSLisandro Dalcin   PetscInt numNodes; /* Size of node array */
227698a79b9SLisandro Dalcin   PetscInt nodes[8]; /* Node array */
228de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
229de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
230de78e4feSLisandro Dalcin } GmshElement;
231de78e4feSLisandro Dalcin 
232698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v22(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes)
233de78e4feSLisandro Dalcin {
234698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
235698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
236de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
237698a79b9SLisandro Dalcin   int            v, num, nid, snum;
238698a79b9SLisandro Dalcin   double         *coordinates;
239de78e4feSLisandro Dalcin   PetscErrorCode ierr;
240de78e4feSLisandro Dalcin 
241de78e4feSLisandro Dalcin   PetscFunctionBegin;
242de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
243698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &num);
244de78e4feSLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
245698a79b9SLisandro Dalcin   ierr = PetscMalloc1(num*3, &coordinates);CHKERRQ(ierr);
246698a79b9SLisandro Dalcin   *numVertices = num;
247698a79b9SLisandro Dalcin   *gmsh_nodes = coordinates;
248698a79b9SLisandro Dalcin   for (v = 0; v < num; ++v) {
249698a79b9SLisandro Dalcin     double *xyz = coordinates + v*3;
25011c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
25111c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
25211c56cb0SLisandro Dalcin     if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift);
25311c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
25411c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
25511c56cb0SLisandro Dalcin   }
25611c56cb0SLisandro Dalcin   PetscFunctionReturn(0);
25711c56cb0SLisandro Dalcin }
25811c56cb0SLisandro Dalcin 
259de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
260de78e4feSLisandro Dalcin    file contents multiple times to figure out the true number of cells and facets
261de78e4feSLisandro Dalcin    in the given mesh. To make this more efficient we read the file contents only
262de78e4feSLisandro Dalcin    once and store them in memory, while determining the true number of cells. */
263698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v22(GmshFile* gmsh, PETSC_UNUSED int shift, PetscInt *numCells, GmshElement **gmsh_elems)
26411c56cb0SLisandro Dalcin {
265698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
266698a79b9SLisandro Dalcin   PetscBool      binary = gmsh->binary;
267698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
268de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
269df032b82SMatthew G. Knepley   GmshElement   *elements;
270698a79b9SLisandro Dalcin   int            i, c, p, num, ibuf[1+4+512], snum;
27111c56cb0SLisandro Dalcin   int            cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
272df032b82SMatthew G. Knepley   PetscErrorCode ierr;
273df032b82SMatthew G. Knepley 
274df032b82SMatthew G. Knepley   PetscFunctionBegin;
275feb237baSPierre Jolivet   ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
276698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &num);
277de78e4feSLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
278698a79b9SLisandro Dalcin   ierr = PetscMalloc1(num, &elements);CHKERRQ(ierr);
279698a79b9SLisandro Dalcin   *numCells = num;
280698a79b9SLisandro Dalcin   *gmsh_elems = elements;
281698a79b9SLisandro Dalcin   for (c = 0; c < num;) {
28211c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
28311c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
284df032b82SMatthew G. Knepley     if (binary) {
285df032b82SMatthew G. Knepley       cellType = ibuf[0];
286df032b82SMatthew G. Knepley       numElem = ibuf[1];
287df032b82SMatthew G. Knepley       numTags = ibuf[2];
288df032b82SMatthew G. Knepley     } else {
289df032b82SMatthew G. Knepley       elements[c].id = ibuf[0];
290df032b82SMatthew G. Knepley       cellType = ibuf[1];
291df032b82SMatthew G. Knepley       numTags = ibuf[2];
292df032b82SMatthew G. Knepley       numElem = 1;
293df032b82SMatthew G. Knepley     }
294bf6ba3a3SMatthew G. Knepley     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
295bf6ba3a3SMatthew G. Knepley     numNodesIgnore = 0;
296df032b82SMatthew G. Knepley     switch (cellType) {
297df032b82SMatthew G. Knepley     case 1: /* 2-node line */
298df032b82SMatthew G. Knepley       dim = 1;
299df032b82SMatthew G. Knepley       numNodes = 2;
300df032b82SMatthew G. Knepley       break;
301df032b82SMatthew G. Knepley     case 2: /* 3-node triangle */
302df032b82SMatthew G. Knepley       dim = 2;
303df032b82SMatthew G. Knepley       numNodes = 3;
304df032b82SMatthew G. Knepley       break;
305df032b82SMatthew G. Knepley     case 3: /* 4-node quadrangle */
306df032b82SMatthew G. Knepley       dim = 2;
307df032b82SMatthew G. Knepley       numNodes = 4;
308df032b82SMatthew G. Knepley       break;
309df032b82SMatthew G. Knepley     case 4: /* 4-node tetrahedron */
310df032b82SMatthew G. Knepley       dim  = 3;
311df032b82SMatthew G. Knepley       numNodes = 4;
312df032b82SMatthew G. Knepley       break;
313df032b82SMatthew G. Knepley     case 5: /* 8-node hexahedron */
314df032b82SMatthew G. Knepley       dim = 3;
315df032b82SMatthew G. Knepley       numNodes = 8;
316df032b82SMatthew G. Knepley       break;
317dea2a3c7SStefano Zampini     case 6: /* 6-node wedge */
318dea2a3c7SStefano Zampini       dim = 3;
319dea2a3c7SStefano Zampini       numNodes = 6;
320dea2a3c7SStefano Zampini       break;
321bf6ba3a3SMatthew G. Knepley     case 8: /* 3-node 2nd order line */
322bf6ba3a3SMatthew G. Knepley       dim = 1;
323bf6ba3a3SMatthew G. Knepley       numNodes = 2;
324bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 1;
325bf6ba3a3SMatthew G. Knepley       break;
326bf6ba3a3SMatthew G. Knepley     case 9: /* 6-node 2nd order triangle */
327bf6ba3a3SMatthew G. Knepley       dim = 2;
328bf6ba3a3SMatthew G. Knepley       numNodes = 3;
329bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 3;
330bf6ba3a3SMatthew G. Knepley       break;
331a1e86c74SStefano Zampini     case 10: /* 9-node 2nd order quadrangle */
332a1e86c74SStefano Zampini       dim = 2;
333a1e86c74SStefano Zampini       numNodes = 4;
334a1e86c74SStefano Zampini       numNodesIgnore = 5;
335a1e86c74SStefano Zampini       break;
336a1e86c74SStefano Zampini     case 11: /* 10-node 2nd order tetrahedron */
337a1e86c74SStefano Zampini       dim  = 3;
338a1e86c74SStefano Zampini       numNodes = 4;
339a1e86c74SStefano Zampini       numNodesIgnore = 6;
340a1e86c74SStefano Zampini       break;
341a1e86c74SStefano Zampini     case 12: /* 27-node 2nd order hexhedron */
342a1e86c74SStefano Zampini       dim  = 3;
343a1e86c74SStefano Zampini       numNodes = 8;
344a1e86c74SStefano Zampini       numNodesIgnore = 19;
345a1e86c74SStefano Zampini       break;
346dea2a3c7SStefano Zampini     case 13: /* 18-node 2nd wedge */
347dea2a3c7SStefano Zampini       dim = 3;
348dea2a3c7SStefano Zampini       numNodes = 6;
349dea2a3c7SStefano Zampini       numNodesIgnore = 12;
350dea2a3c7SStefano Zampini       break;
351df032b82SMatthew G. Knepley     case 15: /* 1-node vertex */
352df032b82SMatthew G. Knepley       dim = 0;
353df032b82SMatthew G. Knepley       numNodes = 1;
354df032b82SMatthew G. Knepley       break;
355bf6ba3a3SMatthew G. Knepley     case 7: /* 5-node pyramid */
356bf6ba3a3SMatthew G. Knepley     case 14: /* 14-node 2nd order pyramid */
357df032b82SMatthew G. Knepley     default:
358df032b82SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
359df032b82SMatthew G. Knepley     }
360df032b82SMatthew G. Knepley     if (binary) {
36111c56cb0SLisandro Dalcin       const int nint = 1 + numTags + numNodes + numNodesIgnore;
36211c56cb0SLisandro Dalcin       /* Loop over element blocks */
363df032b82SMatthew G. Knepley       for (i = 0; i < numElem; ++i, ++c) {
36411c56cb0SLisandro Dalcin         ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
36511c56cb0SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);}
366df032b82SMatthew G. Knepley         elements[c].dim = dim;
367df032b82SMatthew G. Knepley         elements[c].numNodes = numNodes;
368df032b82SMatthew G. Knepley         elements[c].numTags = numTags;
369df032b82SMatthew G. Knepley         elements[c].id = ibuf[0];
370dea2a3c7SStefano Zampini         elements[c].cellType = cellType;
371df032b82SMatthew G. Knepley         for (p = 0; p < numTags;  p++) elements[c].tags[p]  = ibuf[1 + p];
372df032b82SMatthew G. Knepley         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
373df032b82SMatthew G. Knepley       }
374df032b82SMatthew G. Knepley     } else {
375698a79b9SLisandro Dalcin       const int nint = numTags + numNodes + numNodesIgnore;
376df032b82SMatthew G. Knepley       elements[c].dim = dim;
377df032b82SMatthew G. Knepley       elements[c].numNodes = numNodes;
378df032b82SMatthew G. Knepley       elements[c].numTags = numTags;
379dea2a3c7SStefano Zampini       elements[c].cellType = cellType;
380698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
381698a79b9SLisandro Dalcin       for (p = 0; p < numTags;  p++) elements[c].tags[p]  = ibuf[p];
382698a79b9SLisandro Dalcin       for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[numTags + p];
383df032b82SMatthew G. Knepley       c++;
384df032b82SMatthew G. Knepley     }
385df032b82SMatthew G. Knepley   }
386df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
387df032b82SMatthew G. Knepley }
3887d282ae0SMichael Lange 
389de78e4feSLisandro Dalcin /*
390de78e4feSLisandro Dalcin $Entities
391de78e4feSLisandro Dalcin   numPoints(unsigned long) numCurves(unsigned long)
392de78e4feSLisandro Dalcin     numSurfaces(unsigned long) numVolumes(unsigned long)
393de78e4feSLisandro Dalcin   // points
394de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
395de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
396de78e4feSLisandro Dalcin     numPhysicals(unsigned long) phyisicalTag[...](int)
397de78e4feSLisandro Dalcin   ...
398de78e4feSLisandro Dalcin   // curves
399de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
400de78e4feSLisandro Dalcin      boxMaxX(double) boxMaxY(double) boxMaxZ(double)
401de78e4feSLisandro Dalcin      numPhysicals(unsigned long) physicalTag[...](int)
402de78e4feSLisandro Dalcin      numBREPVert(unsigned long) tagBREPVert[...](int)
403de78e4feSLisandro Dalcin   ...
404de78e4feSLisandro Dalcin   // surfaces
405de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
406de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
407de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
408de78e4feSLisandro Dalcin     numBREPCurve(unsigned long) tagBREPCurve[...](int)
409de78e4feSLisandro Dalcin   ...
410de78e4feSLisandro Dalcin   // volumes
411de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
412de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
413de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
414de78e4feSLisandro Dalcin     numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int)
415de78e4feSLisandro Dalcin   ...
416de78e4feSLisandro Dalcin $EndEntities
417de78e4feSLisandro Dalcin */
418698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v40(GmshFile *gmsh, GmshEntities **entities)
419de78e4feSLisandro Dalcin {
420698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
421698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
422698a79b9SLisandro Dalcin   long           index, num, lbuf[4];
423730356f6SLisandro Dalcin   int            dim, eid, numTags, *ibuf, t;
424698a79b9SLisandro Dalcin   PetscInt       count[4], i;
425a5ba3d71SLisandro Dalcin   GmshEntity     *entity = NULL;
426de78e4feSLisandro Dalcin   PetscErrorCode ierr;
427de78e4feSLisandro Dalcin 
428de78e4feSLisandro Dalcin   PetscFunctionBegin;
429698a79b9SLisandro Dalcin   ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr);
430698a79b9SLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);}
431698a79b9SLisandro Dalcin   for (i = 0; i < 4; ++i) count[i] = lbuf[i];
432730356f6SLisandro Dalcin   ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr);
433de78e4feSLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
434730356f6SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
435730356f6SLisandro Dalcin       ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
436730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);}
437730356f6SLisandro Dalcin       ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr);
438de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
439de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);}
440de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
441de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);}
442698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr);
443de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr);
444730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);}
445de78e4feSLisandro Dalcin       entity->numTags = numTags = (int) PetscMin(num, 4);
446de78e4feSLisandro Dalcin       for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t];
447de78e4feSLisandro Dalcin       if (dim == 0) continue;
448de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
449de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);}
450698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr);
451de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr);
452730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);}
453de78e4feSLisandro Dalcin     }
454de78e4feSLisandro Dalcin   }
455de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
456de78e4feSLisandro Dalcin }
457de78e4feSLisandro Dalcin 
458de78e4feSLisandro Dalcin /*
459de78e4feSLisandro Dalcin $Nodes
460de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numNodes(unsigned long)
461de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long)
462de78e4feSLisandro Dalcin     tag(int) x(double) y(double) z(double)
463de78e4feSLisandro Dalcin     ...
464de78e4feSLisandro Dalcin   ...
465de78e4feSLisandro Dalcin $EndNodes
466de78e4feSLisandro Dalcin */
467698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v40(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes)
468de78e4feSLisandro Dalcin {
469698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
470698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
471de78e4feSLisandro Dalcin   long           block, node, v, numEntityBlocks, numTotalNodes, numNodes;
472de78e4feSLisandro Dalcin   int            info[3], nid;
473de78e4feSLisandro Dalcin   double         *coordinates;
474de78e4feSLisandro Dalcin   PetscErrorCode ierr;
475de78e4feSLisandro Dalcin 
476de78e4feSLisandro Dalcin   PetscFunctionBegin;
477de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
478de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);}
479de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
480de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
481de78e4feSLisandro Dalcin   ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr);
482698a79b9SLisandro Dalcin   *numVertices = numTotalNodes;
483de78e4feSLisandro Dalcin   *gmsh_nodes = coordinates;
484de78e4feSLisandro Dalcin   for (v = 0, block = 0; block < numEntityBlocks; ++block) {
485de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
486de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
487de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
488698a79b9SLisandro Dalcin     if (gmsh->binary) {
489277f51e8SBarry Smith       size_t nbytes = sizeof(int) + 3*sizeof(double);
490277f51e8SBarry Smith       char   *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */
491698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr);
492de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr);
493de78e4feSLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++v) {
494de78e4feSLisandro Dalcin         char   *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int);
495de78e4feSLisandro Dalcin         double *xyz = coordinates + v*3;
49630815ce0SLisandro Dalcin         if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr);}
49730815ce0SLisandro Dalcin         if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
498de78e4feSLisandro Dalcin         ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr);
499de78e4feSLisandro Dalcin         ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr);
500de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
501de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
502277f51e8SBarry Smith         if ((long)nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %ld", nid, v+shift);
503de78e4feSLisandro Dalcin       }
504de78e4feSLisandro Dalcin     } else {
505de78e4feSLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++v) {
506de78e4feSLisandro Dalcin         double *xyz = coordinates + v*3;
507de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
508de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
509277f51e8SBarry Smith         if ((long)nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %ld", nid, v+shift);
510de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
511de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
512de78e4feSLisandro Dalcin       }
513de78e4feSLisandro Dalcin     }
514de78e4feSLisandro Dalcin   }
515de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
516de78e4feSLisandro Dalcin }
517de78e4feSLisandro Dalcin 
518de78e4feSLisandro Dalcin /*
519de78e4feSLisandro Dalcin $Elements
520de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numElements(unsigned long)
521de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long)
522de78e4feSLisandro Dalcin     tag(int) numVert[...](int)
523de78e4feSLisandro Dalcin     ...
524de78e4feSLisandro Dalcin   ...
525de78e4feSLisandro Dalcin $EndElements
526de78e4feSLisandro Dalcin */
527698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v40(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems)
528de78e4feSLisandro Dalcin {
529698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
530698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
531de78e4feSLisandro Dalcin   long           c, block, numEntityBlocks, numTotalElements, elem, numElements;
532de78e4feSLisandro Dalcin   int            p, info[3], *ibuf = NULL;
533de78e4feSLisandro Dalcin   int            eid, dim, numTags, *tags, cellType, numNodes;
534a5ba3d71SLisandro Dalcin   GmshEntity     *entity = NULL;
535de78e4feSLisandro Dalcin   GmshElement    *elements;
536de78e4feSLisandro Dalcin   PetscErrorCode ierr;
537de78e4feSLisandro Dalcin 
538de78e4feSLisandro Dalcin   PetscFunctionBegin;
539de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
540de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);}
541de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
542de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);}
543de78e4feSLisandro Dalcin   ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr);
544698a79b9SLisandro Dalcin   *numCells = numTotalElements;
545de78e4feSLisandro Dalcin   *gmsh_elems = elements;
546de78e4feSLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
547de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
548de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);}
549de78e4feSLisandro Dalcin     eid = info[0]; dim = info[1]; cellType = info[2];
550730356f6SLisandro Dalcin     ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr);
551730356f6SLisandro Dalcin     numTags = entity->numTags;
552730356f6SLisandro Dalcin     tags = entity->tags;
553de78e4feSLisandro Dalcin     switch (cellType) {
554de78e4feSLisandro Dalcin     case 1: /* 2-node line */
555de78e4feSLisandro Dalcin       numNodes = 2;
556de78e4feSLisandro Dalcin       break;
557de78e4feSLisandro Dalcin     case 2: /* 3-node triangle */
558698a79b9SLisandro Dalcin       numNodes = 3;
559698a79b9SLisandro Dalcin       break;
560de78e4feSLisandro Dalcin     case 3: /* 4-node quadrangle */
561de78e4feSLisandro Dalcin       numNodes = 4;
562de78e4feSLisandro Dalcin       break;
563de78e4feSLisandro Dalcin     case 4: /* 4-node tetrahedron */
564de78e4feSLisandro Dalcin       numNodes = 4;
565de78e4feSLisandro Dalcin       break;
566de78e4feSLisandro Dalcin     case 5: /* 8-node hexahedron */
567de78e4feSLisandro Dalcin       numNodes = 8;
568de78e4feSLisandro Dalcin       break;
569de78e4feSLisandro Dalcin     case 6: /* 6-node wedge */
570de78e4feSLisandro Dalcin       numNodes = 6;
571de78e4feSLisandro Dalcin       break;
572de78e4feSLisandro Dalcin     case 15: /* 1-node vertex */
573de78e4feSLisandro Dalcin       numNodes = 1;
574de78e4feSLisandro Dalcin       break;
575de78e4feSLisandro Dalcin     default:
576de78e4feSLisandro Dalcin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
577de78e4feSLisandro Dalcin     }
578de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
579de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);}
580698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr);
581de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr);
582de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);}
583de78e4feSLisandro Dalcin     for (elem = 0; elem < numElements; ++elem, ++c) {
584de78e4feSLisandro Dalcin       int *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
585de78e4feSLisandro Dalcin       GmshElement *element = elements + c;
586de78e4feSLisandro Dalcin       element->dim = dim;
587de78e4feSLisandro Dalcin       element->cellType = cellType;
588de78e4feSLisandro Dalcin       element->numNodes = numNodes;
589de78e4feSLisandro Dalcin       element->numTags = numTags;
590de78e4feSLisandro Dalcin       element->id = id[0];
591de78e4feSLisandro Dalcin       for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p];
592de78e4feSLisandro Dalcin       for (p = 0; p < numTags;  p++) element->tags[p]  = tags[p];
593de78e4feSLisandro Dalcin     }
594de78e4feSLisandro Dalcin   }
595698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
596698a79b9SLisandro Dalcin }
597698a79b9SLisandro Dalcin 
598698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v40(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt)
599698a79b9SLisandro Dalcin {
600698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
601698a79b9SLisandro Dalcin   int            fileFormat = gmsh->fileFormat;
602698a79b9SLisandro Dalcin   PetscBool      binary = gmsh->binary;
603698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
604698a79b9SLisandro Dalcin   int            numPeriodic, snum, i;
605698a79b9SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
606698a79b9SLisandro Dalcin   PetscErrorCode ierr;
607698a79b9SLisandro Dalcin 
608698a79b9SLisandro Dalcin   PetscFunctionBegin;
609698a79b9SLisandro Dalcin   if (fileFormat == 22 || !binary) {
610698a79b9SLisandro Dalcin     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
611698a79b9SLisandro Dalcin     snum = sscanf(line, "%d", &numPeriodic);
612698a79b9SLisandro Dalcin     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
613698a79b9SLisandro Dalcin   } else {
614698a79b9SLisandro Dalcin     ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
615698a79b9SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);}
616698a79b9SLisandro Dalcin   }
617698a79b9SLisandro Dalcin   for (i = 0; i < numPeriodic; i++) {
618698a79b9SLisandro Dalcin     int    ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode;
619698a79b9SLisandro Dalcin     long   j, nNodes;
620698a79b9SLisandro Dalcin     double affine[16];
621698a79b9SLisandro Dalcin 
622698a79b9SLisandro Dalcin     if (fileFormat == 22 || !binary) {
623698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
624698a79b9SLisandro Dalcin       snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag);
625698a79b9SLisandro Dalcin       if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
626698a79b9SLisandro Dalcin     } else {
627698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
628698a79b9SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
629698a79b9SLisandro Dalcin       slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2];
630698a79b9SLisandro Dalcin     }
631698a79b9SLisandro Dalcin     (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */
632698a79b9SLisandro Dalcin 
633698a79b9SLisandro Dalcin     if (fileFormat == 22 || !binary) {
634698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
635698a79b9SLisandro Dalcin       snum = sscanf(line, "%ld", &nNodes);
6364c500f23SPierre Jolivet       if (snum != 1) { /* discard transformation and try again */
637698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
638698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
639698a79b9SLisandro Dalcin         snum = sscanf(line, "%ld", &nNodes);
640698a79b9SLisandro Dalcin         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
641698a79b9SLisandro Dalcin       }
642698a79b9SLisandro Dalcin     } else {
643698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
644698a79b9SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
6454c500f23SPierre Jolivet       if (nNodes == -1) { /* discard transformation and try again */
646698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
647698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
648698a79b9SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
649698a79b9SLisandro Dalcin       }
650698a79b9SLisandro Dalcin     }
651698a79b9SLisandro Dalcin 
652698a79b9SLisandro Dalcin     for (j = 0; j < nNodes; j++) {
653698a79b9SLisandro Dalcin       if (fileFormat == 22 || !binary) {
654698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
655698a79b9SLisandro Dalcin         snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
656698a79b9SLisandro Dalcin         if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
657698a79b9SLisandro Dalcin       } else {
658698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr);
659698a79b9SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);}
660698a79b9SLisandro Dalcin         slaveNode = ibuf[0]; masterNode = ibuf[1];
661698a79b9SLisandro Dalcin       }
662698a79b9SLisandro Dalcin       slaveMap[slaveNode - shift] = masterNode - shift;
663698a79b9SLisandro Dalcin       ierr = PetscBTSet(bt, slaveNode - shift);CHKERRQ(ierr);
664698a79b9SLisandro Dalcin       ierr = PetscBTSet(bt, masterNode - shift);CHKERRQ(ierr);
665698a79b9SLisandro Dalcin     }
666698a79b9SLisandro Dalcin   }
667698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
668698a79b9SLisandro Dalcin }
669698a79b9SLisandro Dalcin 
670698a79b9SLisandro Dalcin /*
671698a79b9SLisandro Dalcin $Entities
672698a79b9SLisandro Dalcin   numPoints(size_t) numCurves(size_t)
673698a79b9SLisandro Dalcin     numSurfaces(size_t) numVolumes(size_t)
674698a79b9SLisandro Dalcin   pointTag(int) X(double) Y(double) Z(double)
675698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
676698a79b9SLisandro Dalcin   ...
677698a79b9SLisandro Dalcin   curveTag(int) minX(double) minY(double) minZ(double)
678698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
679698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
680698a79b9SLisandro Dalcin     numBoundingPoints(size_t) pointTag(int) ...
681698a79b9SLisandro Dalcin   ...
682698a79b9SLisandro Dalcin   surfaceTag(int) minX(double) minY(double) minZ(double)
683698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
684698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
685698a79b9SLisandro Dalcin     numBoundingCurves(size_t) curveTag(int) ...
686698a79b9SLisandro Dalcin   ...
687698a79b9SLisandro Dalcin   volumeTag(int) minX(double) minY(double) minZ(double)
688698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
689698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
690698a79b9SLisandro Dalcin     numBoundngSurfaces(size_t) surfaceTag(int) ...
691698a79b9SLisandro Dalcin   ...
692698a79b9SLisandro Dalcin $EndEntities
693698a79b9SLisandro Dalcin */
694698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v41(GmshFile *gmsh, GmshEntities **entities)
695698a79b9SLisandro Dalcin {
696698a79b9SLisandro Dalcin   PetscInt       count[4], index, numTags, i;
697698a79b9SLisandro Dalcin   int            dim, eid, *tags = NULL;
698698a79b9SLisandro Dalcin   GmshEntity     *entity = NULL;
699698a79b9SLisandro Dalcin   PetscErrorCode ierr;
700698a79b9SLisandro Dalcin 
701698a79b9SLisandro Dalcin   PetscFunctionBegin;
702698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr);
703698a79b9SLisandro Dalcin   ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr);
704698a79b9SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
705698a79b9SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
706698a79b9SLisandro Dalcin       ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr);
707698a79b9SLisandro Dalcin       ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr);
708698a79b9SLisandro Dalcin       ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr);
709698a79b9SLisandro Dalcin       ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr);
710698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr);
711698a79b9SLisandro Dalcin       ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr);
712698a79b9SLisandro Dalcin       entity->numTags = PetscMin(numTags, 4);
713698a79b9SLisandro Dalcin       for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i];
714698a79b9SLisandro Dalcin       if (dim == 0) continue;
715698a79b9SLisandro Dalcin       ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr);
716698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr);
717698a79b9SLisandro Dalcin       ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr);
718698a79b9SLisandro Dalcin     }
719698a79b9SLisandro Dalcin   }
720698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
721698a79b9SLisandro Dalcin }
722698a79b9SLisandro Dalcin 
723*33a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
724698a79b9SLisandro Dalcin $Nodes
725698a79b9SLisandro Dalcin   numEntityBlocks(size_t) numNodes(size_t)
726698a79b9SLisandro Dalcin     minNodeTag(size_t) maxNodeTag(size_t)
727698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t)
728698a79b9SLisandro Dalcin     nodeTag(size_t)
729698a79b9SLisandro Dalcin     ...
730698a79b9SLisandro Dalcin     x(double) y(double) z(double)
731698a79b9SLisandro Dalcin        < u(double; if parametric and entityDim = 1 or entityDim = 2) >
732698a79b9SLisandro Dalcin        < v(double; if parametric and entityDim = 2) >
733698a79b9SLisandro Dalcin     ...
734698a79b9SLisandro Dalcin   ...
735698a79b9SLisandro Dalcin $EndNodes
736*33a088b6SMatthew G. Knepley 
737*33a088b6SMatthew G. Knepley We need to allow holes in the node numbering, so we must have a map between the file numbering and a contiguous numbering.
738698a79b9SLisandro Dalcin */
739698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v41(GmshFile *gmsh, int shift, PetscInt *numVertices, double **gmsh_nodes)
740698a79b9SLisandro Dalcin {
741698a79b9SLisandro Dalcin   int            info[3];
742*33a088b6SMatthew G. Knepley   PetscInt       sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, block, *nodeTag = NULL, node, i;
743698a79b9SLisandro Dalcin   double         *coordinates;
744698a79b9SLisandro Dalcin   PetscErrorCode ierr;
745698a79b9SLisandro Dalcin 
746698a79b9SLisandro Dalcin   PetscFunctionBegin;
747698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr);
748*33a088b6SMatthew G. Knepley   numEntityBlocks = sizes[0]; numNodes = sizes[1]; gmsh->nodeTagMin = sizes[2]; gmsh->nodeTagMax = sizes[3]+1;
749698a79b9SLisandro Dalcin   ierr = PetscMalloc1(numNodes*3, &coordinates);CHKERRQ(ierr);
750*33a088b6SMatthew G. Knepley   ierr = PetscMalloc1(gmsh->nodeTagMax-gmsh->nodeTagMin, &gmsh->nodeMap);CHKERRQ(ierr);
751*33a088b6SMatthew G. Knepley   for (i = 0; i < gmsh->nodeTagMax-gmsh->nodeTagMin; ++i) gmsh->nodeMap[i] = -1;
752698a79b9SLisandro Dalcin   *numVertices = numNodes;
753698a79b9SLisandro Dalcin   *gmsh_nodes = coordinates;
754698a79b9SLisandro Dalcin   for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) {
755698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr);
756698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr);
757698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, numNodesBlock, sizeof(PetscInt), &nodeTag);CHKERRQ(ierr);
758698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, nodeTag, numNodesBlock);CHKERRQ(ierr);
759*33a088b6SMatthew G. Knepley     for (i = 0; i < numNodesBlock; ++i) {
760*33a088b6SMatthew G. Knepley       const PetscInt off = nodeTag[i] - gmsh->nodeTagMin;
761*33a088b6SMatthew G. Knepley 
762*33a088b6SMatthew G. Knepley       if (gmsh->nodeMap[off] != -1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Repeated node tag %D", nodeTag[i]);
763*33a088b6SMatthew G. Knepley       gmsh->nodeMap[off] = node+i+shift;
764*33a088b6SMatthew G. Knepley     }
765698a79b9SLisandro Dalcin     if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported");
766698a79b9SLisandro Dalcin     ierr = GmshReadDouble(gmsh, coordinates+node*3, numNodesBlock*3);CHKERRQ(ierr);
767698a79b9SLisandro Dalcin   }
768698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
769698a79b9SLisandro Dalcin }
770698a79b9SLisandro Dalcin 
771*33a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
772698a79b9SLisandro Dalcin $Elements
773698a79b9SLisandro Dalcin   numEntityBlocks(size_t) numElements(size_t)
774698a79b9SLisandro Dalcin     minElementTag(size_t) maxElementTag(size_t)
775698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t)
776698a79b9SLisandro Dalcin     elementTag(size_t) nodeTag(size_t) ...
777698a79b9SLisandro Dalcin     ...
778698a79b9SLisandro Dalcin   ...
779698a79b9SLisandro Dalcin $EndElements
780698a79b9SLisandro Dalcin */
781698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v41(GmshFile *gmsh, PETSC_UNUSED int shift, GmshEntities *entities, PetscInt *numCells, GmshElement **gmsh_elems)
782698a79b9SLisandro Dalcin {
783698a79b9SLisandro Dalcin   int            info[3], eid, dim, cellType, *tags;
784*33a088b6SMatthew G. Knepley   PetscInt       nodeTagMin = gmsh->nodeTagMin, *nodeMap = gmsh->nodeMap, numNodes;
785*33a088b6SMatthew G. Knepley   PetscInt       sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numTags, block, elem, c, p;
786698a79b9SLisandro Dalcin   GmshEntity     *entity = NULL;
787698a79b9SLisandro Dalcin   GmshElement    *elements;
788698a79b9SLisandro Dalcin   PetscErrorCode ierr;
789698a79b9SLisandro Dalcin 
790698a79b9SLisandro Dalcin   PetscFunctionBegin;
791698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr);
792698a79b9SLisandro Dalcin   numEntityBlocks = sizes[0]; numElements = sizes[1];
793698a79b9SLisandro Dalcin   ierr = PetscCalloc1(numElements, &elements);CHKERRQ(ierr);
794698a79b9SLisandro Dalcin   *numCells = numElements;
795698a79b9SLisandro Dalcin   *gmsh_elems = elements;
796698a79b9SLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
797698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr);
798698a79b9SLisandro Dalcin     dim = info[0]; eid = info[1]; cellType = info[2];
799698a79b9SLisandro Dalcin     ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr);
800698a79b9SLisandro Dalcin     numTags = entity->numTags;
801698a79b9SLisandro Dalcin     tags = entity->tags;
802698a79b9SLisandro Dalcin     switch (cellType) {
803698a79b9SLisandro Dalcin     case 1: /* 2-node line */
804698a79b9SLisandro Dalcin       numNodes = 2;
805698a79b9SLisandro Dalcin       break;
806698a79b9SLisandro Dalcin     case 2: /* 3-node triangle */
807698a79b9SLisandro Dalcin       numNodes = 3;
808698a79b9SLisandro Dalcin       break;
809698a79b9SLisandro Dalcin     case 3: /* 4-node quadrangle */
810698a79b9SLisandro Dalcin       numNodes = 4;
811698a79b9SLisandro Dalcin       break;
812698a79b9SLisandro Dalcin     case 4: /* 4-node tetrahedron */
813698a79b9SLisandro Dalcin       numNodes = 4;
814698a79b9SLisandro Dalcin       break;
815698a79b9SLisandro Dalcin     case 5: /* 8-node hexahedron */
816698a79b9SLisandro Dalcin       numNodes = 8;
817698a79b9SLisandro Dalcin       break;
818698a79b9SLisandro Dalcin     case 6: /* 6-node wedge */
819698a79b9SLisandro Dalcin       numNodes = 6;
820698a79b9SLisandro Dalcin       break;
821698a79b9SLisandro Dalcin     case 15: /* 1-node vertex */
822698a79b9SLisandro Dalcin       numNodes = 1;
823698a79b9SLisandro Dalcin       break;
824698a79b9SLisandro Dalcin     default:
825698a79b9SLisandro Dalcin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
826698a79b9SLisandro Dalcin     }
827698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr);
828698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr);
829698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr);
830698a79b9SLisandro Dalcin     for (elem = 0; elem < numBlockElements; ++elem, ++c) {
831698a79b9SLisandro Dalcin       GmshElement *element = elements + c;
832698a79b9SLisandro Dalcin       PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
833698a79b9SLisandro Dalcin       element->id       = id[0];
834698a79b9SLisandro Dalcin       element->dim      = dim;
835698a79b9SLisandro Dalcin       element->cellType = cellType;
836698a79b9SLisandro Dalcin       element->numNodes = numNodes;
837698a79b9SLisandro Dalcin       element->numTags  = numTags;
838*33a088b6SMatthew G. Knepley       for (p = 0; p < numNodes; p++) element->nodes[p] = nodeMap[nodes[p]-nodeTagMin];
839698a79b9SLisandro Dalcin       for (p = 0; p < numTags;  p++) element->tags[p]  = tags[p];
840698a79b9SLisandro Dalcin     }
841698a79b9SLisandro Dalcin   }
842698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
843698a79b9SLisandro Dalcin }
844698a79b9SLisandro Dalcin 
845698a79b9SLisandro Dalcin /*
846698a79b9SLisandro Dalcin $Periodic
847698a79b9SLisandro Dalcin   numPeriodicLinks(size_t)
848698a79b9SLisandro Dalcin  entityDim(int) entityTag(int) entityTagMaster(int)
849698a79b9SLisandro Dalcin   numAffine(size_t) value(double) ...
850698a79b9SLisandro Dalcin   numCorrespondingNodes(size_t)
851698a79b9SLisandro Dalcin     nodeTag(size_t) nodeTagMaster(size_t)
852698a79b9SLisandro Dalcin     ...
853698a79b9SLisandro Dalcin   ...
854698a79b9SLisandro Dalcin $EndPeriodic
855698a79b9SLisandro Dalcin */
856698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPeriodic_v41(GmshFile *gmsh, int shift, PetscInt slaveMap[], PetscBT bt)
857698a79b9SLisandro Dalcin {
858698a79b9SLisandro Dalcin   int            info[3];
859698a79b9SLisandro Dalcin   PetscInt       numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node;
860698a79b9SLisandro Dalcin   double         dbuf[16];
861698a79b9SLisandro Dalcin   PetscErrorCode ierr;
862698a79b9SLisandro Dalcin 
863698a79b9SLisandro Dalcin   PetscFunctionBegin;
864698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr);
865698a79b9SLisandro Dalcin   for (link = 0; link < numPeriodicLinks; ++link) {
866698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr);
867698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr);
868698a79b9SLisandro Dalcin     ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr);
869698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr);
870698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr);
871698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr);
872698a79b9SLisandro Dalcin     for (node = 0; node < numCorrespondingNodes; ++node) {
873698a79b9SLisandro Dalcin       PetscInt slaveNode  = nodeTags[node*2+0] - shift;
874698a79b9SLisandro Dalcin       PetscInt masterNode = nodeTags[node*2+1] - shift;
875698a79b9SLisandro Dalcin       slaveMap[slaveNode] = masterNode;
876698a79b9SLisandro Dalcin       ierr = PetscBTSet(bt, slaveNode);CHKERRQ(ierr);
877698a79b9SLisandro Dalcin       ierr = PetscBTSet(bt, masterNode);CHKERRQ(ierr);
878698a79b9SLisandro Dalcin     }
879698a79b9SLisandro Dalcin   }
880698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
881698a79b9SLisandro Dalcin }
882698a79b9SLisandro Dalcin 
883d6689ca9SLisandro Dalcin /*
884d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2
885d6689ca9SLisandro Dalcin   version(ASCII double; currently 4.1)
886d6689ca9SLisandro Dalcin   file-type(ASCII int; 0 for ASCII mode, 1 for binary mode)
887d6689ca9SLisandro Dalcin   data-size(ASCII int; sizeof(size_t))
888d6689ca9SLisandro Dalcin   < int with value one; only in binary mode, to detect endianness >
889d6689ca9SLisandro Dalcin $EndMeshFormat
890d6689ca9SLisandro Dalcin */
891698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadMeshFormat(GmshFile *gmsh)
892698a79b9SLisandro Dalcin {
893698a79b9SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
894698a79b9SLisandro Dalcin   int            snum, fileType, fileFormat, dataSize, checkEndian;
895698a79b9SLisandro Dalcin   float          version;
896698a79b9SLisandro Dalcin   PetscErrorCode ierr;
897698a79b9SLisandro Dalcin 
898698a79b9SLisandro Dalcin   PetscFunctionBegin;
899698a79b9SLisandro Dalcin   ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr);
900698a79b9SLisandro Dalcin   snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
901698a79b9SLisandro Dalcin   if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
902698a79b9SLisandro 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);
903698a79b9SLisandro Dalcin   if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version);
904698a79b9SLisandro 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);
905698a79b9SLisandro Dalcin   if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is binary but Gmsh file is ASCII");
906698a79b9SLisandro Dalcin   if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Viewer is ASCII but Gmsh file is binary");
907698a79b9SLisandro Dalcin   fileFormat = (int)(version*10); /* XXX Should use (int)roundf(version*10) ? */
908698a79b9SLisandro 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);
909698a79b9SLisandro 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);
910698a79b9SLisandro Dalcin   gmsh->fileFormat = fileFormat;
911698a79b9SLisandro Dalcin   gmsh->dataSize = dataSize;
912698a79b9SLisandro Dalcin   gmsh->byteSwap = PETSC_FALSE;
913698a79b9SLisandro Dalcin   if (gmsh->binary) {
914698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr);
915698a79b9SLisandro Dalcin     if (checkEndian != 1) {
916698a79b9SLisandro Dalcin       ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr);
917698a79b9SLisandro Dalcin       if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to detect endianness in Gmsh file header: %s", line);
918698a79b9SLisandro Dalcin       gmsh->byteSwap = PETSC_TRUE;
919698a79b9SLisandro Dalcin     }
920698a79b9SLisandro Dalcin   }
921698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
922698a79b9SLisandro Dalcin }
923698a79b9SLisandro Dalcin 
9241f643af3SLisandro Dalcin /*
9251f643af3SLisandro Dalcin PhysicalNames
9261f643af3SLisandro Dalcin   numPhysicalNames(ASCII int)
9271f643af3SLisandro Dalcin   dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max)
9281f643af3SLisandro Dalcin   ...
9291f643af3SLisandro Dalcin $EndPhysicalNames
9301f643af3SLisandro Dalcin */
931698a79b9SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadPhysicalNames(GmshFile *gmsh)
932698a79b9SLisandro Dalcin {
9331f643af3SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q;
9341f643af3SLisandro Dalcin   int            snum, numRegions, region, dim, tag;
935698a79b9SLisandro Dalcin   PetscErrorCode ierr;
936698a79b9SLisandro Dalcin 
937698a79b9SLisandro Dalcin   PetscFunctionBegin;
938698a79b9SLisandro Dalcin   ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
939698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &numRegions);
940698a79b9SLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
941698a79b9SLisandro Dalcin   for (region = 0; region < numRegions; ++region) {
9421f643af3SLisandro Dalcin     ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr);
9431f643af3SLisandro Dalcin     snum = sscanf(line, "%d %d", &dim, &tag);
9441f643af3SLisandro Dalcin     if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
9451f643af3SLisandro Dalcin     ierr = GmshReadString(gmsh, line, -(PetscInt)sizeof(line));CHKERRQ(ierr);
9461f643af3SLisandro Dalcin     ierr = PetscStrchr(line, '"', &p);CHKERRQ(ierr);
9471f643af3SLisandro Dalcin     if (!p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
9481f643af3SLisandro Dalcin     ierr = PetscStrrchr(line, '"', &q);CHKERRQ(ierr);
9491f643af3SLisandro Dalcin     if (q == p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
9501f643af3SLisandro Dalcin     ierr = PetscStrncpy(name, p+1, (size_t)(q-p-1));CHKERRQ(ierr);
951698a79b9SLisandro Dalcin   }
952de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
953de78e4feSLisandro Dalcin }
954de78e4feSLisandro Dalcin 
95596ca5757SLisandro Dalcin 
95696ca5757SLisandro Dalcin PETSC_STATIC_INLINE DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt ctGmsh)
95796ca5757SLisandro Dalcin {
95896ca5757SLisandro Dalcin   switch (ctGmsh) {
95996ca5757SLisandro Dalcin     case 15:
96096ca5757SLisandro Dalcin       return DM_POLYTOPE_POINT;
96196ca5757SLisandro Dalcin     case 1:
96296ca5757SLisandro Dalcin     case 8:
96396ca5757SLisandro Dalcin       return DM_POLYTOPE_SEGMENT;
96496ca5757SLisandro Dalcin     case 2:
96596ca5757SLisandro Dalcin     case 9:
96696ca5757SLisandro Dalcin       return DM_POLYTOPE_TRIANGLE;
96796ca5757SLisandro Dalcin     case 3:
96896ca5757SLisandro Dalcin     case 10:
96996ca5757SLisandro Dalcin       return DM_POLYTOPE_QUADRILATERAL;
97096ca5757SLisandro Dalcin     case 4:
97196ca5757SLisandro Dalcin     case 11:
97296ca5757SLisandro Dalcin       return DM_POLYTOPE_TETRAHEDRON;
97396ca5757SLisandro Dalcin     case 5:
97496ca5757SLisandro Dalcin     case 12:
97596ca5757SLisandro Dalcin       return DM_POLYTOPE_HEXAHEDRON;
97696ca5757SLisandro Dalcin     case 6:
97796ca5757SLisandro Dalcin     case 13:
97896ca5757SLisandro Dalcin       return DM_POLYTOPE_TRI_PRISM;
97996ca5757SLisandro Dalcin     case 7:
98096ca5757SLisandro Dalcin     case 14:
98196ca5757SLisandro Dalcin       return DM_POLYTOPE_UNKNOWN; /* Pyramid */
98296ca5757SLisandro Dalcin     default:
98396ca5757SLisandro Dalcin       return DM_POLYTOPE_UNKNOWN;
98496ca5757SLisandro Dalcin   }
98596ca5757SLisandro Dalcin }
98696ca5757SLisandro Dalcin 
987d6689ca9SLisandro Dalcin /*@C
988d6689ca9SLisandro Dalcin   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
989d6689ca9SLisandro Dalcin 
990d6689ca9SLisandro Dalcin + comm        - The MPI communicator
991d6689ca9SLisandro Dalcin . filename    - Name of the Gmsh file
992d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh
993d6689ca9SLisandro Dalcin 
994d6689ca9SLisandro Dalcin   Output Parameter:
995d6689ca9SLisandro Dalcin . dm  - The DM object representing the mesh
996d6689ca9SLisandro Dalcin 
997d6689ca9SLisandro Dalcin   Level: beginner
998d6689ca9SLisandro Dalcin 
999d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
1000d6689ca9SLisandro Dalcin @*/
1001d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1002d6689ca9SLisandro Dalcin {
1003d6689ca9SLisandro Dalcin   PetscViewer     viewer;
1004d6689ca9SLisandro Dalcin   PetscMPIInt     rank;
1005d6689ca9SLisandro Dalcin   int             fileType;
1006d6689ca9SLisandro Dalcin   PetscViewerType vtype;
1007d6689ca9SLisandro Dalcin   PetscErrorCode  ierr;
1008d6689ca9SLisandro Dalcin 
1009d6689ca9SLisandro Dalcin   PetscFunctionBegin;
1010d6689ca9SLisandro Dalcin   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1011d6689ca9SLisandro Dalcin 
1012d6689ca9SLisandro Dalcin   /* Determine Gmsh file type (ASCII or binary) from file header */
1013d6689ca9SLisandro Dalcin   if (!rank) {
1014d6689ca9SLisandro Dalcin     GmshFile    gmsh_, *gmsh = &gmsh_;
1015d6689ca9SLisandro Dalcin     char        line[PETSC_MAX_PATH_LEN];
1016d6689ca9SLisandro Dalcin     int         snum;
1017d6689ca9SLisandro Dalcin     float       version;
1018d6689ca9SLisandro Dalcin 
1019580bdb30SBarry Smith     ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr);
1020d6689ca9SLisandro Dalcin     ierr = PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer);CHKERRQ(ierr);
1021d6689ca9SLisandro Dalcin     ierr = PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
1022d6689ca9SLisandro Dalcin     ierr = PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ);CHKERRQ(ierr);
1023d6689ca9SLisandro Dalcin     ierr = PetscViewerFileSetName(gmsh->viewer, filename);CHKERRQ(ierr);
1024d6689ca9SLisandro Dalcin     /* Read only the first two lines of the Gmsh file */
1025d6689ca9SLisandro Dalcin     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1026d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr);
1027d6689ca9SLisandro Dalcin     ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr);
1028d6689ca9SLisandro Dalcin     snum = sscanf(line, "%f %d", &version, &fileType);
1029d6689ca9SLisandro Dalcin     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
1030d6689ca9SLisandro 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);
1031d6689ca9SLisandro Dalcin     if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version %3.1f not supported", (double)version);
1032d6689ca9SLisandro 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);
1033d6689ca9SLisandro Dalcin     ierr = PetscViewerDestroy(&gmsh->viewer);CHKERRQ(ierr);
1034d6689ca9SLisandro Dalcin   }
1035d6689ca9SLisandro Dalcin   ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
1036d6689ca9SLisandro Dalcin   vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY;
1037d6689ca9SLisandro Dalcin 
1038d6689ca9SLisandro Dalcin   /* Create appropriate viewer and build plex */
1039d6689ca9SLisandro Dalcin   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
1040d6689ca9SLisandro Dalcin   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
1041d6689ca9SLisandro Dalcin   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
1042d6689ca9SLisandro Dalcin   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
1043d6689ca9SLisandro Dalcin   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
1044d6689ca9SLisandro Dalcin   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1045d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
1046d6689ca9SLisandro Dalcin }
1047d6689ca9SLisandro Dalcin 
1048331830f3SMatthew G. Knepley /*@
10497d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
1050331830f3SMatthew G. Knepley 
1051d083f849SBarry Smith   Collective
1052331830f3SMatthew G. Knepley 
1053331830f3SMatthew G. Knepley   Input Parameters:
1054331830f3SMatthew G. Knepley + comm  - The MPI communicator
1055331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
1056331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
1057331830f3SMatthew G. Knepley 
1058331830f3SMatthew G. Knepley   Output Parameter:
1059331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
1060331830f3SMatthew G. Knepley 
1061698a79b9SLisandro Dalcin   Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format
1062331830f3SMatthew G. Knepley 
1063331830f3SMatthew G. Knepley   Level: beginner
1064331830f3SMatthew G. Knepley 
1065331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
1066331830f3SMatthew G. Knepley @*/
1067331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
1068331830f3SMatthew G. Knepley {
106911c56cb0SLisandro Dalcin   PetscViewer    parentviewer = NULL;
107011c56cb0SLisandro Dalcin   double        *coordsIn = NULL;
1071730356f6SLisandro Dalcin   GmshEntities  *entities = NULL;
10723c67d5adSSatish Balay   GmshElement   *gmsh_elem = NULL;
1073331830f3SMatthew G. Knepley   PetscSection   coordSection;
1074331830f3SMatthew G. Knepley   Vec            coordinates;
10756fbe17bfSStefano Zampini   PetscBT        periodicV = NULL, periodicC = NULL;
107684572febSToby Isaac   PetscScalar   *coords;
1077412e9a14SMatthew G. Knepley   PetscInt       dim = 0, embedDim = -1, coordSize, c, v, d, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL;
1078698a79b9SLisandro Dalcin   PetscInt       numVertices = 0, numCells = 0, trueNumCells = 0;
1079698a79b9SLisandro Dalcin   int            i, shift = 1;
108011c56cb0SLisandro Dalcin   PetscMPIInt    rank;
1081ef12879bSLisandro Dalcin   PetscBool      binary, zerobase = PETSC_FALSE, usemarker = PETSC_FALSE;
10820db3fc9eSLisandro Dalcin   PetscBool      enable_hybrid = interpolate, periodic = PETSC_TRUE;
108396ca5757SLisandro Dalcin   PetscBool      hasTetra = PETSC_FALSE;
1084331830f3SMatthew G. Knepley   PetscErrorCode ierr;
1085331830f3SMatthew G. Knepley 
1086331830f3SMatthew G. Knepley   PetscFunctionBegin;
1087331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1088ef12879bSLisandro Dalcin   ierr = PetscObjectOptionsBegin((PetscObject)viewer);CHKERRQ(ierr);
1089ef12879bSLisandro Dalcin   ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Gmsh options");CHKERRQ(ierr);
1090ef12879bSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", enable_hybrid, &enable_hybrid, NULL);CHKERRQ(ierr);
1091ef12879bSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL);CHKERRQ(ierr);
1092ef12879bSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL);CHKERRQ(ierr);
1093ef12879bSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_zero_base", "Read Gmsh file with zero base indices", "DMPlexCreateGmsh", zerobase, &zerobase, NULL);CHKERRQ(ierr);
10945a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", embedDim, &embedDim, NULL,PETSC_DECIDE);CHKERRQ(ierr);
1095ef12879bSLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
1096ef12879bSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
109711c56cb0SLisandro Dalcin   if (zerobase) shift = 0;
109811c56cb0SLisandro Dalcin 
1099331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1100331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
11013b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
110211c56cb0SLisandro Dalcin 
110311c56cb0SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
110411c56cb0SLisandro Dalcin 
110511c56cb0SLisandro Dalcin   /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
11063b46f529SLisandro Dalcin   if (binary) {
110711c56cb0SLisandro Dalcin     parentviewer = viewer;
110811c56cb0SLisandro Dalcin     ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
110911c56cb0SLisandro Dalcin   }
111011c56cb0SLisandro Dalcin 
11113b46f529SLisandro Dalcin   if (!rank) {
1112698a79b9SLisandro Dalcin     GmshFile  gmsh_, *gmsh = &gmsh_;
1113698a79b9SLisandro Dalcin     char      line[PETSC_MAX_PATH_LEN];
1114698a79b9SLisandro Dalcin     PetscBool match;
1115331830f3SMatthew G. Knepley 
1116580bdb30SBarry Smith     ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr);
1117698a79b9SLisandro Dalcin     gmsh->viewer = viewer;
1118698a79b9SLisandro Dalcin     gmsh->binary = binary;
1119698a79b9SLisandro Dalcin 
1120698a79b9SLisandro Dalcin     /* Read mesh format */
1121d6689ca9SLisandro Dalcin     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1122d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr);
1123698a79b9SLisandro Dalcin     ierr = DMPlexCreateGmsh_ReadMeshFormat(gmsh);CHKERRQ(ierr);
1124d6689ca9SLisandro Dalcin     ierr = GmshReadEndSection(gmsh, "$EndMeshFormat", line);CHKERRQ(ierr);
112511c56cb0SLisandro Dalcin 
1126dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
1127d6689ca9SLisandro Dalcin     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1128d6689ca9SLisandro Dalcin     ierr = GmshMatch(gmsh,"$PhysicalNames", line, &match);CHKERRQ(ierr);
1129dc0ede3bSMatthew G. Knepley     if (match) {
1130698a79b9SLisandro Dalcin       ierr = DMPlexCreateGmsh_ReadPhysicalNames(gmsh);CHKERRQ(ierr);
1131d6689ca9SLisandro Dalcin       ierr = GmshReadEndSection(gmsh, "$EndPhysicalNames", line);CHKERRQ(ierr);
1132698a79b9SLisandro Dalcin       /* Initial read for entity section */
1133d6689ca9SLisandro Dalcin       ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1134dc0ede3bSMatthew G. Knepley     }
113511c56cb0SLisandro Dalcin 
1136de78e4feSLisandro Dalcin     /* Read entities */
1137698a79b9SLisandro Dalcin     if (gmsh->fileFormat >= 40) {
1138d6689ca9SLisandro Dalcin       ierr = GmshExpect(gmsh, "$Entities", line);CHKERRQ(ierr);
1139698a79b9SLisandro Dalcin       switch (gmsh->fileFormat) {
1140698a79b9SLisandro Dalcin       case 41: ierr = DMPlexCreateGmsh_ReadEntities_v41(gmsh, &entities);CHKERRQ(ierr); break;
1141698a79b9SLisandro Dalcin       default: ierr = DMPlexCreateGmsh_ReadEntities_v40(gmsh, &entities);CHKERRQ(ierr); break;
1142698a79b9SLisandro Dalcin       }
1143d6689ca9SLisandro Dalcin       ierr = GmshReadEndSection(gmsh, "$EndEntities", line);CHKERRQ(ierr);
1144698a79b9SLisandro Dalcin       /* Initial read for nodes section */
1145d6689ca9SLisandro Dalcin       ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1146de78e4feSLisandro Dalcin     }
1147de78e4feSLisandro Dalcin 
1148698a79b9SLisandro Dalcin     /* Read nodes */
1149d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$Nodes", line);CHKERRQ(ierr);
1150698a79b9SLisandro Dalcin     switch (gmsh->fileFormat) {
1151698a79b9SLisandro Dalcin     case 41: ierr = DMPlexCreateGmsh_ReadNodes_v41(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break;
1152698a79b9SLisandro Dalcin     case 40: ierr = DMPlexCreateGmsh_ReadNodes_v40(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break;
1153698a79b9SLisandro Dalcin     default: ierr = DMPlexCreateGmsh_ReadNodes_v22(gmsh, shift, &numVertices, &coordsIn);CHKERRQ(ierr); break;
1154de78e4feSLisandro Dalcin     }
1155d6689ca9SLisandro Dalcin     ierr = GmshReadEndSection(gmsh, "$EndNodes", line);CHKERRQ(ierr);
115611c56cb0SLisandro Dalcin 
1157698a79b9SLisandro Dalcin     /* Read elements */
1158feb237baSPierre Jolivet     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1159d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$Elements", line);CHKERRQ(ierr);
1160698a79b9SLisandro Dalcin     switch (gmsh->fileFormat) {
1161698a79b9SLisandro Dalcin     case 41: ierr = DMPlexCreateGmsh_ReadElements_v41(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break;
1162698a79b9SLisandro Dalcin     case 40: ierr = DMPlexCreateGmsh_ReadElements_v40(gmsh, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); break;
1163d6689ca9SLisandro Dalcin     default: ierr = DMPlexCreateGmsh_ReadElements_v22(gmsh, shift, &numCells, &gmsh_elem);CHKERRQ(ierr); break;
1164de78e4feSLisandro Dalcin     }
1165d6689ca9SLisandro Dalcin     ierr = GmshReadEndSection(gmsh, "$EndElements", line);CHKERRQ(ierr);
1166de78e4feSLisandro Dalcin 
1167698a79b9SLisandro Dalcin     /* OPTIONAL Read periodic section */
1168698a79b9SLisandro Dalcin     if (periodic) {
1169ef12879bSLisandro Dalcin       ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1170ef12879bSLisandro Dalcin       ierr = GmshMatch(gmsh, "$Periodic", line, &periodic);CHKERRQ(ierr);
1171ef12879bSLisandro Dalcin     }
1172ef12879bSLisandro Dalcin     if (periodic) {
1173698a79b9SLisandro Dalcin       PetscInt pVert, *periodicMapT, *aux;
1174de78e4feSLisandro Dalcin 
1175698a79b9SLisandro Dalcin       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
1176698a79b9SLisandro Dalcin       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
1177698a79b9SLisandro Dalcin       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
1178698a79b9SLisandro Dalcin 
1179d6689ca9SLisandro Dalcin       ierr = GmshExpect(gmsh, "$Periodic", line);CHKERRQ(ierr);
1180698a79b9SLisandro Dalcin       switch (gmsh->fileFormat) {
1181698a79b9SLisandro Dalcin       case 41: ierr = DMPlexCreateGmsh_ReadPeriodic_v41(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break;
1182698a79b9SLisandro Dalcin       default: ierr = DMPlexCreateGmsh_ReadPeriodic_v40(gmsh, shift, periodicMapT, periodicV);CHKERRQ(ierr); break;
1183698a79b9SLisandro Dalcin       }
1184d6689ca9SLisandro Dalcin       ierr = GmshReadEndSection(gmsh, "$EndPeriodic", line);CHKERRQ(ierr);
1185698a79b9SLisandro Dalcin 
1186698a79b9SLisandro Dalcin       /* we may have slaves of slaves */
1187698a79b9SLisandro Dalcin       for (i = 0; i < numVertices; i++) {
1188698a79b9SLisandro Dalcin         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
1189698a79b9SLisandro Dalcin           periodicMapT[i] = periodicMapT[periodicMapT[i]];
1190698a79b9SLisandro Dalcin         }
1191698a79b9SLisandro Dalcin       }
1192698a79b9SLisandro Dalcin       /* periodicMap : from old to new numbering (periodic vertices excluded)
1193698a79b9SLisandro Dalcin          periodicMapI: from new to old numbering */
1194698a79b9SLisandro Dalcin       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
1195698a79b9SLisandro Dalcin       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
1196698a79b9SLisandro Dalcin       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
1197698a79b9SLisandro Dalcin       for (i = 0, pVert = 0; i < numVertices; i++) {
1198698a79b9SLisandro Dalcin         if (periodicMapT[i] != i) {
1199698a79b9SLisandro Dalcin           pVert++;
1200698a79b9SLisandro Dalcin         } else {
1201698a79b9SLisandro Dalcin           aux[i] = i - pVert;
1202698a79b9SLisandro Dalcin           periodicMapI[i - pVert] = i;
1203698a79b9SLisandro Dalcin         }
1204698a79b9SLisandro Dalcin       }
1205698a79b9SLisandro Dalcin       for (i = 0 ; i < numVertices; i++) {
1206698a79b9SLisandro Dalcin         periodicMap[i] = aux[periodicMapT[i]];
1207698a79b9SLisandro Dalcin       }
1208698a79b9SLisandro Dalcin       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
1209698a79b9SLisandro Dalcin       ierr = PetscFree(aux);CHKERRQ(ierr);
1210698a79b9SLisandro Dalcin       /* remove periodic vertices */
1211698a79b9SLisandro Dalcin       numVertices = numVertices - pVert;
1212698a79b9SLisandro Dalcin     }
1213698a79b9SLisandro Dalcin 
1214698a79b9SLisandro Dalcin     ierr = GmshEntitiesDestroy(&entities);CHKERRQ(ierr);
1215698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr);
1216698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr);
1217*33a088b6SMatthew G. Knepley     ierr = PetscFree(gmsh->nodeMap);CHKERRQ(ierr);
1218698a79b9SLisandro Dalcin   }
1219698a79b9SLisandro Dalcin 
1220698a79b9SLisandro Dalcin   if (parentviewer) {
1221698a79b9SLisandro Dalcin     ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
1222698a79b9SLisandro Dalcin   }
1223698a79b9SLisandro Dalcin 
1224698a79b9SLisandro Dalcin   if (!rank) {
1225698a79b9SLisandro Dalcin     PetscBool hybrid   = PETSC_FALSE;
12260db3fc9eSLisandro Dalcin     PetscInt  cellType = -1;
1227698a79b9SLisandro Dalcin 
1228a4bb7517SMichael Lange     for (trueNumCells = 0, c = 0; c < numCells; ++c) {
12290db3fc9eSLisandro Dalcin       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0; hybrid = PETSC_FALSE; cellType = -1;}
12300db3fc9eSLisandro Dalcin       if (gmsh_elem[c].dim < dim) continue;
12310db3fc9eSLisandro Dalcin       if (cellType == -1) cellType = gmsh_elem[c].cellType;
12320db3fc9eSLisandro Dalcin       /* different cell type indicate an hybrid mesh in PLEX */
12330db3fc9eSLisandro Dalcin       if (cellType != gmsh_elem[c].cellType) hybrid = PETSC_TRUE;
1234dea2a3c7SStefano Zampini       /* wedges always indicate an hybrid mesh in PLEX */
12350db3fc9eSLisandro Dalcin       if (cellType == 6 || cellType == 13) hybrid = PETSC_TRUE;
123696ca5757SLisandro Dalcin       if (cellType == 4 || cellType == 11) hasTetra = PETSC_TRUE;
12370db3fc9eSLisandro Dalcin       trueNumCells++;
1238db397164SMichael Lange     }
1239dea2a3c7SStefano Zampini     /* Renumber cells for hybrid grids */
1240dea2a3c7SStefano Zampini     if (hybrid && enable_hybrid) {
1241dea2a3c7SStefano Zampini       PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL;
1242dea2a3c7SStefano Zampini       PetscInt cell, tn, *tp;
1243dea2a3c7SStefano Zampini       int      n1 = 0,n2 = 0;
1244dea2a3c7SStefano Zampini 
1245dea2a3c7SStefano Zampini       ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr);
1246dea2a3c7SStefano Zampini       ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr);
1247dea2a3c7SStefano Zampini       for (cell = 0, c = 0; c < numCells; ++c) {
1248dea2a3c7SStefano Zampini         if (gmsh_elem[c].dim == dim) {
1249dea2a3c7SStefano Zampini           if (!n1) n1 = gmsh_elem[c].cellType;
1250dea2a3c7SStefano Zampini           else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType;
1251dea2a3c7SStefano Zampini 
1252dea2a3c7SStefano Zampini           if      (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell;
1253dea2a3c7SStefano Zampini           else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell;
1254dea2a3c7SStefano Zampini           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types");
1255dea2a3c7SStefano Zampini           cell++;
1256dea2a3c7SStefano Zampini         }
1257dea2a3c7SStefano Zampini       }
1258dea2a3c7SStefano Zampini 
1259dea2a3c7SStefano Zampini       switch (n1) {
1260dea2a3c7SStefano Zampini       case 2: /* triangles */
1261dea2a3c7SStefano Zampini       case 9:
1262dea2a3c7SStefano Zampini         switch (n2) {
1263dea2a3c7SStefano Zampini         case 0: /* single type mesh */
1264dea2a3c7SStefano Zampini         case 3: /* quads */
1265dea2a3c7SStefano Zampini         case 10:
1266dea2a3c7SStefano Zampini           break;
1267dea2a3c7SStefano Zampini         default:
1268dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
1269dea2a3c7SStefano Zampini         }
1270dea2a3c7SStefano Zampini         break;
1271a5b208b6SMatthew G. Knepley       case 3: /* quadrilateral */
1272dea2a3c7SStefano Zampini       case 10:
1273dea2a3c7SStefano Zampini         switch (n2) {
1274dea2a3c7SStefano Zampini         case 0: /* single type mesh */
1275dea2a3c7SStefano Zampini         case 2: /* swap since we list simplices first */
1276dea2a3c7SStefano Zampini         case 9:
1277dea2a3c7SStefano Zampini           tn  = hc1;
1278dea2a3c7SStefano Zampini           hc1 = hc2;
1279dea2a3c7SStefano Zampini           hc2 = tn;
1280dea2a3c7SStefano Zampini 
1281dea2a3c7SStefano Zampini           tp           = hybridCells1;
1282dea2a3c7SStefano Zampini           hybridCells1 = hybridCells2;
1283dea2a3c7SStefano Zampini           hybridCells2 = tp;
1284dea2a3c7SStefano Zampini           break;
1285dea2a3c7SStefano Zampini         default:
1286dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
1287dea2a3c7SStefano Zampini         }
1288dea2a3c7SStefano Zampini         break;
1289dea2a3c7SStefano Zampini       case 4: /* tetrahedra */
1290dea2a3c7SStefano Zampini       case 11:
1291dea2a3c7SStefano Zampini         switch (n2) {
1292dea2a3c7SStefano Zampini         case 0: /* single type mesh */
1293dea2a3c7SStefano Zampini         case 6: /* wedges */
1294dea2a3c7SStefano Zampini         case 13:
1295dea2a3c7SStefano Zampini           break;
1296dea2a3c7SStefano Zampini         default:
1297dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
1298dea2a3c7SStefano Zampini         }
1299dea2a3c7SStefano Zampini         break;
1300a5b208b6SMatthew G. Knepley       case 5: /* hexahedra */
1301a5b208b6SMatthew G. Knepley       case 12:
1302a5b208b6SMatthew G. Knepley         switch (n2) {
1303a5b208b6SMatthew G. Knepley         case 0: /* single type mesh */
1304a5b208b6SMatthew G. Knepley         case 6: /* wedges */
1305a5b208b6SMatthew G. Knepley         case 13:
1306a5b208b6SMatthew G. Knepley           break;
1307a5b208b6SMatthew G. Knepley         default:
1308a5b208b6SMatthew G. Knepley           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
1309a5b208b6SMatthew G. Knepley         }
1310a5b208b6SMatthew G. Knepley         break;
1311a5b208b6SMatthew G. Knepley       case 6: /* wedge */
1312dea2a3c7SStefano Zampini       case 13:
1313dea2a3c7SStefano Zampini         switch (n2) {
1314dea2a3c7SStefano Zampini         case 0: /* single type mesh */
1315a5b208b6SMatthew G. Knepley         case 4: /* tetrahedra: swap since we list simplices first */
1316dea2a3c7SStefano Zampini         case 11:
1317a5b208b6SMatthew G. Knepley         case 5: /* hexahedra */
1318a5b208b6SMatthew G. Knepley         case 12:
1319dea2a3c7SStefano Zampini           tn  = hc1;
1320dea2a3c7SStefano Zampini           hc1 = hc2;
1321dea2a3c7SStefano Zampini           hc2 = tn;
1322dea2a3c7SStefano Zampini 
1323dea2a3c7SStefano Zampini           tp           = hybridCells1;
1324dea2a3c7SStefano Zampini           hybridCells1 = hybridCells2;
1325dea2a3c7SStefano Zampini           hybridCells2 = tp;
1326dea2a3c7SStefano Zampini           break;
1327dea2a3c7SStefano Zampini         default:
1328dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
1329dea2a3c7SStefano Zampini         }
1330dea2a3c7SStefano Zampini         break;
1331dea2a3c7SStefano Zampini       default:
1332dea2a3c7SStefano Zampini         SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
1333dea2a3c7SStefano Zampini       }
1334dea2a3c7SStefano Zampini       ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr);
1335dea2a3c7SStefano Zampini       for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell;
1336dea2a3c7SStefano Zampini       for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1;
1337dea2a3c7SStefano Zampini       ierr = PetscFree(hybridCells1);CHKERRQ(ierr);
1338dea2a3c7SStefano Zampini       ierr = PetscFree(hybridCells2);CHKERRQ(ierr);
1339dea2a3c7SStefano Zampini     }
1340dea2a3c7SStefano Zampini 
134111c56cb0SLisandro Dalcin   }
134211c56cb0SLisandro Dalcin 
1343a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
1344412e9a14SMatthew G. Knepley   /*   We do not want this label automatically computed, instead we compute it here */
1345412e9a14SMatthew G. Knepley   ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr);
1346db397164SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
1347a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
1348a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
134996ca5757SLisandro Dalcin       PetscInt       ucell = hybridMap ? hybridMap[cell] : cell;
135096ca5757SLisandro Dalcin       DMPolytopeType ctype = DMPolytopeTypeFromGmsh(gmsh_elem[c].cellType);
135196ca5757SLisandro Dalcin       if (hybridMap && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR;
135296ca5757SLisandro Dalcin       ierr = DMPlexSetConeSize(*dm, ucell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
135396ca5757SLisandro Dalcin       ierr = DMPlexSetCellType(*dm, ucell, ctype);CHKERRQ(ierr);
1354a4bb7517SMichael Lange       cell++;
1355db397164SMichael Lange     }
1356331830f3SMatthew G. Knepley   }
135796ca5757SLisandro Dalcin   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
135896ca5757SLisandro Dalcin     ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr);
135996ca5757SLisandro Dalcin   }
1360331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
136196ca5757SLisandro Dalcin 
1362a4bb7517SMichael Lange   /* Add cell-vertex connections */
1363a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
1364a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
136511c56cb0SLisandro Dalcin       PetscInt pcone[8], corner;
136696ca5757SLisandro Dalcin       PetscInt ucell = hybridMap ? hybridMap[cell] : cell;
1367a4bb7517SMichael Lange       for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
1368dd169d64SMatthew G. Knepley         const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
1369917266a4SLisandro Dalcin         pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells;
1370db397164SMichael Lange       }
137196ca5757SLisandro Dalcin       ierr = DMPlexReorderCell(*dm, ucell, pcone);CHKERRQ(ierr);
137296ca5757SLisandro Dalcin       ierr = DMPlexSetCone(*dm, ucell, pcone);CHKERRQ(ierr);
1373a4bb7517SMichael Lange       cell++;
1374331830f3SMatthew G. Knepley     }
1375a4bb7517SMichael Lange   }
137696ca5757SLisandro Dalcin 
13776228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1378c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
1379331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
1380331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
1381331830f3SMatthew G. Knepley   if (interpolate) {
13825fd9971aSMatthew G. Knepley     DM idm;
1383331830f3SMatthew G. Knepley 
1384331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1385331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
1386331830f3SMatthew G. Knepley     *dm  = idm;
1387331830f3SMatthew G. Knepley   }
1388917266a4SLisandro Dalcin 
1389f6c8a31fSLisandro Dalcin   if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation");
1390917266a4SLisandro Dalcin   if (!rank && usemarker) {
1391d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
1392d3f73514SStefano Zampini 
1393d3f73514SStefano Zampini     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
1394d3f73514SStefano Zampini     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1395d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
1396d3f73514SStefano Zampini       PetscInt suppSize;
1397d3f73514SStefano Zampini 
1398d3f73514SStefano Zampini       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
1399d3f73514SStefano Zampini       if (suppSize == 1) {
1400d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
1401d3f73514SStefano Zampini 
1402d3f73514SStefano Zampini         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
1403d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
1404d3f73514SStefano Zampini           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
1405d3f73514SStefano Zampini         }
1406d3f73514SStefano Zampini         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
1407d3f73514SStefano Zampini       }
1408d3f73514SStefano Zampini     }
1409d3f73514SStefano Zampini   }
141016ad7e67SMichael Lange 
141116ad7e67SMichael Lange   if (!rank) {
141211c56cb0SLisandro Dalcin     PetscInt vStart, vEnd;
1413d1a54cefSMatthew G. Knepley 
141416ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
141511c56cb0SLisandro Dalcin     for (cell = 0, c = 0; c < numCells; ++c) {
141611c56cb0SLisandro Dalcin 
141711c56cb0SLisandro Dalcin       /* Create face sets */
14185964b3dcSLisandro Dalcin       if (interpolate && gmsh_elem[c].dim == dim-1) {
141916ad7e67SMichael Lange         const PetscInt *join;
142011c56cb0SLisandro Dalcin         PetscInt        joinSize, pcone[8], corner;
142111c56cb0SLisandro Dalcin         /* Find the relevant facet with vertex joins */
1422a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
1423dd169d64SMatthew G. Knepley           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
1424917266a4SLisandro Dalcin           pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart;
142516ad7e67SMichael Lange         }
142611c56cb0SLisandro Dalcin         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr);
1427f10f1c67SMatthew 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);
1428c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
1429a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
143016ad7e67SMichael Lange       }
14316e1dd89aSLawrence Mitchell 
14326e1dd89aSLawrence Mitchell       /* Create cell sets */
14336e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
14346e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
1435dea2a3c7SStefano Zampini           ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
14366e1dd89aSLawrence Mitchell         }
1437917266a4SLisandro Dalcin         cell++;
14386e1dd89aSLawrence Mitchell       }
14390c070f12SSander Arens 
14400c070f12SSander Arens       /* Create vertex sets */
14410c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
14420c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
1443917266a4SLisandro Dalcin           const PetscInt cc = gmsh_elem[c].nodes[0] - shift;
144411c56cb0SLisandro Dalcin           const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart;
1445d08df55aSStefano Zampini           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
14460c070f12SSander Arens         }
14470c070f12SSander Arens       }
14480c070f12SSander Arens     }
144916ad7e67SMichael Lange   }
145016ad7e67SMichael Lange 
145111c56cb0SLisandro Dalcin   /* Create coordinates */
1452dea2a3c7SStefano Zampini   if (embedDim < 0) embedDim = dim;
1453dea2a3c7SStefano Zampini   ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr);
1454979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1455331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
14561d64f2ccSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
1457f45c9589SStefano Zampini   if (periodicMap) { /* we need to localize coordinates on cells */
1458f45c9589SStefano Zampini     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
1459f45c9589SStefano Zampini   } else {
146075b5763bSMichael Lange     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
1461f45c9589SStefano Zampini   }
146275b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
14631d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
14641d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
1465331830f3SMatthew G. Knepley   }
146611c56cb0SLisandro Dalcin   if (periodicMap) {
1467437bdd5bSStefano Zampini     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
1468f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
1469f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
147096ca5757SLisandro Dalcin         PetscInt  corner, pc = PETSC_FALSE;
1471437bdd5bSStefano Zampini         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
1472917266a4SLisandro Dalcin           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift));
1473437bdd5bSStefano Zampini         }
1474437bdd5bSStefano Zampini         if (pc) {
147511c56cb0SLisandro Dalcin           PetscInt dof = gmsh_elem[c].numNodes*embedDim;
1476dea2a3c7SStefano Zampini           PetscInt ucell = hybridMap ? hybridMap[cell] : cell;
1477dea2a3c7SStefano Zampini           ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr);
1478dea2a3c7SStefano Zampini           ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr);
1479dea2a3c7SStefano Zampini           ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr);
14806fbe17bfSStefano Zampini         }
1481f45c9589SStefano Zampini         cell++;
1482f45c9589SStefano Zampini       }
1483f45c9589SStefano Zampini     }
1484f45c9589SStefano Zampini   }
1485331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1486331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
14878b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
1488331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1489331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
14901d64f2ccSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
1491331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
1492331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1493f45c9589SStefano Zampini   if (periodicMap) {
1494f45c9589SStefano Zampini     PetscInt off;
1495f45c9589SStefano Zampini 
1496f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
1497f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
149896ca5757SLisandro Dalcin         PetscInt pcone[8], corner;
1499dea2a3c7SStefano Zampini         PetscInt ucell = hybridMap ? hybridMap[cell] : cell;
1500dea2a3c7SStefano Zampini         if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) {
1501f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
1502917266a4SLisandro Dalcin             pcone[corner] = gmsh_elem[c].nodes[corner] - shift;
1503f45c9589SStefano Zampini           }
150496ca5757SLisandro Dalcin           ierr = DMPlexReorderCell(*dm, ucell, pcone);CHKERRQ(ierr);
1505dea2a3c7SStefano Zampini           ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr);
150696ca5757SLisandro Dalcin           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner)
150796ca5757SLisandro Dalcin             for (v = pcone[corner], d = 0; d < embedDim; ++d)
150811c56cb0SLisandro Dalcin               coords[off++] = (PetscReal) coordsIn[v*3+d];
1509f45c9589SStefano Zampini         }
1510f45c9589SStefano Zampini         cell++;
1511f45c9589SStefano Zampini       }
1512f45c9589SStefano Zampini     }
1513f45c9589SStefano Zampini     for (v = 0; v < numVertices; ++v) {
1514f45c9589SStefano Zampini       ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
1515dd169d64SMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
151611c56cb0SLisandro Dalcin         coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d];
1517f45c9589SStefano Zampini       }
1518f45c9589SStefano Zampini     }
1519f45c9589SStefano Zampini   } else {
1520331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
15211d64f2ccSMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
152211c56cb0SLisandro Dalcin         coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d];
1523331830f3SMatthew G. Knepley       }
1524331830f3SMatthew G. Knepley     }
1525331830f3SMatthew G. Knepley   }
1526331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1527331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
1528ef12879bSLisandro Dalcin 
1529ef12879bSLisandro Dalcin   periodic = periodicMap ? PETSC_TRUE : PETSC_FALSE;
1530ef12879bSLisandro Dalcin   ierr = MPI_Bcast(&periodic, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
153111c56cb0SLisandro Dalcin   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
153211c56cb0SLisandro Dalcin 
153311c56cb0SLisandro Dalcin   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
153411c56cb0SLisandro Dalcin   ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
1535dea2a3c7SStefano Zampini   ierr = PetscFree(hybridMap);CHKERRQ(ierr);
1536d08df55aSStefano Zampini   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
1537fcd9ca0aSStefano Zampini   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
15386fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
15396fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
154011c56cb0SLisandro Dalcin   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
154111c56cb0SLisandro Dalcin 
15423b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
1543331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
1544331830f3SMatthew G. Knepley }
1545