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