xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision f6c8a31f3a5641ae13cb90805eced0e729fa0106)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3331830f3SMatthew G. Knepley 
47d282ae0SMichael Lange /*@C
57d282ae0SMichael Lange   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
67d282ae0SMichael Lange 
77d282ae0SMichael Lange + comm        - The MPI communicator
87d282ae0SMichael Lange . filename    - Name of the Gmsh file
97d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh
107d282ae0SMichael Lange 
117d282ae0SMichael Lange   Output Parameter:
127d282ae0SMichael Lange . dm  - The DM object representing the mesh
137d282ae0SMichael Lange 
147d282ae0SMichael Lange   Level: beginner
157d282ae0SMichael Lange 
167d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
177d282ae0SMichael Lange @*/
187d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
197d282ae0SMichael Lange {
2011c56cb0SLisandro Dalcin   PetscViewer     viewer;
21abc86ac4SMichael Lange   PetscMPIInt     rank;
2211c56cb0SLisandro Dalcin   int             fileType;
237d282ae0SMichael Lange   PetscViewerType vtype;
247d282ae0SMichael Lange   PetscErrorCode  ierr;
257d282ae0SMichael Lange 
267d282ae0SMichael Lange   PetscFunctionBegin;
27abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2811c56cb0SLisandro Dalcin 
297d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
3011c56cb0SLisandro Dalcin   if (!rank) {
3111c56cb0SLisandro Dalcin     PetscViewer vheader;
3211c56cb0SLisandro Dalcin     char        line[PETSC_MAX_PATH_LEN];
3311c56cb0SLisandro Dalcin     PetscBool   match;
3411c56cb0SLisandro Dalcin     int         snum;
3511c56cb0SLisandro Dalcin     float       version;
3611c56cb0SLisandro Dalcin 
3711c56cb0SLisandro Dalcin     ierr = PetscViewerCreate(PETSC_COMM_SELF, &vheader);CHKERRQ(ierr);
387d282ae0SMichael Lange     ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
397d282ae0SMichael Lange     ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
407d282ae0SMichael Lange     ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
417d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
42060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
4311c56cb0SLisandro Dalcin     ierr = PetscStrncmp(line, "$MeshFormat", sizeof(line), &match);CHKERRQ(ierr);
447d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
45060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
4611c56cb0SLisandro Dalcin     snum = sscanf(line, "%f %d", &version, &fileType);
47f6ac7a6aSMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
48f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
4911c56cb0SLisandro Dalcin     ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
50abc86ac4SMichael Lange   }
5111c56cb0SLisandro Dalcin   ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
5211c56cb0SLisandro Dalcin   vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY;
5311c56cb0SLisandro Dalcin 
547d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
557d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
577d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
587d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
597d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
607d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
617d282ae0SMichael Lange   PetscFunctionReturn(0);
627d282ae0SMichael Lange }
637d282ae0SMichael Lange 
6411c56cb0SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int numVertices, double **coordinates)
65df032b82SMatthew G. Knepley {
6611c56cb0SLisandro Dalcin   int            v,nid;
6711c56cb0SLisandro Dalcin   PetscErrorCode ierr;
6811c56cb0SLisandro Dalcin 
6911c56cb0SLisandro Dalcin   PetscFunctionBegin;
7011c56cb0SLisandro Dalcin   ierr = PetscMalloc1(numVertices*3, coordinates);CHKERRQ(ierr);
7111c56cb0SLisandro Dalcin   for (v = 0; v < numVertices; ++v) {
7211c56cb0SLisandro Dalcin     double *xyz = *coordinates + v*3;
7311c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
7411c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
7511c56cb0SLisandro Dalcin     if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift);
7611c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
7711c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
7811c56cb0SLisandro Dalcin   }
7911c56cb0SLisandro Dalcin   PetscFunctionReturn(0);
8011c56cb0SLisandro Dalcin }
8111c56cb0SLisandro Dalcin 
8211c56cb0SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, int numCells, GmshElement **gmsh_elems)
8311c56cb0SLisandro Dalcin {
84df032b82SMatthew G. Knepley   GmshElement   *elements;
8511c56cb0SLisandro Dalcin   int            i, c, p, ibuf[1+4+512];
8611c56cb0SLisandro Dalcin   int            cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
87df032b82SMatthew G. Knepley   PetscErrorCode ierr;
88df032b82SMatthew G. Knepley 
89df032b82SMatthew G. Knepley   PetscFunctionBegin;
90df032b82SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
91df032b82SMatthew G. Knepley   for (c = 0; c < numCells;) {
9211c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
9311c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
94df032b82SMatthew G. Knepley     if (binary) {
95df032b82SMatthew G. Knepley       cellType = ibuf[0];
96df032b82SMatthew G. Knepley       numElem = ibuf[1];
97df032b82SMatthew G. Knepley       numTags = ibuf[2];
98df032b82SMatthew G. Knepley     } else {
99df032b82SMatthew G. Knepley       elements[c].id = ibuf[0];
100df032b82SMatthew G. Knepley       cellType = ibuf[1];
101df032b82SMatthew G. Knepley       numTags = ibuf[2];
102df032b82SMatthew G. Knepley       numElem = 1;
103df032b82SMatthew G. Knepley     }
104bf6ba3a3SMatthew G. Knepley     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
105bf6ba3a3SMatthew G. Knepley     numNodesIgnore = 0;
106df032b82SMatthew G. Knepley     switch (cellType) {
107df032b82SMatthew G. Knepley     case 1: /* 2-node line */
108df032b82SMatthew G. Knepley       dim = 1;
109df032b82SMatthew G. Knepley       numNodes = 2;
110df032b82SMatthew G. Knepley       break;
111df032b82SMatthew G. Knepley     case 2: /* 3-node triangle */
112df032b82SMatthew G. Knepley       dim = 2;
113df032b82SMatthew G. Knepley       numNodes = 3;
114df032b82SMatthew G. Knepley       break;
115df032b82SMatthew G. Knepley     case 3: /* 4-node quadrangle */
116df032b82SMatthew G. Knepley       dim = 2;
117df032b82SMatthew G. Knepley       numNodes = 4;
118df032b82SMatthew G. Knepley       break;
119df032b82SMatthew G. Knepley     case 4: /* 4-node tetrahedron */
120df032b82SMatthew G. Knepley       dim  = 3;
121df032b82SMatthew G. Knepley       numNodes = 4;
122df032b82SMatthew G. Knepley       break;
123df032b82SMatthew G. Knepley     case 5: /* 8-node hexahedron */
124df032b82SMatthew G. Knepley       dim = 3;
125df032b82SMatthew G. Knepley       numNodes = 8;
126df032b82SMatthew G. Knepley       break;
127bf6ba3a3SMatthew G. Knepley     case 8: /* 3-node 2nd order line */
128bf6ba3a3SMatthew G. Knepley       dim = 1;
129bf6ba3a3SMatthew G. Knepley       numNodes = 2;
130bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 1;
131bf6ba3a3SMatthew G. Knepley       break;
132bf6ba3a3SMatthew G. Knepley     case 9: /* 6-node 2nd order triangle */
133bf6ba3a3SMatthew G. Knepley       dim = 2;
134bf6ba3a3SMatthew G. Knepley       numNodes = 3;
135bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 3;
136bf6ba3a3SMatthew G. Knepley       break;
137df032b82SMatthew G. Knepley     case 15: /* 1-node vertex */
138df032b82SMatthew G. Knepley       dim = 0;
139df032b82SMatthew G. Knepley       numNodes = 1;
140df032b82SMatthew G. Knepley       break;
141bf6ba3a3SMatthew G. Knepley     case 6: /* 6-node prism */
142bf6ba3a3SMatthew G. Knepley     case 7: /* 5-node pyramid */
143bf6ba3a3SMatthew G. Knepley     case 10: /* 9-node 2nd order quadrangle */
144bf6ba3a3SMatthew G. Knepley     case 11: /* 10-node 2nd order tetrahedron */
145bf6ba3a3SMatthew G. Knepley     case 12: /* 27-node 2nd order hexhedron */
146bf6ba3a3SMatthew G. Knepley     case 13: /* 19-node 2nd order prism */
147bf6ba3a3SMatthew G. Knepley     case 14: /* 14-node 2nd order pyramid */
148df032b82SMatthew G. Knepley     default:
149df032b82SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
150df032b82SMatthew G. Knepley     }
151df032b82SMatthew G. Knepley     if (binary) {
15211c56cb0SLisandro Dalcin       const int nint = 1 + numTags + numNodes + numNodesIgnore;
15311c56cb0SLisandro Dalcin       /* Loop over element blocks */
154df032b82SMatthew G. Knepley       for (i = 0; i < numElem; ++i, ++c) {
15511c56cb0SLisandro Dalcin         ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
15611c56cb0SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);}
157df032b82SMatthew G. Knepley         elements[c].dim = dim;
158df032b82SMatthew G. Knepley         elements[c].numNodes = numNodes;
159df032b82SMatthew G. Knepley         elements[c].numTags = numTags;
160df032b82SMatthew G. Knepley         elements[c].id = ibuf[0];
161df032b82SMatthew G. Knepley         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
162df032b82SMatthew G. Knepley         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
163df032b82SMatthew G. Knepley       }
164df032b82SMatthew G. Knepley     } else {
165df032b82SMatthew G. Knepley       elements[c].dim = dim;
166df032b82SMatthew G. Knepley       elements[c].numNodes = numNodes;
167df032b82SMatthew G. Knepley       elements[c].numTags = numTags;
168df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
169df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
17011c56cb0SLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr);
171df032b82SMatthew G. Knepley       c++;
172df032b82SMatthew G. Knepley     }
173df032b82SMatthew G. Knepley   }
174df032b82SMatthew G. Knepley   *gmsh_elems = elements;
175df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
176df032b82SMatthew G. Knepley }
1777d282ae0SMichael Lange 
178331830f3SMatthew G. Knepley /*@
1797d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
180331830f3SMatthew G. Knepley 
181331830f3SMatthew G. Knepley   Collective on comm
182331830f3SMatthew G. Knepley 
183331830f3SMatthew G. Knepley   Input Parameters:
184331830f3SMatthew G. Knepley + comm  - The MPI communicator
185331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
186331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
187331830f3SMatthew G. Knepley 
188331830f3SMatthew G. Knepley   Output Parameter:
189331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
190331830f3SMatthew G. Knepley 
191331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
1923d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
193331830f3SMatthew G. Knepley 
194331830f3SMatthew G. Knepley   Level: beginner
195331830f3SMatthew G. Knepley 
196331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
197331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
198331830f3SMatthew G. Knepley @*/
199331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
200331830f3SMatthew G. Knepley {
20111c56cb0SLisandro Dalcin   PetscViewer    parentviewer = NULL;
20211c56cb0SLisandro Dalcin   double        *coordsIn = NULL;
2033c67d5adSSatish Balay   GmshElement   *gmsh_elem = NULL;
204331830f3SMatthew G. Knepley   PetscSection   coordSection;
205331830f3SMatthew G. Knepley   Vec            coordinates;
2066fbe17bfSStefano Zampini   PetscBT        periodicV = NULL, periodicC = NULL;
20784572febSToby Isaac   PetscScalar   *coords;
208dd169d64SMatthew G. Knepley   PetscInt       dim = 0, embedDim, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL;
2091d64f2ccSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1;
21011c56cb0SLisandro Dalcin   PetscMPIInt    rank;
211331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
21211c56cb0SLisandro Dalcin   PetscBool      binary, mpiio, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, isbd = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE;
213331830f3SMatthew G. Knepley   PetscErrorCode ierr;
214331830f3SMatthew G. Knepley 
215331830f3SMatthew G. Knepley   PetscFunctionBegin;
216331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
21711c56cb0SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr);
2185964b3dcSLisandro Dalcin   ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr);
21911c56cb0SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr);
22011c56cb0SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_bd", &isbd, NULL);CHKERRQ(ierr);
22111c56cb0SLisandro Dalcin   if (zerobase) shift = 0;
22211c56cb0SLisandro Dalcin 
223331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
224331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2253b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
22611c56cb0SLisandro Dalcin 
22711c56cb0SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
22811c56cb0SLisandro Dalcin   ierr = PetscViewerBinaryGetUseMPIIO(viewer, &mpiio);CHKERRQ(ierr);
22911c56cb0SLisandro Dalcin 
23011c56cb0SLisandro Dalcin   /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
23111c56cb0SLisandro Dalcin   if (binary && !mpiio) {
23211c56cb0SLisandro Dalcin     parentviewer = viewer;
23311c56cb0SLisandro Dalcin     ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
23411c56cb0SLisandro Dalcin   }
23511c56cb0SLisandro Dalcin 
23611c56cb0SLisandro Dalcin   if (!rank || mpiio) {
237331830f3SMatthew G. Knepley     PetscBool match;
238331830f3SMatthew G. Knepley     int       fileType, dataSize;
239f6ac7a6aSMichael Lange     float     version;
240331830f3SMatthew G. Knepley 
241331830f3SMatthew G. Knepley     /* Read header */
242060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
24304d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
244331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
245060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
246f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
247f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
248f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
249331830f3SMatthew G. Knepley     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
25004d1ad83SMichael Lange     if (binary) {
251b9eae255SMichael Lange       int checkInt;
252060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
25304d1ad83SMichael Lange       if (checkInt != 1) {
254b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
25511c56cb0SLisandro Dalcin         if (checkInt == 1) byteSwap = PETSC_TRUE;
25604d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
25704d1ad83SMichael Lange       }
2580877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
259060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
26004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
261331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
26211c56cb0SLisandro Dalcin 
263dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
264060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
265dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
266dc0ede3bSMatthew G. Knepley     if (match) {
267dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
268dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
269dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
27011c56cb0SLisandro Dalcin       for (r = 0; r < numRegions; ++r) {
27111c56cb0SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
27211c56cb0SLisandro Dalcin       }
273dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
274dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
275dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
276dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
277dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
278dc0ede3bSMatthew G. Knepley     }
27911c56cb0SLisandro Dalcin 
280dc0ede3bSMatthew G. Knepley     /* Read vertices */
28104d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
282331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
283060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
2840877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
2850877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
28611c56cb0SLisandro Dalcin     ierr = DMPlexCreateGmsh_ReadNodes(viewer, binary, byteSwap, shift, numVertices, &coordsIn);CHKERRQ(ierr);
287060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
28804d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
289331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
29011c56cb0SLisandro Dalcin 
291331830f3SMatthew G. Knepley     /* Read cells */
292060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
29304d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
294331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
295060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
2960877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
2970877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
298a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
299a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
300a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
301a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
30211c56cb0SLisandro Dalcin     ierr = DMPlexCreateGmsh_ReadElements(viewer, binary, byteSwap, shift, numCells, &gmsh_elem);CHKERRQ(ierr);
303a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
304a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
305a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
306db397164SMichael Lange     }
307d08df55aSStefano Zampini     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
3081b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
309331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
31011c56cb0SLisandro Dalcin 
31111c56cb0SLisandro Dalcin     /* OPTIONAL Read periodic section */
312d08df55aSStefano Zampini     if (periodic) {
313f45c9589SStefano Zampini       PetscInt pVert, *periodicMapT, *aux;
314d08df55aSStefano Zampini       int      numPeriodic;
315d08df55aSStefano Zampini 
316d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
317d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
318d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
319f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
3206fbe17bfSStefano Zampini       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
321f45c9589SStefano Zampini       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
322d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
323d08df55aSStefano Zampini       snum = sscanf(line, "%d", &numPeriodic);
324d08df55aSStefano Zampini       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
325d08df55aSStefano Zampini       for (i = 0; i < numPeriodic; i++) {
326da83f57bSLisandro Dalcin         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes, slaveNode, masterNode;
327d08df55aSStefano Zampini 
328d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
329d08df55aSStefano Zampini         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
330d08df55aSStefano Zampini         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
331da83f57bSLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
332da83f57bSLisandro Dalcin         snum = sscanf(line, "%d", &nNodes);
333da83f57bSLisandro Dalcin         if (snum != 1) { /* discard tranformation and try again */
33472ffbcc9SStefano Zampini           ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
335d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
336d08df55aSStefano Zampini           snum = sscanf(line, "%d", &nNodes);
337d08df55aSStefano Zampini           if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
338da83f57bSLisandro Dalcin         }
339d08df55aSStefano Zampini         for (j = 0; j < nNodes; j++) {
340d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
34111c56cb0SLisandro Dalcin           snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
342d08df55aSStefano Zampini           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
343917266a4SLisandro Dalcin           periodicMapT[slaveNode - shift] = masterNode - shift;
344917266a4SLisandro Dalcin           ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr);
345917266a4SLisandro Dalcin           ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr);
346d08df55aSStefano Zampini         }
347d08df55aSStefano Zampini       }
348d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
349d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
350d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
351d08df55aSStefano Zampini       /* we may have slaves of slaves */
352d08df55aSStefano Zampini       for (i = 0; i < numVertices; i++) {
353f45c9589SStefano Zampini         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
354f45c9589SStefano Zampini           periodicMapT[i] = periodicMapT[periodicMapT[i]];
355d08df55aSStefano Zampini         }
356d08df55aSStefano Zampini       }
357f45c9589SStefano Zampini       /* periodicMap : from old to new numbering (periodic vertices excluded)
358f45c9589SStefano Zampini          periodicMapI: from new to old numbering */
359f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
360f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
361f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
362f45c9589SStefano Zampini       for (i = 0, pVert = 0; i < numVertices; i++) {
363f45c9589SStefano Zampini         if (periodicMapT[i] != i) {
364f45c9589SStefano Zampini           pVert++;
365f45c9589SStefano Zampini         } else {
366f45c9589SStefano Zampini           aux[i] = i - pVert;
367f45c9589SStefano Zampini           periodicMapI[i - pVert] = i;
368f45c9589SStefano Zampini         }
369f45c9589SStefano Zampini       }
370f45c9589SStefano Zampini       for (i = 0 ; i < numVertices; i++) {
371f45c9589SStefano Zampini         periodicMap[i] = aux[periodicMapT[i]];
372f45c9589SStefano Zampini       }
373f45c9589SStefano Zampini       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
374f45c9589SStefano Zampini       ierr = PetscFree(aux);CHKERRQ(ierr);
375f45c9589SStefano Zampini       /* remove periodic vertices */
376f45c9589SStefano Zampini       numVertices = numVertices - pVert;
377d08df55aSStefano Zampini     }
378331830f3SMatthew G. Knepley   }
37911c56cb0SLisandro Dalcin 
38011c56cb0SLisandro Dalcin   /* For (binary) MPI I/O we read on all ranks, but only build the plex on rank 0 */
38111c56cb0SLisandro Dalcin   if (mpiio && rank) {
38211c56cb0SLisandro Dalcin     numVertices = 0; numCells = 0; trueNumCells = 0;
38311c56cb0SLisandro Dalcin     ierr = PetscFree(periodicMap);CHKERRQ(ierr);
38411c56cb0SLisandro Dalcin     ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
38511c56cb0SLisandro Dalcin   }
38611c56cb0SLisandro Dalcin 
38711c56cb0SLisandro Dalcin   if (parentviewer) {
38811c56cb0SLisandro Dalcin     ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
38911c56cb0SLisandro Dalcin   }
39011c56cb0SLisandro Dalcin 
391a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
392331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
393a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
394a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
395a4bb7517SMichael Lange       ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
396a4bb7517SMichael Lange       cell++;
397331830f3SMatthew G. Knepley     }
398331830f3SMatthew G. Knepley   }
399331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
400a4bb7517SMichael Lange   /* Add cell-vertex connections */
401a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
402a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
40311c56cb0SLisandro Dalcin       PetscInt pcone[8], corner;
404a4bb7517SMichael Lange       for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
405dd169d64SMatthew G. Knepley         const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
406917266a4SLisandro Dalcin         pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells;
407331830f3SMatthew G. Knepley       }
40897ed6be6Ssarens       if (dim == 3) {
40997ed6be6Ssarens         /* Tetrahedra are inverted */
41097ed6be6Ssarens         if (gmsh_elem[c].numNodes == 4) {
41197ed6be6Ssarens           PetscInt tmp = pcone[0];
41297ed6be6Ssarens           pcone[0] = pcone[1];
41397ed6be6Ssarens           pcone[1] = tmp;
41497ed6be6Ssarens         }
41597ed6be6Ssarens         /* Hexahedra are inverted */
41697ed6be6Ssarens         if (gmsh_elem[c].numNodes == 8) {
41797ed6be6Ssarens           PetscInt tmp = pcone[1];
41897ed6be6Ssarens           pcone[1] = pcone[3];
41997ed6be6Ssarens           pcone[3] = tmp;
42097ed6be6Ssarens         }
42197ed6be6Ssarens       }
422a4bb7517SMichael Lange       ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
423a4bb7517SMichael Lange       cell++;
424331830f3SMatthew G. Knepley     }
425a4bb7517SMichael Lange   }
4266228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
427c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
428331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
429331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
430331830f3SMatthew G. Knepley   if (interpolate) {
431e56d480eSMatthew G. Knepley     DM idm = NULL;
432331830f3SMatthew G. Knepley 
433331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
434331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
435331830f3SMatthew G. Knepley     *dm  = idm;
436331830f3SMatthew G. Knepley   }
437917266a4SLisandro Dalcin 
438*f6c8a31fSLisandro Dalcin   if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation");
439917266a4SLisandro Dalcin   if (!rank && usemarker) {
440d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
441d3f73514SStefano Zampini 
442d3f73514SStefano Zampini     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
443d3f73514SStefano Zampini     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
444d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
445d3f73514SStefano Zampini       PetscInt suppSize;
446d3f73514SStefano Zampini 
447d3f73514SStefano Zampini       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
448d3f73514SStefano Zampini       if (suppSize == 1) {
449d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
450d3f73514SStefano Zampini 
451d3f73514SStefano Zampini         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
452d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
453d3f73514SStefano Zampini           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
454d3f73514SStefano Zampini         }
455d3f73514SStefano Zampini         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
456d3f73514SStefano Zampini       }
457d3f73514SStefano Zampini     }
458d3f73514SStefano Zampini   }
45916ad7e67SMichael Lange 
46016ad7e67SMichael Lange   if (!rank) {
46111c56cb0SLisandro Dalcin     PetscInt vStart, vEnd;
462d1a54cefSMatthew G. Knepley 
46316ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
46411c56cb0SLisandro Dalcin     for (cell = 0, c = 0; c < numCells; ++c) {
46511c56cb0SLisandro Dalcin 
46611c56cb0SLisandro Dalcin       /* Create face sets */
4675964b3dcSLisandro Dalcin       if (interpolate && gmsh_elem[c].dim == dim-1) {
46816ad7e67SMichael Lange         const PetscInt *join;
46911c56cb0SLisandro Dalcin         PetscInt        joinSize, pcone[8], corner;
47011c56cb0SLisandro Dalcin         /* Find the relevant facet with vertex joins */
471a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
472dd169d64SMatthew G. Knepley           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
473917266a4SLisandro Dalcin           pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart;
47416ad7e67SMichael Lange         }
47511c56cb0SLisandro Dalcin         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr);
476a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
477c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
478a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
47916ad7e67SMichael Lange       }
4806e1dd89aSLawrence Mitchell 
4816e1dd89aSLawrence Mitchell       /* Create cell sets */
4826e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
4836e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
4846e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4856e1dd89aSLawrence Mitchell         }
486917266a4SLisandro Dalcin         cell++;
4876e1dd89aSLawrence Mitchell       }
4880c070f12SSander Arens 
4890c070f12SSander Arens       /* Create vertex sets */
4900c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
4910c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
492917266a4SLisandro Dalcin           const PetscInt cc = gmsh_elem[c].nodes[0] - shift;
49311c56cb0SLisandro Dalcin           const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart;
494d08df55aSStefano Zampini           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4950c070f12SSander Arens         }
4960c070f12SSander Arens       }
4970c070f12SSander Arens     }
49816ad7e67SMichael Lange   }
49916ad7e67SMichael Lange 
50011c56cb0SLisandro Dalcin   /* Create coordinates */
50111c56cb0SLisandro Dalcin   embedDim = isbd ? dim+1 : dim;
502979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
503331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
5041d64f2ccSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
505f45c9589SStefano Zampini   if (periodicMap) { /* we need to localize coordinates on cells */
506f45c9589SStefano Zampini     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
507f45c9589SStefano Zampini   } else {
50875b5763bSMichael Lange     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
509f45c9589SStefano Zampini   }
51075b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
5111d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
5121d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
513331830f3SMatthew G. Knepley   }
51411c56cb0SLisandro Dalcin   if (periodicMap) {
515437bdd5bSStefano Zampini     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
516f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
517f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
518437bdd5bSStefano Zampini         PetscInt  corner;
51911c56cb0SLisandro Dalcin         PetscBool pc = PETSC_FALSE;
520437bdd5bSStefano Zampini         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
521917266a4SLisandro Dalcin           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift));
522437bdd5bSStefano Zampini         }
523437bdd5bSStefano Zampini         if (pc) {
52411c56cb0SLisandro Dalcin           PetscInt dof = gmsh_elem[c].numNodes*embedDim;
525437bdd5bSStefano Zampini           ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr);
52611c56cb0SLisandro Dalcin           ierr = PetscSectionSetDof(coordSection, cell, dof);CHKERRQ(ierr);
52711c56cb0SLisandro Dalcin           ierr = PetscSectionSetFieldDof(coordSection, cell, 0, dof);CHKERRQ(ierr);
5286fbe17bfSStefano Zampini         }
529f45c9589SStefano Zampini         cell++;
530f45c9589SStefano Zampini       }
531f45c9589SStefano Zampini     }
532f45c9589SStefano Zampini   }
533331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
534331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
5358b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
536331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
537331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
5381d64f2ccSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
539331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
540331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
541f45c9589SStefano Zampini   if (periodicMap) {
542f45c9589SStefano Zampini     PetscInt off;
543f45c9589SStefano Zampini 
544f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
545f45c9589SStefano Zampini       PetscInt pcone[8], corner;
546f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
5476fbe17bfSStefano Zampini         if (PetscUnlikely(PetscBTLookup(periodicC, cell))) {
548f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
549917266a4SLisandro Dalcin             pcone[corner] = gmsh_elem[c].nodes[corner] - shift;
550f45c9589SStefano Zampini           }
551f45c9589SStefano Zampini           if (dim == 3) {
552f45c9589SStefano Zampini             /* Tetrahedra are inverted */
553f45c9589SStefano Zampini             if (gmsh_elem[c].numNodes == 4) {
554f45c9589SStefano Zampini               PetscInt tmp = pcone[0];
555f45c9589SStefano Zampini               pcone[0] = pcone[1];
556f45c9589SStefano Zampini               pcone[1] = tmp;
557f45c9589SStefano Zampini             }
558f45c9589SStefano Zampini             /* Hexahedra are inverted */
559f45c9589SStefano Zampini             if (gmsh_elem[c].numNodes == 8) {
560f45c9589SStefano Zampini               PetscInt tmp = pcone[1];
561f45c9589SStefano Zampini               pcone[1] = pcone[3];
562f45c9589SStefano Zampini               pcone[3] = tmp;
563f45c9589SStefano Zampini             }
564f45c9589SStefano Zampini           }
565f45c9589SStefano Zampini           ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
566f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
56711c56cb0SLisandro Dalcin             v = pcone[corner];
568dd169d64SMatthew G. Knepley             for (d = 0; d < embedDim; ++d) {
56911c56cb0SLisandro Dalcin               coords[off++] = (PetscReal) coordsIn[v*3+d];
570f45c9589SStefano Zampini             }
571f45c9589SStefano Zampini           }
5726fbe17bfSStefano Zampini         }
573f45c9589SStefano Zampini         cell++;
574f45c9589SStefano Zampini       }
575f45c9589SStefano Zampini     }
576f45c9589SStefano Zampini     for (v = 0; v < numVertices; ++v) {
577f45c9589SStefano Zampini       ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
578dd169d64SMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
57911c56cb0SLisandro Dalcin         coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d];
580f45c9589SStefano Zampini       }
581f45c9589SStefano Zampini     }
582f45c9589SStefano Zampini   } else {
583331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
5841d64f2ccSMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
58511c56cb0SLisandro Dalcin         coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d];
586331830f3SMatthew G. Knepley       }
587331830f3SMatthew G. Knepley     }
588331830f3SMatthew G. Knepley   }
589331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
590331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
59111c56cb0SLisandro Dalcin   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
59211c56cb0SLisandro Dalcin 
59311c56cb0SLisandro Dalcin   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
59411c56cb0SLisandro Dalcin   ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
595d08df55aSStefano Zampini   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
596fcd9ca0aSStefano Zampini   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
5976fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
5986fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
59911c56cb0SLisandro Dalcin   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
60011c56cb0SLisandro Dalcin 
6013b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
602331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
603331830f3SMatthew G. Knepley }
604