xref: /petsc/src/dm/impls/plex/plexcreate.c (revision a6dfd86ebdf75fa024070af26d51b62fd16650a3)
1552f7358SJed Brown #define PETSCDM_DLL
2552f7358SJed Brown #include <petsc-private/pleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <petscdmda.h>
4552f7358SJed Brown 
5552f7358SJed Brown #undef __FUNCT__
6552f7358SJed Brown #define __FUNCT__ "DMSetFromOptions_Plex"
7552f7358SJed Brown PetscErrorCode  DMSetFromOptions_Plex(DM dm)
8552f7358SJed Brown {
9552f7358SJed Brown   DM_Plex    *mesh = (DM_Plex *) dm->data;
10552f7358SJed Brown   PetscErrorCode ierr;
11552f7358SJed Brown 
12552f7358SJed Brown   PetscFunctionBegin;
13552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14552f7358SJed Brown   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
15552f7358SJed Brown     /* Handle DMPlex refinement */
16552f7358SJed Brown     /* Handle associated vectors */
17552f7358SJed Brown     /* Handle viewing */
18552f7358SJed Brown     ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, PETSC_NULL);CHKERRQ(ierr);
19552f7358SJed Brown     ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, PETSC_NULL);CHKERRQ(ierr);
20552f7358SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
21552f7358SJed Brown   PetscFunctionReturn(0);
22552f7358SJed Brown }
23552f7358SJed Brown 
24552f7358SJed Brown #undef __FUNCT__
25552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary"
26552f7358SJed Brown /*
27552f7358SJed Brown  Simple square boundary:
28552f7358SJed Brown 
29552f7358SJed Brown  18--5-17--4--16
30552f7358SJed Brown   |     |     |
31552f7358SJed Brown   6    10     3
32552f7358SJed Brown   |     |     |
33552f7358SJed Brown  19-11-20--9--15
34552f7358SJed Brown   |     |     |
35552f7358SJed Brown   7     8     2
36552f7358SJed Brown   |     |     |
37552f7358SJed Brown  12--0-13--1--14
38552f7358SJed Brown */
39552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
40552f7358SJed Brown {
41552f7358SJed Brown   PetscInt       numVertices = (edges[0]+1)*(edges[1]+1);
42552f7358SJed Brown   PetscInt       numEdges    = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
43552f7358SJed Brown   PetscInt       markerTop      = 1;
44552f7358SJed Brown   PetscInt       markerBottom   = 1;
45552f7358SJed Brown   PetscInt       markerRight    = 1;
46552f7358SJed Brown   PetscInt       markerLeft     = 1;
47552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
48552f7358SJed Brown   Vec            coordinates;
49552f7358SJed Brown   PetscSection   coordSection;
50552f7358SJed Brown   PetscScalar   *coords;
51552f7358SJed Brown   PetscInt       coordSize;
52552f7358SJed Brown   PetscMPIInt    rank;
53552f7358SJed Brown   PetscInt       v, vx, vy;
54552f7358SJed Brown   PetscErrorCode ierr;
55552f7358SJed Brown 
56552f7358SJed Brown   PetscFunctionBegin;
57552f7358SJed Brown   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, PETSC_NULL);CHKERRQ(ierr);
58552f7358SJed Brown   if (markerSeparate) {
59552f7358SJed Brown     markerTop    = 1;
60552f7358SJed Brown     markerBottom = 0;
61552f7358SJed Brown     markerRight  = 0;
62552f7358SJed Brown     markerLeft   = 0;
63552f7358SJed Brown   }
64552f7358SJed Brown   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
65552f7358SJed Brown   if (!rank) {
66552f7358SJed Brown     PetscInt e, ex, ey;
67552f7358SJed Brown 
68552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
69552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
70552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
71552f7358SJed Brown     }
72552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
73552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
74552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
75552f7358SJed Brown         PetscInt edge    = vx*edges[1] + ey + edges[0]*(edges[1]+1);
76552f7358SJed Brown         PetscInt vertex  = ey*(edges[0]+1) + vx + numEdges;
77552f7358SJed Brown         PetscInt cone[2] = {vertex, vertex+edges[0]+1};
78552f7358SJed Brown 
79552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
80552f7358SJed Brown         if (vx == edges[0]) {
81552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
82552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
83552f7358SJed Brown           if (ey == edges[1]-1) {
84552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
85552f7358SJed Brown           }
86552f7358SJed Brown         } else if (vx == 0) {
87552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
88552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
89552f7358SJed Brown           if (ey == edges[1]-1) {
90552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
91552f7358SJed Brown           }
92552f7358SJed Brown         }
93552f7358SJed Brown       }
94552f7358SJed Brown     }
95552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
96552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
97552f7358SJed Brown         PetscInt edge    = vy*edges[0]     + ex;
98552f7358SJed Brown         PetscInt vertex  = vy*(edges[0]+1) + ex + numEdges;
99552f7358SJed Brown         PetscInt cone[2] = {vertex, vertex+1};
100552f7358SJed Brown 
101552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
102552f7358SJed Brown         if (vy == edges[1]) {
103552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
104552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
105552f7358SJed Brown           if (ex == edges[0]-1) {
106552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
107552f7358SJed Brown           }
108552f7358SJed Brown         } else if (vy == 0) {
109552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
110552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
111552f7358SJed Brown           if (ex == edges[0]-1) {
112552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
113552f7358SJed Brown           }
114552f7358SJed Brown         }
115552f7358SJed Brown       }
116552f7358SJed Brown     }
117552f7358SJed Brown   }
118552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
119552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
120552f7358SJed Brown   /* Build coordinates */
121552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
122552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
123552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
124552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
125552f7358SJed Brown   }
126552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
127552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
128552f7358SJed Brown   ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
129552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
130552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
131552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
132552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
133552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
134552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
135552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
136552f7358SJed Brown     }
137552f7358SJed Brown   }
138552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
139552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
140552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
141552f7358SJed Brown   PetscFunctionReturn(0);
142552f7358SJed Brown }
143552f7358SJed Brown 
144552f7358SJed Brown #undef __FUNCT__
145552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
146552f7358SJed Brown /*
147552f7358SJed Brown  Simple cubic boundary:
148552f7358SJed Brown 
149552f7358SJed Brown      2-------3
150552f7358SJed Brown     /|      /|
151552f7358SJed Brown    6-------7 |
152552f7358SJed Brown    | |     | |
153552f7358SJed Brown    | 0-----|-1
154552f7358SJed Brown    |/      |/
155552f7358SJed Brown    4-------5
156552f7358SJed Brown */
157552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
158552f7358SJed Brown {
159552f7358SJed Brown   PetscInt       numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
160552f7358SJed Brown   PetscInt       numFaces    = 6;
161552f7358SJed Brown   Vec            coordinates;
162552f7358SJed Brown   PetscSection   coordSection;
163552f7358SJed Brown   PetscScalar   *coords;
164552f7358SJed Brown   PetscInt       coordSize;
165552f7358SJed Brown   PetscMPIInt    rank;
166552f7358SJed Brown   PetscInt       v, vx, vy, vz;
167552f7358SJed Brown   PetscErrorCode ierr;
168552f7358SJed Brown 
169552f7358SJed Brown   PetscFunctionBegin;
170552f7358SJed Brown   if ((faces[0] < 1) || (faces[1] < 1) || (faces[2] < 1)) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Must have at least 1 face per side");
171552f7358SJed Brown   if ((faces[0] > 1) || (faces[1] > 1) || (faces[2] > 1)) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Currently can't handle more than 1 face per side");
172552f7358SJed Brown   ierr = PetscMalloc(numVertices*2 * sizeof(PetscReal), &coords);CHKERRQ(ierr);
173552f7358SJed Brown   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
174552f7358SJed Brown   if (!rank) {
175552f7358SJed Brown     PetscInt f;
176552f7358SJed Brown 
177552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
178552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
179552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
180552f7358SJed Brown     }
181552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
182552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
183552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
184552f7358SJed Brown     }
185552f7358SJed Brown     { /* Side 0 (Front) */
186552f7358SJed Brown       PetscInt cone[4] = {numFaces+4, numFaces+5, numFaces+7, numFaces+6};
187552f7358SJed Brown       ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr);
188552f7358SJed Brown     }
189552f7358SJed Brown     { /* Side 1 (Back) */
190552f7358SJed Brown       PetscInt cone[4] = {numFaces+1, numFaces+0, numFaces+2, numFaces+3};
191552f7358SJed Brown       ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr);
192552f7358SJed Brown     }
193552f7358SJed Brown     { /* Side 2 (Bottom) */
194552f7358SJed Brown       PetscInt cone[4] = {numFaces+0, numFaces+1, numFaces+5, numFaces+4};
195552f7358SJed Brown       ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr);
196552f7358SJed Brown     }
197552f7358SJed Brown     { /* Side 3 (Top) */
198552f7358SJed Brown       PetscInt cone[4] = {numFaces+6, numFaces+7, numFaces+3, numFaces+2};
199552f7358SJed Brown       ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr);
200552f7358SJed Brown     }
201552f7358SJed Brown     { /* Side 4 (Left) */
202552f7358SJed Brown       PetscInt cone[4] = {numFaces+0, numFaces+4, numFaces+6, numFaces+2};
203552f7358SJed Brown       ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr);
204552f7358SJed Brown     }
205552f7358SJed Brown     { /* Side 5 (Right) */
206552f7358SJed Brown       PetscInt cone[4] = {numFaces+5, numFaces+1, numFaces+3, numFaces+7};
207552f7358SJed Brown       ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr);
208552f7358SJed Brown     }
209552f7358SJed Brown   }
210552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
211552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
212552f7358SJed Brown   /* Build coordinates */
213552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
214552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
215552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
216552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
217552f7358SJed Brown   }
218552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
219552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
220552f7358SJed Brown   ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
221552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
222552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
223552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
224552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
225552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
226552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
227552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
228552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
229552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
230552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
231552f7358SJed Brown       }
232552f7358SJed Brown     }
233552f7358SJed Brown   }
234552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
235552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
236552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
237552f7358SJed Brown   PetscFunctionReturn(0);
238552f7358SJed Brown }
239552f7358SJed Brown 
240552f7358SJed Brown #undef __FUNCT__
241552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
242552f7358SJed Brown /*
243552f7358SJed Brown  Simple square mesh:
244552f7358SJed Brown 
245552f7358SJed Brown  22--8-23--9--24
246552f7358SJed Brown   |     |     |
247552f7358SJed Brown  13  2 14  3  15
248552f7358SJed Brown   |     |     |
249552f7358SJed Brown  19--6-20--7--21
250552f7358SJed Brown   |     |     |
251552f7358SJed Brown  10  0 11  1 12
252552f7358SJed Brown   |     |     |
253552f7358SJed Brown  16--4-17--5--18
254552f7358SJed Brown */
255552f7358SJed Brown PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
256552f7358SJed Brown {
257552f7358SJed Brown   PetscInt       markerTop      = 1;
258552f7358SJed Brown   PetscInt       markerBottom   = 1;
259552f7358SJed Brown   PetscInt       markerRight    = 1;
260552f7358SJed Brown   PetscInt       markerLeft     = 1;
261552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
262552f7358SJed Brown   PetscMPIInt    rank;
263552f7358SJed Brown   PetscErrorCode ierr;
264552f7358SJed Brown 
265552f7358SJed Brown   PetscFunctionBegin;
266552f7358SJed Brown   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
267552f7358SJed Brown   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, PETSC_NULL);CHKERRQ(ierr);
268552f7358SJed Brown   if (markerSeparate) {
269552f7358SJed Brown     markerTop    = 3;
270552f7358SJed Brown     markerBottom = 1;
271552f7358SJed Brown     markerRight  = 2;
272552f7358SJed Brown     markerLeft   = 4;
273552f7358SJed Brown   }
274552f7358SJed Brown   {
275552f7358SJed Brown     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
276552f7358SJed Brown     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
277552f7358SJed Brown     const PetscInt numXVertices = !rank ? edges[0]+1 : 0;
278552f7358SJed Brown     const PetscInt numYVertices = !rank ? edges[1]+1 : 0;
279552f7358SJed Brown     const PetscInt numTotXEdges = numXEdges*numYVertices;
280552f7358SJed Brown     const PetscInt numTotYEdges = numYEdges*numXVertices;
281552f7358SJed Brown     const PetscInt numVertices  = numXVertices*numYVertices;
282552f7358SJed Brown     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
283552f7358SJed Brown     const PetscInt numFaces     = numXEdges*numYEdges;
284552f7358SJed Brown     const PetscInt firstVertex  = numFaces;
285552f7358SJed Brown     const PetscInt firstXEdge   = numFaces + numVertices;
286552f7358SJed Brown     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
287552f7358SJed Brown     Vec            coordinates;
288552f7358SJed Brown     PetscSection   coordSection;
289552f7358SJed Brown     PetscScalar   *coords;
290552f7358SJed Brown     PetscInt       coordSize;
291552f7358SJed Brown     PetscInt       v, vx, vy;
292552f7358SJed Brown     PetscInt       f, fx, fy, e, ex, ey;
293552f7358SJed Brown 
294552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
295552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
296552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
297552f7358SJed Brown     }
298552f7358SJed Brown     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
299552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
300552f7358SJed Brown     }
301552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
302552f7358SJed Brown     /* Build faces */
303552f7358SJed Brown     for (fy = 0; fy < numYEdges; fy++) {
304552f7358SJed Brown       for (fx = 0; fx < numXEdges; fx++) {
305552f7358SJed Brown         const PetscInt face    = fy*numXEdges + fx;
306552f7358SJed Brown         const PetscInt edgeL   = firstYEdge + fx*numYEdges + fy;
307552f7358SJed Brown         const PetscInt edgeB   = firstXEdge + fy*numXEdges + fx;
308552f7358SJed Brown         const PetscInt cone[4] = {edgeB, edgeL+numYEdges, edgeB+numXEdges, edgeL};
309552f7358SJed Brown         const PetscInt ornt[4] = {0,     0,               -2,              -2};
310552f7358SJed Brown 
311552f7358SJed Brown         ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
312552f7358SJed Brown         ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
313552f7358SJed Brown       }
314552f7358SJed Brown     }
315552f7358SJed Brown     /* Build Y edges*/
316552f7358SJed Brown     for (vx = 0; vx < numXVertices; vx++) {
317552f7358SJed Brown       for (ey = 0; ey < numYEdges; ey++) {
318552f7358SJed Brown         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
319552f7358SJed Brown         const PetscInt vertex  = firstVertex + ey*numXVertices + vx;
320552f7358SJed Brown         const PetscInt cone[2] = {vertex, vertex+numXVertices};
321552f7358SJed Brown 
322552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
323552f7358SJed Brown         if (vx == numXVertices-1) {
324552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
325552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
326552f7358SJed Brown           if (ey == numYEdges-1) {
327552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
328552f7358SJed Brown           }
329552f7358SJed Brown         } else if (vx == 0) {
330552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
331552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
332552f7358SJed Brown           if (ey == numYEdges-1) {
333552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
334552f7358SJed Brown           }
335552f7358SJed Brown         }
336552f7358SJed Brown       }
337552f7358SJed Brown     }
338552f7358SJed Brown     /* Build X edges*/
339552f7358SJed Brown     for (vy = 0; vy < numYVertices; vy++) {
340552f7358SJed Brown       for (ex = 0; ex < numXEdges; ex++) {
341552f7358SJed Brown         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
342552f7358SJed Brown         const PetscInt vertex  = firstVertex + vy*numXVertices + ex;
343552f7358SJed Brown         const PetscInt cone[2] = {vertex, vertex+1};
344552f7358SJed Brown 
345552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
346552f7358SJed Brown         if (vy == numYVertices-1) {
347552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
348552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
349552f7358SJed Brown           if (ex == numXEdges-1) {
350552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
351552f7358SJed Brown           }
352552f7358SJed Brown         } else if (vy == 0) {
353552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
354552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
355552f7358SJed Brown           if (ex == numXEdges-1) {
356552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
357552f7358SJed Brown           }
358552f7358SJed Brown         }
359552f7358SJed Brown       }
360552f7358SJed Brown     }
361552f7358SJed Brown     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
362552f7358SJed Brown     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
363552f7358SJed Brown     /* Build coordinates */
364552f7358SJed Brown     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
365552f7358SJed Brown     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
366552f7358SJed Brown     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
367552f7358SJed Brown       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
368552f7358SJed Brown     }
369552f7358SJed Brown     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
370552f7358SJed Brown     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
371552f7358SJed Brown     ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
372552f7358SJed Brown     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
373552f7358SJed Brown     ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
374552f7358SJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
375552f7358SJed Brown     for (vy = 0; vy < numYVertices; ++vy) {
376552f7358SJed Brown       for (vx = 0; vx < numXVertices; ++vx) {
377552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
378552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
379552f7358SJed Brown       }
380552f7358SJed Brown     }
381552f7358SJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
382552f7358SJed Brown     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
383552f7358SJed Brown     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
384552f7358SJed Brown   }
385552f7358SJed Brown   PetscFunctionReturn(0);
386552f7358SJed Brown }
387552f7358SJed Brown 
388552f7358SJed Brown #undef __FUNCT__
389552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
3900adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
3910adebc6cSBarry Smith {
392552f7358SJed Brown   DM             boundary;
393552f7358SJed Brown   PetscErrorCode ierr;
394552f7358SJed Brown 
395552f7358SJed Brown   PetscFunctionBegin;
396552f7358SJed Brown   PetscValidPointer(dm, 4);
397552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
398552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
399552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
400552f7358SJed Brown   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
401552f7358SJed Brown   switch(dim) {
402552f7358SJed Brown   case 2:
403552f7358SJed Brown   {
404552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
405552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
406552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
407552f7358SJed Brown 
408552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
409552f7358SJed Brown     break;
410552f7358SJed Brown   }
411552f7358SJed Brown   case 3:
412552f7358SJed Brown   {
413552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
414552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
415552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
416552f7358SJed Brown 
417552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
418552f7358SJed Brown     break;
419552f7358SJed Brown   }
420552f7358SJed Brown   default:
421552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
422552f7358SJed Brown   }
423552f7358SJed Brown   ierr = DMPlexGenerate(boundary, PETSC_NULL, interpolate, dm);CHKERRQ(ierr);
424552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
425552f7358SJed Brown   PetscFunctionReturn(0);
426552f7358SJed Brown }
427552f7358SJed Brown 
428552f7358SJed Brown #undef __FUNCT__
429552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
430*a6dfd86eSKarl Rupp PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DM *dm)
431*a6dfd86eSKarl Rupp {
432552f7358SJed Brown   PetscErrorCode ierr;
433552f7358SJed Brown 
434552f7358SJed Brown   PetscFunctionBegin;
435552f7358SJed Brown   PetscValidPointer(dm, 4);
436552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
437552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
438552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
439552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
440552f7358SJed Brown   switch(dim) {
441552f7358SJed Brown   case 2:
442552f7358SJed Brown   {
443552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
444552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
445552f7358SJed Brown 
446552f7358SJed Brown     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells);CHKERRQ(ierr);
447552f7358SJed Brown     break;
448552f7358SJed Brown   }
449552f7358SJed Brown #if 0
450552f7358SJed Brown   case 3:
451552f7358SJed Brown   {
452552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
453552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
454552f7358SJed Brown 
455552f7358SJed Brown     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells);CHKERRQ(ierr);
456552f7358SJed Brown     break;
457552f7358SJed Brown   }
458552f7358SJed Brown #endif
459552f7358SJed Brown   default:
460552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
461552f7358SJed Brown   }
462552f7358SJed Brown   PetscFunctionReturn(0);
463552f7358SJed Brown }
464552f7358SJed Brown 
465552f7358SJed Brown /* External function declarations here */
466552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
467552f7358SJed Brown extern PetscErrorCode DMCreateMatrix_Plex(DM dm, MatType mtype, Mat *J);
468552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
469552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
470552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
471552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
472552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
473552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
474552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
475552f7358SJed Brown 
476552f7358SJed Brown #undef __FUNCT__
477552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
478552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
479552f7358SJed Brown {
480552f7358SJed Brown   PetscErrorCode ierr;
481552f7358SJed Brown 
482552f7358SJed Brown   PetscFunctionBegin;
483552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
484552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
485552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex);CHKERRQ(ierr);
486552f7358SJed Brown   PetscFunctionReturn(0);
487552f7358SJed Brown }
488552f7358SJed Brown 
489552f7358SJed Brown #undef __FUNCT__
490552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
491552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
492552f7358SJed Brown {
493552f7358SJed Brown   PetscErrorCode ierr;
494552f7358SJed Brown 
495552f7358SJed Brown   PetscFunctionBegin;
496552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
497552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
498552f7358SJed Brown   PetscFunctionReturn(0);
499552f7358SJed Brown }
500552f7358SJed Brown 
501552f7358SJed Brown 
502552f7358SJed Brown #undef __FUNCT__
503552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
504552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
505552f7358SJed Brown {
506552f7358SJed Brown   PetscErrorCode ierr;
507552f7358SJed Brown 
508552f7358SJed Brown   PetscFunctionBegin;
509552f7358SJed Brown   ierr = PetscStrallocpy(VECSTANDARD, (char**)&dm->vectype);CHKERRQ(ierr);
510552f7358SJed Brown   dm->ops->view               = DMView_Plex;
511552f7358SJed Brown   dm->ops->setfromoptions     = DMSetFromOptions_Plex;
512552f7358SJed Brown   dm->ops->setup              = DMSetUp_Plex;
513552f7358SJed Brown   dm->ops->createglobalvector = DMCreateGlobalVector_Plex;
514552f7358SJed Brown   dm->ops->createlocalvector  = DMCreateLocalVector_Plex;
515552f7358SJed Brown   dm->ops->createlocaltoglobalmapping      = PETSC_NULL;
516552f7358SJed Brown   dm->ops->createlocaltoglobalmappingblock = PETSC_NULL;
517552f7358SJed Brown   dm->ops->createfieldis      = PETSC_NULL;
518552f7358SJed Brown   dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex;
519552f7358SJed Brown   dm->ops->getcoloring        = 0;
520552f7358SJed Brown   dm->ops->creatematrix       = DMCreateMatrix_Plex;
521552f7358SJed Brown   dm->ops->createinterpolation= 0;
522552f7358SJed Brown   dm->ops->getaggregates      = 0;
523552f7358SJed Brown   dm->ops->getinjection       = 0;
524552f7358SJed Brown   dm->ops->refine             = DMRefine_Plex;
525552f7358SJed Brown   dm->ops->coarsen            = 0;
526552f7358SJed Brown   dm->ops->refinehierarchy    = 0;
527552f7358SJed Brown   dm->ops->coarsenhierarchy   = 0;
528552f7358SJed Brown   dm->ops->globaltolocalbegin = PETSC_NULL;
529552f7358SJed Brown   dm->ops->globaltolocalend   = PETSC_NULL;
530552f7358SJed Brown   dm->ops->localtoglobalbegin = PETSC_NULL;
531552f7358SJed Brown   dm->ops->localtoglobalend   = PETSC_NULL;
532552f7358SJed Brown   dm->ops->destroy            = DMDestroy_Plex;
533552f7358SJed Brown   dm->ops->createsubdm        = DMCreateSubDM_Plex;
534552f7358SJed Brown   dm->ops->locatepoints       = DMLocatePoints_Plex;
535552f7358SJed Brown   PetscFunctionReturn(0);
536552f7358SJed Brown }
537552f7358SJed Brown 
538552f7358SJed Brown EXTERN_C_BEGIN
539552f7358SJed Brown #undef __FUNCT__
540552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
541552f7358SJed Brown PetscErrorCode DMCreate_Plex(DM dm)
542552f7358SJed Brown {
543552f7358SJed Brown   DM_Plex       *mesh;
544770b213bSMatthew G Knepley   PetscInt       unit, d;
545552f7358SJed Brown   PetscErrorCode ierr;
546552f7358SJed Brown 
547552f7358SJed Brown   PetscFunctionBegin;
548552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
549552f7358SJed Brown   ierr = PetscNewLog(dm, DM_Plex, &mesh);CHKERRQ(ierr);
550552f7358SJed Brown   dm->data = mesh;
551552f7358SJed Brown 
552552f7358SJed Brown   mesh->refct             = 1;
553552f7358SJed Brown   mesh->dim               = 0;
554552f7358SJed Brown   ierr = PetscSectionCreate(((PetscObject) dm)->comm, &mesh->coneSection);CHKERRQ(ierr);
555552f7358SJed Brown   mesh->maxConeSize       = 0;
556552f7358SJed Brown   mesh->cones             = PETSC_NULL;
557552f7358SJed Brown   mesh->coneOrientations  = PETSC_NULL;
558552f7358SJed Brown   ierr = PetscSectionCreate(((PetscObject) dm)->comm, &mesh->supportSection);CHKERRQ(ierr);
559552f7358SJed Brown   mesh->maxSupportSize    = 0;
560552f7358SJed Brown   mesh->supports          = PETSC_NULL;
561552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
562552f7358SJed Brown   mesh->refinementLimit   = -1.0;
563552f7358SJed Brown 
564552f7358SJed Brown   mesh->facesTmp         = PETSC_NULL;
565552f7358SJed Brown 
566552f7358SJed Brown   mesh->subpointMap      = PETSC_NULL;
567552f7358SJed Brown 
568552f7358SJed Brown   for(unit = 0; unit < NUM_PETSC_UNITS; ++unit) {
569552f7358SJed Brown     mesh->scale[unit]    = 1.0;
570552f7358SJed Brown   }
571552f7358SJed Brown 
572552f7358SJed Brown   mesh->labels               = PETSC_NULL;
573552f7358SJed Brown   mesh->globalVertexNumbers  = PETSC_NULL;
574552f7358SJed Brown   mesh->globalCellNumbers    = PETSC_NULL;
575770b213bSMatthew G Knepley   for(d = 0; d < 8; ++d) {
576770b213bSMatthew G Knepley     mesh->hybridPointMax[d]  = PETSC_DETERMINE;
577770b213bSMatthew G Knepley   }
578552f7358SJed Brown   mesh->vtkCellHeight        = 0;
579552f7358SJed Brown   mesh->preallocCenterDim    = -1;
580552f7358SJed Brown 
581552f7358SJed Brown   mesh->integrateResidualFEM       = PETSC_NULL;
582552f7358SJed Brown   mesh->integrateJacobianActionFEM = PETSC_NULL;
583552f7358SJed Brown   mesh->integrateJacobianFEM       = PETSC_NULL;
584552f7358SJed Brown 
585552f7358SJed Brown   mesh->printSetValues       = PETSC_FALSE;
586552f7358SJed Brown   mesh->printFEM             = 0;
587552f7358SJed Brown 
588552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
589552f7358SJed Brown   PetscFunctionReturn(0);
590552f7358SJed Brown }
591552f7358SJed Brown EXTERN_C_END
592552f7358SJed Brown 
593552f7358SJed Brown #undef __FUNCT__
594552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
595552f7358SJed Brown /*@
596552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
597552f7358SJed Brown 
598552f7358SJed Brown   Collective on MPI_Comm
599552f7358SJed Brown 
600552f7358SJed Brown   Input Parameter:
601552f7358SJed Brown . comm - The communicator for the DMPlex object
602552f7358SJed Brown 
603552f7358SJed Brown   Output Parameter:
604552f7358SJed Brown . mesh  - The DMPlex object
605552f7358SJed Brown 
606552f7358SJed Brown   Level: beginner
607552f7358SJed Brown 
608552f7358SJed Brown .keywords: DMPlex, create
609552f7358SJed Brown @*/
610552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
611552f7358SJed Brown {
612552f7358SJed Brown   PetscErrorCode ierr;
613552f7358SJed Brown 
614552f7358SJed Brown   PetscFunctionBegin;
615552f7358SJed Brown   PetscValidPointer(mesh,2);
616552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
617552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
618552f7358SJed Brown   PetscFunctionReturn(0);
619552f7358SJed Brown }
620552f7358SJed Brown 
621552f7358SJed Brown #undef __FUNCT__
622552f7358SJed Brown #define __FUNCT__ "DMPlexClone"
623552f7358SJed Brown /*@
624552f7358SJed Brown   DMPlexClone - Creates a DMPlex object with the same mesh as the original.
625552f7358SJed Brown 
626552f7358SJed Brown   Collective on MPI_Comm
627552f7358SJed Brown 
628552f7358SJed Brown   Input Parameter:
629552f7358SJed Brown . dm - The original DMPlex object
630552f7358SJed Brown 
631552f7358SJed Brown   Output Parameter:
632552f7358SJed Brown . newdm  - The new DMPlex object
633552f7358SJed Brown 
634552f7358SJed Brown   Level: beginner
635552f7358SJed Brown 
636552f7358SJed Brown .keywords: DMPlex, create
637552f7358SJed Brown @*/
638552f7358SJed Brown PetscErrorCode DMPlexClone(DM dm, DM *newdm)
639552f7358SJed Brown {
640552f7358SJed Brown   DM_Plex    *mesh;
641552f7358SJed Brown   void          *ctx;
642552f7358SJed Brown   PetscErrorCode ierr;
643552f7358SJed Brown 
644552f7358SJed Brown   PetscFunctionBegin;
645552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
646552f7358SJed Brown   PetscValidPointer(newdm,2);
647552f7358SJed Brown   ierr = DMCreate(((PetscObject) dm)->comm, newdm);CHKERRQ(ierr);
648552f7358SJed Brown   ierr = PetscSFDestroy(&(*newdm)->sf);CHKERRQ(ierr);
649552f7358SJed Brown   ierr = PetscObjectReference((PetscObject) dm->sf);CHKERRQ(ierr);
650552f7358SJed Brown   (*newdm)->sf = dm->sf;
651552f7358SJed Brown   mesh = (DM_Plex *) dm->data;
652552f7358SJed Brown   mesh->refct++;
653552f7358SJed Brown   (*newdm)->data = mesh;
654552f7358SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
655552f7358SJed Brown   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
656552f7358SJed Brown   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
657552f7358SJed Brown   ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr);
658552f7358SJed Brown   PetscFunctionReturn(0);
659552f7358SJed Brown }
660