xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 96509d17ab28f09a5113ddd6b708eeb68d54db76)
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 = 0;
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   PetscBT        periodicV = NULL, periodicC = NULL;
190   PetscScalar   *coords;
191   PetscReal     *coordsIn = NULL;
192   PetscInt       dim = 0, embedDim, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL;
193   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1;
194   PetscMPIInt    num_proc, rank;
195   char           line[PETSC_MAX_PATH_LEN];
196   PetscBool      zerobase = PETSC_FALSE, isbd = PETSC_FALSE, match, binary, bswap = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE;
197   PetscErrorCode ierr;
198 
199   PetscFunctionBegin;
200   ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr);
201   if (zerobase) shift = 0;
202   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
203   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
204   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
205   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
206   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
207   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
208   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
209   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_periodic",&periodic,NULL);CHKERRQ(ierr);
210   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_usemarker",&usemarker,NULL);CHKERRQ(ierr);
211   if (!rank || binary) {
212     PetscBool match;
213     int       fileType, dataSize;
214     float     version;
215 
216     /* Read header */
217     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
218     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
219     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
220     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
221     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
222     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
223     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
224     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
225     if (binary) {
226       int checkInt;
227       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
228       if (checkInt != 1) {
229         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
230         if (checkInt == 1) bswap = PETSC_TRUE;
231         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
232       }
233     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
234     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
235     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
236     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
237     /* OPTIONAL Read physical names */
238     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
239     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
240     if (match) {
241       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
242       snum = sscanf(line, "%d", &numRegions);
243       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
244       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
245       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
246       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
247       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
248       /* Initial read for vertex section */
249       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
250     }
251     /* Read vertices */
252     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
253     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
254     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
255     snum = sscanf(line, "%d", &numVertices);
256     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
257     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
258     if (binary) {
259       int      tInt;
260       double   tDouble[3];
261 
262       for (v = 0; v < numVertices; ++v) {
263         ierr = PetscViewerRead(viewer, &tInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
264         ierr = PetscViewerRead(viewer, tDouble, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
265         if (bswap) {ierr = PetscByteSwap(tDouble, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
266         coordsIn[v*3+0] = (PetscReal) tDouble[0];
267         coordsIn[v*3+1] = (PetscReal) tDouble[1];
268         coordsIn[v*3+2] = (PetscReal) tDouble[2];
269       }
270     } else {
271       for (v = 0; v < numVertices; ++v) {
272         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
273         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
274         if (i != (int) v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+shift);
275       }
276     }
277     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
278     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
279     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
280     /* Read cells */
281     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
282     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
283     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
284     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
285     snum = sscanf(line, "%d", &numCells);
286     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
287   }
288 
289   if (!rank || binary) {
290     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
291        file contents multiple times to figure out the true number of cells and facets
292        in the given mesh. To make this more efficient we read the file contents only
293        once and store them in memory, while determining the true number of cells. */
294     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
295     for (trueNumCells=0, c = 0; c < numCells; ++c) {
296       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
297       if (gmsh_elem[c].dim == dim) trueNumCells++;
298     }
299     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
300     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
301     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
302     if (periodic) {
303       PetscInt pVert, *periodicMapT, *aux;
304       int      numPeriodic;
305 
306       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
307       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
308       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
309       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
310       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
311       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
312       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
313       snum = sscanf(line, "%d", &numPeriodic);
314       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
315       for (i = 0; i < numPeriodic; i++) {
316         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes;
317 
318         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
319         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
320         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
321         /* Type of tranformation, discarded */
322         ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
323         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
324         snum = sscanf(line, "%d", &nNodes);
325         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
326         for (j = 0; j < nNodes; j++) {
327           PetscInt slaveNode, masterNode;
328 
329           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
330           snum = sscanf(line, "%" PetscInt_FMT " %" PetscInt_FMT, &slaveNode, &masterNode);
331           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
332           periodicMapT[slaveNode - shift] = masterNode - shift;
333           ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr);
334           ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr);
335         }
336       }
337       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
338       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
339       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
340       /* we may have slaves of slaves */
341       for (i = 0; i < numVertices; i++) {
342         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
343           periodicMapT[i] = periodicMapT[periodicMapT[i]];
344         }
345       }
346       /* periodicMap : from old to new numbering (periodic vertices excluded)
347          periodicMapI: from new to old numbering */
348       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
349       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
350       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
351       for (i = 0, pVert = 0; i < numVertices; i++) {
352         if (periodicMapT[i] != i) {
353           pVert++;
354         } else {
355           aux[i] = i - pVert;
356           periodicMapI[i - pVert] = i;
357         }
358       }
359       for (i = 0 ; i < numVertices; i++) {
360         periodicMap[i] = aux[periodicMapT[i]];
361       }
362       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
363       ierr = PetscFree(aux);CHKERRQ(ierr);
364       /* remove periodic vertices */
365       numVertices = numVertices - pVert;
366     }
367   }
368   /* For binary we read on all ranks, but only build the plex on rank 0 */
369   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
370   /* Allocate the cell-vertex mesh */
371   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
372   if (!rank) {
373     for (cell = 0, c = 0; c < numCells; ++c) {
374       if (gmsh_elem[c].dim == dim) {
375         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
376         cell++;
377       }
378     }
379   }
380   ierr = DMSetUp(*dm);CHKERRQ(ierr);
381   /* Add cell-vertex connections */
382   if (!rank) {
383     PetscInt pcone[8], corner;
384 
385     for (cell = 0, c = 0; c < numCells; ++c) {
386       if (gmsh_elem[c].dim == dim) {
387         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
388           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
389           pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells;
390         }
391         if (dim == 3) {
392           /* Tetrahedra are inverted */
393           if (gmsh_elem[c].numNodes == 4) {
394             PetscInt tmp = pcone[0];
395             pcone[0] = pcone[1];
396             pcone[1] = tmp;
397           }
398           /* Hexahedra are inverted */
399           if (gmsh_elem[c].numNodes == 8) {
400             PetscInt tmp = pcone[1];
401             pcone[1] = pcone[3];
402             pcone[3] = tmp;
403           }
404         }
405         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
406         cell++;
407       }
408     }
409   }
410   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
411   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
412   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
413   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
414   if (interpolate) {
415     DM idm = NULL;
416 
417     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
418     ierr = DMDestroy(dm);CHKERRQ(ierr);
419     *dm  = idm;
420   }
421 
422   if (!rank && usemarker) {
423     PetscInt f, fStart, fEnd;
424 
425     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
426     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
427     for (f = fStart; f < fEnd; ++f) {
428       PetscInt suppSize;
429 
430       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
431       if (suppSize == 1) {
432         PetscInt *cone = NULL, coneSize, p;
433 
434         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
435         for (p = 0; p < coneSize; p += 2) {
436           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
437         }
438         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
439       }
440     }
441   }
442 
443   if (!rank) {
444     /* Apply boundary IDs by finding the relevant facets with vertex joins */
445     PetscInt pcone[8], corner, vStart, vEnd;
446 
447     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
448     for (c = 0; c < numCells; ++c) {
449       if (gmsh_elem[c].dim == dim-1) {
450         const PetscInt *join;
451         PetscInt        joinSize;
452 
453         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
454           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
455           pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart;
456         }
457         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
458         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
459         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
460         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
461       }
462     }
463 
464     /* Create cell sets */
465     for (cell = 0, c = 0; c < numCells; ++c) {
466       if (gmsh_elem[c].dim == dim) {
467         if (gmsh_elem[c].numTags > 0) {
468           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
469         }
470         cell++;
471       }
472     }
473 
474     /* Create vertex sets */
475     for (c = 0; c < numCells; ++c) {
476       if (gmsh_elem[c].dim == 0) {
477         if (gmsh_elem[c].numTags > 0) {
478           const PetscInt cc = gmsh_elem[c].nodes[0] - shift;
479           PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart;
480           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
481         }
482       }
483     }
484   }
485 
486   /* Read coordinates */
487   embedDim = dim;
488   ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_bd", &isbd, NULL);CHKERRQ(ierr);
489   if (isbd) embedDim = dim+1;
490   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
491   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
492   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
493   if (periodicMap) { /* we need to localize coordinates on cells */
494     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
495   } else {
496     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
497   }
498   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
499     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
500     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
501   }
502   if (!rank && periodicMap) {
503     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
504     for (cell = 0, c = 0; c < numCells; ++c) {
505       if (gmsh_elem[c].dim == dim) {
506         PetscBool pc = PETSC_FALSE;
507         PetscInt  corner;
508         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
509           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift));
510         }
511         if (pc) {
512           ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr);
513           ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
514           ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
515         }
516         cell++;
517       }
518     }
519   }
520   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
521   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
522   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
523   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
524   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
525   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
526   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
527   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
528   if (!rank) {
529 
530     if (periodicMap) {
531       PetscInt off;
532 
533       for (cell = 0, c = 0; c < numCells; ++c) {
534         PetscInt pcone[8], corner;
535         if (gmsh_elem[c].dim == dim) {
536           if (PetscUnlikely(PetscBTLookup(periodicC,cell))) {
537             for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
538               pcone[corner] = gmsh_elem[c].nodes[corner] - shift;
539             }
540             if (dim == 3) {
541               /* Tetrahedra are inverted */
542               if (gmsh_elem[c].numNodes == 4) {
543                 PetscInt tmp = pcone[0];
544                 pcone[0] = pcone[1];
545                 pcone[1] = tmp;
546               }
547               /* Hexahedra are inverted */
548               if (gmsh_elem[c].numNodes == 8) {
549                 PetscInt tmp = pcone[1];
550                 pcone[1] = pcone[3];
551                 pcone[3] = tmp;
552               }
553             }
554             ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
555             for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
556               PetscInt v = pcone[corner];
557               for (d = 0; d < embedDim; ++d) {
558                 coords[off++] = coordsIn[v*3+d];
559               }
560             }
561           }
562           cell++;
563         }
564       }
565       for (v = 0; v < numVertices; ++v) {
566         ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
567         for (d = 0; d < embedDim; ++d) {
568           coords[off+d] = coordsIn[periodicMapI[v]*3+d];
569         }
570       }
571     } else {
572       for (v = 0; v < numVertices; ++v) {
573         for (d = 0; d < embedDim; ++d) {
574           coords[v*embedDim+d] = coordsIn[v*3+d];
575         }
576       }
577     }
578   }
579   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
580   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
581   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
582   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
583   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
584   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
585   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
586   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
587   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
588   /* Clean up intermediate storage */
589   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
590   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593