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