xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 3b46f529e54e44e362952d229e16c6ede533652a)
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];
212*3b46f529SLisandro Dalcin   PetscBool      binary, 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 
22911c56cb0SLisandro Dalcin   /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
230*3b46f529SLisandro Dalcin   if (binary) {
23111c56cb0SLisandro Dalcin     parentviewer = viewer;
23211c56cb0SLisandro Dalcin     ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
23311c56cb0SLisandro Dalcin   }
23411c56cb0SLisandro Dalcin 
235*3b46f529SLisandro Dalcin   if (!rank) {
236331830f3SMatthew G. Knepley     PetscBool match;
237331830f3SMatthew G. Knepley     int       fileType, dataSize;
238f6ac7a6aSMichael Lange     float     version;
239331830f3SMatthew G. Knepley 
240331830f3SMatthew G. Knepley     /* Read header */
241060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
24204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
243331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
244060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
245f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
246f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
247f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
248331830f3SMatthew 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);
24904d1ad83SMichael Lange     if (binary) {
250b9eae255SMichael Lange       int checkInt;
251060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
25204d1ad83SMichael Lange       if (checkInt != 1) {
253b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
25411c56cb0SLisandro Dalcin         if (checkInt == 1) byteSwap = PETSC_TRUE;
25504d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
25604d1ad83SMichael Lange       }
2570877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
258060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
25904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
260331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
26111c56cb0SLisandro Dalcin 
262dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
263060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
264dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
265dc0ede3bSMatthew G. Knepley     if (match) {
266dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
267dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
268dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
26911c56cb0SLisandro Dalcin       for (r = 0; r < numRegions; ++r) {
27011c56cb0SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
27111c56cb0SLisandro Dalcin       }
272dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
273dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
274dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
275dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
276dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
277dc0ede3bSMatthew G. Knepley     }
27811c56cb0SLisandro Dalcin 
279dc0ede3bSMatthew G. Knepley     /* Read vertices */
28004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
281331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
282060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
2830877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
2840877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
28511c56cb0SLisandro Dalcin     ierr = DMPlexCreateGmsh_ReadNodes(viewer, binary, byteSwap, shift, numVertices, &coordsIn);CHKERRQ(ierr);
286060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
28704d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
288331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
28911c56cb0SLisandro Dalcin 
290331830f3SMatthew G. Knepley     /* Read cells */
291060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
29204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
293331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
294060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
2950877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
2960877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
297a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
298a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
299a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
300a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
30111c56cb0SLisandro Dalcin     ierr = DMPlexCreateGmsh_ReadElements(viewer, binary, byteSwap, shift, numCells, &gmsh_elem);CHKERRQ(ierr);
302a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
303a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
304a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
305db397164SMichael Lange     }
306d08df55aSStefano Zampini     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
3071b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
308331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
30911c56cb0SLisandro Dalcin 
31011c56cb0SLisandro Dalcin     /* OPTIONAL Read periodic section */
311d08df55aSStefano Zampini     if (periodic) {
312f45c9589SStefano Zampini       PetscInt pVert, *periodicMapT, *aux;
313d08df55aSStefano Zampini       int      numPeriodic;
314d08df55aSStefano Zampini 
315d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
316d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
317d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
318f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
3196fbe17bfSStefano Zampini       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
320f45c9589SStefano Zampini       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
321d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
322d08df55aSStefano Zampini       snum = sscanf(line, "%d", &numPeriodic);
323d08df55aSStefano Zampini       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
324d08df55aSStefano Zampini       for (i = 0; i < numPeriodic; i++) {
325da83f57bSLisandro Dalcin         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes, slaveNode, masterNode;
326d08df55aSStefano Zampini 
327d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
328d08df55aSStefano Zampini         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
329d08df55aSStefano Zampini         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
330da83f57bSLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
331da83f57bSLisandro Dalcin         snum = sscanf(line, "%d", &nNodes);
332da83f57bSLisandro Dalcin         if (snum != 1) { /* discard tranformation and try again */
33372ffbcc9SStefano Zampini           ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
334d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
335d08df55aSStefano Zampini           snum = sscanf(line, "%d", &nNodes);
336d08df55aSStefano Zampini           if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
337da83f57bSLisandro Dalcin         }
338d08df55aSStefano Zampini         for (j = 0; j < nNodes; j++) {
339d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
34011c56cb0SLisandro Dalcin           snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
341d08df55aSStefano Zampini           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
342917266a4SLisandro Dalcin           periodicMapT[slaveNode - shift] = masterNode - shift;
343917266a4SLisandro Dalcin           ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr);
344917266a4SLisandro Dalcin           ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr);
345d08df55aSStefano Zampini         }
346d08df55aSStefano Zampini       }
347d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
348d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
349d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
350d08df55aSStefano Zampini       /* we may have slaves of slaves */
351d08df55aSStefano Zampini       for (i = 0; i < numVertices; i++) {
352f45c9589SStefano Zampini         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
353f45c9589SStefano Zampini           periodicMapT[i] = periodicMapT[periodicMapT[i]];
354d08df55aSStefano Zampini         }
355d08df55aSStefano Zampini       }
356f45c9589SStefano Zampini       /* periodicMap : from old to new numbering (periodic vertices excluded)
357f45c9589SStefano Zampini          periodicMapI: from new to old numbering */
358f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
359f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
360f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
361f45c9589SStefano Zampini       for (i = 0, pVert = 0; i < numVertices; i++) {
362f45c9589SStefano Zampini         if (periodicMapT[i] != i) {
363f45c9589SStefano Zampini           pVert++;
364f45c9589SStefano Zampini         } else {
365f45c9589SStefano Zampini           aux[i] = i - pVert;
366f45c9589SStefano Zampini           periodicMapI[i - pVert] = i;
367f45c9589SStefano Zampini         }
368f45c9589SStefano Zampini       }
369f45c9589SStefano Zampini       for (i = 0 ; i < numVertices; i++) {
370f45c9589SStefano Zampini         periodicMap[i] = aux[periodicMapT[i]];
371f45c9589SStefano Zampini       }
372f45c9589SStefano Zampini       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
373f45c9589SStefano Zampini       ierr = PetscFree(aux);CHKERRQ(ierr);
374f45c9589SStefano Zampini       /* remove periodic vertices */
375f45c9589SStefano Zampini       numVertices = numVertices - pVert;
376d08df55aSStefano Zampini     }
377331830f3SMatthew G. Knepley   }
37811c56cb0SLisandro Dalcin 
37911c56cb0SLisandro Dalcin   if (parentviewer) {
38011c56cb0SLisandro Dalcin     ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
38111c56cb0SLisandro Dalcin   }
38211c56cb0SLisandro Dalcin 
383a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
384331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
385a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
386a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
387a4bb7517SMichael Lange       ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
388a4bb7517SMichael Lange       cell++;
389331830f3SMatthew G. Knepley     }
390331830f3SMatthew G. Knepley   }
391331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
392a4bb7517SMichael Lange   /* Add cell-vertex connections */
393a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
394a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
39511c56cb0SLisandro Dalcin       PetscInt pcone[8], corner;
396a4bb7517SMichael Lange       for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
397dd169d64SMatthew G. Knepley         const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
398917266a4SLisandro Dalcin         pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells;
399331830f3SMatthew G. Knepley       }
40097ed6be6Ssarens       if (dim == 3) {
40197ed6be6Ssarens         /* Tetrahedra are inverted */
40297ed6be6Ssarens         if (gmsh_elem[c].numNodes == 4) {
40397ed6be6Ssarens           PetscInt tmp = pcone[0];
40497ed6be6Ssarens           pcone[0] = pcone[1];
40597ed6be6Ssarens           pcone[1] = tmp;
40697ed6be6Ssarens         }
40797ed6be6Ssarens         /* Hexahedra are inverted */
40897ed6be6Ssarens         if (gmsh_elem[c].numNodes == 8) {
40997ed6be6Ssarens           PetscInt tmp = pcone[1];
41097ed6be6Ssarens           pcone[1] = pcone[3];
41197ed6be6Ssarens           pcone[3] = tmp;
41297ed6be6Ssarens         }
41397ed6be6Ssarens       }
414a4bb7517SMichael Lange       ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
415a4bb7517SMichael Lange       cell++;
416331830f3SMatthew G. Knepley     }
417a4bb7517SMichael Lange   }
4186228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
419c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
420331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
421331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
422331830f3SMatthew G. Knepley   if (interpolate) {
4235fd9971aSMatthew G. Knepley     DM idm;
424331830f3SMatthew G. Knepley 
425331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
426331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
427331830f3SMatthew G. Knepley     *dm  = idm;
428331830f3SMatthew G. Knepley   }
429917266a4SLisandro Dalcin 
430f6c8a31fSLisandro Dalcin   if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation");
431917266a4SLisandro Dalcin   if (!rank && usemarker) {
432d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
433d3f73514SStefano Zampini 
434d3f73514SStefano Zampini     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
435d3f73514SStefano Zampini     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
436d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
437d3f73514SStefano Zampini       PetscInt suppSize;
438d3f73514SStefano Zampini 
439d3f73514SStefano Zampini       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
440d3f73514SStefano Zampini       if (suppSize == 1) {
441d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
442d3f73514SStefano Zampini 
443d3f73514SStefano Zampini         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
444d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
445d3f73514SStefano Zampini           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
446d3f73514SStefano Zampini         }
447d3f73514SStefano Zampini         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
448d3f73514SStefano Zampini       }
449d3f73514SStefano Zampini     }
450d3f73514SStefano Zampini   }
45116ad7e67SMichael Lange 
45216ad7e67SMichael Lange   if (!rank) {
45311c56cb0SLisandro Dalcin     PetscInt vStart, vEnd;
454d1a54cefSMatthew G. Knepley 
45516ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
45611c56cb0SLisandro Dalcin     for (cell = 0, c = 0; c < numCells; ++c) {
45711c56cb0SLisandro Dalcin 
45811c56cb0SLisandro Dalcin       /* Create face sets */
4595964b3dcSLisandro Dalcin       if (interpolate && gmsh_elem[c].dim == dim-1) {
46016ad7e67SMichael Lange         const PetscInt *join;
46111c56cb0SLisandro Dalcin         PetscInt        joinSize, pcone[8], corner;
46211c56cb0SLisandro Dalcin         /* Find the relevant facet with vertex joins */
463a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
464dd169d64SMatthew G. Knepley           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
465917266a4SLisandro Dalcin           pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart;
46616ad7e67SMichael Lange         }
46711c56cb0SLisandro Dalcin         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr);
468a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
469c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
470a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
47116ad7e67SMichael Lange       }
4726e1dd89aSLawrence Mitchell 
4736e1dd89aSLawrence Mitchell       /* Create cell sets */
4746e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
4756e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
4766e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4776e1dd89aSLawrence Mitchell         }
478917266a4SLisandro Dalcin         cell++;
4796e1dd89aSLawrence Mitchell       }
4800c070f12SSander Arens 
4810c070f12SSander Arens       /* Create vertex sets */
4820c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
4830c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
484917266a4SLisandro Dalcin           const PetscInt cc = gmsh_elem[c].nodes[0] - shift;
48511c56cb0SLisandro Dalcin           const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart;
486d08df55aSStefano Zampini           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4870c070f12SSander Arens         }
4880c070f12SSander Arens       }
4890c070f12SSander Arens     }
49016ad7e67SMichael Lange   }
49116ad7e67SMichael Lange 
49211c56cb0SLisandro Dalcin   /* Create coordinates */
49311c56cb0SLisandro Dalcin   embedDim = isbd ? dim+1 : dim;
494979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
495331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
4961d64f2ccSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
497f45c9589SStefano Zampini   if (periodicMap) { /* we need to localize coordinates on cells */
498f45c9589SStefano Zampini     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
499f45c9589SStefano Zampini   } else {
50075b5763bSMichael Lange     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
501f45c9589SStefano Zampini   }
50275b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
5031d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
5041d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
505331830f3SMatthew G. Knepley   }
50611c56cb0SLisandro Dalcin   if (periodicMap) {
507437bdd5bSStefano Zampini     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
508f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
509f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
510437bdd5bSStefano Zampini         PetscInt  corner;
51111c56cb0SLisandro Dalcin         PetscBool pc = PETSC_FALSE;
512437bdd5bSStefano Zampini         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
513917266a4SLisandro Dalcin           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift));
514437bdd5bSStefano Zampini         }
515437bdd5bSStefano Zampini         if (pc) {
51611c56cb0SLisandro Dalcin           PetscInt dof = gmsh_elem[c].numNodes*embedDim;
517437bdd5bSStefano Zampini           ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr);
51811c56cb0SLisandro Dalcin           ierr = PetscSectionSetDof(coordSection, cell, dof);CHKERRQ(ierr);
51911c56cb0SLisandro Dalcin           ierr = PetscSectionSetFieldDof(coordSection, cell, 0, dof);CHKERRQ(ierr);
5206fbe17bfSStefano Zampini         }
521f45c9589SStefano Zampini         cell++;
522f45c9589SStefano Zampini       }
523f45c9589SStefano Zampini     }
524f45c9589SStefano Zampini   }
525331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
526331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
5278b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
528331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
529331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
5301d64f2ccSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
531331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
532331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
533f45c9589SStefano Zampini   if (periodicMap) {
534f45c9589SStefano Zampini     PetscInt off;
535f45c9589SStefano Zampini 
536f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
537f45c9589SStefano Zampini       PetscInt pcone[8], corner;
538f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
5396fbe17bfSStefano Zampini         if (PetscUnlikely(PetscBTLookup(periodicC, cell))) {
540f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
541917266a4SLisandro Dalcin             pcone[corner] = gmsh_elem[c].nodes[corner] - shift;
542f45c9589SStefano Zampini           }
543f45c9589SStefano Zampini           if (dim == 3) {
544f45c9589SStefano Zampini             /* Tetrahedra are inverted */
545f45c9589SStefano Zampini             if (gmsh_elem[c].numNodes == 4) {
546f45c9589SStefano Zampini               PetscInt tmp = pcone[0];
547f45c9589SStefano Zampini               pcone[0] = pcone[1];
548f45c9589SStefano Zampini               pcone[1] = tmp;
549f45c9589SStefano Zampini             }
550f45c9589SStefano Zampini             /* Hexahedra are inverted */
551f45c9589SStefano Zampini             if (gmsh_elem[c].numNodes == 8) {
552f45c9589SStefano Zampini               PetscInt tmp = pcone[1];
553f45c9589SStefano Zampini               pcone[1] = pcone[3];
554f45c9589SStefano Zampini               pcone[3] = tmp;
555f45c9589SStefano Zampini             }
556f45c9589SStefano Zampini           }
557f45c9589SStefano Zampini           ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
558f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
55911c56cb0SLisandro Dalcin             v = pcone[corner];
560dd169d64SMatthew G. Knepley             for (d = 0; d < embedDim; ++d) {
56111c56cb0SLisandro Dalcin               coords[off++] = (PetscReal) coordsIn[v*3+d];
562f45c9589SStefano Zampini             }
563f45c9589SStefano Zampini           }
5646fbe17bfSStefano Zampini         }
565f45c9589SStefano Zampini         cell++;
566f45c9589SStefano Zampini       }
567f45c9589SStefano Zampini     }
568f45c9589SStefano Zampini     for (v = 0; v < numVertices; ++v) {
569f45c9589SStefano Zampini       ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
570dd169d64SMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
57111c56cb0SLisandro Dalcin         coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d];
572f45c9589SStefano Zampini       }
573f45c9589SStefano Zampini     }
574f45c9589SStefano Zampini   } else {
575331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
5761d64f2ccSMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
57711c56cb0SLisandro Dalcin         coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d];
578331830f3SMatthew G. Knepley       }
579331830f3SMatthew G. Knepley     }
580331830f3SMatthew G. Knepley   }
581331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
582331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
58311c56cb0SLisandro Dalcin   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
58411c56cb0SLisandro Dalcin 
58511c56cb0SLisandro Dalcin   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
58611c56cb0SLisandro Dalcin   ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
587d08df55aSStefano Zampini   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
588fcd9ca0aSStefano Zampini   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
5896fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
5906fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
59111c56cb0SLisandro Dalcin   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
59211c56cb0SLisandro Dalcin 
5933b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
594331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
595331830f3SMatthew G. Knepley }
596