xref: /petsc/src/dm/impls/plex/plexcreate.c (revision d79735219aee4f8cd85b3ffe9a0d82178a16a9a9)
1 #define PETSCDM_DLL
2 #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 #include <petscdmda.h>
4 #include <petscsf.h>
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "DMPlexCreateDoublet"
8 /*@
9   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
10 
11   Collective on MPI_Comm
12 
13   Input Parameters:
14 + comm - The communicator for the DM object
15 . dim - The spatial dimension
16 . simplex - Flag for simplicial cells, otherwise they are tensor product cells
17 . interpolate - Flag to create intermediate mesh pieces (edges, faces)
18 . refinementUniform - Flag for uniform parallel refinement
19 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
20 
21   Output Parameter:
22 . dm  - The DM object
23 
24   Level: beginner
25 
26 .keywords: DM, create
27 .seealso: DMSetType(), DMCreate()
28 @*/
29 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
30 {
31   DM             dm;
32   PetscInt       p;
33   PetscMPIInt    rank;
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
38   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
39   ierr = DMPlexSetDimension(dm, dim);CHKERRQ(ierr);
40   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
41   switch (dim) {
42   case 2:
43     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
44     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
45     break;
46   case 3:
47     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
48     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
49     break;
50   default:
51     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
52   }
53   if (rank) {
54     PetscInt numPoints[2] = {0, 0};
55     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
56   } else {
57     switch (dim) {
58     case 2:
59       if (simplex) {
60         PetscInt    numPoints[2]        = {4, 2};
61         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
62         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
63         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
64         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
65         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
66 
67         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
68         for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
69       } else {
70         PetscInt    numPoints[2]        = {6, 2};
71         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
72         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
73         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
74         PetscScalar vertexCoords[12]    = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  -1.0, 0.5,  1.0, -0.5,  1.0, 0.5};
75 
76         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
77       }
78       break;
79     case 3:
80       if (simplex) {
81         PetscInt    numPoints[2]        = {5, 2};
82         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
83         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
84         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
85         PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  1.0, 0.0, 0.0};
86         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
87 
88         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
89         for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
90       } else {
91         PetscInt    numPoints[2]         = {12, 2};
92         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
93         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
94         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
95         PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5,  -1.0,  0.5, -0.5,  0.0,  0.5, -0.5,   0.0, -0.5, -0.5,
96                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
97                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
98 
99         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
100       }
101       break;
102     default:
103       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
104     }
105   }
106   *newdm = dm;
107   if (refinementLimit > 0.0) {
108     DM rdm;
109     const char *name;
110 
111     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
112     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
113     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
114     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
115     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
116     ierr = DMDestroy(newdm);CHKERRQ(ierr);
117     *newdm = rdm;
118   }
119   if (interpolate) {
120     DM idm = NULL;
121     const char *name;
122 
123     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
124     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
125     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
126     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
127     ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr);
128     ierr = DMDestroy(newdm);CHKERRQ(ierr);
129     *newdm = idm;
130   }
131   {
132     DM refinedMesh     = NULL;
133     DM distributedMesh = NULL;
134 
135     /* Distribute mesh over processes */
136     ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr);
137     if (distributedMesh) {
138       ierr = DMDestroy(newdm);CHKERRQ(ierr);
139       *newdm = distributedMesh;
140     }
141     if (refinementUniform) {
142       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
143       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
144       if (refinedMesh) {
145         ierr = DMDestroy(newdm);CHKERRQ(ierr);
146         *newdm = refinedMesh;
147       }
148     }
149   }
150   PetscFunctionReturn(0);
151 }
152 
153 #undef __FUNCT__
154 #define __FUNCT__ "DMPlexCreateSquareBoundary"
155 /*@
156   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
157 
158   Collective on MPI_Comm
159 
160   Input Parameters:
161 + comm  - The communicator for the DM object
162 . lower - The lower left corner coordinates
163 . upper - The upper right corner coordinates
164 - edges - The number of cells in each direction
165 
166   Output Parameter:
167 . dm  - The DM object
168 
169   Note: Here is the numbering returned for 2 cells in each direction:
170 $ 18--5-17--4--16
171 $  |     |     |
172 $  6    10     3
173 $  |     |     |
174 $ 19-11-20--9--15
175 $  |     |     |
176 $  7     8     2
177 $  |     |     |
178 $ 12--0-13--1--14
179 
180   Level: beginner
181 
182 .keywords: DM, create
183 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
184 @*/
185 PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
186 {
187   PetscInt       numVertices    = (edges[0]+1)*(edges[1]+1);
188   PetscInt       numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
189   PetscInt       markerTop      = 1;
190   PetscInt       markerBottom   = 1;
191   PetscInt       markerRight    = 1;
192   PetscInt       markerLeft     = 1;
193   PetscBool      markerSeparate = PETSC_FALSE;
194   Vec            coordinates;
195   PetscSection   coordSection;
196   PetscScalar    *coords;
197   PetscInt       coordSize;
198   PetscMPIInt    rank;
199   PetscInt       v, vx, vy;
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
204   if (markerSeparate) {
205     markerTop    = 1;
206     markerBottom = 0;
207     markerRight  = 0;
208     markerLeft   = 0;
209   }
210   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
211   if (!rank) {
212     PetscInt e, ex, ey;
213 
214     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
215     for (e = 0; e < numEdges; ++e) {
216       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
217     }
218     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
219     for (vx = 0; vx <= edges[0]; vx++) {
220       for (ey = 0; ey < edges[1]; ey++) {
221         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
222         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
223         PetscInt cone[2];
224 
225         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
226         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
227         if (vx == edges[0]) {
228           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
229           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
230           if (ey == edges[1]-1) {
231             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
232           }
233         } else if (vx == 0) {
234           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
235           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
236           if (ey == edges[1]-1) {
237             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
238           }
239         }
240       }
241     }
242     for (vy = 0; vy <= edges[1]; vy++) {
243       for (ex = 0; ex < edges[0]; ex++) {
244         PetscInt edge   = vy*edges[0]     + ex;
245         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
246         PetscInt cone[2];
247 
248         cone[0] = vertex; cone[1] = vertex+1;
249         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
250         if (vy == edges[1]) {
251           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
252           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
253           if (ex == edges[0]-1) {
254             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
255           }
256         } else if (vy == 0) {
257           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
258           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
259           if (ex == edges[0]-1) {
260             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
261           }
262         }
263       }
264     }
265   }
266   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
267   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
268   /* Build coordinates */
269   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
270   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
271   for (v = numEdges; v < numEdges+numVertices; ++v) {
272     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
273   }
274   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
275   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
276   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
277   ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
278   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);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 #undef __FUNCT__
294 #define __FUNCT__ "DMPlexCreateCubeBoundary"
295 /*@
296   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
297 
298   Collective on MPI_Comm
299 
300   Input Parameters:
301 + comm  - The communicator for the DM object
302 . lower - The lower left front corner coordinates
303 . upper - The upper right back corner coordinates
304 - edges - The number of cells in each direction
305 
306   Output Parameter:
307 . dm  - The DM object
308 
309   Level: beginner
310 
311 .keywords: DM, create
312 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
313 @*/
314 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
315 {
316   PetscInt       vertices[3], numVertices;
317   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
318   Vec            coordinates;
319   PetscSection   coordSection;
320   PetscScalar    *coords;
321   PetscInt       coordSize;
322   PetscMPIInt    rank;
323   PetscInt       v, vx, vy, vz;
324   PetscInt       voffset, iface=0, cone[4];
325   PetscErrorCode ierr;
326 
327   PetscFunctionBegin;
328   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");
329   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
330   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
331   numVertices = vertices[0]*vertices[1]*vertices[2];
332   if (!rank) {
333     PetscInt f;
334 
335     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
336     for (f = 0; f < numFaces; ++f) {
337       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
338     }
339     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
340     for (v = 0; v < numFaces+numVertices; ++v) {
341       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
342     }
343 
344     /* Side 0 (Top) */
345     for (vy = 0; vy < faces[1]; vy++) {
346       for (vx = 0; vx < faces[0]; vx++) {
347         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
348         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
349         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
350         iface++;
351       }
352     }
353 
354     /* Side 1 (Bottom) */
355     for (vy = 0; vy < faces[1]; vy++) {
356       for (vx = 0; vx < faces[0]; vx++) {
357         voffset = numFaces + vy*(faces[0]+1) + vx;
358         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
359         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
360         iface++;
361       }
362     }
363 
364     /* Side 2 (Front) */
365     for (vz = 0; vz < faces[2]; vz++) {
366       for (vx = 0; vx < faces[0]; vx++) {
367         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
368         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
369         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
370         iface++;
371       }
372     }
373 
374     /* Side 3 (Back) */
375     for (vz = 0; vz < faces[2]; vz++) {
376       for (vx = 0; vx < faces[0]; vx++) {
377         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
378         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
379         cone[2] = voffset+1; cone[3] = voffset;
380         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
381         iface++;
382       }
383     }
384 
385     /* Side 4 (Left) */
386     for (vz = 0; vz < faces[2]; vz++) {
387       for (vy = 0; vy < faces[1]; vy++) {
388         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
389         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
390         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
391         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
392         iface++;
393       }
394     }
395 
396     /* Side 5 (Right) */
397     for (vz = 0; vz < faces[2]; vz++) {
398       for (vy = 0; vy < faces[1]; vy++) {
399         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx;
400         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
401         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
402         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
403         iface++;
404       }
405     }
406   }
407   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
408   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
409   /* Build coordinates */
410   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
411   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
412   for (v = numFaces; v < numFaces+numVertices; ++v) {
413     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
414   }
415   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
416   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
417   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
418   ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr);
419   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
420   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
421   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
422   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
423   for (vz = 0; vz <= faces[2]; ++vz) {
424     for (vy = 0; vy <= faces[1]; ++vy) {
425       for (vx = 0; vx <= faces[0]; ++vx) {
426         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
427         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
428         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
429       }
430     }
431   }
432   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
433   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
434   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
435   PetscFunctionReturn(0);
436 }
437 
438 #undef __FUNCT__
439 #define __FUNCT__ "DMPlexCreateSquareMesh"
440 /*@
441   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
442 
443   Collective on MPI_Comm
444 
445   Input Parameters:
446 + comm  - The communicator for the DM object
447 . lower - The lower left corner coordinates
448 . upper - The upper right corner coordinates
449 . edges - The number of cells in each direction
450 . bdX   - The boundary type for the X direction
451 - bdY   - The boundary type for the Y direction
452 
453   Output Parameter:
454 . dm  - The DM object
455 
456   Note: Here is the numbering returned for 2 cells in each direction:
457 $ 22--8-23--9--24
458 $  |     |     |
459 $ 13  2 14  3  15
460 $  |     |     |
461 $ 19--6-20--7--21
462 $  |     |     |
463 $ 10  0 11  1 12
464 $  |     |     |
465 $ 16--4-17--5--18
466 
467   Level: beginner
468 
469 .keywords: DM, create
470 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
471 @*/
472 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
473 {
474   PetscInt       markerTop      = 1;
475   PetscInt       markerBottom   = 1;
476   PetscInt       markerRight    = 1;
477   PetscInt       markerLeft     = 1;
478   PetscBool      markerSeparate = PETSC_FALSE;
479   PetscMPIInt    rank;
480   PetscErrorCode ierr;
481 
482   PetscFunctionBegin;
483   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
484   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
485   if (markerSeparate) {
486     markerTop    = 3;
487     markerBottom = 1;
488     markerRight  = 2;
489     markerLeft   = 4;
490   }
491   {
492     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
493     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
494     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
495     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
496     const PetscInt numTotXEdges = numXEdges*numYVertices;
497     const PetscInt numTotYEdges = numYEdges*numXVertices;
498     const PetscInt numVertices  = numXVertices*numYVertices;
499     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
500     const PetscInt numFaces     = numXEdges*numYEdges;
501     const PetscInt firstVertex  = numFaces;
502     const PetscInt firstXEdge   = numFaces + numVertices;
503     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
504     Vec            coordinates;
505     PetscSection   coordSection;
506     PetscScalar   *coords;
507     PetscInt       coordSize;
508     PetscInt       v, vx, vy;
509     PetscInt       f, fx, fy, e, ex, ey;
510 
511     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
512     for (f = 0; f < numFaces; ++f) {
513       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
514     }
515     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
516       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
517     }
518     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
519     /* Build faces */
520     for (fy = 0; fy < numYEdges; ++fy) {
521       for (fx = 0; fx < numXEdges; ++fx) {
522         PetscInt face    = fy*numXEdges + fx;
523         PetscInt edgeL   = firstYEdge +   fx                 *numYEdges + fy;
524         PetscInt edgeR   = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy;
525         PetscInt edgeB   = firstXEdge +   fy                 *numXEdges + fx;
526         PetscInt edgeT   = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx;
527         PetscInt ornt[4] = {0, 0, -2, -2};
528         PetscInt cone[4];
529 
530         if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
531         if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
532         cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
533         ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
534         ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
535       }
536     }
537     /* Build Y edges*/
538     for (vx = 0; vx < numXVertices; vx++) {
539       for (ey = 0; ey < numYEdges; ey++) {
540         const PetscInt nextv   = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx;
541         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
542         const PetscInt vertexB = firstVertex + ey*numXVertices + vx;
543         const PetscInt vertexT = firstVertex + nextv;
544         PetscInt       cone[2];
545 
546         cone[0] = vertexB; cone[1] = vertexT;
547         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
548         if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
549           if (vx == numXVertices-1) {
550             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
551             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
552             if (ey == numYEdges-1) {
553               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
554             }
555           } else if (vx == 0) {
556             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
557             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
558             if (ey == numYEdges-1) {
559               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
560             }
561           }
562         }
563       }
564     }
565     /* Build X edges*/
566     for (vy = 0; vy < numYVertices; vy++) {
567       for (ex = 0; ex < numXEdges; ex++) {
568         const PetscInt nextv   = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices;
569         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
570         const PetscInt vertexL = firstVertex + vy*numXVertices + ex;
571         const PetscInt vertexR = firstVertex + nextv;
572         PetscInt       cone[2];
573 
574         cone[0] = vertexL; cone[1] = vertexR;
575         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
576         if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
577           if (vy == numYVertices-1) {
578             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
579             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
580             if (ex == numXEdges-1) {
581               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
582             }
583           } else if (vy == 0) {
584             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
585             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
586             if (ex == numXEdges-1) {
587               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
588             }
589           }
590         }
591       }
592     }
593     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
594     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
595     /* Build coordinates */
596     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
597     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
598     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
599       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
600     }
601     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
602     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
603     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
604     ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
605     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
606     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
607     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
608     for (vy = 0; vy < numYVertices; ++vy) {
609       for (vx = 0; vx < numXVertices; ++vx) {
610         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
611         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
612       }
613     }
614     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
615     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
616     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
617   }
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "DMPlexCreateBoxMesh"
623 /*@
624   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
625 
626   Collective on MPI_Comm
627 
628   Input Parameters:
629 + comm - The communicator for the DM object
630 . dim - The spatial dimension
631 - interpolate - Flag to create intermediate mesh pieces (edges, faces)
632 
633   Output Parameter:
634 . dm  - The DM object
635 
636   Level: beginner
637 
638 .keywords: DM, create
639 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
640 @*/
641 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
642 {
643   DM             boundary;
644   PetscErrorCode ierr;
645 
646   PetscFunctionBegin;
647   PetscValidPointer(dm, 4);
648   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
649   PetscValidLogicalCollectiveInt(boundary,dim,2);
650   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
651   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
652   switch (dim) {
653   case 2:
654   {
655     PetscReal lower[2] = {0.0, 0.0};
656     PetscReal upper[2] = {1.0, 1.0};
657     PetscInt  edges[2] = {2, 2};
658 
659     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
660     break;
661   }
662   case 3:
663   {
664     PetscReal lower[3] = {0.0, 0.0, 0.0};
665     PetscReal upper[3] = {1.0, 1.0, 1.0};
666     PetscInt  faces[3] = {1, 1, 1};
667 
668     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
669     break;
670   }
671   default:
672     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
673   }
674   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
675   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
676   PetscFunctionReturn(0);
677 }
678 
679 #undef __FUNCT__
680 #define __FUNCT__ "DMPlexCreateHexBoxMesh"
681 /*@
682   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
683 
684   Collective on MPI_Comm
685 
686   Input Parameters:
687 + comm  - The communicator for the DM object
688 . dim   - The spatial dimension
689 . periodicX - The boundary type for the X direction
690 . periodicY - The boundary type for the Y direction
691 . periodicZ - The boundary type for the Z direction
692 - cells - The number of cells in each direction
693 
694   Output Parameter:
695 . dm  - The DM object
696 
697   Level: beginner
698 
699 .keywords: DM, create
700 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
701 @*/
702 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
703 {
704   PetscErrorCode ierr;
705 
706   PetscFunctionBegin;
707   PetscValidPointer(dm, 4);
708   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
709   PetscValidLogicalCollectiveInt(*dm,dim,2);
710   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
711   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
712   switch (dim) {
713   case 2:
714   {
715     PetscReal lower[2] = {0.0, 0.0};
716     PetscReal upper[2] = {1.0, 1.0};
717 
718     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
719     break;
720   }
721 #if 0
722   case 3:
723   {
724     PetscReal lower[3] = {0.0, 0.0, 0.0};
725     PetscReal upper[3] = {1.0, 1.0, 1.0};
726 
727     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
728     break;
729   }
730 #endif
731   default:
732     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
733   }
734   PetscFunctionReturn(0);
735 }
736 
737 /* External function declarations here */
738 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
739 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx);
740 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
741 extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
742 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
743 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
744 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
745 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
746 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
747 extern PetscErrorCode DMSetUp_Plex(DM dm);
748 extern PetscErrorCode DMDestroy_Plex(DM dm);
749 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
750 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
751 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
752 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "DMPlexReplace_Static"
756 /* Replace dm with the contents of dmNew
757    - Share the DM_Plex structure
758    - Share the coordinates
759    - Share the SF
760 */
761 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
762 {
763   PetscSF        sf;
764   PetscSection   coordSection;
765   Vec            coords;
766   PetscErrorCode ierr;
767 
768   PetscFunctionBegin;
769   ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr);
770   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
771   ierr = DMGetCoordinateSection(dmNew, &coordSection);CHKERRQ(ierr);
772   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
773   ierr = DMSetCoordinateSection(dm, coordSection);CHKERRQ(ierr);
774   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
775   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
776   dm->data = dmNew->data;
777   ((DM_Plex *) dmNew->data)->refct++;
778   PetscFunctionReturn(0);
779 }
780 
781 #undef __FUNCT__
782 #define __FUNCT__ "DMPlexSwap_Static"
783 /* Swap dm with the contents of dmNew
784    - Swap the DM_Plex structure
785    - Swap the coordinates
786 */
787 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
788 {
789   DM             coordDMA, coordDMB;
790   Vec            coordsA,  coordsB;
791   void          *tmp;
792   PetscErrorCode ierr;
793 
794   PetscFunctionBegin;
795   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
796   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
797   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
798   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
799   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
800   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
801 
802   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
803   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
804   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
805   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
806   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
807   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
808   tmp       = dmA->data;
809   dmA->data = dmB->data;
810   dmB->data = tmp;
811   PetscFunctionReturn(0);
812 }
813 
814 #undef __FUNCT__
815 #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex"
816 PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(DM dm)
817 {
818   DM_Plex       *mesh = (DM_Plex*) dm->data;
819   DMBoundary     b;
820   PetscErrorCode ierr;
821 
822   PetscFunctionBegin;
823   /* Handle boundary conditions */
824   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr);
825   for (b = mesh->boundary; b; b = b->next) {
826     char      optname[1024];
827     PetscInt  ids[1024], len = 1024, i;
828     PetscBool flg;
829 
830     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
831     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
832     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
833     if (flg) {
834       DMLabel label;
835 
836       ierr = DMPlexGetLabel(dm, b->labelname, &label);CHKERRQ(ierr);
837       for (i = 0; i < len; ++i) {
838         PetscBool has;
839 
840         ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr);
841         if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name);
842       }
843       b->numids = len;
844       ierr = PetscFree(b->ids);CHKERRQ(ierr);
845       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
846       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
847     }
848   }
849   ierr = PetscOptionsEnd();CHKERRQ(ierr);
850   /* Handle viewing */
851   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
852   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
853   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr);
854   PetscFunctionReturn(0);
855 }
856 
857 #undef __FUNCT__
858 #define __FUNCT__ "DMSetFromOptions_Plex"
859 PetscErrorCode  DMSetFromOptions_Plex(DM dm)
860 {
861   PetscInt       refine = 0, r;
862   PetscBool      isHierarchy;
863   PetscErrorCode ierr;
864 
865   PetscFunctionBegin;
866   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
867   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
868   /* Handle DMPlex refinement */
869   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
870   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
871   ierr = DMPlexSetRefinementUniform(dm, refine ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
872   if (refine && isHierarchy) {
873     DM *dms;
874 
875     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
876     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
877     /* Total hack since we do not pass in a pointer */
878     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
879     if (refine == 1) {
880       ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
881     } else {
882       ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
883       ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
884     }
885     /* Free DMs */
886     for (r = 0; r < refine; ++r) {
887       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
888       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
889     }
890     ierr = PetscFree(dms);CHKERRQ(ierr);
891   } else {
892     for (r = 0; r < refine; ++r) {
893       DM refinedMesh;
894 
895       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
896       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
897       /* Total hack since we do not pass in a pointer */
898       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
899       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
900     }
901   }
902   ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
903   ierr = PetscOptionsTail();CHKERRQ(ierr);
904   PetscFunctionReturn(0);
905 }
906 
907 #undef __FUNCT__
908 #define __FUNCT__ "DMCreateGlobalVector_Plex"
909 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
910 {
911   PetscErrorCode ierr;
912 
913   PetscFunctionBegin;
914   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
915   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
916   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr);
917   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr);
918   PetscFunctionReturn(0);
919 }
920 
921 #undef __FUNCT__
922 #define __FUNCT__ "DMCreateLocalVector_Plex"
923 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
924 {
925   PetscErrorCode ierr;
926 
927   PetscFunctionBegin;
928   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
929   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
930   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr);
931   PetscFunctionReturn(0);
932 }
933 
934 #undef __FUNCT__
935 #define __FUNCT__ "DMInitialize_Plex"
936 PetscErrorCode DMInitialize_Plex(DM dm)
937 {
938   PetscFunctionBegin;
939   dm->ops->view                            = DMView_Plex;
940   dm->ops->load                            = DMLoad_Plex;
941   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
942   dm->ops->clone                           = DMClone_Plex;
943   dm->ops->setup                           = DMSetUp_Plex;
944   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
945   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
946   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
947   dm->ops->getlocaltoglobalmapping         = NULL;
948   dm->ops->getlocaltoglobalmappingblock    = NULL;
949   dm->ops->createfieldis                   = NULL;
950   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
951   dm->ops->getcoloring                     = NULL;
952   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
953   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
954   dm->ops->getaggregates                   = NULL;
955   dm->ops->getinjection                    = DMCreateInjection_Plex;
956   dm->ops->refine                          = DMRefine_Plex;
957   dm->ops->coarsen                         = DMCoarsen_Plex;
958   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
959   dm->ops->coarsenhierarchy                = NULL;
960   dm->ops->globaltolocalbegin              = NULL;
961   dm->ops->globaltolocalend                = NULL;
962   dm->ops->localtoglobalbegin              = NULL;
963   dm->ops->localtoglobalend                = NULL;
964   dm->ops->destroy                         = DMDestroy_Plex;
965   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
966   dm->ops->locatepoints                    = DMLocatePoints_Plex;
967   PetscFunctionReturn(0);
968 }
969 
970 #undef __FUNCT__
971 #define __FUNCT__ "DMClone_Plex"
972 PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
973 {
974   DM_Plex        *mesh = (DM_Plex *) dm->data;
975   PetscErrorCode ierr;
976 
977   PetscFunctionBegin;
978   mesh->refct++;
979   (*newdm)->data = mesh;
980   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
981   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
982   PetscFunctionReturn(0);
983 }
984 
985 /*MC
986   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
987                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
988                     specified by a PetscSection object. Ownership in the global representation is determined by
989                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
990 
991   Level: intermediate
992 
993 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
994 M*/
995 
996 #undef __FUNCT__
997 #define __FUNCT__ "DMCreate_Plex"
998 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
999 {
1000   DM_Plex        *mesh;
1001   PetscInt       unit, d;
1002   PetscErrorCode ierr;
1003 
1004   PetscFunctionBegin;
1005   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1006   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
1007   dm->data = mesh;
1008 
1009   mesh->refct             = 1;
1010   mesh->dim               = 0;
1011   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
1012   mesh->maxConeSize       = 0;
1013   mesh->cones             = NULL;
1014   mesh->coneOrientations  = NULL;
1015   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
1016   mesh->maxSupportSize    = 0;
1017   mesh->supports          = NULL;
1018   mesh->refinementUniform = PETSC_TRUE;
1019   mesh->refinementLimit   = -1.0;
1020 
1021   mesh->facesTmp = NULL;
1022 
1023   mesh->subpointMap = NULL;
1024 
1025   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
1026 
1027   mesh->labels              = NULL;
1028   mesh->depthLabel          = NULL;
1029   mesh->globalVertexNumbers = NULL;
1030   mesh->globalCellNumbers   = NULL;
1031   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1032   mesh->vtkCellHeight       = 0;
1033   mesh->useCone             = PETSC_FALSE;
1034   mesh->useClosure          = PETSC_TRUE;
1035 
1036   mesh->printSetValues = PETSC_FALSE;
1037   mesh->printFEM       = 0;
1038   mesh->printTol       = 1.0e-10;
1039 
1040   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 #undef __FUNCT__
1045 #define __FUNCT__ "DMPlexCreate"
1046 /*@
1047   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1048 
1049   Collective on MPI_Comm
1050 
1051   Input Parameter:
1052 . comm - The communicator for the DMPlex object
1053 
1054   Output Parameter:
1055 . mesh  - The DMPlex object
1056 
1057   Level: beginner
1058 
1059 .keywords: DMPlex, create
1060 @*/
1061 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1062 {
1063   PetscErrorCode ierr;
1064 
1065   PetscFunctionBegin;
1066   PetscValidPointer(mesh,2);
1067   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1068   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1069   PetscFunctionReturn(0);
1070 }
1071 
1072 #undef __FUNCT__
1073 #define __FUNCT__ "DMPlexBuildFromCellList_Private"
1074 /*
1075   This takes as input the common mesh generator output, a list of the vertices for each cell
1076 */
1077 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
1078 {
1079   PetscInt      *cone, c, p;
1080   PetscErrorCode ierr;
1081 
1082   PetscFunctionBegin;
1083   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
1084   for (c = 0; c < numCells; ++c) {
1085     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
1086   }
1087   ierr = DMSetUp(dm);CHKERRQ(ierr);
1088   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1089   for (c = 0; c < numCells; ++c) {
1090     for (p = 0; p < numCorners; ++p) {
1091       cone[p] = cells[c*numCorners+p]+numCells;
1092     }
1093     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
1094   }
1095   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1096   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1097   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 #undef __FUNCT__
1102 #define __FUNCT__ "DMPlexBuildCoordinates_Private"
1103 /*
1104   This takes as input the coordinates for each vertex
1105 */
1106 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
1107 {
1108   PetscSection   coordSection;
1109   Vec            coordinates;
1110   PetscScalar   *coords;
1111   PetscInt       coordSize, v, d;
1112   PetscErrorCode ierr;
1113 
1114   PetscFunctionBegin;
1115   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1116   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1117   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
1118   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
1119   for (v = numCells; v < numCells+numVertices; ++v) {
1120     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
1121     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
1122   }
1123   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1124   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1125   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
1126   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
1127   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1128   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1129   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1130   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1131   for (v = 0; v < numVertices; ++v) {
1132     for (d = 0; d < spaceDim; ++d) {
1133       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
1134     }
1135   }
1136   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1137   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1138   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "DMPlexCreateFromCellList"
1144 /*@C
1145   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
1146 
1147   Input Parameters:
1148 + comm - The communicator
1149 . dim - The topological dimension of the mesh
1150 . numCells - The number of cells
1151 . numVertices - The number of vertices
1152 . numCorners - The number of vertices for each cell
1153 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1154 . cells - An array of numCells*numCorners numbers, the vertices for each cell
1155 . spaceDim - The spatial dimension used for coordinates
1156 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
1157 
1158   Output Parameter:
1159 . dm - The DM
1160 
1161   Note: Two triangles sharing a face
1162 $
1163 $        2
1164 $      / | \
1165 $     /  |  \
1166 $    /   |   \
1167 $   0  0 | 1  3
1168 $    \   |   /
1169 $     \  |  /
1170 $      \ | /
1171 $        1
1172 would have input
1173 $  numCells = 2, numVertices = 4
1174 $  cells = [0 1 2  1 3 2]
1175 $
1176 which would result in the DMPlex
1177 $
1178 $        4
1179 $      / | \
1180 $     /  |  \
1181 $    /   |   \
1182 $   2  0 | 1  5
1183 $    \   |   /
1184 $     \  |  /
1185 $      \ | /
1186 $        3
1187 
1188   Level: beginner
1189 
1190 .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
1191 @*/
1192 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)
1193 {
1194   PetscErrorCode ierr;
1195 
1196   PetscFunctionBegin;
1197   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1198   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1199   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
1200   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
1201   if (interpolate) {
1202     DM idm = NULL;
1203 
1204     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1205     ierr = DMDestroy(dm);CHKERRQ(ierr);
1206     *dm  = idm;
1207   }
1208   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 #undef __FUNCT__
1213 #define __FUNCT__ "DMPlexCreateFromDAG"
1214 /*@
1215   DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM
1216 
1217   Input Parameters:
1218 + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension()
1219 . depth - The depth of the DAG
1220 . numPoints - The number of points at each depth
1221 . coneSize - The cone size of each point
1222 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1223 . coneOrientations - The orientation of each cone point
1224 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1225 
1226   Output Parameter:
1227 . dm - The DM
1228 
1229   Note: Two triangles sharing a face would have input
1230 $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1231 $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1232 $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1233 $
1234 which would result in the DMPlex
1235 $
1236 $        4
1237 $      / | \
1238 $     /  |  \
1239 $    /   |   \
1240 $   2  0 | 1  5
1241 $    \   |   /
1242 $     \  |  /
1243 $      \ | /
1244 $        3
1245 $
1246 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1247 
1248   Level: advanced
1249 
1250 .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1251 @*/
1252 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
1253 {
1254   Vec            coordinates;
1255   PetscSection   coordSection;
1256   PetscScalar    *coords;
1257   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
1258   PetscErrorCode ierr;
1259 
1260   PetscFunctionBegin;
1261   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1262   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
1263   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
1264   for (p = pStart; p < pEnd; ++p) {
1265     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
1266   }
1267   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
1268   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
1269     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
1270     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
1271   }
1272   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1273   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1274   /* Build coordinates */
1275   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1276   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1277   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
1278   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
1279   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
1280     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
1281     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
1282   }
1283   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1284   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1285   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
1286   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
1287   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1288   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1289   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1290   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1291   for (v = 0; v < numPoints[0]; ++v) {
1292     PetscInt off;
1293 
1294     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
1295     for (d = 0; d < dim; ++d) {
1296       coords[off+d] = vertexCoords[v*dim+d];
1297     }
1298   }
1299   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1300   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1301   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1302   PetscFunctionReturn(0);
1303 }
1304