xref: /petsc/src/dm/impls/plex/plexcreate.c (revision 7b59f5a92658f217c01d5163df9c7e3326a23d6c)
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;
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 = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
278   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
279   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
280   for (vy = 0; vy <= edges[1]; ++vy) {
281     for (vx = 0; vx <= edges[0]; ++vx) {
282       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
283       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
284     }
285   }
286   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
287   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
288   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
289   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "DMPlexCreateCubeBoundary"
294 /*@
295   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
296 
297   Collective on MPI_Comm
298 
299   Input Parameters:
300 + comm  - The communicator for the DM object
301 . lower - The lower left front corner coordinates
302 . upper - The upper right back corner coordinates
303 - edges - The number of cells in each direction
304 
305   Output Parameter:
306 . dm  - The DM object
307 
308   Level: beginner
309 
310 .keywords: DM, create
311 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
312 @*/
313 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
314 {
315   PetscInt       vertices[3] = {faces[0]+1, faces[1]+1, faces[2]+1};
316   PetscInt       numVertices = vertices[0]*vertices[1]*vertices[2];
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   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 = DMPlexSetLabelValue(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(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
416   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
417   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
418   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
419   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
420   for (vz = 0; vz <= faces[2]; ++vz) {
421     for (vy = 0; vy <= faces[1]; ++vy) {
422       for (vx = 0; vx <= faces[0]; ++vx) {
423         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
424         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
425         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
426       }
427     }
428   }
429   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
430   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
431   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "DMPlexCreateSquareMesh"
437 /*@
438   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
439 
440   Collective on MPI_Comm
441 
442   Input Parameters:
443 + comm  - The communicator for the DM object
444 . lower - The lower left corner coordinates
445 . upper - The upper right corner coordinates
446 . edges - The number of cells in each direction
447 . bdX   - The boundary type for the X direction
448 - bdY   - The boundary type for the Y direction
449 
450   Output Parameter:
451 . dm  - The DM object
452 
453   Note: Here is the numbering returned for 2 cells in each direction:
454 $ 22--8-23--9--24
455 $  |     |     |
456 $ 13  2 14  3  15
457 $  |     |     |
458 $ 19--6-20--7--21
459 $  |     |     |
460 $ 10  0 11  1 12
461 $  |     |     |
462 $ 16--4-17--5--18
463 
464   Level: beginner
465 
466 .keywords: DM, create
467 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
468 @*/
469 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
470 {
471   PetscInt       markerTop      = 1;
472   PetscInt       markerBottom   = 1;
473   PetscInt       markerRight    = 1;
474   PetscInt       markerLeft     = 1;
475   PetscBool      markerSeparate = PETSC_FALSE;
476   PetscMPIInt    rank;
477   PetscErrorCode ierr;
478 
479   PetscFunctionBegin;
480   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
481   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
482   if (markerSeparate) {
483     markerTop    = 3;
484     markerBottom = 1;
485     markerRight  = 2;
486     markerLeft   = 4;
487   }
488   {
489     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
490     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
491     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
492     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
493     const PetscInt numTotXEdges = numXEdges*numYVertices;
494     const PetscInt numTotYEdges = numYEdges*numXVertices;
495     const PetscInt numVertices  = numXVertices*numYVertices;
496     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
497     const PetscInt numFaces     = numXEdges*numYEdges;
498     const PetscInt firstVertex  = numFaces;
499     const PetscInt firstXEdge   = numFaces + numVertices;
500     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
501     Vec            coordinates;
502     PetscSection   coordSection;
503     PetscScalar   *coords;
504     PetscInt       coordSize;
505     PetscInt       v, vx, vy;
506     PetscInt       f, fx, fy, e, ex, ey;
507 
508     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
509     for (f = 0; f < numFaces; ++f) {
510       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
511     }
512     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
513       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
514     }
515     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
516     /* Build faces */
517     for (fy = 0; fy < numYEdges; fy++) {
518       for (fx = 0; fx < numXEdges; fx++) {
519         const PetscInt face    = fy*numXEdges + fx;
520         const PetscInt edgeL   = firstYEdge +   fx                 *numYEdges + fy;
521         const PetscInt edgeR   = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy;
522         const PetscInt edgeB   = firstXEdge +   fy                 *numXEdges + fx;
523         const PetscInt edgeT   = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx;
524         const PetscInt ornt[4] = {0, 0, -2, -2};
525         PetscInt       cone[4];
526 
527         cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
528         ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
529         ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
530       }
531     }
532     /* Build Y edges*/
533     for (vx = 0; vx < numXVertices; vx++) {
534       for (ey = 0; ey < numYEdges; ey++) {
535         const PetscInt nextv   = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx;
536         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
537         const PetscInt vertexB = firstVertex + ey*numXVertices + vx;
538         const PetscInt vertexT = firstVertex + nextv;
539         PetscInt       cone[2];
540 
541         cone[0] = vertexB; cone[1] = vertexT;
542         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
543         if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
544           if (vx == numXVertices-1) {
545             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
546             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
547             if (ey == numYEdges-1) {
548               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
549             }
550           } else if (vx == 0) {
551             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
552             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
553             if (ey == numYEdges-1) {
554               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
555             }
556           }
557         }
558       }
559     }
560     /* Build X edges*/
561     for (vy = 0; vy < numYVertices; vy++) {
562       for (ex = 0; ex < numXEdges; ex++) {
563         const PetscInt nextv   = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices;
564         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
565         const PetscInt vertexL = firstVertex + vy*numXVertices + ex;
566         const PetscInt vertexR = firstVertex + nextv;
567         PetscInt       cone[2];
568 
569         cone[0] = vertexL; cone[1] = vertexR;
570         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
571         if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
572           if (vy == numYVertices-1) {
573             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
574             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
575             if (ex == numXEdges-1) {
576               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
577             }
578           } else if (vy == 0) {
579             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
580             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
581             if (ex == numXEdges-1) {
582               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
583             }
584           }
585         }
586       }
587     }
588     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
589     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
590     /* Build coordinates */
591     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
592     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
593     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
594       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
595     }
596     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
597     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
598     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
599     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
600     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
601     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
602     for (vy = 0; vy < numYVertices; ++vy) {
603       for (vx = 0; vx < numXVertices; ++vx) {
604         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
605         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
606       }
607     }
608     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
609     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
610     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
611   }
612   PetscFunctionReturn(0);
613 }
614 
615 #undef __FUNCT__
616 #define __FUNCT__ "DMPlexCreateBoxMesh"
617 /*@
618   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
619 
620   Collective on MPI_Comm
621 
622   Input Parameters:
623 + comm - The communicator for the DM object
624 . dim - The spatial dimension
625 - interpolate - Flag to create intermediate mesh pieces (edges, faces)
626 
627   Output Parameter:
628 . dm  - The DM object
629 
630   Level: beginner
631 
632 .keywords: DM, create
633 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
634 @*/
635 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
636 {
637   DM             boundary;
638   PetscErrorCode ierr;
639 
640   PetscFunctionBegin;
641   PetscValidPointer(dm, 4);
642   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
643   PetscValidLogicalCollectiveInt(boundary,dim,2);
644   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
645   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
646   switch (dim) {
647   case 2:
648   {
649     PetscReal lower[2] = {0.0, 0.0};
650     PetscReal upper[2] = {1.0, 1.0};
651     PetscInt  edges[2] = {2, 2};
652 
653     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
654     break;
655   }
656   case 3:
657   {
658     PetscReal lower[3] = {0.0, 0.0, 0.0};
659     PetscReal upper[3] = {1.0, 1.0, 1.0};
660     PetscInt  faces[3] = {1, 1, 1};
661 
662     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
663     break;
664   }
665   default:
666     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
667   }
668   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
669   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
670   PetscFunctionReturn(0);
671 }
672 
673 #undef __FUNCT__
674 #define __FUNCT__ "DMPlexCreateHexBoxMesh"
675 /*@
676   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
677 
678   Collective on MPI_Comm
679 
680   Input Parameters:
681 + comm  - The communicator for the DM object
682 . dim   - The spatial dimension
683 . periodicX - The boundary type for the X direction
684 . periodicY - The boundary type for the Y direction
685 . periodicZ - The boundary type for the Z direction
686 - cells - The number of cells in each direction
687 
688   Output Parameter:
689 . dm  - The DM object
690 
691   Level: beginner
692 
693 .keywords: DM, create
694 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
695 @*/
696 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
697 {
698   PetscErrorCode ierr;
699 
700   PetscFunctionBegin;
701   PetscValidPointer(dm, 4);
702   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
703   PetscValidLogicalCollectiveInt(*dm,dim,2);
704   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
705   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
706   switch (dim) {
707   case 2:
708   {
709     PetscReal lower[2] = {0.0, 0.0};
710     PetscReal upper[2] = {1.0, 1.0};
711 
712     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
713     break;
714   }
715 #if 0
716   case 3:
717   {
718     PetscReal lower[3] = {0.0, 0.0, 0.0};
719     PetscReal upper[3] = {1.0, 1.0, 1.0};
720 
721     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
722     break;
723   }
724 #endif
725   default:
726     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
727   }
728   PetscFunctionReturn(0);
729 }
730 
731 /* External function declarations here */
732 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
733 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx);
734 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
735 extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
736 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
737 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
738 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
739 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
740 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
741 extern PetscErrorCode DMSetUp_Plex(DM dm);
742 extern PetscErrorCode DMDestroy_Plex(DM dm);
743 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
744 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
745 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "DMPlexReplace_Static"
749 /* Replace dm with the contents of dmNew
750    - Share the DM_Plex structure
751    - Share the coordinates
752 */
753 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
754 {
755   PetscSection   coordSection;
756   Vec            coords;
757   PetscErrorCode ierr;
758 
759   PetscFunctionBegin;
760   ierr = DMGetCoordinateSection(dmNew, &coordSection);CHKERRQ(ierr);
761   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
762   ierr = DMSetCoordinateSection(dm, coordSection);CHKERRQ(ierr);
763   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
764   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
765   dm->data = dmNew->data;
766   ((DM_Plex *) dmNew->data)->refct++;
767   PetscFunctionReturn(0);
768 }
769 
770 #undef __FUNCT__
771 #define __FUNCT__ "DMPlexSwap_Static"
772 /* Swap dm with the contents of dmNew
773    - Swap the DM_Plex structure
774    - Swap the coordinates
775 */
776 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
777 {
778   DM             coordDMA, coordDMB;
779   Vec            coordsA,  coordsB;
780   void          *tmp;
781   PetscErrorCode ierr;
782 
783   PetscFunctionBegin;
784   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
785   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
786   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
787   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
788   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
789   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
790 
791   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
792   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
793   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
794   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
795   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
796   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
797   tmp       = dmA->data;
798   dmA->data = dmB->data;
799   dmB->data = tmp;
800   PetscFunctionReturn(0);
801 }
802 
803 #undef __FUNCT__
804 #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex"
805 PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(DM dm)
806 {
807   DM_Plex       *mesh = (DM_Plex*) dm->data;
808   DMBoundary     b;
809   PetscErrorCode ierr;
810 
811   PetscFunctionBegin;
812   /* Handle boundary conditions */
813   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr);
814   for (b = mesh->boundary; b; b = b->next) {
815     char      optname[1024];
816     PetscInt  ids[1024], len = 1024, i;
817     PetscBool flg;
818 
819     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
820     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
821     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
822     if (flg) {
823       DMLabel label;
824 
825       ierr = DMPlexGetLabel(dm, b->name, &label);CHKERRQ(ierr);
826       for (i = 0; i < len; ++i) {
827         PetscBool has;
828 
829         ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr);
830         if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name);
831       }
832       b->numids = len;
833       ierr = PetscFree(b->ids);CHKERRQ(ierr);
834       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
835       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
836     }
837   }
838   ierr = PetscOptionsEnd();CHKERRQ(ierr);
839   /* Handle viewing */
840   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
841   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
842   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr);
843   PetscFunctionReturn(0);
844 }
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "DMSetFromOptions_Plex"
848 PetscErrorCode  DMSetFromOptions_Plex(DM dm)
849 {
850   PetscInt       refine = 0, r;
851   PetscBool      isHierarchy;
852   PetscErrorCode ierr;
853 
854   PetscFunctionBegin;
855   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
856   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
857   /* Handle DMPlex refinement */
858   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
859   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
860   ierr = DMPlexSetRefinementUniform(dm, refine ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
861   if (refine && isHierarchy) {
862     DM *dms;
863 
864     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
865     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
866     /* Total hack since we do not pass in a pointer */
867     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
868     if (refine == 1) {
869       ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
870     } else {
871       ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
872       ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
873     }
874     /* Free DMs */
875     for (r = 0; r < refine; ++r) {
876       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
877       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
878     }
879     ierr = PetscFree(dms);CHKERRQ(ierr);
880   } else {
881     for (r = 0; r < refine; ++r) {
882       DM refinedMesh;
883 
884       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
885       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
886       /* Total hack since we do not pass in a pointer */
887       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
888       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
889     }
890   }
891   ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
892   ierr = PetscOptionsTail();CHKERRQ(ierr);
893   PetscFunctionReturn(0);
894 }
895 
896 #undef __FUNCT__
897 #define __FUNCT__ "DMCreateGlobalVector_Plex"
898 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
899 {
900   PetscErrorCode ierr;
901 
902   PetscFunctionBegin;
903   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
904   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
905   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr);
906   PetscFunctionReturn(0);
907 }
908 
909 #undef __FUNCT__
910 #define __FUNCT__ "DMCreateLocalVector_Plex"
911 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
912 {
913   PetscErrorCode ierr;
914 
915   PetscFunctionBegin;
916   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
917   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
918   PetscFunctionReturn(0);
919 }
920 
921 #undef __FUNCT__
922 #define __FUNCT__ "DMInitialize_Plex"
923 PetscErrorCode DMInitialize_Plex(DM dm)
924 {
925   PetscFunctionBegin;
926   dm->ops->view                            = DMView_Plex;
927   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
928   dm->ops->clone                           = DMClone_Plex;
929   dm->ops->setup                           = DMSetUp_Plex;
930   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
931   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
932   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
933   dm->ops->getlocaltoglobalmapping         = NULL;
934   dm->ops->getlocaltoglobalmappingblock    = NULL;
935   dm->ops->createfieldis                   = NULL;
936   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
937   dm->ops->getcoloring                     = NULL;
938   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
939   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
940   dm->ops->getaggregates                   = NULL;
941   dm->ops->getinjection                    = DMCreateInjection_Plex;
942   dm->ops->refine                          = DMRefine_Plex;
943   dm->ops->coarsen                         = DMCoarsen_Plex;
944   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
945   dm->ops->coarsenhierarchy                = NULL;
946   dm->ops->globaltolocalbegin              = NULL;
947   dm->ops->globaltolocalend                = NULL;
948   dm->ops->localtoglobalbegin              = NULL;
949   dm->ops->localtoglobalend                = NULL;
950   dm->ops->destroy                         = DMDestroy_Plex;
951   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
952   dm->ops->locatepoints                    = DMLocatePoints_Plex;
953   PetscFunctionReturn(0);
954 }
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "DMClone_Plex"
958 PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
959 {
960   DM_Plex        *mesh = (DM_Plex *) dm->data;
961   PetscErrorCode ierr;
962 
963   PetscFunctionBegin;
964   mesh->refct++;
965   (*newdm)->data = mesh;
966   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
967   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
968   PetscFunctionReturn(0);
969 }
970 
971 /*MC
972   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
973                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
974                     specified by a PetscSection object. Ownership in the global representation is determined by
975                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
976 
977   Level: intermediate
978 
979 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
980 M*/
981 
982 #undef __FUNCT__
983 #define __FUNCT__ "DMCreate_Plex"
984 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
985 {
986   DM_Plex        *mesh;
987   PetscInt       unit, d;
988   PetscErrorCode ierr;
989 
990   PetscFunctionBegin;
991   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
992   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
993   dm->data = mesh;
994 
995   mesh->refct             = 1;
996   mesh->dim               = 0;
997   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
998   mesh->maxConeSize       = 0;
999   mesh->cones             = NULL;
1000   mesh->coneOrientations  = NULL;
1001   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
1002   mesh->maxSupportSize    = 0;
1003   mesh->supports          = NULL;
1004   mesh->refinementUniform = PETSC_TRUE;
1005   mesh->refinementLimit   = -1.0;
1006 
1007   mesh->facesTmp = NULL;
1008 
1009   mesh->subpointMap = NULL;
1010 
1011   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
1012 
1013   mesh->labels              = NULL;
1014   mesh->depthLabel          = NULL;
1015   mesh->globalVertexNumbers = NULL;
1016   mesh->globalCellNumbers   = NULL;
1017   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1018   mesh->vtkCellHeight     = 0;
1019   mesh->preallocCenterDim = -1;
1020 
1021   mesh->printSetValues = PETSC_FALSE;
1022   mesh->printFEM       = 0;
1023   mesh->printTol       = 1.0e-10;
1024 
1025   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 #undef __FUNCT__
1030 #define __FUNCT__ "DMPlexCreate"
1031 /*@
1032   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1033 
1034   Collective on MPI_Comm
1035 
1036   Input Parameter:
1037 . comm - The communicator for the DMPlex object
1038 
1039   Output Parameter:
1040 . mesh  - The DMPlex object
1041 
1042   Level: beginner
1043 
1044 .keywords: DMPlex, create
1045 @*/
1046 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1047 {
1048   PetscErrorCode ierr;
1049 
1050   PetscFunctionBegin;
1051   PetscValidPointer(mesh,2);
1052   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1053   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 #undef __FUNCT__
1058 #define __FUNCT__ "DMPlexBuildFromCellList_Private"
1059 /*
1060   This takes as input the common mesh generator output, a list of the vertices for each cell
1061 */
1062 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
1063 {
1064   PetscInt      *cone, c, p;
1065   PetscErrorCode ierr;
1066 
1067   PetscFunctionBegin;
1068   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
1069   for (c = 0; c < numCells; ++c) {
1070     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
1071   }
1072   ierr = DMSetUp(dm);CHKERRQ(ierr);
1073   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1074   for (c = 0; c < numCells; ++c) {
1075     for (p = 0; p < numCorners; ++p) {
1076       cone[p] = cells[c*numCorners+p]+numCells;
1077     }
1078     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
1079   }
1080   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1081   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1082   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 #undef __FUNCT__
1087 #define __FUNCT__ "DMPlexBuildCoordinates_Private"
1088 /*
1089   This takes as input the coordinates for each vertex
1090 */
1091 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
1092 {
1093   PetscSection   coordSection;
1094   Vec            coordinates;
1095   PetscScalar   *coords;
1096   PetscInt       coordSize, v, d;
1097   PetscErrorCode ierr;
1098 
1099   PetscFunctionBegin;
1100   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1101   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1102   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
1103   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
1104   for (v = numCells; v < numCells+numVertices; ++v) {
1105     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
1106     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
1107   }
1108   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1109   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1110   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
1111   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1112   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1113   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1114   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1115   for (v = 0; v < numVertices; ++v) {
1116     for (d = 0; d < spaceDim; ++d) {
1117       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
1118     }
1119   }
1120   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1121   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1122   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "DMPlexCreateFromCellList"
1128 /*@C
1129   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
1130 
1131   Input Parameters:
1132 + comm - The communicator
1133 . dim - The topological dimension of the mesh
1134 . numCells - The number of cells
1135 . numVertices - The number of vertices
1136 . numCorners - The number of vertices for each cell
1137 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1138 . cells - An array of numCells*numCorners numbers, the vertices for each cell
1139 . spaceDim - The spatial dimension used for coordinates
1140 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
1141 
1142   Output Parameter:
1143 . dm - The DM
1144 
1145   Note: Two triangles sharing a face
1146 $
1147 $        2
1148 $      / | \
1149 $     /  |  \
1150 $    /   |   \
1151 $   0  0 | 1  3
1152 $    \   |   /
1153 $     \  |  /
1154 $      \ | /
1155 $        1
1156 would have input
1157 $  numCells = 2, numVertices = 4
1158 $  cells = [0 1 2  1 3 2]
1159 $
1160 which would result in the DMPlex
1161 $
1162 $        4
1163 $      / | \
1164 $     /  |  \
1165 $    /   |   \
1166 $   2  0 | 1  5
1167 $    \   |   /
1168 $     \  |  /
1169 $      \ | /
1170 $        3
1171 
1172   Level: beginner
1173 
1174 .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
1175 @*/
1176 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)
1177 {
1178   PetscErrorCode ierr;
1179 
1180   PetscFunctionBegin;
1181   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1182   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1183   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
1184   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
1185   if (interpolate) {
1186     DM idm;
1187 
1188     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1189     ierr = DMDestroy(dm);CHKERRQ(ierr);
1190     *dm  = idm;
1191   }
1192   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 #undef __FUNCT__
1197 #define __FUNCT__ "DMPlexCreateFromDAG"
1198 /*@
1199   DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM
1200 
1201   Input Parameters:
1202 + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension()
1203 . depth - The depth of the DAG
1204 . numPoints - The number of points at each depth
1205 . coneSize - The cone size of each point
1206 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1207 . coneOrientations - The orientation of each cone point
1208 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1209 
1210   Output Parameter:
1211 . dm - The DM
1212 
1213   Note: Two triangles sharing a face would have input
1214 $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1215 $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1216 $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1217 $
1218 which would result in the DMPlex
1219 $
1220 $        4
1221 $      / | \
1222 $     /  |  \
1223 $    /   |   \
1224 $   2  0 | 1  5
1225 $    \   |   /
1226 $     \  |  /
1227 $      \ | /
1228 $        3
1229 $
1230 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1231 
1232   Level: advanced
1233 
1234 .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1235 @*/
1236 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
1237 {
1238   Vec            coordinates;
1239   PetscSection   coordSection;
1240   PetscScalar    *coords;
1241   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
1242   PetscErrorCode ierr;
1243 
1244   PetscFunctionBegin;
1245   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1246   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
1247   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
1248   for (p = pStart; p < pEnd; ++p) {
1249     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
1250   }
1251   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
1252   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
1253     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
1254     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
1255   }
1256   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1257   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1258   /* Build coordinates */
1259   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1260   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1261   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
1262   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
1263   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
1264     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
1265     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
1266   }
1267   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1268   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1269   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
1270   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1271   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1272   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1273   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1274   for (v = 0; v < numPoints[0]; ++v) {
1275     PetscInt off;
1276 
1277     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
1278     for (d = 0; d < dim; ++d) {
1279       coords[off+d] = vertexCoords[v*dim+d];
1280     }
1281   }
1282   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1283   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1284   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1285   PetscFunctionReturn(0);
1286 }
1287