xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 275fc1aefd25bc5577e6bb3daa1fca383bc51f02)
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       PetscInt elementSize;
260       char     *buffer;
261       double   *baseptr;
262 
263       elementSize = sizeof(int) + 3*sizeof(double);
264       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
265       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
266       for (v = 0; v < numVertices; ++v) {
267         baseptr = ((double*)(buffer+v*elementSize+sizeof(int)));
268         if (bswap) {ierr = PetscByteSwap(baseptr, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
269         coordsIn[v*3+0] = (PetscReal) baseptr[0];
270         coordsIn[v*3+1] = (PetscReal) baseptr[1];
271         coordsIn[v*3+2] = (PetscReal) baseptr[2];
272       }
273       ierr = PetscFree(buffer);CHKERRQ(ierr);
274     } else {
275       for (v = 0; v < numVertices; ++v) {
276         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
277         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
278         if (i != (int) v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+shift);
279       }
280     }
281     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
282     ierr = PetscStrncmp(line, "$EndNodes", 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     /* Read cells */
285     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
286     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
287     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
288     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
289     snum = sscanf(line, "%d", &numCells);
290     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
291   }
292 
293   if (!rank || binary) {
294     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
295        file contents multiple times to figure out the true number of cells and facets
296        in the given mesh. To make this more efficient we read the file contents only
297        once and store them in memory, while determining the true number of cells. */
298     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
299     for (trueNumCells=0, c = 0; c < numCells; ++c) {
300       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
301       if (gmsh_elem[c].dim == dim) trueNumCells++;
302     }
303     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
304     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
305     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
306     if (periodic) {
307       PetscInt pVert, *periodicMapT, *aux;
308       int      numPeriodic;
309 
310       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
311       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
312       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
313       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
314       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
315       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
316       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
317       snum = sscanf(line, "%d", &numPeriodic);
318       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
319       for (i = 0; i < numPeriodic; i++) {
320         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes;
321 
322         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
323         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
324         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
325         /* Type of tranformation, discarded */
326         ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
327         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
328         snum = sscanf(line, "%d", &nNodes);
329         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
330         for (j = 0; j < nNodes; j++) {
331           PetscInt slaveNode, masterNode;
332 
333           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
334           snum = sscanf(line, "%" PetscInt_FMT " %" PetscInt_FMT, &slaveNode, &masterNode);
335           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
336           periodicMapT[slaveNode - shift] = masterNode - shift;
337           ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr);
338           ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr);
339         }
340       }
341       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
342       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
343       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
344       /* we may have slaves of slaves */
345       for (i = 0; i < numVertices; i++) {
346         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
347           periodicMapT[i] = periodicMapT[periodicMapT[i]];
348         }
349       }
350       /* periodicMap : from old to new numbering (periodic vertices excluded)
351          periodicMapI: from new to old numbering */
352       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
353       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
354       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
355       for (i = 0, pVert = 0; i < numVertices; i++) {
356         if (periodicMapT[i] != i) {
357           pVert++;
358         } else {
359           aux[i] = i - pVert;
360           periodicMapI[i - pVert] = i;
361         }
362       }
363       for (i = 0 ; i < numVertices; i++) {
364         periodicMap[i] = aux[periodicMapT[i]];
365       }
366       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
367       ierr = PetscFree(aux);CHKERRQ(ierr);
368       /* remove periodic vertices */
369       numVertices = numVertices - pVert;
370     }
371   }
372   /* For binary we read on all ranks, but only build the plex on rank 0 */
373   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
374   /* Allocate the cell-vertex mesh */
375   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
376   if (!rank) {
377     for (cell = 0, c = 0; c < numCells; ++c) {
378       if (gmsh_elem[c].dim == dim) {
379         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
380         cell++;
381       }
382     }
383   }
384   ierr = DMSetUp(*dm);CHKERRQ(ierr);
385   /* Add cell-vertex connections */
386   if (!rank) {
387     PetscInt pcone[8], corner;
388 
389     for (cell = 0, c = 0; c < numCells; ++c) {
390       if (gmsh_elem[c].dim == dim) {
391         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
392           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
393           pcone[corner] = (periodicMap ? periodicMap[cc] : 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 
426   if (!rank && usemarker) {
427     PetscInt f, fStart, fEnd;
428 
429     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
430     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
431     for (f = fStart; f < fEnd; ++f) {
432       PetscInt suppSize;
433 
434       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
435       if (suppSize == 1) {
436         PetscInt *cone = NULL, coneSize, p;
437 
438         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
439         for (p = 0; p < coneSize; p += 2) {
440           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
441         }
442         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
443       }
444     }
445   }
446 
447   if (!rank) {
448     /* Apply boundary IDs by finding the relevant facets with vertex joins */
449     PetscInt pcone[8], corner, vStart, vEnd;
450 
451     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
452     for (c = 0; c < numCells; ++c) {
453       if (gmsh_elem[c].dim == dim-1) {
454         const PetscInt *join;
455         PetscInt        joinSize;
456 
457         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
458           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
459           pcone[corner] = (periodicMap ? periodicMap[cc] : 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         }
474         cell++;
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           const PetscInt cc = gmsh_elem[c].nodes[0] - shift;
483           PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart;
484           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
485         }
486       }
487     }
488   }
489 
490   /* Read coordinates */
491   embedDim = dim;
492   ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_bd", &isbd, NULL);CHKERRQ(ierr);
493   if (isbd) embedDim = dim+1;
494   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
495   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
496   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
497   if (periodicMap) { /* we need to localize coordinates on cells */
498     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
499   } else {
500     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
501   }
502   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
503     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
504     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
505   }
506   if (periodicMap) {
507     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
508     for (cell = 0, c = 0; c < numCells; ++c) {
509       if (gmsh_elem[c].dim == dim) {
510         PetscBool pc = PETSC_FALSE;
511         PetscInt  corner;
512         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
513           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift));
514         }
515         if (pc) {
516           ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr);
517           ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
518           ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
519         }
520         cell++;
521       }
522     }
523   }
524   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
525   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
526   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
527   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
528   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
529   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
530   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
531   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
532   if (!rank) {
533 
534     if (periodicMap) {
535       PetscInt off;
536 
537       for (cell = 0, c = 0; c < numCells; ++c) {
538         PetscInt pcone[8], corner;
539         if (gmsh_elem[c].dim == dim) {
540           if (PetscUnlikely(PetscBTLookup(periodicC,cell))) {
541             for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
542               pcone[corner] = gmsh_elem[c].nodes[corner] - shift;
543             }
544             if (dim == 3) {
545               /* Tetrahedra are inverted */
546               if (gmsh_elem[c].numNodes == 4) {
547                 PetscInt tmp = pcone[0];
548                 pcone[0] = pcone[1];
549                 pcone[1] = tmp;
550               }
551               /* Hexahedra are inverted */
552               if (gmsh_elem[c].numNodes == 8) {
553                 PetscInt tmp = pcone[1];
554                 pcone[1] = pcone[3];
555                 pcone[3] = tmp;
556               }
557             }
558             ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
559             for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
560               PetscInt v = pcone[corner];
561               for (d = 0; d < embedDim; ++d) {
562                 coords[off++] = coordsIn[v*3+d];
563               }
564             }
565           }
566           cell++;
567         }
568       }
569       for (v = 0; v < numVertices; ++v) {
570         ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
571         for (d = 0; d < embedDim; ++d) {
572           coords[off+d] = coordsIn[periodicMapI[v]*3+d];
573         }
574       }
575     } else {
576       for (v = 0; v < numVertices; ++v) {
577         for (d = 0; d < embedDim; ++d) {
578           coords[v*embedDim+d] = coordsIn[v*3+d];
579         }
580       }
581     }
582   }
583   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
584   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
585   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
586   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
587   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
588   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
589   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
590   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
591   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
592   /* Clean up intermediate storage */
593   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
594   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597