xref: /petsc/src/dm/impls/plex/plexcreate.c (revision c8c68bd81ffedb5843f812e67c4542d57b420fc5)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 #include <petscdmda.h>
4 #include <petscsf.h>
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "DMPlexCreateDoublet"
8 /*@
9   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
10 
11   Collective on MPI_Comm
12 
13   Input Parameters:
14 + comm - The communicator for the DM object
15 . dim - The spatial dimension
16 . simplex - Flag for simplicial cells, otherwise they are tensor product cells
17 . interpolate - Flag to create intermediate mesh pieces (edges, faces)
18 . refinementUniform - Flag for uniform parallel refinement
19 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
20 
21   Output Parameter:
22 . dm  - The DM object
23 
24   Level: beginner
25 
26 .keywords: DM, create
27 .seealso: DMSetType(), DMCreate()
28 @*/
29 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
30 {
31   DM             dm;
32   PetscInt       p;
33   PetscMPIInt    rank;
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
38   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
39   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
40   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
41   switch (dim) {
42   case 2:
43     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
44     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
45     break;
46   case 3:
47     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
48     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
49     break;
50   default:
51     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
52   }
53   if (rank) {
54     PetscInt numPoints[2] = {0, 0};
55     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
56   } else {
57     switch (dim) {
58     case 2:
59       if (simplex) {
60         PetscInt    numPoints[2]        = {4, 2};
61         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
62         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
63         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
64         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
65         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
66 
67         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
68         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
69       } else {
70         PetscInt    numPoints[2]        = {6, 2};
71         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
72         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
73         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
74         PetscScalar vertexCoords[12]    = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  -1.0, 0.5,  1.0, -0.5,  1.0, 0.5};
75 
76         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
77       }
78       break;
79     case 3:
80       if (simplex) {
81         PetscInt    numPoints[2]        = {5, 2};
82         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
83         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
84         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
85         PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  1.0, 0.0, 0.0};
86         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
87 
88         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
89         for (p = 0; p < 5; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
90       } else {
91         PetscInt    numPoints[2]         = {12, 2};
92         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
93         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
94         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
95         PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5,  -1.0,  0.5, -0.5,  0.0,  0.5, -0.5,   0.0, -0.5, -0.5,
96                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
97                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
98 
99         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
100       }
101       break;
102     default:
103       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
104     }
105   }
106   *newdm = dm;
107   if (refinementLimit > 0.0) {
108     DM rdm;
109     const char *name;
110 
111     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
112     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
113     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
114     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
115     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
116     ierr = DMDestroy(newdm);CHKERRQ(ierr);
117     *newdm = rdm;
118   }
119   if (interpolate) {
120     DM idm = NULL;
121     const char *name;
122 
123     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
124     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
125     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
126     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
127     ierr = DMCopyLabels(*newdm, idm);CHKERRQ(ierr);
128     ierr = DMDestroy(newdm);CHKERRQ(ierr);
129     *newdm = idm;
130   }
131   {
132     DM refinedMesh     = NULL;
133     DM distributedMesh = NULL;
134 
135     /* Distribute mesh over processes */
136     ierr = DMPlexDistribute(*newdm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
137     if (distributedMesh) {
138       ierr = DMDestroy(newdm);CHKERRQ(ierr);
139       *newdm = distributedMesh;
140     }
141     if (refinementUniform) {
142       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
143       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
144       if (refinedMesh) {
145         ierr = DMDestroy(newdm);CHKERRQ(ierr);
146         *newdm = refinedMesh;
147       }
148     }
149   }
150   PetscFunctionReturn(0);
151 }
152 
153 #undef __FUNCT__
154 #define __FUNCT__ "DMPlexCreateSquareBoundary"
155 /*@
156   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
157 
158   Collective on MPI_Comm
159 
160   Input Parameters:
161 + comm  - The communicator for the DM object
162 . lower - The lower left corner coordinates
163 . upper - The upper right corner coordinates
164 - edges - The number of cells in each direction
165 
166   Output Parameter:
167 . dm  - The DM object
168 
169   Note: Here is the numbering returned for 2 cells in each direction:
170 $ 18--5-17--4--16
171 $  |     |     |
172 $  6    10     3
173 $  |     |     |
174 $ 19-11-20--9--15
175 $  |     |     |
176 $  7     8     2
177 $  |     |     |
178 $ 12--0-13--1--14
179 
180   Level: beginner
181 
182 .keywords: DM, create
183 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
184 @*/
185 PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
186 {
187   const PetscInt numVertices    = (edges[0]+1)*(edges[1]+1);
188   const PetscInt numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
189   PetscInt       markerTop      = 1;
190   PetscInt       markerBottom   = 1;
191   PetscInt       markerRight    = 1;
192   PetscInt       markerLeft     = 1;
193   PetscBool      markerSeparate = PETSC_FALSE;
194   Vec            coordinates;
195   PetscSection   coordSection;
196   PetscScalar    *coords;
197   PetscInt       coordSize;
198   PetscMPIInt    rank;
199   PetscInt       v, vx, vy;
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
204   if (markerSeparate) {
205     markerTop    = 3;
206     markerBottom = 1;
207     markerRight  = 2;
208     markerLeft   = 4;
209   }
210   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
211   if (!rank) {
212     PetscInt e, ex, ey;
213 
214     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
215     for (e = 0; e < numEdges; ++e) {
216       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
217     }
218     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
219     for (vx = 0; vx <= edges[0]; vx++) {
220       for (ey = 0; ey < edges[1]; ey++) {
221         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
222         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
223         PetscInt cone[2];
224 
225         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
226         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
227         if (vx == edges[0]) {
228           ierr = DMSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
229           ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
230           if (ey == edges[1]-1) {
231             ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
232           }
233         } else if (vx == 0) {
234           ierr = DMSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
235           ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
236           if (ey == edges[1]-1) {
237             ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
238           }
239         }
240       }
241     }
242     for (vy = 0; vy <= edges[1]; vy++) {
243       for (ex = 0; ex < edges[0]; ex++) {
244         PetscInt edge   = vy*edges[0]     + ex;
245         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
246         PetscInt cone[2];
247 
248         cone[0] = vertex; cone[1] = vertex+1;
249         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
250         if (vy == edges[1]) {
251           ierr = DMSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
252           ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
253           if (ex == edges[0]-1) {
254             ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
255           }
256         } else if (vy == 0) {
257           ierr = DMSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
258           ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
259           if (ex == edges[0]-1) {
260             ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
261           }
262         }
263       }
264     }
265   }
266   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
267   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
268   /* Build coordinates */
269   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
270   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
271   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
272   ierr = PetscSectionSetFieldComponents(coordSection, 0, 2);CHKERRQ(ierr);
273   for (v = numEdges; v < numEdges+numVertices; ++v) {
274     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
275     ierr = PetscSectionSetFieldDof(coordSection, v, 0, 2);CHKERRQ(ierr);
276   }
277   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
278   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
279   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
280   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
281   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
282   ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
283   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
284   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
285   for (vy = 0; vy <= edges[1]; ++vy) {
286     for (vx = 0; vx <= edges[0]; ++vx) {
287       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
288       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
289     }
290   }
291   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
292   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
293   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "DMPlexCreateCubeBoundary"
299 /*@
300   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
301 
302   Collective on MPI_Comm
303 
304   Input Parameters:
305 + comm  - The communicator for the DM object
306 . lower - The lower left front corner coordinates
307 . upper - The upper right back corner coordinates
308 - edges - The number of cells in each direction
309 
310   Output Parameter:
311 . dm  - The DM object
312 
313   Level: beginner
314 
315 .keywords: DM, create
316 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
317 @*/
318 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
319 {
320   PetscInt       vertices[3], numVertices;
321   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
322   Vec            coordinates;
323   PetscSection   coordSection;
324   PetscScalar    *coords;
325   PetscInt       coordSize;
326   PetscMPIInt    rank;
327   PetscInt       v, vx, vy, vz;
328   PetscInt       voffset, iface=0, cone[4];
329   PetscErrorCode ierr;
330 
331   PetscFunctionBegin;
332   if ((faces[0] < 1) || (faces[1] < 1) || (faces[2] < 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side");
333   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
334   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
335   numVertices = vertices[0]*vertices[1]*vertices[2];
336   if (!rank) {
337     PetscInt f;
338 
339     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
340     for (f = 0; f < numFaces; ++f) {
341       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
342     }
343     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
344     for (v = 0; v < numFaces+numVertices; ++v) {
345       ierr = DMSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
346     }
347 
348     /* Side 0 (Top) */
349     for (vy = 0; vy < faces[1]; vy++) {
350       for (vx = 0; vx < faces[0]; vx++) {
351         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
352         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
353         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
354         iface++;
355       }
356     }
357 
358     /* Side 1 (Bottom) */
359     for (vy = 0; vy < faces[1]; vy++) {
360       for (vx = 0; vx < faces[0]; vx++) {
361         voffset = numFaces + vy*(faces[0]+1) + vx;
362         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
363         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
364         iface++;
365       }
366     }
367 
368     /* Side 2 (Front) */
369     for (vz = 0; vz < faces[2]; vz++) {
370       for (vx = 0; vx < faces[0]; vx++) {
371         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
372         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
373         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
374         iface++;
375       }
376     }
377 
378     /* Side 3 (Back) */
379     for (vz = 0; vz < faces[2]; vz++) {
380       for (vx = 0; vx < faces[0]; vx++) {
381         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
382         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
383         cone[2] = voffset+1; cone[3] = voffset;
384         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
385         iface++;
386       }
387     }
388 
389     /* Side 4 (Left) */
390     for (vz = 0; vz < faces[2]; vz++) {
391       for (vy = 0; vy < faces[1]; vy++) {
392         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
393         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
394         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
395         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
396         iface++;
397       }
398     }
399 
400     /* Side 5 (Right) */
401     for (vz = 0; vz < faces[2]; vz++) {
402       for (vy = 0; vy < faces[1]; vy++) {
403         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx;
404         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
405         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
406         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
407         iface++;
408       }
409     }
410   }
411   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
412   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
413   /* Build coordinates */
414   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
415   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
416   for (v = numFaces; v < numFaces+numVertices; ++v) {
417     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
418   }
419   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
420   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
421   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
422   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
423   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
424   ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr);
425   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
426   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
427   for (vz = 0; vz <= faces[2]; ++vz) {
428     for (vy = 0; vy <= faces[1]; ++vy) {
429       for (vx = 0; vx <= faces[0]; ++vx) {
430         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
431         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
432         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
433       }
434     }
435   }
436   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
437   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
438   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "DMPlexCreateCubeMesh_Internal"
444 static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
445 {
446   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
447   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
448   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
449   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
450   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
451   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
452   PetscInt       dim;
453   PetscBool      markerSeparate = PETSC_FALSE;
454   PetscMPIInt    rank;
455   PetscErrorCode ierr;
456 
457   PetscFunctionBegin;
458   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
459   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
460   ierr = DMCreateLabel(dm,"marker");CHKERRQ(ierr);
461   ierr = DMCreateLabel(dm,"Face Sets");CHKERRQ(ierr);
462   switch (dim) {
463   case 2:
464     faceMarkerTop    = 3;
465     faceMarkerBottom = 1;
466     faceMarkerRight  = 2;
467     faceMarkerLeft   = 4;
468     break;
469   case 3:
470     faceMarkerBottom = 1;
471     faceMarkerTop    = 2;
472     faceMarkerFront  = 3;
473     faceMarkerBack   = 4;
474     faceMarkerRight  = 5;
475     faceMarkerLeft   = 6;
476     break;
477   default:
478     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim);
479     break;
480   }
481   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
482   if (markerSeparate) {
483     markerBottom = faceMarkerBottom;
484     markerTop    = faceMarkerTop;
485     markerFront  = faceMarkerFront;
486     markerBack   = faceMarkerBack;
487     markerRight  = faceMarkerRight;
488     markerLeft   = faceMarkerLeft;
489   }
490   {
491     const PetscInt numXEdges    = !rank ? edges[0] : 0;
492     const PetscInt numYEdges    = !rank ? edges[1] : 0;
493     const PetscInt numZEdges    = !rank ? edges[2] : 0;
494     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
495     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
496     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
497     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
498     const PetscInt numXFaces    = numYEdges*numZEdges;
499     const PetscInt numYFaces    = numXEdges*numZEdges;
500     const PetscInt numZFaces    = numXEdges*numYEdges;
501     const PetscInt numTotXFaces = numXVertices*numXFaces;
502     const PetscInt numTotYFaces = numYVertices*numYFaces;
503     const PetscInt numTotZFaces = numZVertices*numZFaces;
504     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
505     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
506     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
507     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
508     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
509     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
510     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
511     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
512     const PetscInt firstYFace   = firstXFace + numTotXFaces;
513     const PetscInt firstZFace   = firstYFace + numTotYFaces;
514     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
515     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
516     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
517     Vec            coordinates;
518     PetscSection   coordSection;
519     PetscScalar   *coords;
520     PetscInt       coordSize;
521     PetscInt       v, vx, vy, vz;
522     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;
523 
524     ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr);
525     for (c = 0; c < numCells; c++) {
526       ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr);
527     }
528     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
529       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
530     }
531     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
532       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
533     }
534     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
535     /* Build cells */
536     for (fz = 0; fz < numZEdges; ++fz) {
537       for (fy = 0; fy < numYEdges; ++fy) {
538         for (fx = 0; fx < numXEdges; ++fx) {
539           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
540           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
541           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
542           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
543           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
544           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
545           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
546                             /* B,  T,  F,  K,  R,  L */
547           PetscInt ornt[6] = {-4,  0,  0, -1,  0, -4}; /* ??? */
548           PetscInt cone[6];
549 
550           /* no boundary twisting in 3D */
551           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
552           ierr    = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr);
553           ierr    = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr);
554         }
555       }
556     }
557     /* Build x faces */
558     for (fz = 0; fz < numZEdges; ++fz) {
559       for (fy = 0; fy < numYEdges; ++fy) {
560         for (fx = 0; fx < numXVertices; ++fx) {
561           PetscInt face    = firstXFace + (fz*numYEdges+fy)*numXVertices + fx;
562           PetscInt edgeL   = firstZEdge + (  fy*                 numXVertices+fx)*numZEdges + fz;
563           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
564           PetscInt edgeB   = firstYEdge + (  fz*                 numXVertices+fx)*numYEdges + fy;
565           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
566           PetscInt ornt[4] = {0, 0, -2, -2};
567           PetscInt cone[4];
568 
569           if (dim == 3) {
570             /* markers */
571             if (bdX != DM_BOUNDARY_PERIODIC) {
572               if (fx == numXVertices-1) {
573                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);CHKERRQ(ierr);
574                 ierr = DMSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr);
575               }
576               else if (fx == 0) {
577                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);CHKERRQ(ierr);
578                 ierr = DMSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr);
579               }
580             }
581           }
582           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
583           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
584           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
585         }
586       }
587     }
588     /* Build y faces */
589     for (fz = 0; fz < numZEdges; ++fz) {
590       for (fx = 0; fx < numXEdges; ++fx) {
591         for (fy = 0; fy < numYVertices; ++fy) {
592           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
593           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx                 )*numZEdges + fz;
594           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
595           PetscInt edgeB   = firstXEdge + (  fz                 *numYVertices+fy)*numXEdges + fx;
596           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
597           PetscInt ornt[4] = {0, 0, -2, -2};
598           PetscInt cone[4];
599 
600           if (dim == 3) {
601             /* markers */
602             if (bdY != DM_BOUNDARY_PERIODIC) {
603               if (fy == numYVertices-1) {
604                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);CHKERRQ(ierr);
605                 ierr = DMSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr);
606               }
607               else if (fy == 0) {
608                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);CHKERRQ(ierr);
609                 ierr = DMSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr);
610               }
611             }
612           }
613           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
614           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
615           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
616         }
617       }
618     }
619     /* Build z faces */
620     for (fy = 0; fy < numYEdges; ++fy) {
621       for (fx = 0; fx < numXEdges; ++fx) {
622         for (fz = 0; fz < numZVertices; fz++) {
623           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
624           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx                 )*numYEdges + fy;
625           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
626           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy                 )*numXEdges + fx;
627           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
628           PetscInt ornt[4] = {0, 0, -2, -2};
629           PetscInt cone[4];
630 
631           if (dim == 2) {
632             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
633             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
634           }
635           else {
636             /* markers */
637             if (bdZ != DM_BOUNDARY_PERIODIC) {
638               if (fz == numZVertices-1) {
639                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);CHKERRQ(ierr);
640                 ierr = DMSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr);
641               }
642               else if (fz == 0) {
643                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);CHKERRQ(ierr);
644                 ierr = DMSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr);
645               }
646             }
647           }
648           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
649           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
650           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
651         }
652       }
653     }
654     /* Build Z edges*/
655     for (vy = 0; vy < numYVertices; vy++) {
656       for (vx = 0; vx < numXVertices; vx++) {
657         for (ez = 0; ez < numZEdges; ez++) {
658           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
659           const PetscInt vertexB = firstVertex + (  ez                 *numYVertices+vy)*numXVertices + vx;
660           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
661           PetscInt       cone[2];
662 
663           if (dim == 3) {
664             if (bdX != DM_BOUNDARY_PERIODIC) {
665               if (vx == numXVertices-1) {
666                 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
667               }
668               else if (vx == 0) {
669                 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
670               }
671             }
672             if (bdY != DM_BOUNDARY_PERIODIC) {
673               if (vy == numYVertices-1) {
674                 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
675               }
676               else if (vy == 0) {
677                 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
678               }
679             }
680           }
681           cone[0] = vertexB; cone[1] = vertexT;
682           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
683         }
684       }
685     }
686     /* Build Y edges*/
687     for (vz = 0; vz < numZVertices; vz++) {
688       for (vx = 0; vx < numXVertices; vx++) {
689         for (ey = 0; ey < numYEdges; ey++) {
690           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
691           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
692           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
693           const PetscInt vertexK = firstVertex + nextv;
694           PetscInt       cone[2];
695 
696           cone[0] = vertexF; cone[1] = vertexK;
697           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
698           if (dim == 2) {
699             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
700               if (vx == numXVertices-1) {
701                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);CHKERRQ(ierr);
702                 ierr = DMSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
703                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
704                 if (ey == numYEdges-1) {
705                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
706                 }
707               }
708               else if (vx == 0) {
709                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);CHKERRQ(ierr);
710                 ierr = DMSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
711                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
712                 if (ey == numYEdges-1) {
713                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
714                 }
715               }
716             }
717           }
718           else {
719             if (bdX != DM_BOUNDARY_PERIODIC) {
720               if (vx == numXVertices-1) {
721                 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
722               }
723               else if (vx == 0) {
724                 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
725               }
726             }
727             if (bdZ != DM_BOUNDARY_PERIODIC) {
728               if (vz == numZVertices-1) {
729                 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
730               }
731               else if (vz == 0) {
732                 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
733               }
734             }
735           }
736         }
737       }
738     }
739     /* Build X edges*/
740     for (vz = 0; vz < numZVertices; vz++) {
741       for (vy = 0; vy < numYVertices; vy++) {
742         for (ex = 0; ex < numXEdges; ex++) {
743           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
744           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
745           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
746           const PetscInt vertexR = firstVertex + nextv;
747           PetscInt       cone[2];
748 
749           cone[0] = vertexL; cone[1] = vertexR;
750           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
751           if (dim == 2) {
752             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
753               if (vy == numYVertices-1) {
754                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);CHKERRQ(ierr);
755                 ierr = DMSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
756                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
757                 if (ex == numXEdges-1) {
758                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
759                 }
760               }
761               else if (vy == 0) {
762                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);CHKERRQ(ierr);
763                 ierr = DMSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
764                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
765                 if (ex == numXEdges-1) {
766                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
767                 }
768               }
769             }
770           }
771           else {
772             if (bdY != DM_BOUNDARY_PERIODIC) {
773               if (vy == numYVertices-1) {
774                 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
775               }
776               else if (vy == 0) {
777                 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
778               }
779             }
780             if (bdZ != DM_BOUNDARY_PERIODIC) {
781               if (vz == numZVertices-1) {
782                 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
783               }
784               else if (vz == 0) {
785                 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
786               }
787             }
788           }
789         }
790       }
791     }
792     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
793     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
794     /* Build coordinates */
795     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
796     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
797     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
798     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
799     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
800       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
801       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
802     }
803     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
804     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
805     ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
806     ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
807     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
808     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
809     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
810     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
811     for (vz = 0; vz < numZVertices; ++vz) {
812       for (vy = 0; vy < numYVertices; ++vy) {
813         for (vx = 0; vx < numXVertices; ++vx) {
814           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
815           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
816           if (dim == 3) {
817             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
818           }
819         }
820       }
821     }
822     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
823     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
824     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
825   }
826   PetscFunctionReturn(0);
827 }
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "DMPlexCreateSquareMesh"
831 /*@
832   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
833 
834   Collective on MPI_Comm
835 
836   Input Parameters:
837 + comm  - The communicator for the DM object
838 . lower - The lower left corner coordinates
839 . upper - The upper right corner coordinates
840 . edges - The number of cells in each direction
841 . bdX   - The boundary type for the X direction
842 - bdY   - The boundary type for the Y direction
843 
844   Output Parameter:
845 . dm  - The DM object
846 
847   Note: Here is the numbering returned for 2 cells in each direction:
848 $ 22--8-23--9--24
849 $  |     |     |
850 $ 13  2 14  3  15
851 $  |     |     |
852 $ 19--6-20--7--21
853 $  |     |     |
854 $ 10  0 11  1 12
855 $  |     |     |
856 $ 16--4-17--5--18
857 
858   Level: beginner
859 
860 .keywords: DM, create
861 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
862 @*/
863 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
864 {
865   PetscReal      lower3[3], upper3[3];
866   PetscInt       edges3[3];
867   PetscErrorCode ierr;
868 
869   PetscFunctionBegin;
870   lower3[0] = lower[0]; lower3[1] = lower[1]; lower3[2] = 0.;
871   upper3[0] = upper[0]; upper3[1] = upper[1]; upper3[2] = 0.;
872   edges3[0] = edges[0]; edges3[1] = edges[1]; edges3[2] = 0;
873   ierr = DMPlexCreateCubeMesh_Internal(dm, lower3, upper3, edges3, bdX, bdY, DM_BOUNDARY_NONE);CHKERRQ(ierr);
874   PetscFunctionReturn(0);
875 }
876 
877 #undef __FUNCT__
878 #define __FUNCT__ "DMPlexCreateBoxMesh"
879 /*@
880   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
881 
882   Collective on MPI_Comm
883 
884   Input Parameters:
885 + comm - The communicator for the DM object
886 . dim - The spatial dimension
887 . numFaces - Number of faces per dimension
888 - interpolate - Flag to create intermediate mesh pieces (edges, faces)
889 
890   Output Parameter:
891 . dm  - The DM object
892 
893   Level: beginner
894 
895 .keywords: DM, create
896 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
897 @*/
898 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm)
899 {
900   DM             boundary;
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   PetscValidPointer(dm, 4);
905   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
906   PetscValidLogicalCollectiveInt(boundary,dim,2);
907   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
908   ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr);
909   ierr = DMSetCoordinateDim(boundary, dim);CHKERRQ(ierr);
910   switch (dim) {
911   case 2:
912   {
913     PetscReal lower[2] = {0.0, 0.0};
914     PetscReal upper[2] = {1.0, 1.0};
915     PetscInt  edges[2];
916 
917     edges[0] = numFaces; edges[1] = numFaces;
918     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
919     break;
920   }
921   case 3:
922   {
923     PetscReal lower[3] = {0.0, 0.0, 0.0};
924     PetscReal upper[3] = {1.0, 1.0, 1.0};
925     PetscInt  faces[3];
926 
927     faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces;
928     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
929     break;
930   }
931   default:
932     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
933   }
934   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
935   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
936   PetscFunctionReturn(0);
937 }
938 
939 #undef __FUNCT__
940 #define __FUNCT__ "DMPlexCreateHexBoxMesh"
941 /*@
942   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
943 
944   Collective on MPI_Comm
945 
946   Input Parameters:
947 + comm  - The communicator for the DM object
948 . dim   - The spatial dimension
949 . periodicX - The boundary type for the X direction
950 . periodicY - The boundary type for the Y direction
951 . periodicZ - The boundary type for the Z direction
952 - cells - The number of cells in each direction
953 
954   Output Parameter:
955 . dm  - The DM object
956 
957   Level: beginner
958 
959 .keywords: DM, create
960 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
961 @*/
962 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
963 {
964   PetscInt       i;
965   PetscErrorCode ierr;
966 
967   PetscFunctionBegin;
968   PetscValidPointer(dm, 7);
969   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
970   PetscValidLogicalCollectiveInt(*dm,dim,2);
971   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
972   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
973   switch (dim) {
974   case 2:
975   {
976     PetscReal lower[2] = {0.0, 0.0};
977     PetscReal upper[2] = {1.0, 1.0};
978 
979     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
980     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
981         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST) {
982       PetscReal      L[2];
983       PetscReal      maxCell[2];
984       DMBoundaryType bdType[2] = {periodicX, periodicY};
985 
986       for (i = 0; i < dim; i++) {
987         L[i]       = upper[i] - lower[i];
988         maxCell[i] = 1.1 * (L[i] / cells[i]);
989       }
990 
991       ierr = DMSetPeriodicity(*dm,maxCell,L,bdType);CHKERRQ(ierr);
992     }
993     break;
994   }
995   case 3:
996   {
997     PetscReal lower[3] = {0.0, 0.0, 0.0};
998     PetscReal upper[3] = {1.0, 1.0, 1.0};
999 
1000     ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
1001     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
1002         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST ||
1003         periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1004       PetscReal      L[3];
1005       PetscReal      maxCell[3];
1006       DMBoundaryType bdType[3] = {periodicX, periodicY, periodicZ};
1007 
1008       for (i = 0; i < dim; i++) {
1009         L[i]       = upper[i] - lower[i];
1010         maxCell[i] = 1.1 * (L[i] / cells[i]);
1011       }
1012 
1013       ierr = DMSetPeriodicity(*dm,maxCell,L,bdType);CHKERRQ(ierr);
1014     }
1015     break;
1016   }
1017   default:
1018     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
1019   }
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 /* External function declarations here */
1024 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
1025 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
1026 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
1027 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
1028 extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
1029 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
1030 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
1031 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
1032 extern PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmRefined);
1033 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
1034 extern PetscErrorCode DMCoarsenHierarchy_Plex(DM dm, PetscInt nlevels, DM dmCoarsened[]);
1035 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
1036 extern PetscErrorCode DMSetUp_Plex(DM dm);
1037 extern PetscErrorCode DMDestroy_Plex(DM dm);
1038 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
1039 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
1040 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
1041 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF);
1042 extern PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
1043 extern PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
1044 extern PetscErrorCode DMProjectFieldLocal_Plex(DM,PetscReal,Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscScalar[]),InsertMode,Vec);
1045 extern PetscErrorCode DMProjectFieldLabelLocal_Plex(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscScalar[]),InsertMode,Vec);
1046 extern PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
1047 extern PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
1048 extern PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
1049 
1050 #undef __FUNCT__
1051 #define __FUNCT__ "DMPlexReplace_Static"
1052 /* Replace dm with the contents of dmNew
1053    - Share the DM_Plex structure
1054    - Share the coordinates
1055    - Share the SF
1056 */
1057 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
1058 {
1059   PetscSF          sf;
1060   DM               coordDM, coarseDM;
1061   Vec              coords;
1062   const PetscReal *maxCell, *L;
1063   const DMBoundaryType *bd;
1064   PetscErrorCode   ierr;
1065 
1066   PetscFunctionBegin;
1067   ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr);
1068   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
1069   ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr);
1070   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
1071   ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr);
1072   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
1073   ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr);
1074   if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L, bd);CHKERRQ(ierr);}
1075   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
1076   dm->data = dmNew->data;
1077   ((DM_Plex *) dmNew->data)->refct++;
1078   dmNew->labels->refct++;
1079   if (!--(dm->labels->refct)) {
1080     DMLabelLink next = dm->labels->next;
1081 
1082     /* destroy the labels */
1083     while (next) {
1084       DMLabelLink tmp = next->next;
1085 
1086       ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr);
1087       ierr = PetscFree(next);CHKERRQ(ierr);
1088       next = tmp;
1089     }
1090     ierr = PetscFree(dm->labels);CHKERRQ(ierr);
1091   }
1092   dm->labels = dmNew->labels;
1093   dm->depthLabel = dmNew->depthLabel;
1094   ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr);
1095   ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr);
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 #undef __FUNCT__
1100 #define __FUNCT__ "DMPlexSwap_Static"
1101 /* Swap dm with the contents of dmNew
1102    - Swap the DM_Plex structure
1103    - Swap the coordinates
1104    - Swap the point PetscSF
1105 */
1106 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
1107 {
1108   DM              coordDMA, coordDMB;
1109   Vec             coordsA,  coordsB;
1110   PetscSF         sfA,      sfB;
1111   void            *tmp;
1112   DMLabelLinkList listTmp;
1113   DMLabel         depthTmp;
1114   PetscErrorCode  ierr;
1115 
1116   PetscFunctionBegin;
1117   ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr);
1118   ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr);
1119   ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr);
1120   ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr);
1121   ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr);
1122   ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr);
1123 
1124   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
1125   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
1126   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
1127   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
1128   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
1129   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
1130 
1131   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
1132   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
1133   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
1134   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
1135   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
1136   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
1137 
1138   tmp       = dmA->data;
1139   dmA->data = dmB->data;
1140   dmB->data = tmp;
1141   listTmp   = dmA->labels;
1142   dmA->labels = dmB->labels;
1143   dmB->labels = listTmp;
1144   depthTmp  = dmA->depthLabel;
1145   dmA->depthLabel = dmB->depthLabel;
1146   dmB->depthLabel = depthTmp;
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 #undef __FUNCT__
1151 #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex"
1152 PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
1153 {
1154   DM_Plex       *mesh = (DM_Plex*) dm->data;
1155   PetscErrorCode ierr;
1156 
1157   PetscFunctionBegin;
1158   /* Handle viewing */
1159   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
1160   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
1161   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr);
1162   /* Point Location */
1163   ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr);
1164   /* Generation and remeshing */
1165   ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr);
1166   /* Projection behavior */
1167   ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr);
1168   ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr);
1169   PetscFunctionReturn(0);
1170 }
1171 
1172 #undef __FUNCT__
1173 #define __FUNCT__ "DMSetFromOptions_Plex"
1174 PetscErrorCode  DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
1175 {
1176   PetscInt       refine = 0, coarsen = 0, r;
1177   PetscBool      isHierarchy;
1178   PetscErrorCode ierr;
1179 
1180   PetscFunctionBegin;
1181   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1182   ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr);
1183   /* Handle DMPlex refinement */
1184   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
1185   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
1186   if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);}
1187   if (refine && isHierarchy) {
1188     DM *dms, coarseDM;
1189 
1190     ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr);
1191     ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
1192     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
1193     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
1194     /* Total hack since we do not pass in a pointer */
1195     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
1196     if (refine == 1) {
1197       ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
1198       ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr);
1199     } else {
1200       ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
1201       ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr);
1202       ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
1203       ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr);
1204     }
1205     ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr);
1206     ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr);
1207     /* Free DMs */
1208     for (r = 0; r < refine; ++r) {
1209       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr);
1210       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
1211     }
1212     ierr = PetscFree(dms);CHKERRQ(ierr);
1213   } else {
1214     for (r = 0; r < refine; ++r) {
1215       DM refinedMesh;
1216 
1217       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
1218       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
1219       /* Total hack since we do not pass in a pointer */
1220       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
1221       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
1222       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
1223     }
1224   }
1225   /* Handle DMPlex coarsening */
1226   ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr);
1227   ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr);
1228   if (coarsen && isHierarchy) {
1229     DM *dms;
1230 
1231     ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr);
1232     ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr);
1233     /* Free DMs */
1234     for (r = 0; r < coarsen; ++r) {
1235       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr);
1236       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
1237     }
1238     ierr = PetscFree(dms);CHKERRQ(ierr);
1239   } else {
1240     for (r = 0; r < coarsen; ++r) {
1241       DM coarseMesh;
1242 
1243       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
1244       ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr);
1245       /* Total hack since we do not pass in a pointer */
1246       ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr);
1247       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
1248       ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr);
1249     }
1250   }
1251   /* Handle */
1252   ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
1253   ierr = PetscOptionsTail();CHKERRQ(ierr);
1254   PetscFunctionReturn(0);
1255 }
1256 
1257 #undef __FUNCT__
1258 #define __FUNCT__ "DMCreateGlobalVector_Plex"
1259 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
1260 {
1261   PetscErrorCode ierr;
1262 
1263   PetscFunctionBegin;
1264   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
1265   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
1266   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr);
1267   ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr);
1268   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr);
1269   ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr);
1270   PetscFunctionReturn(0);
1271 }
1272 
1273 #undef __FUNCT__
1274 #define __FUNCT__ "DMCreateLocalVector_Plex"
1275 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
1276 {
1277   PetscErrorCode ierr;
1278 
1279   PetscFunctionBegin;
1280   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
1281   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
1282   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr);
1283   PetscFunctionReturn(0);
1284 }
1285 
1286 #undef __FUNCT__
1287 #define __FUNCT__ "DMGetDimPoints_Plex"
1288 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
1289 {
1290   PetscInt       depth, d;
1291   PetscErrorCode ierr;
1292 
1293   PetscFunctionBegin;
1294   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1295   if (depth == 1) {
1296     ierr = DMGetDimension(dm, &d);CHKERRQ(ierr);
1297     if (dim == 0)      {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);}
1298     else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);}
1299     else               {*pStart = 0; *pEnd = 0;}
1300   } else {
1301     ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);
1302   }
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 #undef __FUNCT__
1307 #define __FUNCT__ "DMGetNeighors_Plex"
1308 static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
1309 {
1310   PetscSF        sf;
1311   PetscErrorCode ierr;
1312 
1313   PetscFunctionBegin;
1314   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
1315   ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr);
1316   PetscFunctionReturn(0);
1317 }
1318 
1319 #undef __FUNCT__
1320 #define __FUNCT__ "DMInitialize_Plex"
1321 PetscErrorCode DMInitialize_Plex(DM dm)
1322 {
1323   PetscErrorCode ierr;
1324 
1325   PetscFunctionBegin;
1326   dm->ops->view                            = DMView_Plex;
1327   dm->ops->load                            = DMLoad_Plex;
1328   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
1329   dm->ops->clone                           = DMClone_Plex;
1330   dm->ops->setup                           = DMSetUp_Plex;
1331   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
1332   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
1333   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
1334   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
1335   dm->ops->getlocaltoglobalmapping         = NULL;
1336   dm->ops->createfieldis                   = NULL;
1337   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
1338   dm->ops->getcoloring                     = NULL;
1339   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
1340   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
1341   dm->ops->getaggregates                   = NULL;
1342   dm->ops->getinjection                    = DMCreateInjection_Plex;
1343   dm->ops->refine                          = DMRefine_Plex;
1344   dm->ops->coarsen                         = DMCoarsen_Plex;
1345   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
1346   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
1347   dm->ops->globaltolocalbegin              = NULL;
1348   dm->ops->globaltolocalend                = NULL;
1349   dm->ops->localtoglobalbegin              = NULL;
1350   dm->ops->localtoglobalend                = NULL;
1351   dm->ops->destroy                         = DMDestroy_Plex;
1352   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
1353   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
1354   dm->ops->locatepoints                    = DMLocatePoints_Plex;
1355   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
1356   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
1357   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
1358   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
1359   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
1360   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
1361   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
1362   dm->ops->getneighbors                    = DMGetNeighors_Plex;
1363   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Plex);CHKERRQ(ierr);
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 #undef __FUNCT__
1368 #define __FUNCT__ "DMClone_Plex"
1369 PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
1370 {
1371   DM_Plex        *mesh = (DM_Plex *) dm->data;
1372   PetscErrorCode ierr;
1373 
1374   PetscFunctionBegin;
1375   mesh->refct++;
1376   (*newdm)->data = mesh;
1377   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
1378   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
1379   PetscFunctionReturn(0);
1380 }
1381 
1382 /*MC
1383   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
1384                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
1385                     specified by a PetscSection object. Ownership in the global representation is determined by
1386                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
1387 
1388   Level: intermediate
1389 
1390 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
1391 M*/
1392 
1393 #undef __FUNCT__
1394 #define __FUNCT__ "DMCreate_Plex"
1395 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1396 {
1397   DM_Plex        *mesh;
1398   PetscInt       unit, d;
1399   PetscErrorCode ierr;
1400 
1401   PetscFunctionBegin;
1402   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1403   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
1404   dm->dim  = 0;
1405   dm->data = mesh;
1406 
1407   mesh->refct             = 1;
1408   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
1409   mesh->maxConeSize       = 0;
1410   mesh->cones             = NULL;
1411   mesh->coneOrientations  = NULL;
1412   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
1413   mesh->maxSupportSize    = 0;
1414   mesh->supports          = NULL;
1415   mesh->refinementUniform = PETSC_TRUE;
1416   mesh->refinementLimit   = -1.0;
1417 
1418   mesh->facesTmp = NULL;
1419 
1420   mesh->tetgenOpts   = NULL;
1421   mesh->triangleOpts = NULL;
1422   ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr);
1423   ierr = PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);CHKERRQ(ierr);
1424   mesh->remeshBd     = PETSC_FALSE;
1425 
1426   mesh->subpointMap = NULL;
1427 
1428   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
1429 
1430   mesh->regularRefinement   = PETSC_FALSE;
1431   mesh->depthState          = -1;
1432   mesh->globalVertexNumbers = NULL;
1433   mesh->globalCellNumbers   = NULL;
1434   mesh->anchorSection       = NULL;
1435   mesh->anchorIS            = NULL;
1436   mesh->createanchors       = NULL;
1437   mesh->computeanchormatrix = NULL;
1438   mesh->parentSection       = NULL;
1439   mesh->parents             = NULL;
1440   mesh->childIDs            = NULL;
1441   mesh->childSection        = NULL;
1442   mesh->children            = NULL;
1443   mesh->referenceTree       = NULL;
1444   mesh->getchildsymmetry    = NULL;
1445   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1446   mesh->vtkCellHeight       = 0;
1447   mesh->useCone             = PETSC_FALSE;
1448   mesh->useClosure          = PETSC_TRUE;
1449   mesh->useAnchors          = PETSC_FALSE;
1450 
1451   mesh->maxProjectionHeight = 0;
1452 
1453   mesh->printSetValues = PETSC_FALSE;
1454   mesh->printFEM       = 0;
1455   mesh->printTol       = 1.0e-10;
1456 
1457   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
1458   PetscFunctionReturn(0);
1459 }
1460 
1461 #undef __FUNCT__
1462 #define __FUNCT__ "DMPlexCreate"
1463 /*@
1464   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1465 
1466   Collective on MPI_Comm
1467 
1468   Input Parameter:
1469 . comm - The communicator for the DMPlex object
1470 
1471   Output Parameter:
1472 . mesh  - The DMPlex object
1473 
1474   Level: beginner
1475 
1476 .keywords: DMPlex, create
1477 @*/
1478 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1479 {
1480   PetscErrorCode ierr;
1481 
1482   PetscFunctionBegin;
1483   PetscValidPointer(mesh,2);
1484   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1485   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1486   PetscFunctionReturn(0);
1487 }
1488 
1489 #undef __FUNCT__
1490 #define __FUNCT__ "DMPlexBuildFromCellList_Parallel_Private"
1491 /*
1492   This takes as input the common mesh generator output, a list of the vertices for each cell, but vertex numbers are global and an SF is built for them
1493 */
1494 static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert)
1495 {
1496   PetscSF         sfPoint;
1497   PetscLayout     vLayout;
1498   PetscHashI      vhash;
1499   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
1500   const PetscInt *vrange;
1501   PetscInt        numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
1502   PETSC_UNUSED PetscHashIIter ret, iter;
1503   PetscMPIInt     rank, numProcs;
1504   PetscErrorCode  ierr;
1505 
1506   PetscFunctionBegin;
1507   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
1508   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
1509   /* Partition vertices */
1510   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr);
1511   ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr);
1512   ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
1513   ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
1514   ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr);
1515   /* Count vertices and map them to procs */
1516   PetscHashICreate(vhash);
1517   for (c = 0; c < numCells; ++c) {
1518     for (p = 0; p < numCorners; ++p) {
1519       PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter);
1520     }
1521   }
1522   PetscHashISize(vhash, numVerticesAdj);
1523   ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr);
1524   off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr);
1525   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
1526   ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr);
1527   for (v = 0; v < numVerticesAdj; ++v) {
1528     const PetscInt gv = verticesAdj[v];
1529     PetscInt       vrank;
1530 
1531     ierr = PetscFindInt(gv, numProcs+1, vrange, &vrank);CHKERRQ(ierr);
1532     vrank = vrank < 0 ? -(vrank+2) : vrank;
1533     remoteVerticesAdj[v].index = gv - vrange[vrank];
1534     remoteVerticesAdj[v].rank  = vrank;
1535   }
1536   /* Create cones */
1537   ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr);
1538   for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);}
1539   ierr = DMSetUp(dm);CHKERRQ(ierr);
1540   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1541   for (c = 0; c < numCells; ++c) {
1542     for (p = 0; p < numCorners; ++p) {
1543       const PetscInt gv = cells[c*numCorners+p];
1544       PetscInt       lv;
1545 
1546       ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr);
1547       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
1548       cone[p] = lv+numCells;
1549     }
1550     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
1551   }
1552   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1553   /* Create SF for vertices */
1554   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr);
1555   ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr);
1556   ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr);
1557   ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr);
1558   ierr = PetscFree(verticesAdj);CHKERRQ(ierr);
1559   /* Build pointSF */
1560   ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr);
1561   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
1562   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
1563   ierr = PetscSFReduceBegin(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr);
1564   ierr = PetscSFReduceEnd(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr);
1565   for (v = 0; v < numVertices;    ++v) if (vertexOwner[v].rank < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global vertex %d on rank %d was unclaimed", v + vrange[rank], rank);
1566   ierr = PetscSFBcastBegin(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr);
1567   ierr = PetscSFBcastEnd(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr);
1568   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
1569   ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr);
1570   ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr);
1571   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
1572     if (vertexLocal[v].rank != rank) {
1573       localVertex[g]        = v+numCells;
1574       remoteVertex[g].index = vertexLocal[v].index;
1575       remoteVertex[g].rank  = vertexLocal[v].rank;
1576       ++g;
1577     }
1578   }
1579   ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr);
1580   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
1581   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1582   ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr);
1583   ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr);
1584   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
1585   PetscHashIDestroy(vhash);
1586   /* Fill in the rest of the topology structure */
1587   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1588   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1589   PetscFunctionReturn(0);
1590 }
1591 
1592 #undef __FUNCT__
1593 #define __FUNCT__ "DMPlexBuildCoordinates_Parallel_Private"
1594 /*
1595   This takes as input the coordinates for each owned vertex
1596 */
1597 static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscSF sfVert, const double vertexCoords[])
1598 {
1599   PetscSection   coordSection;
1600   Vec            coordinates;
1601   PetscScalar   *coords;
1602   PetscInt       numVertices, numVerticesAdj, coordSize, v;
1603   PetscErrorCode ierr;
1604 
1605   PetscFunctionBegin;
1606   ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr);
1607   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1608   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1609   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
1610   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr);
1611   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
1612     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
1613     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
1614   }
1615   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1616   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1617   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
1618   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
1619   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1620   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1621   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1622   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1623   {
1624     MPI_Datatype coordtype;
1625 
1626     /* Need a temp buffer for coords if we have complex/single */
1627     ierr = MPI_Type_contiguous(spaceDim, MPI_DOUBLE, &coordtype);CHKERRQ(ierr);
1628     ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr);
1629     ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr);
1630     ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr);
1631     ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr);
1632   }
1633   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1634   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1635   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1636   PetscFunctionReturn(0);
1637 }
1638 
1639 #undef __FUNCT__
1640 #define __FUNCT__ "DMPlexCreateFromCellListParallel"
1641 /*@C
1642   DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
1643 
1644   Input Parameters:
1645 + comm - The communicator
1646 . dim - The topological dimension of the mesh
1647 . numCells - The number of cells owned by this process
1648 . numVertices - The number of vertices owned by this process
1649 . numCorners - The number of vertices for each cell
1650 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1651 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
1652 . spaceDim - The spatial dimension used for coordinates
1653 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
1654 
1655   Output Parameter:
1656 . dm - The DM
1657 
1658   Note: Two triangles sharing a face
1659 $
1660 $        2
1661 $      / | \
1662 $     /  |  \
1663 $    /   |   \
1664 $   0  0 | 1  3
1665 $    \   |   /
1666 $     \  |  /
1667 $      \ | /
1668 $        1
1669 would have input
1670 $  numCells = 2, numVertices = 4
1671 $  cells = [0 1 2  1 3 2]
1672 $
1673 which would result in the DMPlex
1674 $
1675 $        4
1676 $      / | \
1677 $     /  |  \
1678 $    /   |   \
1679 $   2  0 | 1  5
1680 $    \   |   /
1681 $     \  |  /
1682 $      \ | /
1683 $        3
1684 
1685   Level: beginner
1686 
1687 .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
1688 @*/
1689 PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm)
1690 {
1691   PetscSF        sfVert;
1692   PetscErrorCode ierr;
1693 
1694   PetscFunctionBegin;
1695   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1696   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1697   PetscValidLogicalCollectiveInt(*dm, dim, 2);
1698   PetscValidLogicalCollectiveInt(*dm, spaceDim, 8);
1699   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
1700   ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr);
1701   if (interpolate) {
1702     DM idm = NULL;
1703 
1704     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1705     ierr = DMDestroy(dm);CHKERRQ(ierr);
1706     *dm  = idm;
1707   }
1708   ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, sfVert, vertexCoords);CHKERRQ(ierr);
1709   ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);
1710   PetscFunctionReturn(0);
1711 }
1712 
1713 #undef __FUNCT__
1714 #define __FUNCT__ "DMPlexBuildFromCellList_Private"
1715 /*
1716   This takes as input the common mesh generator output, a list of the vertices for each cell
1717 */
1718 static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
1719 {
1720   PetscInt      *cone, c, p;
1721   PetscErrorCode ierr;
1722 
1723   PetscFunctionBegin;
1724   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
1725   for (c = 0; c < numCells; ++c) {
1726     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
1727   }
1728   ierr = DMSetUp(dm);CHKERRQ(ierr);
1729   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1730   for (c = 0; c < numCells; ++c) {
1731     for (p = 0; p < numCorners; ++p) {
1732       cone[p] = cells[c*numCorners+p]+numCells;
1733     }
1734     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
1735   }
1736   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1737   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1738   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1739   PetscFunctionReturn(0);
1740 }
1741 
1742 #undef __FUNCT__
1743 #define __FUNCT__ "DMPlexBuildCoordinates_Private"
1744 /*
1745   This takes as input the coordinates for each vertex
1746 */
1747 static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
1748 {
1749   PetscSection   coordSection;
1750   Vec            coordinates;
1751   DM             cdm;
1752   PetscScalar   *coords;
1753   PetscInt       v, d;
1754   PetscErrorCode ierr;
1755 
1756   PetscFunctionBegin;
1757   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1758   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1759   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
1760   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
1761   for (v = numCells; v < numCells+numVertices; ++v) {
1762     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
1763     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
1764   }
1765   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1766 
1767   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1768   ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr);
1769   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
1770   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1771   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1772   for (v = 0; v < numVertices; ++v) {
1773     for (d = 0; d < spaceDim; ++d) {
1774       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
1775     }
1776   }
1777   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1778   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1779   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1780   PetscFunctionReturn(0);
1781 }
1782 
1783 #undef __FUNCT__
1784 #define __FUNCT__ "DMPlexCreateFromCellList"
1785 /*@C
1786   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
1787 
1788   Input Parameters:
1789 + comm - The communicator
1790 . dim - The topological dimension of the mesh
1791 . numCells - The number of cells
1792 . numVertices - The number of vertices
1793 . numCorners - The number of vertices for each cell
1794 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1795 . cells - An array of numCells*numCorners numbers, the vertices for each cell
1796 . spaceDim - The spatial dimension used for coordinates
1797 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
1798 
1799   Output Parameter:
1800 . dm - The DM
1801 
1802   Note: Two triangles sharing a face
1803 $
1804 $        2
1805 $      / | \
1806 $     /  |  \
1807 $    /   |   \
1808 $   0  0 | 1  3
1809 $    \   |   /
1810 $     \  |  /
1811 $      \ | /
1812 $        1
1813 would have input
1814 $  numCells = 2, numVertices = 4
1815 $  cells = [0 1 2  1 3 2]
1816 $
1817 which would result in the DMPlex
1818 $
1819 $        4
1820 $      / | \
1821 $     /  |  \
1822 $    /   |   \
1823 $   2  0 | 1  5
1824 $    \   |   /
1825 $     \  |  /
1826 $      \ | /
1827 $        3
1828 
1829   Level: beginner
1830 
1831 .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
1832 @*/
1833 PetscErrorCode DMPlexCreateFromCellList(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm)
1834 {
1835   PetscErrorCode ierr;
1836 
1837   PetscFunctionBegin;
1838   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1839   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1840   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
1841   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
1842   if (interpolate) {
1843     DM idm = NULL;
1844 
1845     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1846     ierr = DMDestroy(dm);CHKERRQ(ierr);
1847     *dm  = idm;
1848   }
1849   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
1850   PetscFunctionReturn(0);
1851 }
1852 
1853 #undef __FUNCT__
1854 #define __FUNCT__ "DMPlexCreateFromDAG"
1855 /*@
1856   DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM
1857 
1858   Input Parameters:
1859 + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
1860 . depth - The depth of the DAG
1861 . numPoints - The number of points at each depth
1862 . coneSize - The cone size of each point
1863 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1864 . coneOrientations - The orientation of each cone point
1865 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1866 
1867   Output Parameter:
1868 . dm - The DM
1869 
1870   Note: Two triangles sharing a face would have input
1871 $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1872 $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1873 $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1874 $
1875 which would result in the DMPlex
1876 $
1877 $        4
1878 $      / | \
1879 $     /  |  \
1880 $    /   |   \
1881 $   2  0 | 1  5
1882 $    \   |   /
1883 $     \  |  /
1884 $      \ | /
1885 $        3
1886 $
1887 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1888 
1889   Level: advanced
1890 
1891 .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1892 @*/
1893 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
1894 {
1895   Vec            coordinates;
1896   PetscSection   coordSection;
1897   PetscScalar    *coords;
1898   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;
1899   PetscErrorCode ierr;
1900 
1901   PetscFunctionBegin;
1902   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1903   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1904   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
1905   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
1906   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
1907   for (p = pStart; p < pEnd; ++p) {
1908     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
1909     if (firstVertex < 0 && !coneSize[p - pStart]) {
1910       firstVertex = p - pStart;
1911     }
1912   }
1913   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
1914   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
1915   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
1916     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
1917     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
1918   }
1919   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1920   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1921   /* Build coordinates */
1922   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1923   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1924   ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr);
1925   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
1926   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
1927     ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr);
1928     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr);
1929   }
1930   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1931   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1932   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
1933   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1934   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1935   ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr);
1936   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1937   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1938   for (v = 0; v < numPoints[0]; ++v) {
1939     PetscInt off;
1940 
1941     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
1942     for (d = 0; d < dimEmbed; ++d) {
1943       coords[off+d] = vertexCoords[v*dimEmbed+d];
1944     }
1945   }
1946   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1947   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1948   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1949   PetscFunctionReturn(0);
1950 }
1951 
1952 #undef __FUNCT__
1953 #define __FUNCT__ "DMPlexCreateFromFile"
1954 /*@C
1955   DMPlexCreateFromFile - This takes a filename and produces a DM
1956 
1957   Input Parameters:
1958 + comm - The communicator
1959 . filename - A file name
1960 - interpolate - Flag to create intermediate mesh pieces (edges, faces)
1961 
1962   Output Parameter:
1963 . dm - The DM
1964 
1965   Level: beginner
1966 
1967 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
1968 @*/
1969 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1970 {
1971   const char    *extGmsh   = ".msh";
1972   const char    *extCGNS   = ".cgns";
1973   const char    *extExodus = ".exo";
1974   const char    *extFluent = ".cas";
1975   const char    *extHDF5   = ".h5";
1976   size_t         len;
1977   PetscBool      isGmsh, isCGNS, isExodus, isFluent, isHDF5;
1978   PetscMPIInt    rank;
1979   PetscErrorCode ierr;
1980 
1981   PetscFunctionBegin;
1982   PetscValidPointer(filename, 2);
1983   PetscValidPointer(dm, 4);
1984   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1985   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
1986   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
1987   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,   4, &isGmsh);CHKERRQ(ierr);
1988   ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,   5, &isCGNS);CHKERRQ(ierr);
1989   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr);
1990   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr);
1991   ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,   3, &isHDF5);CHKERRQ(ierr);
1992   if (isGmsh) {
1993     ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
1994   } else if (isCGNS) {
1995     ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
1996   } else if (isExodus) {
1997     ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
1998   } else if (isFluent) {
1999     ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2000   } else if (isHDF5) {
2001     PetscViewer viewer;
2002 
2003     ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
2004     ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr);
2005     ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
2006     ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
2007     ierr = DMCreate(comm, dm);CHKERRQ(ierr);
2008     ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2009     ierr = DMLoad(*dm, viewer);CHKERRQ(ierr);
2010     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2011   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
2012   PetscFunctionReturn(0);
2013 }
2014 
2015 #undef __FUNCT__
2016 #define __FUNCT__ "DMPlexCreateReferenceCell"
2017 /*@
2018   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
2019 
2020   Collective on comm
2021 
2022   Input Parameters:
2023 + comm    - The communicator
2024 . dim     - The spatial dimension
2025 - simplex - Flag for simplex, otherwise use a tensor-product cell
2026 
2027   Output Parameter:
2028 . refdm - The reference cell
2029 
2030   Level: intermediate
2031 
2032 .keywords: reference cell
2033 .seealso:
2034 @*/
2035 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
2036 {
2037   DM             rdm;
2038   Vec            coords;
2039   PetscErrorCode ierr;
2040 
2041   PetscFunctionBeginUser;
2042   ierr = DMCreate(comm, &rdm);CHKERRQ(ierr);
2043   ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr);
2044   ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr);
2045   switch (dim) {
2046   case 0:
2047   {
2048     PetscInt    numPoints[1]        = {1};
2049     PetscInt    coneSize[1]         = {0};
2050     PetscInt    cones[1]            = {0};
2051     PetscInt    coneOrientations[1] = {0};
2052     PetscScalar vertexCoords[1]     = {0.0};
2053 
2054     ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
2055   }
2056   break;
2057   case 1:
2058   {
2059     PetscInt    numPoints[2]        = {2, 1};
2060     PetscInt    coneSize[3]         = {2, 0, 0};
2061     PetscInt    cones[2]            = {1, 2};
2062     PetscInt    coneOrientations[2] = {0, 0};
2063     PetscScalar vertexCoords[2]     = {-1.0,  1.0};
2064 
2065     ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
2066   }
2067   break;
2068   case 2:
2069     if (simplex) {
2070       PetscInt    numPoints[2]        = {3, 1};
2071       PetscInt    coneSize[4]         = {3, 0, 0, 0};
2072       PetscInt    cones[3]            = {1, 2, 3};
2073       PetscInt    coneOrientations[3] = {0, 0, 0};
2074       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};
2075 
2076       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
2077     } else {
2078       PetscInt    numPoints[2]        = {4, 1};
2079       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2080       PetscInt    cones[4]            = {1, 2, 3, 4};
2081       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2082       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};
2083 
2084       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
2085     }
2086   break;
2087   case 3:
2088     if (simplex) {
2089       PetscInt    numPoints[2]        = {4, 1};
2090       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2091       PetscInt    cones[4]            = {1, 3, 2, 4};
2092       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2093       PetscScalar vertexCoords[12]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,  -1.0, -1.0, 1.0};
2094 
2095       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
2096     } else {
2097       PetscInt    numPoints[2]        = {8, 1};
2098       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
2099       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
2100       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
2101       PetscScalar vertexCoords[24]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  1.0, 1.0, -1.0,  -1.0, 1.0, -1.0,
2102                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};
2103 
2104       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
2105     }
2106   break;
2107   default:
2108     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
2109   }
2110   *refdm = NULL;
2111   ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr);
2112   if (rdm->coordinateDM) {
2113     DM           ncdm;
2114     PetscSection cs;
2115     PetscInt     pEnd = -1;
2116 
2117     ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr);
2118     if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);}
2119     if (pEnd >= 0) {
2120       ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr);
2121       ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr);
2122       ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr);
2123       ierr = DMDestroy(&ncdm);CHKERRQ(ierr);
2124     }
2125   }
2126   ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr);
2127   if (coords) {
2128     ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr);
2129   } else {
2130     ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr);
2131     if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);}
2132   }
2133   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
2134   PetscFunctionReturn(0);
2135 }
2136