xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision d3f7351471e23353a9a945e50e663e7cbbe770dd)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 
4 /*@C
5   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
6 
7 + comm        - The MPI communicator
8 . filename    - Name of the Gmsh file
9 - interpolate - Create faces and edges in the mesh
10 
11   Output Parameter:
12 . dm  - The DM object representing the mesh
13 
14   Level: beginner
15 
16 .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
17 @*/
18 PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
19 {
20   PetscViewer     viewer, vheader;
21   PetscMPIInt     rank;
22   PetscViewerType vtype;
23   char            line[PETSC_MAX_PATH_LEN];
24   int             snum;
25   PetscBool       match;
26   int             fT;
27   PetscInt        fileType;
28   float           version;
29   PetscErrorCode  ierr;
30 
31   PetscFunctionBegin;
32   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
33   /* Determine Gmsh file type (ASCII or binary) from file header */
34   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
35   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
36   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
37   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
38   if (!rank) {
39     /* Read only the first two lines of the Gmsh file */
40     ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
41     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
42     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
43     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
44     snum = sscanf(line, "%f %d", &version, &fT);
45     fileType = (PetscInt) fT;
46     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
47     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
48   }
49   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
50   /* Create appropriate viewer and build plex */
51   if (fileType == 0) vtype = PETSCVIEWERASCII;
52   else vtype = PETSCVIEWERBINARY;
53   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
54   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
55   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
56   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
57   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
58   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
59   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
60   PetscFunctionReturn(0);
61 }
62 
63 static PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
64 {
65   PetscInt       c, p;
66   GmshElement   *elements;
67   int            i, cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
68   PetscInt       pibuf[64];
69   int            ibuf[16];
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
74   for (c = 0; c < numCells;) {
75     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
76     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
77     if (binary) {
78       cellType = ibuf[0];
79       numElem = ibuf[1];
80       numTags = ibuf[2];
81     } else {
82       elements[c].id = ibuf[0];
83       cellType = ibuf[1];
84       numTags = ibuf[2];
85       numElem = 1;
86     }
87     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
88     numNodesIgnore = 0;
89     switch (cellType) {
90     case 1: /* 2-node line */
91       dim = 1;
92       numNodes = 2;
93       break;
94     case 2: /* 3-node triangle */
95       dim = 2;
96       numNodes = 3;
97       break;
98     case 3: /* 4-node quadrangle */
99       dim = 2;
100       numNodes = 4;
101       break;
102     case 4: /* 4-node tetrahedron */
103       dim  = 3;
104       numNodes = 4;
105       break;
106     case 5: /* 8-node hexahedron */
107       dim = 3;
108       numNodes = 8;
109       break;
110     case 8: /* 3-node 2nd order line */
111       dim = 1;
112       numNodes = 2;
113       numNodesIgnore = 1;
114       break;
115     case 9: /* 6-node 2nd order triangle */
116       dim = 2;
117       numNodes = 3;
118       numNodesIgnore = 3;
119       break;
120     case 15: /* 1-node vertex */
121       dim = 0;
122       numNodes = 1;
123       break;
124     case 6: /* 6-node prism */
125     case 7: /* 5-node pyramid */
126     case 10: /* 9-node 2nd order quadrangle */
127     case 11: /* 10-node 2nd order tetrahedron */
128     case 12: /* 27-node 2nd order hexhedron */
129     case 13: /* 19-node 2nd order prism */
130     case 14: /* 14-node 2nd order pyramid */
131     default:
132       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
133     }
134     if (binary) {
135       const PetscInt nint = numNodes + numTags + 1 + numNodesIgnore;
136       for (i = 0; i < numElem; ++i, ++c) {
137         /* Loop over inner binary element block */
138         elements[c].dim = dim;
139         elements[c].numNodes = numNodes;
140         elements[c].numTags = numTags;
141 
142         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
143         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
144         elements[c].id = ibuf[0];
145         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
146         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
147       }
148     } else {
149       elements[c].dim = dim;
150       elements[c].numNodes = numNodes;
151       elements[c].numTags = numTags;
152       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
153       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
154       ierr = PetscViewerRead(viewer, pibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr);
155       c++;
156     }
157   }
158   *gmsh_elems = elements;
159   PetscFunctionReturn(0);
160 }
161 
162 /*@
163   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
164 
165   Collective on comm
166 
167   Input Parameters:
168 + comm  - The MPI communicator
169 . viewer - The Viewer associated with a Gmsh file
170 - interpolate - Create faces and edges in the mesh
171 
172   Output Parameter:
173 . dm  - The DM object representing the mesh
174 
175   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
176   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
177 
178   Level: beginner
179 
180 .keywords: mesh,Gmsh
181 .seealso: DMPLEX, DMCreate()
182 @*/
183 PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
184 {
185   PetscViewerType vtype;
186   GmshElement   *gmsh_elem = NULL;
187   PetscSection   coordSection;
188   Vec            coordinates;
189   PetscScalar   *coords, *coordsIn = NULL;
190   PetscInt       dim = 0, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL;
191   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum;
192   PetscMPIInt    num_proc, rank;
193   char           line[PETSC_MAX_PATH_LEN];
194   PetscBool      match, binary, bswap = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE;
195   PetscErrorCode ierr;
196 
197   PetscFunctionBegin;
198   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
199   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
200   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
201   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
202   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
203   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
204   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
205   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_periodic",&periodic,NULL);CHKERRQ(ierr);
206   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_usemarker",&usemarker,NULL);CHKERRQ(ierr);
207   if (!rank || binary) {
208     PetscBool match;
209     int       fileType, dataSize;
210     float     version;
211 
212     /* Read header */
213     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
214     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
215     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
216     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
217     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
218     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
219     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
220     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
221     if (binary) {
222       int checkInt;
223       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
224       if (checkInt != 1) {
225         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
226         if (checkInt == 1) bswap = PETSC_TRUE;
227         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
228       }
229     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
230     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
231     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
232     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
233     /* OPTIONAL Read physical names */
234     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
235     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
236     if (match) {
237       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
238       snum = sscanf(line, "%d", &numRegions);
239       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
240       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
241       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
242       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
243       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
244       /* Initial read for vertex section */
245       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
246     }
247     /* Read vertices */
248     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
249     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
250     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
251     snum = sscanf(line, "%d", &numVertices);
252     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
253     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
254     if (binary) {
255       size_t doubleSize, intSize;
256       PetscInt elementSize;
257       char *buffer;
258       PetscScalar *baseptr;
259       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
260       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
261       elementSize = (intSize + 3*doubleSize);
262       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
263       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
264       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
265       for (v = 0; v < numVertices; ++v) {
266         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
267         coordsIn[v*3+0] = baseptr[0];
268         coordsIn[v*3+1] = baseptr[1];
269         coordsIn[v*3+2] = baseptr[2];
270       }
271       ierr = PetscFree(buffer);CHKERRQ(ierr);
272     } else {
273       for (v = 0; v < numVertices; ++v) {
274         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
275         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
276         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
277       }
278     }
279     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
280     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
281     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
282     /* Read cells */
283     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
284     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
285     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
286     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
287     snum = sscanf(line, "%d", &numCells);
288     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
289   }
290 
291   if (!rank || binary) {
292     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
293        file contents multiple times to figure out the true number of cells and facets
294        in the given mesh. To make this more efficient we read the file contents only
295        once and store them in memory, while determining the true number of cells. */
296     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
297     for (trueNumCells=0, c = 0; c < numCells; ++c) {
298       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
299       if (gmsh_elem[c].dim == dim) trueNumCells++;
300     }
301     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
302     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
303     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
304     if (periodic) {
305       PetscInt pVert, *periodicMapT, *aux;
306       int      numPeriodic;
307 
308       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
309       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
310       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
311       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
312       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
313       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
314       snum = sscanf(line, "%d", &numPeriodic);
315       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
316       for (i = 0; i < numPeriodic; i++) {
317         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes;
318 
319         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
320         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
321         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
322         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
323         snum = sscanf(line, "%d", &nNodes);
324         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
325         for (j = 0; j < nNodes; j++) {
326           PetscInt slaveNode, masterNode;
327 
328           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
329           snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
330           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
331           periodicMapT[slaveNode - 1] = masterNode - 1;
332         }
333       }
334       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
335       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
336       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
337       /* we may have slaves of slaves */
338       for (i = 0; i < numVertices; i++) {
339         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
340           periodicMapT[i] = periodicMapT[periodicMapT[i]];
341         }
342       }
343       /* periodicMap : from old to new numbering (periodic vertices excluded)
344          periodicMapI: from new to old numbering */
345       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
346       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
347       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
348       for (i = 0, pVert = 0; i < numVertices; i++) {
349         if (periodicMapT[i] != i) {
350           pVert++;
351         } else {
352           aux[i] = i - pVert;
353           periodicMapI[i - pVert] = i;
354         }
355       }
356       for (i = 0 ; i < numVertices; i++) {
357         periodicMap[i] = aux[periodicMapT[i]];
358       }
359       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
360       ierr = PetscFree(aux);CHKERRQ(ierr);
361       /* remove periodic vertices */
362       numVertices = numVertices - pVert;
363     }
364   }
365   /* For binary we read on all ranks, but only build the plex on rank 0 */
366   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
367   /* Allocate the cell-vertex mesh */
368   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
369   if (!rank) {
370     for (cell = 0, c = 0; c < numCells; ++c) {
371       if (gmsh_elem[c].dim == dim) {
372         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
373         cell++;
374       }
375     }
376   }
377   ierr = DMSetUp(*dm);CHKERRQ(ierr);
378   /* Add cell-vertex connections */
379   if (!rank) {
380     PetscInt pcone[8], corner;
381     for (cell = 0, c = 0; c < numCells; ++c) {
382       if (gmsh_elem[c].dim == dim) {
383         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
384           pcone[corner] = periodicMap ? periodicMap[gmsh_elem[c].nodes[corner] - 1] + trueNumCells
385                                       : gmsh_elem[c].nodes[corner] -1 + trueNumCells;
386         }
387         if (dim == 3) {
388           /* Tetrahedra are inverted */
389           if (gmsh_elem[c].numNodes == 4) {
390             PetscInt tmp = pcone[0];
391             pcone[0] = pcone[1];
392             pcone[1] = tmp;
393           }
394           /* Hexahedra are inverted */
395           if (gmsh_elem[c].numNodes == 8) {
396             PetscInt tmp = pcone[1];
397             pcone[1] = pcone[3];
398             pcone[3] = tmp;
399           }
400         }
401         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
402         cell++;
403       }
404     }
405   }
406   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
407   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
408   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
409   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
410   if (interpolate) {
411     DM idm = NULL;
412 
413     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
414     ierr = DMDestroy(dm);CHKERRQ(ierr);
415     *dm  = idm;
416   }
417   if (usemarker) {
418     PetscInt f, fStart, fEnd;
419 
420     ierr = DMCreateLabel(*dm,"marker");CHKERRQ(ierr);
421     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
422     for (f = fStart; f < fEnd; ++f) {
423       PetscInt suppSize;
424 
425       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
426       if (suppSize == 1) {
427         PetscInt *cone = NULL, coneSize, p;
428 
429         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
430         for (p = 0; p < coneSize; p += 2) {
431           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
432         }
433         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
434       }
435     }
436   }
437 
438   if (!rank) {
439     /* Apply boundary IDs by finding the relevant facets with vertex joins */
440     PetscInt pcone[8], corner, vStart, vEnd;
441 
442     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
443     for (c = 0; c < numCells; ++c) {
444       if (gmsh_elem[c].dim == dim-1) {
445         PetscInt joinSize;
446         const PetscInt *join;
447         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
448           pcone[corner] = periodicMap ? periodicMap[gmsh_elem[c].nodes[corner] - 1] + vStart
449                                       : gmsh_elem[c].nodes[corner] - 1 + vStart;
450         }
451         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
452         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
453         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
454         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
455       }
456     }
457 
458     /* Create cell sets */
459     for (cell = 0, c = 0; c < numCells; ++c) {
460       if (gmsh_elem[c].dim == dim) {
461         if (gmsh_elem[c].numTags > 0) {
462           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
463           cell++;
464         }
465       }
466     }
467 
468     /* Create vertex sets */
469     for (c = 0; c < numCells; ++c) {
470       if (gmsh_elem[c].dim == 0) {
471         if (gmsh_elem[c].numTags > 0) {
472           PetscInt vid =  periodicMap ? periodicMap[gmsh_elem[c].nodes[0] - 1] + vStart
473                                       : gmsh_elem[c].nodes[0] - 1 + vStart;
474           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
475         }
476       }
477     }
478   }
479 
480   /* Read coordinates */
481   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
482   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
483   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
484   if (periodicMap) { /* we need to localize coordinates on cells */
485     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
486   } else {
487     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
488   }
489   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
490     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
491     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
492   }
493   if (periodicMap) {
494     for (cell = 0, c = 0; c < numCells; ++c) {
495       if (gmsh_elem[c].dim == dim) {
496         ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
497         ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
498         cell++;
499       }
500     }
501   }
502   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
503   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
504   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
505   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
506   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
507   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
508   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
509   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
510   if (!rank) {
511 
512     if (periodicMap) {
513       PetscInt off;
514 
515       for (cell = 0, c = 0; c < numCells; ++c) {
516         PetscInt pcone[8], corner;
517         if (gmsh_elem[c].dim == dim) {
518           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
519             pcone[corner] = gmsh_elem[c].nodes[corner] - 1;
520           }
521           if (dim == 3) {
522             /* Tetrahedra are inverted */
523             if (gmsh_elem[c].numNodes == 4) {
524               PetscInt tmp = pcone[0];
525               pcone[0] = pcone[1];
526               pcone[1] = tmp;
527             }
528             /* Hexahedra are inverted */
529             if (gmsh_elem[c].numNodes == 8) {
530               PetscInt tmp = pcone[1];
531               pcone[1] = pcone[3];
532               pcone[3] = tmp;
533             }
534           }
535           ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
536           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
537             PetscInt v = pcone[corner];
538             for (d = 0; d < dim; ++d) {
539               coords[off++] = coordsIn[v*3+d];
540             }
541           }
542           cell++;
543         }
544       }
545       for (v = 0; v < numVertices; ++v) {
546         ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
547         for (d = 0; d < dim; ++d) {
548           coords[off+d] = coordsIn[periodicMapI[v]*3+d];
549         }
550       }
551     } else {
552       for (v = 0; v < numVertices; ++v) {
553         for (d = 0; d < dim; ++d) {
554           coords[v*dim+d] = coordsIn[v*3+d];
555         }
556       }
557     }
558   }
559   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
560   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
561   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
562   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
563   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
564   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
565   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
566   /* Clean up intermediate storage */
567   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
568   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
569   PetscFunctionReturn(0);
570 }
571