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