xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision fcd9ca0a7dbc0ad76a157a10a362aad97c3cc452)
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 {
207d282ae0SMichael Lange   PetscViewer     viewer, vheader;
21abc86ac4SMichael Lange   PetscMPIInt     rank;
227d282ae0SMichael Lange   PetscViewerType vtype;
237d282ae0SMichael Lange   char            line[PETSC_MAX_PATH_LEN];
247d282ae0SMichael Lange   int             snum;
257d282ae0SMichael Lange   PetscBool       match;
26c1c22fd2SMatthew G. Knepley   int             fT;
27c1c22fd2SMatthew G. Knepley   PetscInt        fileType;
28f6ac7a6aSMichael Lange   float           version;
297d282ae0SMichael Lange   PetscErrorCode  ierr;
307d282ae0SMichael Lange 
317d282ae0SMichael Lange   PetscFunctionBegin;
32abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
337d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
347d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
357d282ae0SMichael Lange   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
367d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
377d282ae0SMichael Lange   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
38abc86ac4SMichael Lange   if (!rank) {
397d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
40060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
417d282ae0SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
427d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
43060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
44c1c22fd2SMatthew G. Knepley     snum = sscanf(line, "%f %d", &version, &fT);
45c1c22fd2SMatthew G. Knepley     fileType = (PetscInt) fT;
46f6ac7a6aSMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
47f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
48abc86ac4SMichael Lange   }
49abc86ac4SMichael Lange   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
507d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
517d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
527d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
537d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
547d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
557d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
577d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
587d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
597d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
607d282ae0SMichael Lange   PetscFunctionReturn(0);
617d282ae0SMichael Lange }
627d282ae0SMichael Lange 
63df032b82SMatthew G. Knepley static PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
64df032b82SMatthew G. Knepley {
65df032b82SMatthew G. Knepley   PetscInt       c, p;
66df032b82SMatthew G. Knepley   GmshElement   *elements;
67bf6ba3a3SMatthew G. Knepley   int            i, cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
68bf6ba3a3SMatthew G. Knepley   PetscInt       pibuf[64];
69df032b82SMatthew G. Knepley   int            ibuf[16];
70df032b82SMatthew G. Knepley   PetscErrorCode ierr;
71df032b82SMatthew G. Knepley 
72df032b82SMatthew G. Knepley   PetscFunctionBegin;
73df032b82SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
74df032b82SMatthew G. Knepley   for (c = 0; c < numCells;) {
75df032b82SMatthew G. Knepley     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
76df032b82SMatthew G. Knepley     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
77df032b82SMatthew G. Knepley     if (binary) {
78df032b82SMatthew G. Knepley       cellType = ibuf[0];
79df032b82SMatthew G. Knepley       numElem = ibuf[1];
80df032b82SMatthew G. Knepley       numTags = ibuf[2];
81df032b82SMatthew G. Knepley     } else {
82df032b82SMatthew G. Knepley       elements[c].id = ibuf[0];
83df032b82SMatthew G. Knepley       cellType = ibuf[1];
84df032b82SMatthew G. Knepley       numTags = ibuf[2];
85df032b82SMatthew G. Knepley       numElem = 1;
86df032b82SMatthew G. Knepley     }
87bf6ba3a3SMatthew G. Knepley     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
88bf6ba3a3SMatthew G. Knepley     numNodesIgnore = 0;
89df032b82SMatthew G. Knepley     switch (cellType) {
90df032b82SMatthew G. Knepley     case 1: /* 2-node line */
91df032b82SMatthew G. Knepley       dim = 1;
92df032b82SMatthew G. Knepley       numNodes = 2;
93df032b82SMatthew G. Knepley       break;
94df032b82SMatthew G. Knepley     case 2: /* 3-node triangle */
95df032b82SMatthew G. Knepley       dim = 2;
96df032b82SMatthew G. Knepley       numNodes = 3;
97df032b82SMatthew G. Knepley       break;
98df032b82SMatthew G. Knepley     case 3: /* 4-node quadrangle */
99df032b82SMatthew G. Knepley       dim = 2;
100df032b82SMatthew G. Knepley       numNodes = 4;
101df032b82SMatthew G. Knepley       break;
102df032b82SMatthew G. Knepley     case 4: /* 4-node tetrahedron */
103df032b82SMatthew G. Knepley       dim  = 3;
104df032b82SMatthew G. Knepley       numNodes = 4;
105df032b82SMatthew G. Knepley       break;
106df032b82SMatthew G. Knepley     case 5: /* 8-node hexahedron */
107df032b82SMatthew G. Knepley       dim = 3;
108df032b82SMatthew G. Knepley       numNodes = 8;
109df032b82SMatthew G. Knepley       break;
110bf6ba3a3SMatthew G. Knepley     case 8: /* 3-node 2nd order line */
111bf6ba3a3SMatthew G. Knepley       dim = 1;
112bf6ba3a3SMatthew G. Knepley       numNodes = 2;
113bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 1;
114bf6ba3a3SMatthew G. Knepley       break;
115bf6ba3a3SMatthew G. Knepley     case 9: /* 6-node 2nd order triangle */
116bf6ba3a3SMatthew G. Knepley       dim = 2;
117bf6ba3a3SMatthew G. Knepley       numNodes = 3;
118bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 3;
119bf6ba3a3SMatthew G. Knepley       break;
120df032b82SMatthew G. Knepley     case 15: /* 1-node vertex */
121df032b82SMatthew G. Knepley       dim = 0;
122df032b82SMatthew G. Knepley       numNodes = 1;
123df032b82SMatthew G. Knepley       break;
124bf6ba3a3SMatthew G. Knepley     case 6: /* 6-node prism */
125bf6ba3a3SMatthew G. Knepley     case 7: /* 5-node pyramid */
126bf6ba3a3SMatthew G. Knepley     case 10: /* 9-node 2nd order quadrangle */
127bf6ba3a3SMatthew G. Knepley     case 11: /* 10-node 2nd order tetrahedron */
128bf6ba3a3SMatthew G. Knepley     case 12: /* 27-node 2nd order hexhedron */
129bf6ba3a3SMatthew G. Knepley     case 13: /* 19-node 2nd order prism */
130bf6ba3a3SMatthew G. Knepley     case 14: /* 14-node 2nd order pyramid */
131df032b82SMatthew G. Knepley     default:
132df032b82SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
133df032b82SMatthew G. Knepley     }
134df032b82SMatthew G. Knepley     if (binary) {
135bf6ba3a3SMatthew G. Knepley       const PetscInt nint = numNodes + numTags + 1 + numNodesIgnore;
136df032b82SMatthew G. Knepley       for (i = 0; i < numElem; ++i, ++c) {
137df032b82SMatthew G. Knepley         /* Loop over inner binary element block */
138df032b82SMatthew G. Knepley         elements[c].dim = dim;
139df032b82SMatthew G. Knepley         elements[c].numNodes = numNodes;
140df032b82SMatthew G. Knepley         elements[c].numTags = numTags;
141df032b82SMatthew G. Knepley 
142df032b82SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
143df032b82SMatthew G. Knepley         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
144df032b82SMatthew G. Knepley         elements[c].id = ibuf[0];
145df032b82SMatthew G. Knepley         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
146df032b82SMatthew G. Knepley         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
147df032b82SMatthew G. Knepley       }
148df032b82SMatthew G. Knepley     } else {
149df032b82SMatthew G. Knepley       elements[c].dim = dim;
150df032b82SMatthew G. Knepley       elements[c].numNodes = numNodes;
151df032b82SMatthew G. Knepley       elements[c].numTags = numTags;
152df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
153df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
154bf6ba3a3SMatthew G. Knepley       ierr = PetscViewerRead(viewer, pibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr);
155df032b82SMatthew G. Knepley       c++;
156df032b82SMatthew G. Knepley     }
157df032b82SMatthew G. Knepley   }
158df032b82SMatthew G. Knepley   *gmsh_elems = elements;
159df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
160df032b82SMatthew G. Knepley }
1617d282ae0SMichael Lange 
162331830f3SMatthew G. Knepley /*@
1637d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
164331830f3SMatthew G. Knepley 
165331830f3SMatthew G. Knepley   Collective on comm
166331830f3SMatthew G. Knepley 
167331830f3SMatthew G. Knepley   Input Parameters:
168331830f3SMatthew G. Knepley + comm  - The MPI communicator
169331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
170331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
171331830f3SMatthew G. Knepley 
172331830f3SMatthew G. Knepley   Output Parameter:
173331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
174331830f3SMatthew G. Knepley 
175331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
1763d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
177331830f3SMatthew G. Knepley 
178331830f3SMatthew G. Knepley   Level: beginner
179331830f3SMatthew G. Knepley 
180331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
181331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
182331830f3SMatthew G. Knepley @*/
183331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
184331830f3SMatthew G. Knepley {
18504d1ad83SMichael Lange   PetscViewerType vtype;
1863c67d5adSSatish Balay   GmshElement   *gmsh_elem = NULL;
187331830f3SMatthew G. Knepley   PetscSection   coordSection;
188331830f3SMatthew G. Knepley   Vec            coordinates;
189331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
190f45c9589SStefano Zampini   PetscInt       dim = 0, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL;
191dc0ede3bSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum;
192331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
193331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
194d08df55aSStefano Zampini   PetscBool      match, binary, bswap = PETSC_FALSE, periodic = PETSC_FALSE;
195331830f3SMatthew G. Knepley   PetscErrorCode ierr;
196331830f3SMatthew G. Knepley 
197331830f3SMatthew G. Knepley   PetscFunctionBegin;
198331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
199331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
200331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
201331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2023b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
20304d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
20404d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
205d08df55aSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_periodic_gmsh",&periodic,NULL);CHKERRQ(ierr);
206abc86ac4SMichael Lange   if (!rank || binary) {
207331830f3SMatthew G. Knepley     PetscBool match;
208331830f3SMatthew G. Knepley     int       fileType, dataSize;
209f6ac7a6aSMichael Lange     float     version;
210331830f3SMatthew G. Knepley 
211331830f3SMatthew G. Knepley     /* Read header */
212060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
21304d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
214331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
215060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
216f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
217f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
218f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
219331830f3SMatthew 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);
22004d1ad83SMichael Lange     if (binary) {
221b9eae255SMichael Lange       int checkInt;
222060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
22304d1ad83SMichael Lange       if (checkInt != 1) {
224b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
22504d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
22604d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
22704d1ad83SMichael Lange       }
2280877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
229060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
23004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
231331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
232dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
233060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
234dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
235dc0ede3bSMatthew G. Knepley     if (match) {
236dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
237dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
238dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
239a8667449SSander Arens       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
240dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
241dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
242dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
243dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
244dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
245dc0ede3bSMatthew G. Knepley     }
246dc0ede3bSMatthew G. Knepley     /* Read vertices */
24704d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
248331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
249060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
2500877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
2510877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
252854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
25313aecfe5SMichael Lange     if (binary) {
25413aecfe5SMichael Lange       size_t doubleSize, intSize;
25513aecfe5SMichael Lange       PetscInt elementSize;
25613aecfe5SMichael Lange       char *buffer;
25713aecfe5SMichael Lange       PetscScalar *baseptr;
25813aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
25913aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
26013aecfe5SMichael Lange       elementSize = (intSize + 3*doubleSize);
26113aecfe5SMichael Lange       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
26213aecfe5SMichael Lange       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
26313aecfe5SMichael Lange       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
264331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
26513aecfe5SMichael Lange         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
26613aecfe5SMichael Lange         coordsIn[v*3+0] = baseptr[0];
26713aecfe5SMichael Lange         coordsIn[v*3+1] = baseptr[1];
26813aecfe5SMichael Lange         coordsIn[v*3+2] = baseptr[2];
26913aecfe5SMichael Lange       }
27013aecfe5SMichael Lange       ierr = PetscFree(buffer);CHKERRQ(ierr);
27113aecfe5SMichael Lange     } else {
27213aecfe5SMichael Lange       for (v = 0; v < numVertices; ++v) {
273060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
274060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
275b9eae255SMichael Lange         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
276331830f3SMatthew G. Knepley       }
27713aecfe5SMichael Lange     }
278060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
27904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
280331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
281331830f3SMatthew G. Knepley     /* Read cells */
282060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
28304d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
284331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
285060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
2860877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
2870877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
288331830f3SMatthew G. Knepley   }
289db397164SMichael Lange 
290abc86ac4SMichael Lange   if (!rank || binary) {
291a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
292a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
293a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
294a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
2953d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
296a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
297a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
298a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
299db397164SMichael Lange     }
300d08df55aSStefano Zampini     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
3011b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
302331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
303d08df55aSStefano Zampini     if (periodic) {
304f45c9589SStefano Zampini       PetscInt pVert, *periodicMapT, *aux;
305d08df55aSStefano Zampini       int      numPeriodic;
306d08df55aSStefano Zampini 
307d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
308d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
309d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
310f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
311f45c9589SStefano Zampini       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
312d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
313d08df55aSStefano Zampini       snum = sscanf(line, "%d", &numPeriodic);
314d08df55aSStefano Zampini       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
315d08df55aSStefano Zampini       for (i = 0; i < numPeriodic; i++) {
316d08df55aSStefano Zampini         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes;
317d08df55aSStefano Zampini 
318d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
319d08df55aSStefano Zampini         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
320d08df55aSStefano Zampini         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
321d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
322d08df55aSStefano Zampini         snum = sscanf(line, "%d", &nNodes);
323d08df55aSStefano Zampini         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
324d08df55aSStefano Zampini         for (j = 0; j < nNodes; j++) {
325d08df55aSStefano Zampini           PetscInt slaveNode, masterNode;
326d08df55aSStefano Zampini 
327d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
328d08df55aSStefano Zampini           snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
329d08df55aSStefano Zampini           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
330f45c9589SStefano Zampini           periodicMapT[slaveNode - 1] = masterNode - 1;
331d08df55aSStefano Zampini         }
332d08df55aSStefano Zampini       }
333d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
334d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
335d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
336d08df55aSStefano Zampini       /* we may have slaves of slaves */
337d08df55aSStefano Zampini       for (i = 0; i < numVertices; i++) {
338f45c9589SStefano Zampini         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
339f45c9589SStefano Zampini           periodicMapT[i] = periodicMapT[periodicMapT[i]];
340d08df55aSStefano Zampini         }
341d08df55aSStefano Zampini       }
342f45c9589SStefano Zampini       /* periodicMap : from old to new numbering (periodic vertices excluded)
343f45c9589SStefano Zampini          periodicMapI: from new to old numbering */
344f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
345f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
346f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
347f45c9589SStefano Zampini       for (i = 0, pVert = 0; i < numVertices; i++) {
348f45c9589SStefano Zampini         if (periodicMapT[i] != i) {
349f45c9589SStefano Zampini           pVert++;
350f45c9589SStefano Zampini         } else {
351f45c9589SStefano Zampini           aux[i] = i - pVert;
352f45c9589SStefano Zampini           periodicMapI[i - pVert] = i;
353f45c9589SStefano Zampini         }
354f45c9589SStefano Zampini       }
355f45c9589SStefano Zampini       for (i = 0 ; i < numVertices; i++) {
356f45c9589SStefano Zampini         periodicMap[i] = aux[periodicMapT[i]];
357f45c9589SStefano Zampini       }
358f45c9589SStefano Zampini       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
359f45c9589SStefano Zampini       ierr = PetscFree(aux);CHKERRQ(ierr);
360f45c9589SStefano Zampini       /* remove periodic vertices */
361f45c9589SStefano Zampini       numVertices = numVertices - pVert;
362d08df55aSStefano Zampini     }
363331830f3SMatthew G. Knepley   }
364abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
365abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
366a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
367331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
368331830f3SMatthew G. Knepley   if (!rank) {
369a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
370a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
371a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
372a4bb7517SMichael Lange         cell++;
373331830f3SMatthew G. Knepley       }
374331830f3SMatthew G. Knepley     }
375331830f3SMatthew G. Knepley   }
376331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
377a4bb7517SMichael Lange   /* Add cell-vertex connections */
378331830f3SMatthew G. Knepley   if (!rank) {
379331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
380a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
381a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
382a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
383d08df55aSStefano Zampini           pcone[corner] = periodicMap ? periodicMap[gmsh_elem[c].nodes[corner] - 1] + trueNumCells
384d08df55aSStefano Zampini                                       : gmsh_elem[c].nodes[corner] -1 + trueNumCells;
385331830f3SMatthew G. Knepley         }
38697ed6be6Ssarens         if (dim == 3) {
38797ed6be6Ssarens           /* Tetrahedra are inverted */
38897ed6be6Ssarens           if (gmsh_elem[c].numNodes == 4) {
38997ed6be6Ssarens             PetscInt tmp = pcone[0];
39097ed6be6Ssarens             pcone[0] = pcone[1];
39197ed6be6Ssarens             pcone[1] = tmp;
39297ed6be6Ssarens           }
39397ed6be6Ssarens           /* Hexahedra are inverted */
39497ed6be6Ssarens           if (gmsh_elem[c].numNodes == 8) {
39597ed6be6Ssarens             PetscInt tmp = pcone[1];
39697ed6be6Ssarens             pcone[1] = pcone[3];
39797ed6be6Ssarens             pcone[3] = tmp;
39897ed6be6Ssarens           }
39997ed6be6Ssarens         }
400a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
401a4bb7517SMichael Lange         cell++;
402331830f3SMatthew G. Knepley       }
403a4bb7517SMichael Lange     }
404331830f3SMatthew G. Knepley   }
4056228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
406c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
407331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
408331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
409331830f3SMatthew G. Knepley   if (interpolate) {
410e56d480eSMatthew G. Knepley     DM idm = NULL;
411331830f3SMatthew G. Knepley 
412331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
413331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
414331830f3SMatthew G. Knepley     *dm  = idm;
415331830f3SMatthew G. Knepley   }
41616ad7e67SMichael Lange 
41716ad7e67SMichael Lange   if (!rank) {
41816ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
41916ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
420d1a54cefSMatthew G. Knepley 
42116ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
42216ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
423a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
42416ad7e67SMichael Lange         PetscInt joinSize;
42516ad7e67SMichael Lange         const PetscInt *join;
426a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
427d08df55aSStefano Zampini           pcone[corner] = periodicMap ? periodicMap[gmsh_elem[c].nodes[corner] - 1] + vStart
428d08df55aSStefano Zampini                                       : gmsh_elem[c].nodes[corner] - 1 + vStart;
42916ad7e67SMichael Lange         }
430a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
431a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
432c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
433a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
43416ad7e67SMichael Lange       }
435a4bb7517SMichael Lange     }
4366e1dd89aSLawrence Mitchell 
4376e1dd89aSLawrence Mitchell     /* Create cell sets */
4386e1dd89aSLawrence Mitchell     for (cell = 0, c = 0; c < numCells; ++c) {
4396e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
4406e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
4416e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4426e1dd89aSLawrence Mitchell           cell++;
4436e1dd89aSLawrence Mitchell         }
4446e1dd89aSLawrence Mitchell       }
4456e1dd89aSLawrence Mitchell     }
4460c070f12SSander Arens 
4470c070f12SSander Arens     /* Create vertex sets */
4480c070f12SSander Arens     for (c = 0; c < numCells; ++c) {
4490c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
4500c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
451d08df55aSStefano Zampini           PetscInt vid =  periodicMap ? periodicMap[gmsh_elem[c].nodes[0] - 1] + vStart
452d08df55aSStefano Zampini                                       : gmsh_elem[c].nodes[0] - 1 + vStart;
453d08df55aSStefano Zampini           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4540c070f12SSander Arens         }
4550c070f12SSander Arens       }
4560c070f12SSander Arens     }
45716ad7e67SMichael Lange   }
45816ad7e67SMichael Lange 
459331830f3SMatthew G. Knepley   /* Read coordinates */
460979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
461331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
462331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
463f45c9589SStefano Zampini   if (periodicMap) { /* we need to localize coordinates on cells */
464f45c9589SStefano Zampini     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
465f45c9589SStefano Zampini   } else {
46675b5763bSMichael Lange     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
467f45c9589SStefano Zampini   }
46875b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
469331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
470331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
471331830f3SMatthew G. Knepley   }
472f45c9589SStefano Zampini   if (periodicMap) {
473f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
474f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
475f45c9589SStefano Zampini         ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
476f45c9589SStefano Zampini         ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
477f45c9589SStefano Zampini         cell++;
478f45c9589SStefano Zampini       }
479f45c9589SStefano Zampini     }
480f45c9589SStefano Zampini   }
481331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
482331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
4838b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
484331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
485331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4867635fff5Ssarens   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
487331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
488331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
489331830f3SMatthew G. Knepley   if (!rank) {
490f45c9589SStefano Zampini 
491f45c9589SStefano Zampini     if (periodicMap) {
492f45c9589SStefano Zampini       PetscInt off;
493f45c9589SStefano Zampini 
494f45c9589SStefano Zampini       for (cell = 0, c = 0; c < numCells; ++c) {
495f45c9589SStefano Zampini         PetscInt pcone[8], corner;
496f45c9589SStefano Zampini         if (gmsh_elem[c].dim == dim) {
497f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
498f45c9589SStefano Zampini             pcone[corner] = gmsh_elem[c].nodes[corner] - 1;
499f45c9589SStefano Zampini           }
500f45c9589SStefano Zampini           if (dim == 3) {
501f45c9589SStefano Zampini             /* Tetrahedra are inverted */
502f45c9589SStefano Zampini             if (gmsh_elem[c].numNodes == 4) {
503f45c9589SStefano Zampini               PetscInt tmp = pcone[0];
504f45c9589SStefano Zampini               pcone[0] = pcone[1];
505f45c9589SStefano Zampini               pcone[1] = tmp;
506f45c9589SStefano Zampini             }
507f45c9589SStefano Zampini             /* Hexahedra are inverted */
508f45c9589SStefano Zampini             if (gmsh_elem[c].numNodes == 8) {
509f45c9589SStefano Zampini               PetscInt tmp = pcone[1];
510f45c9589SStefano Zampini               pcone[1] = pcone[3];
511f45c9589SStefano Zampini               pcone[3] = tmp;
512f45c9589SStefano Zampini             }
513f45c9589SStefano Zampini           }
514f45c9589SStefano Zampini           ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
515f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
516f45c9589SStefano Zampini             PetscInt v = pcone[corner];
517f45c9589SStefano Zampini             for (d = 0; d < dim; ++d) {
518f45c9589SStefano Zampini               coords[off++] = coordsIn[v*3+d];
519f45c9589SStefano Zampini             }
520f45c9589SStefano Zampini           }
521f45c9589SStefano Zampini           cell++;
522f45c9589SStefano Zampini         }
523f45c9589SStefano Zampini       }
524f45c9589SStefano Zampini       for (v = 0; v < numVertices; ++v) {
525f45c9589SStefano Zampini         ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
526f45c9589SStefano Zampini         for (d = 0; d < dim; ++d) {
527f45c9589SStefano Zampini           coords[off+d] = coordsIn[periodicMapI[v]*3+d];
528f45c9589SStefano Zampini         }
529f45c9589SStefano Zampini       }
530f45c9589SStefano Zampini     } else {
531331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
532331830f3SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
533331830f3SMatthew G. Knepley           coords[v*dim+d] = coordsIn[v*3+d];
534331830f3SMatthew G. Knepley         }
535331830f3SMatthew G. Knepley       }
536331830f3SMatthew G. Knepley     }
537f45c9589SStefano Zampini   }
53890b157c4SStefano Zampini   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
539331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
540331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
541331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
542331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
543d08df55aSStefano Zampini   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
544*fcd9ca0aSStefano Zampini   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
545a4bb7517SMichael Lange   /* Clean up intermediate storage */
54629403c87SMichael Lange   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
5473b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
548331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
549331830f3SMatthew G. Knepley }
550