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