xref: /petsc/src/dm/impls/plex/plexcreate.c (revision 552f735842aa6127d09b62f740a903cdd0631614)
1*552f7358SJed Brown #define PETSCDM_DLL
2*552f7358SJed Brown #include <petsc-private/pleximpl.h>    /*I   "petscdmplex.h"   I*/
3*552f7358SJed Brown #include <petscdmda.h>
4*552f7358SJed Brown 
5*552f7358SJed Brown #undef __FUNCT__
6*552f7358SJed Brown #define __FUNCT__ "DMSetFromOptions_Plex"
7*552f7358SJed Brown PetscErrorCode  DMSetFromOptions_Plex(DM dm)
8*552f7358SJed Brown {
9*552f7358SJed Brown   DM_Plex    *mesh = (DM_Plex *) dm->data;
10*552f7358SJed Brown   PetscErrorCode ierr;
11*552f7358SJed Brown 
12*552f7358SJed Brown   PetscFunctionBegin;
13*552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14*552f7358SJed Brown   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
15*552f7358SJed Brown     /* Handle DMPlex refinement */
16*552f7358SJed Brown     /* Handle associated vectors */
17*552f7358SJed Brown     /* Handle viewing */
18*552f7358SJed Brown     ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, PETSC_NULL);CHKERRQ(ierr);
19*552f7358SJed Brown     ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, PETSC_NULL);CHKERRQ(ierr);
20*552f7358SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
21*552f7358SJed Brown   PetscFunctionReturn(0);
22*552f7358SJed Brown }
23*552f7358SJed Brown 
24*552f7358SJed Brown #undef __FUNCT__
25*552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary"
26*552f7358SJed Brown /*
27*552f7358SJed Brown  Simple square boundary:
28*552f7358SJed Brown 
29*552f7358SJed Brown  18--5-17--4--16
30*552f7358SJed Brown   |     |     |
31*552f7358SJed Brown   6    10     3
32*552f7358SJed Brown   |     |     |
33*552f7358SJed Brown  19-11-20--9--15
34*552f7358SJed Brown   |     |     |
35*552f7358SJed Brown   7     8     2
36*552f7358SJed Brown   |     |     |
37*552f7358SJed Brown  12--0-13--1--14
38*552f7358SJed Brown */
39*552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
40*552f7358SJed Brown {
41*552f7358SJed Brown   PetscInt       numVertices = (edges[0]+1)*(edges[1]+1);
42*552f7358SJed Brown   PetscInt       numEdges    = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
43*552f7358SJed Brown   PetscInt       markerTop      = 1;
44*552f7358SJed Brown   PetscInt       markerBottom   = 1;
45*552f7358SJed Brown   PetscInt       markerRight    = 1;
46*552f7358SJed Brown   PetscInt       markerLeft     = 1;
47*552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
48*552f7358SJed Brown   Vec            coordinates;
49*552f7358SJed Brown   PetscSection   coordSection;
50*552f7358SJed Brown   PetscScalar   *coords;
51*552f7358SJed Brown   PetscInt       coordSize;
52*552f7358SJed Brown   PetscMPIInt    rank;
53*552f7358SJed Brown   PetscInt       v, vx, vy;
54*552f7358SJed Brown   PetscErrorCode ierr;
55*552f7358SJed Brown 
56*552f7358SJed Brown   PetscFunctionBegin;
57*552f7358SJed Brown   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, PETSC_NULL);CHKERRQ(ierr);
58*552f7358SJed Brown   if (markerSeparate) {
59*552f7358SJed Brown     markerTop    = 1;
60*552f7358SJed Brown     markerBottom = 0;
61*552f7358SJed Brown     markerRight  = 0;
62*552f7358SJed Brown     markerLeft   = 0;
63*552f7358SJed Brown   }
64*552f7358SJed Brown   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
65*552f7358SJed Brown   if (!rank) {
66*552f7358SJed Brown     PetscInt e, ex, ey;
67*552f7358SJed Brown 
68*552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
69*552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
70*552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
71*552f7358SJed Brown     }
72*552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
73*552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
74*552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
75*552f7358SJed Brown         PetscInt edge    = vx*edges[1] + ey + edges[0]*(edges[1]+1);
76*552f7358SJed Brown         PetscInt vertex  = ey*(edges[0]+1) + vx + numEdges;
77*552f7358SJed Brown         PetscInt cone[2] = {vertex, vertex+edges[0]+1};
78*552f7358SJed Brown 
79*552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
80*552f7358SJed Brown         if (vx == edges[0]) {
81*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
82*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
83*552f7358SJed Brown           if (ey == edges[1]-1) {
84*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
85*552f7358SJed Brown           }
86*552f7358SJed Brown         } else if (vx == 0) {
87*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
88*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
89*552f7358SJed Brown           if (ey == edges[1]-1) {
90*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
91*552f7358SJed Brown           }
92*552f7358SJed Brown         }
93*552f7358SJed Brown       }
94*552f7358SJed Brown     }
95*552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
96*552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
97*552f7358SJed Brown         PetscInt edge    = vy*edges[0]     + ex;
98*552f7358SJed Brown         PetscInt vertex  = vy*(edges[0]+1) + ex + numEdges;
99*552f7358SJed Brown         PetscInt cone[2] = {vertex, vertex+1};
100*552f7358SJed Brown 
101*552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
102*552f7358SJed Brown         if (vy == edges[1]) {
103*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
104*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
105*552f7358SJed Brown           if (ex == edges[0]-1) {
106*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
107*552f7358SJed Brown           }
108*552f7358SJed Brown         } else if (vy == 0) {
109*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
110*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
111*552f7358SJed Brown           if (ex == edges[0]-1) {
112*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
113*552f7358SJed Brown           }
114*552f7358SJed Brown         }
115*552f7358SJed Brown       }
116*552f7358SJed Brown     }
117*552f7358SJed Brown   }
118*552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
119*552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
120*552f7358SJed Brown   /* Build coordinates */
121*552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
122*552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
123*552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
124*552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
125*552f7358SJed Brown   }
126*552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
127*552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
128*552f7358SJed Brown   ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
129*552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
130*552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
131*552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
132*552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
133*552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
134*552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
135*552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
136*552f7358SJed Brown     }
137*552f7358SJed Brown   }
138*552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
139*552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
140*552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
141*552f7358SJed Brown   PetscFunctionReturn(0);
142*552f7358SJed Brown }
143*552f7358SJed Brown 
144*552f7358SJed Brown #undef __FUNCT__
145*552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
146*552f7358SJed Brown /*
147*552f7358SJed Brown  Simple cubic boundary:
148*552f7358SJed Brown 
149*552f7358SJed Brown      2-------3
150*552f7358SJed Brown     /|      /|
151*552f7358SJed Brown    6-------7 |
152*552f7358SJed Brown    | |     | |
153*552f7358SJed Brown    | 0-----|-1
154*552f7358SJed Brown    |/      |/
155*552f7358SJed Brown    4-------5
156*552f7358SJed Brown */
157*552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
158*552f7358SJed Brown {
159*552f7358SJed Brown   PetscInt       numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
160*552f7358SJed Brown   PetscInt       numFaces    = 6;
161*552f7358SJed Brown   Vec            coordinates;
162*552f7358SJed Brown   PetscSection   coordSection;
163*552f7358SJed Brown   PetscScalar   *coords;
164*552f7358SJed Brown   PetscInt       coordSize;
165*552f7358SJed Brown   PetscMPIInt    rank;
166*552f7358SJed Brown   PetscInt       v, vx, vy, vz;
167*552f7358SJed Brown   PetscErrorCode ierr;
168*552f7358SJed Brown 
169*552f7358SJed Brown   PetscFunctionBegin;
170*552f7358SJed 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");
171*552f7358SJed 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");
172*552f7358SJed Brown   ierr = PetscMalloc(numVertices*2 * sizeof(PetscReal), &coords);CHKERRQ(ierr);
173*552f7358SJed Brown   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
174*552f7358SJed Brown   if (!rank) {
175*552f7358SJed Brown     PetscInt f;
176*552f7358SJed Brown 
177*552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
178*552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
179*552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
180*552f7358SJed Brown     }
181*552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
182*552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
183*552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
184*552f7358SJed Brown     }
185*552f7358SJed Brown     { /* Side 0 (Front) */
186*552f7358SJed Brown       PetscInt cone[4] = {numFaces+4, numFaces+5, numFaces+7, numFaces+6};
187*552f7358SJed Brown       ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr);
188*552f7358SJed Brown     }
189*552f7358SJed Brown     { /* Side 1 (Back) */
190*552f7358SJed Brown       PetscInt cone[4] = {numFaces+1, numFaces+0, numFaces+2, numFaces+3};
191*552f7358SJed Brown       ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr);
192*552f7358SJed Brown     }
193*552f7358SJed Brown     { /* Side 2 (Bottom) */
194*552f7358SJed Brown       PetscInt cone[4] = {numFaces+0, numFaces+1, numFaces+5, numFaces+4};
195*552f7358SJed Brown       ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr);
196*552f7358SJed Brown     }
197*552f7358SJed Brown     { /* Side 3 (Top) */
198*552f7358SJed Brown       PetscInt cone[4] = {numFaces+6, numFaces+7, numFaces+3, numFaces+2};
199*552f7358SJed Brown       ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr);
200*552f7358SJed Brown     }
201*552f7358SJed Brown     { /* Side 4 (Left) */
202*552f7358SJed Brown       PetscInt cone[4] = {numFaces+0, numFaces+4, numFaces+6, numFaces+2};
203*552f7358SJed Brown       ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr);
204*552f7358SJed Brown     }
205*552f7358SJed Brown     { /* Side 5 (Right) */
206*552f7358SJed Brown       PetscInt cone[4] = {numFaces+5, numFaces+1, numFaces+3, numFaces+7};
207*552f7358SJed Brown       ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr);
208*552f7358SJed Brown     }
209*552f7358SJed Brown   }
210*552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
211*552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
212*552f7358SJed Brown   /* Build coordinates */
213*552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
214*552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
215*552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
216*552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
217*552f7358SJed Brown   }
218*552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
219*552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
220*552f7358SJed Brown   ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
221*552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
222*552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
223*552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
224*552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
225*552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
226*552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
227*552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
228*552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
229*552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
230*552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
231*552f7358SJed Brown       }
232*552f7358SJed Brown     }
233*552f7358SJed Brown   }
234*552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
235*552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
236*552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
237*552f7358SJed Brown   PetscFunctionReturn(0);
238*552f7358SJed Brown }
239*552f7358SJed Brown 
240*552f7358SJed Brown #undef __FUNCT__
241*552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
242*552f7358SJed Brown /*
243*552f7358SJed Brown  Simple square mesh:
244*552f7358SJed Brown 
245*552f7358SJed Brown  22--8-23--9--24
246*552f7358SJed Brown   |     |     |
247*552f7358SJed Brown  13  2 14  3  15
248*552f7358SJed Brown   |     |     |
249*552f7358SJed Brown  19--6-20--7--21
250*552f7358SJed Brown   |     |     |
251*552f7358SJed Brown  10  0 11  1 12
252*552f7358SJed Brown   |     |     |
253*552f7358SJed Brown  16--4-17--5--18
254*552f7358SJed Brown */
255*552f7358SJed Brown PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
256*552f7358SJed Brown {
257*552f7358SJed Brown   PetscInt       markerTop      = 1;
258*552f7358SJed Brown   PetscInt       markerBottom   = 1;
259*552f7358SJed Brown   PetscInt       markerRight    = 1;
260*552f7358SJed Brown   PetscInt       markerLeft     = 1;
261*552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
262*552f7358SJed Brown   PetscMPIInt    rank;
263*552f7358SJed Brown   PetscErrorCode ierr;
264*552f7358SJed Brown 
265*552f7358SJed Brown   PetscFunctionBegin;
266*552f7358SJed Brown   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
267*552f7358SJed Brown   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, PETSC_NULL);CHKERRQ(ierr);
268*552f7358SJed Brown   if (markerSeparate) {
269*552f7358SJed Brown     markerTop    = 3;
270*552f7358SJed Brown     markerBottom = 1;
271*552f7358SJed Brown     markerRight  = 2;
272*552f7358SJed Brown     markerLeft   = 4;
273*552f7358SJed Brown   }
274*552f7358SJed Brown   {
275*552f7358SJed Brown     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
276*552f7358SJed Brown     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
277*552f7358SJed Brown     const PetscInt numXVertices = !rank ? edges[0]+1 : 0;
278*552f7358SJed Brown     const PetscInt numYVertices = !rank ? edges[1]+1 : 0;
279*552f7358SJed Brown     const PetscInt numTotXEdges = numXEdges*numYVertices;
280*552f7358SJed Brown     const PetscInt numTotYEdges = numYEdges*numXVertices;
281*552f7358SJed Brown     const PetscInt numVertices  = numXVertices*numYVertices;
282*552f7358SJed Brown     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
283*552f7358SJed Brown     const PetscInt numFaces     = numXEdges*numYEdges;
284*552f7358SJed Brown     const PetscInt firstVertex  = numFaces;
285*552f7358SJed Brown     const PetscInt firstXEdge   = numFaces + numVertices;
286*552f7358SJed Brown     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
287*552f7358SJed Brown     Vec            coordinates;
288*552f7358SJed Brown     PetscSection   coordSection;
289*552f7358SJed Brown     PetscScalar   *coords;
290*552f7358SJed Brown     PetscInt       coordSize;
291*552f7358SJed Brown     PetscInt       v, vx, vy;
292*552f7358SJed Brown     PetscInt       f, fx, fy, e, ex, ey;
293*552f7358SJed Brown 
294*552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
295*552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
296*552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
297*552f7358SJed Brown     }
298*552f7358SJed Brown     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
299*552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
300*552f7358SJed Brown     }
301*552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
302*552f7358SJed Brown     /* Build faces */
303*552f7358SJed Brown     for (fy = 0; fy < numYEdges; fy++) {
304*552f7358SJed Brown       for (fx = 0; fx < numXEdges; fx++) {
305*552f7358SJed Brown         const PetscInt face    = fy*numXEdges + fx;
306*552f7358SJed Brown         const PetscInt edgeL   = firstYEdge + fx*numYEdges + fy;
307*552f7358SJed Brown         const PetscInt edgeB   = firstXEdge + fy*numXEdges + fx;
308*552f7358SJed Brown         const PetscInt cone[4] = {edgeB, edgeL+numYEdges, edgeB+numXEdges, edgeL};
309*552f7358SJed Brown         const PetscInt ornt[4] = {0,     0,               -2,              -2};
310*552f7358SJed Brown 
311*552f7358SJed Brown         ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
312*552f7358SJed Brown         ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
313*552f7358SJed Brown       }
314*552f7358SJed Brown     }
315*552f7358SJed Brown     /* Build Y edges*/
316*552f7358SJed Brown     for (vx = 0; vx < numXVertices; vx++) {
317*552f7358SJed Brown       for (ey = 0; ey < numYEdges; ey++) {
318*552f7358SJed Brown         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
319*552f7358SJed Brown         const PetscInt vertex  = firstVertex + ey*numXVertices + vx;
320*552f7358SJed Brown         const PetscInt cone[2] = {vertex, vertex+numXVertices};
321*552f7358SJed Brown 
322*552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
323*552f7358SJed Brown         if (vx == numXVertices-1) {
324*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
325*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
326*552f7358SJed Brown           if (ey == numYEdges-1) {
327*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
328*552f7358SJed Brown           }
329*552f7358SJed Brown         } else if (vx == 0) {
330*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
331*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
332*552f7358SJed Brown           if (ey == numYEdges-1) {
333*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
334*552f7358SJed Brown           }
335*552f7358SJed Brown         }
336*552f7358SJed Brown       }
337*552f7358SJed Brown     }
338*552f7358SJed Brown     /* Build X edges*/
339*552f7358SJed Brown     for (vy = 0; vy < numYVertices; vy++) {
340*552f7358SJed Brown       for (ex = 0; ex < numXEdges; ex++) {
341*552f7358SJed Brown         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
342*552f7358SJed Brown         const PetscInt vertex  = firstVertex + vy*numXVertices + ex;
343*552f7358SJed Brown         const PetscInt cone[2] = {vertex, vertex+1};
344*552f7358SJed Brown 
345*552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
346*552f7358SJed Brown         if (vy == numYVertices-1) {
347*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
348*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
349*552f7358SJed Brown           if (ex == numXEdges-1) {
350*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
351*552f7358SJed Brown           }
352*552f7358SJed Brown         } else if (vy == 0) {
353*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
354*552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
355*552f7358SJed Brown           if (ex == numXEdges-1) {
356*552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
357*552f7358SJed Brown           }
358*552f7358SJed Brown         }
359*552f7358SJed Brown       }
360*552f7358SJed Brown     }
361*552f7358SJed Brown     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
362*552f7358SJed Brown     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
363*552f7358SJed Brown     /* Build coordinates */
364*552f7358SJed Brown     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
365*552f7358SJed Brown     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
366*552f7358SJed Brown     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
367*552f7358SJed Brown       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
368*552f7358SJed Brown     }
369*552f7358SJed Brown     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
370*552f7358SJed Brown     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
371*552f7358SJed Brown     ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
372*552f7358SJed Brown     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
373*552f7358SJed Brown     ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
374*552f7358SJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
375*552f7358SJed Brown     for (vy = 0; vy < numYVertices; ++vy) {
376*552f7358SJed Brown       for (vx = 0; vx < numXVertices; ++vx) {
377*552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
378*552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
379*552f7358SJed Brown       }
380*552f7358SJed Brown     }
381*552f7358SJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
382*552f7358SJed Brown     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
383*552f7358SJed Brown     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
384*552f7358SJed Brown   }
385*552f7358SJed Brown   PetscFunctionReturn(0);
386*552f7358SJed Brown }
387*552f7358SJed Brown 
388*552f7358SJed Brown #undef __FUNCT__
389*552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
390*552f7358SJed Brown PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) {
391*552f7358SJed Brown   DM             boundary;
392*552f7358SJed Brown   PetscErrorCode ierr;
393*552f7358SJed Brown 
394*552f7358SJed Brown   PetscFunctionBegin;
395*552f7358SJed Brown   PetscValidPointer(dm, 4);
396*552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
397*552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
398*552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
399*552f7358SJed Brown   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
400*552f7358SJed Brown   switch(dim) {
401*552f7358SJed Brown   case 2:
402*552f7358SJed Brown   {
403*552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
404*552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
405*552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
406*552f7358SJed Brown 
407*552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
408*552f7358SJed Brown     break;
409*552f7358SJed Brown   }
410*552f7358SJed Brown   case 3:
411*552f7358SJed Brown   {
412*552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
413*552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
414*552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
415*552f7358SJed Brown 
416*552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
417*552f7358SJed Brown     break;
418*552f7358SJed Brown   }
419*552f7358SJed Brown   default:
420*552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
421*552f7358SJed Brown   }
422*552f7358SJed Brown   ierr = DMPlexGenerate(boundary, PETSC_NULL, interpolate, dm);CHKERRQ(ierr);
423*552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
424*552f7358SJed Brown   PetscFunctionReturn(0);
425*552f7358SJed Brown }
426*552f7358SJed Brown 
427*552f7358SJed Brown #undef __FUNCT__
428*552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
429*552f7358SJed Brown PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DM *dm) {
430*552f7358SJed Brown   PetscErrorCode ierr;
431*552f7358SJed Brown 
432*552f7358SJed Brown   PetscFunctionBegin;
433*552f7358SJed Brown   PetscValidPointer(dm, 4);
434*552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
435*552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
436*552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
437*552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
438*552f7358SJed Brown   switch(dim) {
439*552f7358SJed Brown   case 2:
440*552f7358SJed Brown   {
441*552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
442*552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
443*552f7358SJed Brown 
444*552f7358SJed Brown     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells);CHKERRQ(ierr);
445*552f7358SJed Brown     break;
446*552f7358SJed Brown   }
447*552f7358SJed Brown #if 0
448*552f7358SJed Brown   case 3:
449*552f7358SJed Brown   {
450*552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
451*552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
452*552f7358SJed Brown 
453*552f7358SJed Brown     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells);CHKERRQ(ierr);
454*552f7358SJed Brown     break;
455*552f7358SJed Brown   }
456*552f7358SJed Brown #endif
457*552f7358SJed Brown   default:
458*552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
459*552f7358SJed Brown   }
460*552f7358SJed Brown   PetscFunctionReturn(0);
461*552f7358SJed Brown }
462*552f7358SJed Brown 
463*552f7358SJed Brown /* External function declarations here */
464*552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
465*552f7358SJed Brown extern PetscErrorCode DMCreateMatrix_Plex(DM dm, MatType mtype, Mat *J);
466*552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
467*552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
468*552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
469*552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
470*552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
471*552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
472*552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
473*552f7358SJed Brown 
474*552f7358SJed Brown #undef __FUNCT__
475*552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
476*552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
477*552f7358SJed Brown {
478*552f7358SJed Brown   PetscErrorCode ierr;
479*552f7358SJed Brown 
480*552f7358SJed Brown   PetscFunctionBegin;
481*552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
482*552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
483*552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex);CHKERRQ(ierr);
484*552f7358SJed Brown   PetscFunctionReturn(0);
485*552f7358SJed Brown }
486*552f7358SJed Brown 
487*552f7358SJed Brown #undef __FUNCT__
488*552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
489*552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
490*552f7358SJed Brown {
491*552f7358SJed Brown   PetscErrorCode ierr;
492*552f7358SJed Brown 
493*552f7358SJed Brown   PetscFunctionBegin;
494*552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
495*552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
496*552f7358SJed Brown   PetscFunctionReturn(0);
497*552f7358SJed Brown }
498*552f7358SJed Brown 
499*552f7358SJed Brown 
500*552f7358SJed Brown #undef __FUNCT__
501*552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
502*552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
503*552f7358SJed Brown {
504*552f7358SJed Brown   PetscErrorCode ierr;
505*552f7358SJed Brown 
506*552f7358SJed Brown   PetscFunctionBegin;
507*552f7358SJed Brown   ierr = PetscStrallocpy(VECSTANDARD, (char**)&dm->vectype);CHKERRQ(ierr);
508*552f7358SJed Brown   dm->ops->view               = DMView_Plex;
509*552f7358SJed Brown   dm->ops->setfromoptions     = DMSetFromOptions_Plex;
510*552f7358SJed Brown   dm->ops->setup              = DMSetUp_Plex;
511*552f7358SJed Brown   dm->ops->createglobalvector = DMCreateGlobalVector_Plex;
512*552f7358SJed Brown   dm->ops->createlocalvector  = DMCreateLocalVector_Plex;
513*552f7358SJed Brown   dm->ops->createlocaltoglobalmapping      = PETSC_NULL;
514*552f7358SJed Brown   dm->ops->createlocaltoglobalmappingblock = PETSC_NULL;
515*552f7358SJed Brown   dm->ops->createfieldis      = PETSC_NULL;
516*552f7358SJed Brown   dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex;
517*552f7358SJed Brown   dm->ops->getcoloring        = 0;
518*552f7358SJed Brown   dm->ops->creatematrix       = DMCreateMatrix_Plex;
519*552f7358SJed Brown   dm->ops->createinterpolation= 0;
520*552f7358SJed Brown   dm->ops->getaggregates      = 0;
521*552f7358SJed Brown   dm->ops->getinjection       = 0;
522*552f7358SJed Brown   dm->ops->refine             = DMRefine_Plex;
523*552f7358SJed Brown   dm->ops->coarsen            = 0;
524*552f7358SJed Brown   dm->ops->refinehierarchy    = 0;
525*552f7358SJed Brown   dm->ops->coarsenhierarchy   = 0;
526*552f7358SJed Brown   dm->ops->globaltolocalbegin = PETSC_NULL;
527*552f7358SJed Brown   dm->ops->globaltolocalend   = PETSC_NULL;
528*552f7358SJed Brown   dm->ops->localtoglobalbegin = PETSC_NULL;
529*552f7358SJed Brown   dm->ops->localtoglobalend   = PETSC_NULL;
530*552f7358SJed Brown   dm->ops->destroy            = DMDestroy_Plex;
531*552f7358SJed Brown   dm->ops->createsubdm        = DMCreateSubDM_Plex;
532*552f7358SJed Brown   dm->ops->locatepoints       = DMLocatePoints_Plex;
533*552f7358SJed Brown   PetscFunctionReturn(0);
534*552f7358SJed Brown }
535*552f7358SJed Brown 
536*552f7358SJed Brown EXTERN_C_BEGIN
537*552f7358SJed Brown #undef __FUNCT__
538*552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
539*552f7358SJed Brown PetscErrorCode DMCreate_Plex(DM dm)
540*552f7358SJed Brown {
541*552f7358SJed Brown   DM_Plex    *mesh;
542*552f7358SJed Brown   PetscInt       unit;
543*552f7358SJed Brown   PetscErrorCode ierr;
544*552f7358SJed Brown 
545*552f7358SJed Brown   PetscFunctionBegin;
546*552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
547*552f7358SJed Brown   ierr = PetscNewLog(dm, DM_Plex, &mesh);CHKERRQ(ierr);
548*552f7358SJed Brown   dm->data = mesh;
549*552f7358SJed Brown 
550*552f7358SJed Brown   mesh->refct             = 1;
551*552f7358SJed Brown   mesh->dim               = 0;
552*552f7358SJed Brown   ierr = PetscSectionCreate(((PetscObject) dm)->comm, &mesh->coneSection);CHKERRQ(ierr);
553*552f7358SJed Brown   mesh->maxConeSize       = 0;
554*552f7358SJed Brown   mesh->cones             = PETSC_NULL;
555*552f7358SJed Brown   mesh->coneOrientations  = PETSC_NULL;
556*552f7358SJed Brown   ierr = PetscSectionCreate(((PetscObject) dm)->comm, &mesh->supportSection);CHKERRQ(ierr);
557*552f7358SJed Brown   mesh->maxSupportSize    = 0;
558*552f7358SJed Brown   mesh->supports          = PETSC_NULL;
559*552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
560*552f7358SJed Brown   mesh->refinementLimit   = -1.0;
561*552f7358SJed Brown 
562*552f7358SJed Brown   mesh->facesTmp         = PETSC_NULL;
563*552f7358SJed Brown 
564*552f7358SJed Brown   mesh->subpointMap      = PETSC_NULL;
565*552f7358SJed Brown 
566*552f7358SJed Brown   for(unit = 0; unit < NUM_PETSC_UNITS; ++unit) {
567*552f7358SJed Brown     mesh->scale[unit]    = 1.0;
568*552f7358SJed Brown   }
569*552f7358SJed Brown 
570*552f7358SJed Brown   mesh->labels               = PETSC_NULL;
571*552f7358SJed Brown   mesh->globalVertexNumbers  = PETSC_NULL;
572*552f7358SJed Brown   mesh->globalCellNumbers    = PETSC_NULL;
573*552f7358SJed Brown   mesh->vtkCellMax           = PETSC_DETERMINE;
574*552f7358SJed Brown   mesh->vtkVertexMax         = PETSC_DETERMINE;
575*552f7358SJed Brown   mesh->vtkCellHeight        = 0;
576*552f7358SJed Brown   mesh->preallocCenterDim    = -1;
577*552f7358SJed Brown 
578*552f7358SJed Brown   mesh->integrateResidualFEM       = PETSC_NULL;
579*552f7358SJed Brown   mesh->integrateJacobianActionFEM = PETSC_NULL;
580*552f7358SJed Brown   mesh->integrateJacobianFEM       = PETSC_NULL;
581*552f7358SJed Brown 
582*552f7358SJed Brown   mesh->printSetValues       = PETSC_FALSE;
583*552f7358SJed Brown   mesh->printFEM             = 0;
584*552f7358SJed Brown 
585*552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
586*552f7358SJed Brown   PetscFunctionReturn(0);
587*552f7358SJed Brown }
588*552f7358SJed Brown EXTERN_C_END
589*552f7358SJed Brown 
590*552f7358SJed Brown #undef __FUNCT__
591*552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
592*552f7358SJed Brown /*@
593*552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
594*552f7358SJed Brown 
595*552f7358SJed Brown   Collective on MPI_Comm
596*552f7358SJed Brown 
597*552f7358SJed Brown   Input Parameter:
598*552f7358SJed Brown . comm - The communicator for the DMPlex object
599*552f7358SJed Brown 
600*552f7358SJed Brown   Output Parameter:
601*552f7358SJed Brown . mesh  - The DMPlex object
602*552f7358SJed Brown 
603*552f7358SJed Brown   Level: beginner
604*552f7358SJed Brown 
605*552f7358SJed Brown .keywords: DMPlex, create
606*552f7358SJed Brown @*/
607*552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
608*552f7358SJed Brown {
609*552f7358SJed Brown   PetscErrorCode ierr;
610*552f7358SJed Brown 
611*552f7358SJed Brown   PetscFunctionBegin;
612*552f7358SJed Brown   PetscValidPointer(mesh,2);
613*552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
614*552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
615*552f7358SJed Brown   PetscFunctionReturn(0);
616*552f7358SJed Brown }
617*552f7358SJed Brown 
618*552f7358SJed Brown #undef __FUNCT__
619*552f7358SJed Brown #define __FUNCT__ "DMPlexClone"
620*552f7358SJed Brown /*@
621*552f7358SJed Brown   DMPlexClone - Creates a DMPlex object with the same mesh as the original.
622*552f7358SJed Brown 
623*552f7358SJed Brown   Collective on MPI_Comm
624*552f7358SJed Brown 
625*552f7358SJed Brown   Input Parameter:
626*552f7358SJed Brown . dm - The original DMPlex object
627*552f7358SJed Brown 
628*552f7358SJed Brown   Output Parameter:
629*552f7358SJed Brown . newdm  - The new DMPlex object
630*552f7358SJed Brown 
631*552f7358SJed Brown   Level: beginner
632*552f7358SJed Brown 
633*552f7358SJed Brown .keywords: DMPlex, create
634*552f7358SJed Brown @*/
635*552f7358SJed Brown PetscErrorCode DMPlexClone(DM dm, DM *newdm)
636*552f7358SJed Brown {
637*552f7358SJed Brown   DM_Plex    *mesh;
638*552f7358SJed Brown   void          *ctx;
639*552f7358SJed Brown   PetscErrorCode ierr;
640*552f7358SJed Brown 
641*552f7358SJed Brown   PetscFunctionBegin;
642*552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
643*552f7358SJed Brown   PetscValidPointer(newdm,2);
644*552f7358SJed Brown   ierr = DMCreate(((PetscObject) dm)->comm, newdm);CHKERRQ(ierr);
645*552f7358SJed Brown   ierr = PetscSFDestroy(&(*newdm)->sf);CHKERRQ(ierr);
646*552f7358SJed Brown   ierr = PetscObjectReference((PetscObject) dm->sf);CHKERRQ(ierr);
647*552f7358SJed Brown   (*newdm)->sf = dm->sf;
648*552f7358SJed Brown   mesh = (DM_Plex *) dm->data;
649*552f7358SJed Brown   mesh->refct++;
650*552f7358SJed Brown   (*newdm)->data = mesh;
651*552f7358SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
652*552f7358SJed Brown   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
653*552f7358SJed Brown   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
654*552f7358SJed Brown   ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr);
655*552f7358SJed Brown   PetscFunctionReturn(0);
656*552f7358SJed Brown }
657