xref: /petsc/src/dm/impls/plex/plexfluent.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMPlexCreateFluentFromFile"
6 /*@C
7   DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file
8 
9 + comm        - The MPI communicator
10 . filename    - Name of the Fluent mesh file
11 - interpolate - Create faces and edges in the mesh
12 
13   Output Parameter:
14 . dm  - The DM object representing the mesh
15 
16   Level: beginner
17 
18 .seealso: DMPlexCreateFromFile(), DMPlexCreateFluent(), DMPlexCreate()
19 @*/
20 PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
21 {
22   PetscViewer     viewer;
23   PetscErrorCode  ierr;
24 
25   PetscFunctionBegin;
26   /* Create file viewer and build plex */
27   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
28   ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
29   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
30   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
31   ierr = DMPlexCreateFluent(comm, viewer, interpolate, dm);CHKERRQ(ierr);
32   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
33   PetscFunctionReturn(0);
34 }
35 
36 #undef __FUNCT__
37 #define __FUNCT__ "DMPlexCreateFluent_ReadString"
38 PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim)
39 {
40   PetscInt i = 0;
41   PetscErrorCode ierr;
42 
43   PetscFunctionBegin;
44   do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, NULL, PETSC_CHAR);CHKERRQ(ierr);}
45   while (buffer[i-1] != '\0' && buffer[i-1] != delim);
46   buffer[i] = '\0';
47   PetscFunctionReturn(0);
48 }
49 
50 #undef __FUNCT__
51 #define __FUNCT__ "DMPlexCreateFluent_ReadValues"
52 PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary)
53 {
54   int            fdes=0;
55   FILE          *file;
56   PetscInt       i;
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   if (binary) {
61     /* Extract raw file descriptor to read binary block */
62     ierr = PetscViewerASCIIGetPointer(viewer, &file);CHKERRQ(ierr);
63     fflush(file); fdes = fileno(file);
64   }
65 
66   if (!binary && dtype == PETSC_INT) {
67     char     cbuf[256];
68     int      ibuf, snum;
69     /* Parse hexadecimal ascii integers */
70     for (i = 0; i < count; i++) {
71       ierr = PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
72       snum = sscanf(cbuf, "%x", &ibuf);
73       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
74       ((PetscInt*)data)[i] = (PetscInt)ibuf;
75     }
76   } else if (binary && dtype == PETSC_INT) {
77     /* Always read 32-bit ints and cast to PetscInt */
78     int *ibuf;
79     ierr = PetscMalloc1(count, &ibuf);CHKERRQ(ierr);
80     ierr = PetscBinaryRead(fdes, ibuf, count, PETSC_ENUM);CHKERRQ(ierr);
81     ierr = PetscByteSwap(ibuf, PETSC_ENUM, count);CHKERRQ(ierr);
82     for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]);
83     ierr = PetscFree(ibuf);CHKERRQ(ierr);
84 
85  } else if (binary && dtype == PETSC_SCALAR) {
86     float *fbuf;
87     /* Always read 32-bit floats and cast to PetscScalar */
88     ierr = PetscMalloc1(count, &fbuf);CHKERRQ(ierr);
89     ierr = PetscBinaryRead(fdes, fbuf, count, PETSC_FLOAT);CHKERRQ(ierr);
90     ierr = PetscByteSwap(fbuf, PETSC_FLOAT, count);CHKERRQ(ierr);
91     for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]);
92     ierr = PetscFree(fbuf);CHKERRQ(ierr);
93   } else {
94     ierr = PetscViewerASCIIRead(viewer, data, count, NULL, dtype);CHKERRQ(ierr);
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 #undef __FUNCT__
100 #define __FUNCT__ "DMPlexCreateFluent_ReadSection"
101 PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s)
102 {
103   char           buffer[PETSC_MAX_PATH_LEN];
104   int            snum;
105   PetscErrorCode ierr;
106 
107   PetscFunctionBegin;
108   /* Fast-forward to next section and derive its index */
109   ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
110   ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr);
111   snum = sscanf(buffer, "%d", &(s->index));
112   /* If we can't match an index return -1 to signal end-of-file */
113   if (snum < 1) {s->index = -1;   PetscFunctionReturn(0);}
114 
115   if (s->index == 0) {           /* Comment */
116     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
117 
118   } else if (s->index == 2) {    /* Dimension */
119     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
120     snum = sscanf(buffer, "%d", &(s->nd));
121     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
122 
123   } else if (s->index == 10 || s->index == 2010) {   /* Vertices */
124     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
125     snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
126     if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
127     if (s->zoneID > 0) {
128       PetscInt numCoords = s->last - s->first + 1;
129       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
130       ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr);
131       ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
132       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
133     }
134     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
135 
136   } else if (s->index == 12 || s->index == 2012) {   /* Cells */
137     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
138     snum = sscanf(buffer, "(%x", &(s->zoneID));
139     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
140     if (s->zoneID == 0) {  /* Header section */
141       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
142       if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
143     } else {               /* Data section */
144       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
145       if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
146       if (s->nd == 0) {
147         /* Read cell type definitions for mixed cells */
148         PetscInt numCells = s->last - s->first + 1;
149         ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
150         ierr = PetscMalloc1(numCells, (PetscInt**)&s->data);CHKERRQ(ierr);
151         ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
152         ierr = PetscFree(s->data);CHKERRQ(ierr);
153         ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
154       }
155     }
156     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
157 
158   } else if (s->index == 13 || s->index == 2013) {   /* Faces */
159     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
160     snum = sscanf(buffer, "(%x", &(s->zoneID));
161     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
162     if (s->zoneID == 0) {  /* Header section */
163       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
164       if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
165     } else {               /* Data section */
166       PetscInt f, numEntries, numFaces, numFaceVert;
167       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
168       if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
169       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
170       switch (s->nd) {
171       case 0: numEntries = PETSC_DETERMINE; break;
172       case 2: numEntries = 2 + 2; break;  /* linear */
173       case 3: numEntries = 2 + 3; break;  /* triangular */
174       case 4: numEntries = 2 + 4; break;  /* quadrilateral */
175       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file");
176       }
177       numFaces = s->last-s->first + 1;
178       if (numEntries != PETSC_DETERMINE) {
179         /* Allocate space only if we already know the size of the block */
180         ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr);
181       }
182       for (f = 0; f < numFaces; f++) {
183         if (s->nd == 0) {
184           /* Determine the size of the block for "mixed" facets */
185           ierr = DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
186           if (numEntries == PETSC_DETERMINE) {
187             numEntries = numFaceVert + 2;
188             ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr);
189           } else {
190             if (numEntries != numFaceVert + 2) {
191               SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files");
192             }
193           }
194         }
195         ierr = DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
196       }
197       s->nd = numEntries - 2;
198       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
199     }
200     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
201 
202   } else {                       /* Unknown section type */
203     PetscInt depth = 1;
204     do {
205       /* Match parentheses when parsing unknown sections */
206       do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR);CHKERRQ(ierr);}
207       while (buffer[0] != '(' && buffer[0] != ')');
208       if (buffer[0] == '(') depth++;
209       if (buffer[0] == ')') depth--;
210     } while (depth > 0);
211     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr);
212   }
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "DMPlexCreateFluent"
218 /*@C
219   DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file.
220 
221   Collective on comm
222 
223   Input Parameters:
224 + comm  - The MPI communicator
225 . viewer - The Viewer associated with a Fluent mesh file
226 - interpolate - Create faces and edges in the mesh
227 
228   Output Parameter:
229 . dm  - The DM object representing the mesh
230 
231   Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm
232 
233   Level: beginner
234 
235 .keywords: mesh, fluent, case
236 .seealso: DMPLEX, DMCreate()
237 @*/
238 PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
239 {
240   PetscMPIInt    rank;
241   PetscInt       c, f, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE;
242   PetscInt       numFaces = PETSC_DETERMINE, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE;
243   PetscInt      *faces = NULL, *cellVertices, *faceZoneIDs = NULL;
244   PetscInt       d, coordSize;
245   PetscScalar   *coords, *coordsIn = NULL;
246   PetscSection   coordSection;
247   Vec            coordinates;
248   PetscErrorCode ierr;
249 
250   PetscFunctionBegin;
251   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
252 
253   if (!rank) {
254     FluentSection s;
255     do {
256       ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr);
257       if (s.index == 2) {                 /* Dimension */
258         dim = s.nd;
259 
260       } else if (s.index == 10 || s.index == 2010) { /* Vertices */
261         if (s.zoneID == 0) numVertices = s.last;
262         else {
263           if (coordsIn) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files");
264           coordsIn = (PetscScalar *) s.data;
265         }
266 
267       } else if (s.index == 12 || s.index == 2012) { /* Cells */
268         if (s.zoneID == 0) numCells = s.last;
269         else {
270           switch (s.nd) {
271           case 0: numCellVertices = PETSC_DETERMINE; break;
272           case 1: numCellVertices = 3; break;  /* triangular */
273           case 2: numCellVertices = 4; break;  /* tetrahedral */
274           case 3: numCellVertices = 4; break;  /* quadrilateral */
275           case 4: numCellVertices = 8; break;  /* hexahedral */
276           case 5: numCellVertices = 5; break;  /* pyramid */
277           case 6: numCellVertices = 6; break;  /* wedge */
278           default: numCellVertices = PETSC_DETERMINE;
279           }
280         }
281 
282       } else if (s.index == 13 || s.index == 2013) { /* Facets */
283         if (s.zoneID == 0) {  /* Header section */
284           numFaces = s.last - s.first + 1;
285           if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE;
286           else numFaceVertices = s.nd;
287         } else {              /* Data section */
288           if (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices) {
289             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported");
290           }
291           if (numFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file");
292           if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd;
293           numFaceEntries = numFaceVertices + 2;
294           if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);}
295           if (!faceZoneIDs) {ierr = PetscMalloc1(numFaces, &faceZoneIDs);CHKERRQ(ierr);}
296           ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr);
297           /* Record the zoneID for each face set */
298           for (f = s.first -1; f < s.last; f++) faceZoneIDs[f] = s.zoneID;
299           ierr = PetscFree(s.data);CHKERRQ(ierr);
300         }
301       }
302     } while (s.index >= 0);
303   }
304   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
305   if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension");
306 
307   /* Allocate cell-vertex mesh */
308   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
309   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
310   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
311   ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr);
312   if (!rank) {
313     if (numCells < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file");
314     /* If no cell type was given we assume simplices */
315     if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1;
316     for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);}
317   }
318   ierr = DMSetUp(*dm);CHKERRQ(ierr);
319 
320   if (!rank) {
321     /* Derive cell-vertex list from face-vertex and face-cell maps */
322     ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr);
323     for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1;
324     for (f = 0; f < numFaces; f++) {
325       PetscInt *cell;
326       const PetscInt cl = faces[f*numFaceEntries + numFaceVertices];
327       const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1];
328       const PetscInt *face = &(faces[f*numFaceEntries]);
329 
330       if (cl > 0) {
331         cell = &(cellVertices[(cl-1) * numCellVertices]);
332         for (v = 0; v < numFaceVertices; v++) {
333           PetscBool found = PETSC_FALSE;
334           for (c = 0; c < numCellVertices; c++) {
335             if (cell[c] < 0) break;
336             if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;}
337           }
338           if (!found) cell[c] = face[v]-1 + numCells;
339         }
340       }
341       if (cr > 0) {
342         cell = &(cellVertices[(cr-1) * numCellVertices]);
343         for (v = 0; v < numFaceVertices; v++) {
344           PetscBool found = PETSC_FALSE;
345           for (c = 0; c < numCellVertices; c++) {
346             if (cell[c] < 0) break;
347             if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;}
348           }
349           if (!found) cell[c] = face[v]-1 + numCells;
350         }
351       }
352     }
353     for (c = 0; c < numCells; c++) {
354       ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr);
355     }
356   }
357   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
358   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
359   if (interpolate) {
360     DM idm = NULL;
361 
362     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
363     ierr = DMDestroy(dm);CHKERRQ(ierr);
364     *dm  = idm;
365   }
366 
367   if (!rank) {
368     PetscInt fi, joinSize, meetSize, *fverts, cells[2];
369     const PetscInt *join, *meet;
370     ierr = PetscMalloc1(numFaceVertices, &fverts);CHKERRQ(ierr);
371     /* Mark facets by finding the full join of all adjacent vertices */
372     for (f = 0; f < numFaces; f++) {
373       const PetscInt cl = faces[f*numFaceEntries + numFaceVertices] - 1;
374       const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1] - 1;
375       if (cl > 0 && cr > 0) {
376         /* If we know both adjoining cells we can use a single-level meet */
377         cells[0] = cl; cells[1] = cr;
378         ierr = DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet);CHKERRQ(ierr);
379         if (meetSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f);
380         ierr = DMPlexSetLabelValue(*dm, "Face Sets", meet[0], faceZoneIDs[f]);CHKERRQ(ierr);
381         ierr = DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet);CHKERRQ(ierr);
382       } else {
383         for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1;
384         ierr = DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr);
385         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f);
386         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], faceZoneIDs[f]);CHKERRQ(ierr);
387         ierr = DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr);
388       }
389     }
390     ierr = PetscFree(fverts);CHKERRQ(ierr);
391   }
392 
393   /* Read coordinates */
394   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
395   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
396   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
397   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
398   for (v = numCells; v < numCells+numVertices; ++v) {
399     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
400     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
401   }
402   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
403   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
404   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
405   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
406   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
407   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
408   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
409   if (!rank) {
410     for (v = 0; v < numVertices; ++v) {
411       for (d = 0; d < dim; ++d) {
412         coords[v*dim+d] = coordsIn[v*dim+d];
413       }
414     }
415   }
416   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
417   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
418   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
419   if (!rank) {
420     ierr = PetscFree(cellVertices);CHKERRQ(ierr);
421     ierr = PetscFree(faces);CHKERRQ(ierr);
422     ierr = PetscFree(faceZoneIDs);CHKERRQ(ierr);
423     ierr = PetscFree(coordsIn);CHKERRQ(ierr);
424   }
425   PetscFunctionReturn(0);
426 }
427