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