xref: /petsc/src/dm/impls/plex/plexcreate.c (revision 3dfda0b10bbf4cd8677f01069365ca416b9c6d67)
1552f7358SJed Brown #define PETSCDM_DLL
234541f0dSBarry Smith #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <petscdmda.h>
40c312b8eSJed Brown #include <petscsf.h>
5552f7358SJed Brown 
6552f7358SJed Brown #undef __FUNCT__
71df5d5c5SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateDoublet"
81df5d5c5SMatthew G. Knepley /*@
91df5d5c5SMatthew G. Knepley   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
101df5d5c5SMatthew G. Knepley 
111df5d5c5SMatthew G. Knepley   Collective on MPI_Comm
121df5d5c5SMatthew G. Knepley 
131df5d5c5SMatthew G. Knepley   Input Parameters:
141df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object
151df5d5c5SMatthew G. Knepley . dim - The spatial dimension
161df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells
171df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces)
181df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement
191df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
201df5d5c5SMatthew G. Knepley 
211df5d5c5SMatthew G. Knepley   Output Parameter:
221df5d5c5SMatthew G. Knepley . dm  - The DM object
231df5d5c5SMatthew G. Knepley 
241df5d5c5SMatthew G. Knepley   Level: beginner
251df5d5c5SMatthew G. Knepley 
261df5d5c5SMatthew G. Knepley .keywords: DM, create
271df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate()
281df5d5c5SMatthew G. Knepley @*/
291df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
301df5d5c5SMatthew G. Knepley {
311df5d5c5SMatthew G. Knepley   DM             dm;
321df5d5c5SMatthew G. Knepley   PetscInt       p;
331df5d5c5SMatthew G. Knepley   PetscMPIInt    rank;
341df5d5c5SMatthew G. Knepley   PetscErrorCode ierr;
351df5d5c5SMatthew G. Knepley 
361df5d5c5SMatthew G. Knepley   PetscFunctionBegin;
371df5d5c5SMatthew G. Knepley   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
381df5d5c5SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
39c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
401df5d5c5SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
41ce78fa2fSMatthew G. Knepley   switch (dim) {
42ce78fa2fSMatthew G. Knepley   case 2:
43ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
44ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
45ce78fa2fSMatthew G. Knepley     break;
46ce78fa2fSMatthew G. Knepley   case 3:
47ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
48ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
49ce78fa2fSMatthew G. Knepley     break;
50ce78fa2fSMatthew G. Knepley   default:
51ce78fa2fSMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
52ce78fa2fSMatthew G. Knepley   }
531df5d5c5SMatthew G. Knepley   if (rank) {
541df5d5c5SMatthew G. Knepley     PetscInt numPoints[2] = {0, 0};
551df5d5c5SMatthew G. Knepley     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
561df5d5c5SMatthew G. Knepley   } else {
571df5d5c5SMatthew G. Knepley     switch (dim) {
581df5d5c5SMatthew G. Knepley     case 2:
591df5d5c5SMatthew G. Knepley       if (simplex) {
601df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {4, 2};
611df5d5c5SMatthew G. Knepley         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
621df5d5c5SMatthew G. Knepley         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
631df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
641df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
651df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
661df5d5c5SMatthew G. Knepley 
671df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
681df5d5c5SMatthew G. Knepley         for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
691df5d5c5SMatthew G. Knepley       } else {
701df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {6, 2};
711df5d5c5SMatthew G. Knepley         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
721df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
731df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
741df5d5c5SMatthew G. Knepley         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};
751df5d5c5SMatthew G. Knepley 
761df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
771df5d5c5SMatthew G. Knepley       }
781df5d5c5SMatthew G. Knepley       break;
791df5d5c5SMatthew G. Knepley     case 3:
801df5d5c5SMatthew G. Knepley       if (simplex) {
811df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {5, 2};
821df5d5c5SMatthew G. Knepley         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
831df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
841df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
851df5d5c5SMatthew G. Knepley         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};
861df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
871df5d5c5SMatthew G. Knepley 
881df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
891df5d5c5SMatthew G. Knepley         for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
901df5d5c5SMatthew G. Knepley       } else {
911df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]         = {12, 2};
921df5d5c5SMatthew G. Knepley         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
931df5d5c5SMatthew G. Knepley         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
941df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
951df5d5c5SMatthew G. Knepley         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,
961df5d5c5SMatthew G. Knepley                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
971df5d5c5SMatthew G. Knepley                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
981df5d5c5SMatthew G. Knepley 
991df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
1001df5d5c5SMatthew G. Knepley       }
1011df5d5c5SMatthew G. Knepley       break;
1021df5d5c5SMatthew G. Knepley     default:
1031df5d5c5SMatthew G. Knepley       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
1041df5d5c5SMatthew G. Knepley     }
1051df5d5c5SMatthew G. Knepley   }
1061df5d5c5SMatthew G. Knepley   *newdm = dm;
1071df5d5c5SMatthew G. Knepley   if (refinementLimit > 0.0) {
1081df5d5c5SMatthew G. Knepley     DM rdm;
1091df5d5c5SMatthew G. Knepley     const char *name;
1101df5d5c5SMatthew G. Knepley 
1111df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
1121df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
1131df5d5c5SMatthew G. Knepley     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
1141df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1151df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
1161df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1171df5d5c5SMatthew G. Knepley     *newdm = rdm;
1181df5d5c5SMatthew G. Knepley   }
1191df5d5c5SMatthew G. Knepley   if (interpolate) {
120e56d480eSMatthew G. Knepley     DM idm = NULL;
1211df5d5c5SMatthew G. Knepley     const char *name;
1221df5d5c5SMatthew G. Knepley 
1231df5d5c5SMatthew G. Knepley     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
1241df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1251df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
1261df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
1271df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr);
1281df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1291df5d5c5SMatthew G. Knepley     *newdm = idm;
1301df5d5c5SMatthew G. Knepley   }
1311df5d5c5SMatthew G. Knepley   {
1321df5d5c5SMatthew G. Knepley     DM refinedMesh     = NULL;
1331df5d5c5SMatthew G. Knepley     DM distributedMesh = NULL;
1341df5d5c5SMatthew G. Knepley 
1351df5d5c5SMatthew G. Knepley     /* Distribute mesh over processes */
13607104863SJed Brown     ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr);
1371df5d5c5SMatthew G. Knepley     if (distributedMesh) {
1381df5d5c5SMatthew G. Knepley       ierr = DMDestroy(newdm);CHKERRQ(ierr);
1391df5d5c5SMatthew G. Knepley       *newdm = distributedMesh;
1401df5d5c5SMatthew G. Knepley     }
1411df5d5c5SMatthew G. Knepley     if (refinementUniform) {
1421df5d5c5SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
1431df5d5c5SMatthew G. Knepley       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
1441df5d5c5SMatthew G. Knepley       if (refinedMesh) {
1451df5d5c5SMatthew G. Knepley         ierr = DMDestroy(newdm);CHKERRQ(ierr);
1461df5d5c5SMatthew G. Knepley         *newdm = refinedMesh;
1471df5d5c5SMatthew G. Knepley       }
1481df5d5c5SMatthew G. Knepley     }
1491df5d5c5SMatthew G. Knepley   }
1501df5d5c5SMatthew G. Knepley   PetscFunctionReturn(0);
1511df5d5c5SMatthew G. Knepley }
1521df5d5c5SMatthew G. Knepley 
1531df5d5c5SMatthew G. Knepley #undef __FUNCT__
154552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary"
15526492d91SMatthew G. Knepley /*@
15626492d91SMatthew G. Knepley   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
157552f7358SJed Brown 
15826492d91SMatthew G. Knepley   Collective on MPI_Comm
15926492d91SMatthew G. Knepley 
16026492d91SMatthew G. Knepley   Input Parameters:
16126492d91SMatthew G. Knepley + comm  - The communicator for the DM object
16226492d91SMatthew G. Knepley . lower - The lower left corner coordinates
16326492d91SMatthew G. Knepley . upper - The upper right corner coordinates
16426492d91SMatthew G. Knepley - edges - The number of cells in each direction
16526492d91SMatthew G. Knepley 
16626492d91SMatthew G. Knepley   Output Parameter:
16726492d91SMatthew G. Knepley . dm  - The DM object
16826492d91SMatthew G. Knepley 
16926492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
17026492d91SMatthew G. Knepley $ 18--5-17--4--16
17126492d91SMatthew G. Knepley $  |     |     |
17226492d91SMatthew G. Knepley $  6    10     3
17326492d91SMatthew G. Knepley $  |     |     |
17426492d91SMatthew G. Knepley $ 19-11-20--9--15
17526492d91SMatthew G. Knepley $  |     |     |
17626492d91SMatthew G. Knepley $  7     8     2
17726492d91SMatthew G. Knepley $  |     |     |
17826492d91SMatthew G. Knepley $ 12--0-13--1--14
17926492d91SMatthew G. Knepley 
18026492d91SMatthew G. Knepley   Level: beginner
18126492d91SMatthew G. Knepley 
18226492d91SMatthew G. Knepley .keywords: DM, create
18326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
18426492d91SMatthew G. Knepley @*/
185552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
186552f7358SJed Brown {
187552f7358SJed Brown   PetscInt       numVertices    = (edges[0]+1)*(edges[1]+1);
188552f7358SJed Brown   PetscInt       numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
189552f7358SJed Brown   PetscInt       markerTop      = 1;
190552f7358SJed Brown   PetscInt       markerBottom   = 1;
191552f7358SJed Brown   PetscInt       markerRight    = 1;
192552f7358SJed Brown   PetscInt       markerLeft     = 1;
193552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
194552f7358SJed Brown   Vec            coordinates;
195552f7358SJed Brown   PetscSection   coordSection;
196552f7358SJed Brown   PetscScalar    *coords;
197552f7358SJed Brown   PetscInt       coordSize;
198552f7358SJed Brown   PetscMPIInt    rank;
199552f7358SJed Brown   PetscInt       v, vx, vy;
200552f7358SJed Brown   PetscErrorCode ierr;
201552f7358SJed Brown 
202552f7358SJed Brown   PetscFunctionBegin;
2030298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
204552f7358SJed Brown   if (markerSeparate) {
205552f7358SJed Brown     markerTop    = 1;
206552f7358SJed Brown     markerBottom = 0;
207552f7358SJed Brown     markerRight  = 0;
208552f7358SJed Brown     markerLeft   = 0;
209552f7358SJed Brown   }
21082f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
211552f7358SJed Brown   if (!rank) {
212552f7358SJed Brown     PetscInt e, ex, ey;
213552f7358SJed Brown 
214552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
215552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
216552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
217552f7358SJed Brown     }
218552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
219552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
220552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
221552f7358SJed Brown         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
222552f7358SJed Brown         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
223da80777bSKarl Rupp         PetscInt cone[2];
224552f7358SJed Brown 
225da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
226552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
227552f7358SJed Brown         if (vx == edges[0]) {
228552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
229552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
230552f7358SJed Brown           if (ey == edges[1]-1) {
231552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
232552f7358SJed Brown           }
233552f7358SJed Brown         } else if (vx == 0) {
234552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
235552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
236552f7358SJed Brown           if (ey == edges[1]-1) {
237552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
238552f7358SJed Brown           }
239552f7358SJed Brown         }
240552f7358SJed Brown       }
241552f7358SJed Brown     }
242552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
243552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
244552f7358SJed Brown         PetscInt edge   = vy*edges[0]     + ex;
245552f7358SJed Brown         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
246da80777bSKarl Rupp         PetscInt cone[2];
247552f7358SJed Brown 
248da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
249552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
250552f7358SJed Brown         if (vy == edges[1]) {
251552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
252552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
253552f7358SJed Brown           if (ex == edges[0]-1) {
254552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
255552f7358SJed Brown           }
256552f7358SJed Brown         } else if (vy == 0) {
257552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
258552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
259552f7358SJed Brown           if (ex == edges[0]-1) {
260552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
261552f7358SJed Brown           }
262552f7358SJed Brown         }
263552f7358SJed Brown       }
264552f7358SJed Brown     }
265552f7358SJed Brown   }
266552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
267552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
268552f7358SJed Brown   /* Build coordinates */
269c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
270552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
271552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
272552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
273552f7358SJed Brown   }
274552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
275552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
27682f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
2774b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
278552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2792eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
280552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
281552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
282552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
283552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
284552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
285552f7358SJed Brown     }
286552f7358SJed Brown   }
287552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
288552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
289552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
290552f7358SJed Brown   PetscFunctionReturn(0);
291552f7358SJed Brown }
292552f7358SJed Brown 
293552f7358SJed Brown #undef __FUNCT__
294552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
29526492d91SMatthew G. Knepley /*@
29626492d91SMatthew G. Knepley   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
297552f7358SJed Brown 
29826492d91SMatthew G. Knepley   Collective on MPI_Comm
29926492d91SMatthew G. Knepley 
30026492d91SMatthew G. Knepley   Input Parameters:
30126492d91SMatthew G. Knepley + comm  - The communicator for the DM object
30226492d91SMatthew G. Knepley . lower - The lower left front corner coordinates
30326492d91SMatthew G. Knepley . upper - The upper right back corner coordinates
30426492d91SMatthew G. Knepley - edges - The number of cells in each direction
30526492d91SMatthew G. Knepley 
30626492d91SMatthew G. Knepley   Output Parameter:
30726492d91SMatthew G. Knepley . dm  - The DM object
30826492d91SMatthew G. Knepley 
30926492d91SMatthew G. Knepley   Level: beginner
31026492d91SMatthew G. Knepley 
31126492d91SMatthew G. Knepley .keywords: DM, create
31226492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
31326492d91SMatthew G. Knepley @*/
314552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
315552f7358SJed Brown {
3169e8abbc3SMichael Lange   PetscInt       vertices[3], numVertices;
3177b59f5a9SMichael Lange   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
318552f7358SJed Brown   Vec            coordinates;
319552f7358SJed Brown   PetscSection   coordSection;
320552f7358SJed Brown   PetscScalar    *coords;
321552f7358SJed Brown   PetscInt       coordSize;
322552f7358SJed Brown   PetscMPIInt    rank;
323552f7358SJed Brown   PetscInt       v, vx, vy, vz;
3247b59f5a9SMichael Lange   PetscInt       voffset, iface=0, cone[4];
325552f7358SJed Brown   PetscErrorCode ierr;
326552f7358SJed Brown 
327552f7358SJed Brown   PetscFunctionBegin;
32882f516ccSBarry Smith   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");
32982f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
3309e8abbc3SMichael Lange   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
3319e8abbc3SMichael Lange   numVertices = vertices[0]*vertices[1]*vertices[2];
332552f7358SJed Brown   if (!rank) {
333552f7358SJed Brown     PetscInt f;
334552f7358SJed Brown 
335552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
336552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
337552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
338552f7358SJed Brown     }
339552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
340552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
341552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
342552f7358SJed Brown     }
3437b59f5a9SMichael Lange 
3447b59f5a9SMichael Lange     /* Side 0 (Top) */
3457b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3467b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3477b59f5a9SMichael Lange         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
3487b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
3497b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3507b59f5a9SMichael Lange         iface++;
351552f7358SJed Brown       }
352552f7358SJed Brown     }
3537b59f5a9SMichael Lange 
3547b59f5a9SMichael Lange     /* Side 1 (Bottom) */
3557b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3567b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3577b59f5a9SMichael Lange         voffset = numFaces + vy*(faces[0]+1) + vx;
3587b59f5a9SMichael Lange         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
3597b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3607b59f5a9SMichael Lange         iface++;
361552f7358SJed Brown       }
362552f7358SJed Brown     }
3637b59f5a9SMichael Lange 
3647b59f5a9SMichael Lange     /* Side 2 (Front) */
3657b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3667b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3677b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
3687b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
3697b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3707b59f5a9SMichael Lange         iface++;
371552f7358SJed Brown       }
3727b59f5a9SMichael Lange     }
3737b59f5a9SMichael Lange 
3747b59f5a9SMichael Lange     /* Side 3 (Back) */
3757b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3767b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3777b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
3787b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
3797b59f5a9SMichael Lange         cone[2] = voffset+1; cone[3] = voffset;
3807b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3817b59f5a9SMichael Lange         iface++;
3827b59f5a9SMichael Lange       }
3837b59f5a9SMichael Lange     }
3847b59f5a9SMichael Lange 
3857b59f5a9SMichael Lange     /* Side 4 (Left) */
3867b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3877b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3887b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
3897b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
3907b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
3917b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3927b59f5a9SMichael Lange         iface++;
3937b59f5a9SMichael Lange       }
3947b59f5a9SMichael Lange     }
3957b59f5a9SMichael Lange 
3967b59f5a9SMichael Lange     /* Side 5 (Right) */
3977b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3987b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3997b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx;
4007b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
4017b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
4027b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
4037b59f5a9SMichael Lange         iface++;
4047b59f5a9SMichael Lange       }
405552f7358SJed Brown     }
406552f7358SJed Brown   }
407552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
408552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
409552f7358SJed Brown   /* Build coordinates */
410c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
411552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
412552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
413552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
414552f7358SJed Brown   }
415552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
416552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
41782f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
4184b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr);
419552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
420552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4212eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
422552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
423552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
424552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
425552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
426552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
427552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
428552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
429552f7358SJed Brown       }
430552f7358SJed Brown     }
431552f7358SJed Brown   }
432552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
433552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
434552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
435552f7358SJed Brown   PetscFunctionReturn(0);
436552f7358SJed Brown }
437552f7358SJed Brown 
438552f7358SJed Brown #undef __FUNCT__
439*3dfda0b1SToby Isaac #define __FUNCT__ "DMPlexCreateCubeMesh_Internal"
440*3dfda0b1SToby Isaac static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
441*3dfda0b1SToby Isaac {
442*3dfda0b1SToby Isaac   PetscInt       markerTop      = 1;
443*3dfda0b1SToby Isaac   PetscInt       markerBottom   = 1;
444*3dfda0b1SToby Isaac   PetscInt       markerFront    = 1;
445*3dfda0b1SToby Isaac   PetscInt       markerBack     = 1;
446*3dfda0b1SToby Isaac   PetscInt       markerRight    = 1;
447*3dfda0b1SToby Isaac   PetscInt       markerLeft     = 1;
448*3dfda0b1SToby Isaac   PetscInt       dim;
449*3dfda0b1SToby Isaac   PetscBool      markerSeparate = PETSC_FALSE;
450*3dfda0b1SToby Isaac   PetscMPIInt    rank;
451*3dfda0b1SToby Isaac   PetscErrorCode ierr;
452*3dfda0b1SToby Isaac 
453*3dfda0b1SToby Isaac   PetscFunctionBegin;
454*3dfda0b1SToby Isaac   ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr);
455*3dfda0b1SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
456*3dfda0b1SToby Isaac   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
457*3dfda0b1SToby Isaac   if (markerSeparate) {
458*3dfda0b1SToby Isaac     switch (dim) {
459*3dfda0b1SToby Isaac     case 2:
460*3dfda0b1SToby Isaac       markerTop    = 3;
461*3dfda0b1SToby Isaac       markerBottom = 1;
462*3dfda0b1SToby Isaac       markerRight  = 2;
463*3dfda0b1SToby Isaac       markerLeft   = 4;
464*3dfda0b1SToby Isaac       break;
465*3dfda0b1SToby Isaac     case 3:
466*3dfda0b1SToby Isaac       markerBottom = 1;
467*3dfda0b1SToby Isaac       markerTop    = 2;
468*3dfda0b1SToby Isaac       markerFront  = 3;
469*3dfda0b1SToby Isaac       markerBack   = 4;
470*3dfda0b1SToby Isaac       markerRight  = 5;
471*3dfda0b1SToby Isaac       markerLeft   = 6;
472*3dfda0b1SToby Isaac       break;
473*3dfda0b1SToby Isaac     default:
474*3dfda0b1SToby Isaac       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim);
475*3dfda0b1SToby Isaac       break;
476*3dfda0b1SToby Isaac     }
477*3dfda0b1SToby Isaac   }
478*3dfda0b1SToby Isaac   {
479*3dfda0b1SToby Isaac     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
480*3dfda0b1SToby Isaac     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
481*3dfda0b1SToby Isaac     const PetscInt numZEdges    = !rank ? edges[2]   : 0;
482*3dfda0b1SToby Isaac     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
483*3dfda0b1SToby Isaac     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
484*3dfda0b1SToby Isaac     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
485*3dfda0b1SToby Isaac     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
486*3dfda0b1SToby Isaac     const PetscInt numXFaces    = numYEdges*numZEdges;
487*3dfda0b1SToby Isaac     const PetscInt numYFaces    = numXEdges*numZEdges;
488*3dfda0b1SToby Isaac     const PetscInt numZFaces    = numXEdges*numYEdges;
489*3dfda0b1SToby Isaac     const PetscInt numTotXFaces = numXVertices*numXFaces;
490*3dfda0b1SToby Isaac     const PetscInt numTotYFaces = numYVertices*numYFaces;
491*3dfda0b1SToby Isaac     const PetscInt numTotZFaces = numZVertices*numZFaces;
492*3dfda0b1SToby Isaac     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
493*3dfda0b1SToby Isaac     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
494*3dfda0b1SToby Isaac     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
495*3dfda0b1SToby Isaac     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
496*3dfda0b1SToby Isaac     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
497*3dfda0b1SToby Isaac     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
498*3dfda0b1SToby Isaac     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
499*3dfda0b1SToby Isaac     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
500*3dfda0b1SToby Isaac     const PetscInt firstYFace   = firstXFace + numTotXFaces;
501*3dfda0b1SToby Isaac     const PetscInt firstZFace   = firstYFace + numTotYFaces;
502*3dfda0b1SToby Isaac     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
503*3dfda0b1SToby Isaac     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
504*3dfda0b1SToby Isaac     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
505*3dfda0b1SToby Isaac     Vec            coordinates;
506*3dfda0b1SToby Isaac     PetscSection   coordSection;
507*3dfda0b1SToby Isaac     PetscScalar   *coords;
508*3dfda0b1SToby Isaac     PetscInt       coordSize;
509*3dfda0b1SToby Isaac     PetscInt       v, vx, vy, vz;
510*3dfda0b1SToby Isaac     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;
511*3dfda0b1SToby Isaac 
512*3dfda0b1SToby Isaac     ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr);
513*3dfda0b1SToby Isaac     for (c = 0; c < numCells; c++) {
514*3dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr);
515*3dfda0b1SToby Isaac     }
516*3dfda0b1SToby Isaac     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
517*3dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
518*3dfda0b1SToby Isaac     }
519*3dfda0b1SToby Isaac     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
520*3dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
521*3dfda0b1SToby Isaac     }
522*3dfda0b1SToby Isaac     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
523*3dfda0b1SToby Isaac     /* Build cells */
524*3dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
525*3dfda0b1SToby Isaac       for (fy = 0; fy < numYEdges; ++fy) {
526*3dfda0b1SToby Isaac         for (fx = 0; fx < numXEdges; ++fx) {
527*3dfda0b1SToby Isaac           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
528*3dfda0b1SToby Isaac           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
529*3dfda0b1SToby Isaac           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
530*3dfda0b1SToby Isaac           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
531*3dfda0b1SToby Isaac           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
532*3dfda0b1SToby Isaac           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
533*3dfda0b1SToby Isaac           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
534*3dfda0b1SToby Isaac                             /* B,  T,  F,  K,  R,  L */
535*3dfda0b1SToby Isaac           PetscInt ornt[8] = {-4,  0,  0, -1,  0, -4}; /* ??? */
536*3dfda0b1SToby Isaac           PetscInt cone[8];
537*3dfda0b1SToby Isaac 
538*3dfda0b1SToby Isaac           /* no boundary twisting in 3D */
539*3dfda0b1SToby Isaac           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
540*3dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr);
541*3dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr);
542*3dfda0b1SToby Isaac         }
543*3dfda0b1SToby Isaac       }
544*3dfda0b1SToby Isaac     }
545*3dfda0b1SToby Isaac     /* Build x faces */
546*3dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
547*3dfda0b1SToby Isaac       for (fy = 0; fy < numYEdges; ++fy) {
548*3dfda0b1SToby Isaac         for (fx = 0; fx < numXVertices; ++fx) {
549*3dfda0b1SToby Isaac           PetscInt face    = firstXFace + (fz*numYEdges+fy)*numXVertices + fx;
550*3dfda0b1SToby Isaac           PetscInt edgeL   = firstZEdge + (  fy*                 numXVertices+fx)*numZEdges + fz;
551*3dfda0b1SToby Isaac           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
552*3dfda0b1SToby Isaac           PetscInt edgeB   = firstYEdge + (  fz*                 numXVertices+fx)*numYEdges + fy;
553*3dfda0b1SToby Isaac           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
554*3dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
555*3dfda0b1SToby Isaac           PetscInt cone[4];
556*3dfda0b1SToby Isaac 
557*3dfda0b1SToby Isaac           if (dim == 3) {
558*3dfda0b1SToby Isaac             /* markers */
559*3dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
560*3dfda0b1SToby Isaac               if (fx == numXVertices-1) {
561*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr);
562*3dfda0b1SToby Isaac               }
563*3dfda0b1SToby Isaac               else if (fx == 0) {
564*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr);
565*3dfda0b1SToby Isaac               }
566*3dfda0b1SToby Isaac             }
567*3dfda0b1SToby Isaac           }
568*3dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
569*3dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
570*3dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
571*3dfda0b1SToby Isaac         }
572*3dfda0b1SToby Isaac       }
573*3dfda0b1SToby Isaac     }
574*3dfda0b1SToby Isaac     /* Build y faces */
575*3dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
576*3dfda0b1SToby Isaac       for (fx = 0; fx < numYEdges; ++fx) {
577*3dfda0b1SToby Isaac         for (fy = 0; fy < numYVertices; ++fy) {
578*3dfda0b1SToby Isaac           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
579*3dfda0b1SToby Isaac           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx                 )*numZEdges + fz;
580*3dfda0b1SToby Isaac           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
581*3dfda0b1SToby Isaac           PetscInt edgeB   = firstXEdge + (  fz                 *numYVertices+fy)*numXEdges + fx;
582*3dfda0b1SToby Isaac           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
583*3dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
584*3dfda0b1SToby Isaac           PetscInt cone[4];
585*3dfda0b1SToby Isaac 
586*3dfda0b1SToby Isaac           if (dim == 3) {
587*3dfda0b1SToby Isaac             /* markers */
588*3dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
589*3dfda0b1SToby Isaac               if (fy == numYVertices-1) {
590*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr);
591*3dfda0b1SToby Isaac               }
592*3dfda0b1SToby Isaac               else if (fy == 0) {
593*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr);
594*3dfda0b1SToby Isaac               }
595*3dfda0b1SToby Isaac             }
596*3dfda0b1SToby Isaac           }
597*3dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
598*3dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
599*3dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
600*3dfda0b1SToby Isaac         }
601*3dfda0b1SToby Isaac       }
602*3dfda0b1SToby Isaac     }
603*3dfda0b1SToby Isaac     /* Build z faces */
604*3dfda0b1SToby Isaac     for (fy = 0; fy < numYEdges; ++fy) {
605*3dfda0b1SToby Isaac       for (fx = 0; fx < numXEdges; ++fx) {
606*3dfda0b1SToby Isaac         for (fz = 0; fz < numZVertices; fz++) {
607*3dfda0b1SToby Isaac           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
608*3dfda0b1SToby Isaac           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx                 )*numYEdges + fy;
609*3dfda0b1SToby Isaac           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
610*3dfda0b1SToby Isaac           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy                 )*numXEdges + fx;
611*3dfda0b1SToby Isaac           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
612*3dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
613*3dfda0b1SToby Isaac           PetscInt cone[4];
614*3dfda0b1SToby Isaac 
615*3dfda0b1SToby Isaac           if (dim == 2) {
616*3dfda0b1SToby Isaac             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
617*3dfda0b1SToby Isaac             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
618*3dfda0b1SToby Isaac           }
619*3dfda0b1SToby Isaac           else {
620*3dfda0b1SToby Isaac             /* markers */
621*3dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
622*3dfda0b1SToby Isaac               if (fz == numZVertices-1) {
623*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr);
624*3dfda0b1SToby Isaac               }
625*3dfda0b1SToby Isaac               else if (fz == 0) {
626*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr);
627*3dfda0b1SToby Isaac               }
628*3dfda0b1SToby Isaac             }
629*3dfda0b1SToby Isaac           }
630*3dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
631*3dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
632*3dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
633*3dfda0b1SToby Isaac         }
634*3dfda0b1SToby Isaac       }
635*3dfda0b1SToby Isaac     }
636*3dfda0b1SToby Isaac     /* Build Z edges*/
637*3dfda0b1SToby Isaac     for (vy = 0; vy < numYVertices; vy++) {
638*3dfda0b1SToby Isaac       for (vx = 0; vx < numXVertices; vx++) {
639*3dfda0b1SToby Isaac         for (ez = 0; ez < numZEdges; ez++) {
640*3dfda0b1SToby Isaac           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
641*3dfda0b1SToby Isaac           const PetscInt vertexB = firstVertex + (  ez                 *numYVertices+vy)*numXVertices + vx;
642*3dfda0b1SToby Isaac           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
643*3dfda0b1SToby Isaac           PetscInt       cone[2];
644*3dfda0b1SToby Isaac 
645*3dfda0b1SToby Isaac           if (dim == 3) {
646*3dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
647*3dfda0b1SToby Isaac               if (vx == numXVertices-1) {
648*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
649*3dfda0b1SToby Isaac               }
650*3dfda0b1SToby Isaac               else if (vx == 0) {
651*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
652*3dfda0b1SToby Isaac               }
653*3dfda0b1SToby Isaac             }
654*3dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
655*3dfda0b1SToby Isaac               if (vy == numYVertices-1) {
656*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
657*3dfda0b1SToby Isaac               }
658*3dfda0b1SToby Isaac               else if (vy == 0) {
659*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
660*3dfda0b1SToby Isaac               }
661*3dfda0b1SToby Isaac             }
662*3dfda0b1SToby Isaac           }
663*3dfda0b1SToby Isaac           cone[0] = vertexB; cone[1] = vertexT;
664*3dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
665*3dfda0b1SToby Isaac         }
666*3dfda0b1SToby Isaac       }
667*3dfda0b1SToby Isaac     }
668*3dfda0b1SToby Isaac     /* Build Y edges*/
669*3dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; vz++) {
670*3dfda0b1SToby Isaac       for (vx = 0; vx < numXVertices; vx++) {
671*3dfda0b1SToby Isaac         for (ey = 0; ey < numYEdges; ey++) {
672*3dfda0b1SToby Isaac           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
673*3dfda0b1SToby Isaac           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
674*3dfda0b1SToby Isaac           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
675*3dfda0b1SToby Isaac           const PetscInt vertexK = firstVertex + nextv;
676*3dfda0b1SToby Isaac           PetscInt       cone[2];
677*3dfda0b1SToby Isaac 
678*3dfda0b1SToby Isaac           cone[0] = vertexF; cone[1] = vertexK;
679*3dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
680*3dfda0b1SToby Isaac           if (dim == 2) {
681*3dfda0b1SToby Isaac             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
682*3dfda0b1SToby Isaac               if (vx == numXVertices-1) {
683*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
684*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
685*3dfda0b1SToby Isaac                 if (ey == numYEdges-1) {
686*3dfda0b1SToby Isaac                   ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
687*3dfda0b1SToby Isaac                 }
688*3dfda0b1SToby Isaac               }
689*3dfda0b1SToby Isaac               else if (vx == 0) {
690*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
691*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
692*3dfda0b1SToby Isaac                 if (ey == numYEdges-1) {
693*3dfda0b1SToby Isaac                   ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
694*3dfda0b1SToby Isaac                 }
695*3dfda0b1SToby Isaac               }
696*3dfda0b1SToby Isaac             }
697*3dfda0b1SToby Isaac           }
698*3dfda0b1SToby Isaac           else {
699*3dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
700*3dfda0b1SToby Isaac               if (vx == numXVertices-1) {
701*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
702*3dfda0b1SToby Isaac               }
703*3dfda0b1SToby Isaac               else if (vx == 0) {
704*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
705*3dfda0b1SToby Isaac               }
706*3dfda0b1SToby Isaac             }
707*3dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
708*3dfda0b1SToby Isaac               if (vz == numZVertices-1) {
709*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
710*3dfda0b1SToby Isaac               }
711*3dfda0b1SToby Isaac               else if (vz == 0) {
712*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
713*3dfda0b1SToby Isaac               }
714*3dfda0b1SToby Isaac             }
715*3dfda0b1SToby Isaac           }
716*3dfda0b1SToby Isaac         }
717*3dfda0b1SToby Isaac       }
718*3dfda0b1SToby Isaac     }
719*3dfda0b1SToby Isaac     /* Build X edges*/
720*3dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; vz++) {
721*3dfda0b1SToby Isaac       for (vy = 0; vy < numYVertices; vy++) {
722*3dfda0b1SToby Isaac         for (ex = 0; ex < numXEdges; ex++) {
723*3dfda0b1SToby Isaac           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
724*3dfda0b1SToby Isaac           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
725*3dfda0b1SToby Isaac           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
726*3dfda0b1SToby Isaac           const PetscInt vertexR = firstVertex + nextv;
727*3dfda0b1SToby Isaac           PetscInt       cone[2];
728*3dfda0b1SToby Isaac 
729*3dfda0b1SToby Isaac           cone[0] = vertexL; cone[1] = vertexR;
730*3dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
731*3dfda0b1SToby Isaac           if (dim == 2) {
732*3dfda0b1SToby Isaac             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
733*3dfda0b1SToby Isaac               if (vy == numYVertices-1) {
734*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
735*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
736*3dfda0b1SToby Isaac                 if (ex == numXEdges-1) {
737*3dfda0b1SToby Isaac                   ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
738*3dfda0b1SToby Isaac                 }
739*3dfda0b1SToby Isaac               }
740*3dfda0b1SToby Isaac               else if (vy == 0) {
741*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
742*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
743*3dfda0b1SToby Isaac                 if (ex == numXEdges-1) {
744*3dfda0b1SToby Isaac                   ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
745*3dfda0b1SToby Isaac                 }
746*3dfda0b1SToby Isaac               }
747*3dfda0b1SToby Isaac             }
748*3dfda0b1SToby Isaac           }
749*3dfda0b1SToby Isaac           else {
750*3dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
751*3dfda0b1SToby Isaac               if (vy == numYVertices-1) {
752*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
753*3dfda0b1SToby Isaac               }
754*3dfda0b1SToby Isaac               else if (vy == 0) {
755*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
756*3dfda0b1SToby Isaac               }
757*3dfda0b1SToby Isaac             }
758*3dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
759*3dfda0b1SToby Isaac               if (vz == numZVertices-1) {
760*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
761*3dfda0b1SToby Isaac               }
762*3dfda0b1SToby Isaac               else if (vz == 0) {
763*3dfda0b1SToby Isaac                 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
764*3dfda0b1SToby Isaac               }
765*3dfda0b1SToby Isaac             }
766*3dfda0b1SToby Isaac           }
767*3dfda0b1SToby Isaac         }
768*3dfda0b1SToby Isaac       }
769*3dfda0b1SToby Isaac     }
770*3dfda0b1SToby Isaac     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
771*3dfda0b1SToby Isaac     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
772*3dfda0b1SToby Isaac     /* Build coordinates */
773*3dfda0b1SToby Isaac     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
774*3dfda0b1SToby Isaac     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
775*3dfda0b1SToby Isaac     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
776*3dfda0b1SToby Isaac     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
777*3dfda0b1SToby Isaac     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
778*3dfda0b1SToby Isaac       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
779*3dfda0b1SToby Isaac       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
780*3dfda0b1SToby Isaac     }
781*3dfda0b1SToby Isaac     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
782*3dfda0b1SToby Isaac     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
783*3dfda0b1SToby Isaac     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
784*3dfda0b1SToby Isaac     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
785*3dfda0b1SToby Isaac     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
786*3dfda0b1SToby Isaac     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
787*3dfda0b1SToby Isaac     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
788*3dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; ++vz) {
789*3dfda0b1SToby Isaac       for (vy = 0; vy < numYVertices; ++vy) {
790*3dfda0b1SToby Isaac         for (vx = 0; vx < numXVertices; ++vx) {
791*3dfda0b1SToby Isaac           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
792*3dfda0b1SToby Isaac           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
793*3dfda0b1SToby Isaac           if (dim == 3) {
794*3dfda0b1SToby Isaac             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
795*3dfda0b1SToby Isaac           }
796*3dfda0b1SToby Isaac         }
797*3dfda0b1SToby Isaac       }
798*3dfda0b1SToby Isaac     }
799*3dfda0b1SToby Isaac     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
800*3dfda0b1SToby Isaac     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
801*3dfda0b1SToby Isaac     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
802*3dfda0b1SToby Isaac   }
803*3dfda0b1SToby Isaac   PetscFunctionReturn(0);
804*3dfda0b1SToby Isaac }
805*3dfda0b1SToby Isaac 
806*3dfda0b1SToby Isaac #undef __FUNCT__
807552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
80826492d91SMatthew G. Knepley /*@
80926492d91SMatthew G. Knepley   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
810552f7358SJed Brown 
81126492d91SMatthew G. Knepley   Collective on MPI_Comm
81226492d91SMatthew G. Knepley 
81326492d91SMatthew G. Knepley   Input Parameters:
81426492d91SMatthew G. Knepley + comm  - The communicator for the DM object
81526492d91SMatthew G. Knepley . lower - The lower left corner coordinates
81626492d91SMatthew G. Knepley . upper - The upper right corner coordinates
817fbf5b483SMatthew G. Knepley . edges - The number of cells in each direction
818fbf5b483SMatthew G. Knepley . bdX   - The boundary type for the X direction
819fbf5b483SMatthew G. Knepley - bdY   - The boundary type for the Y direction
82026492d91SMatthew G. Knepley 
82126492d91SMatthew G. Knepley   Output Parameter:
82226492d91SMatthew G. Knepley . dm  - The DM object
82326492d91SMatthew G. Knepley 
82426492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
82526492d91SMatthew G. Knepley $ 22--8-23--9--24
82626492d91SMatthew G. Knepley $  |     |     |
82726492d91SMatthew G. Knepley $ 13  2 14  3  15
82826492d91SMatthew G. Knepley $  |     |     |
82926492d91SMatthew G. Knepley $ 19--6-20--7--21
83026492d91SMatthew G. Knepley $  |     |     |
83126492d91SMatthew G. Knepley $ 10  0 11  1 12
83226492d91SMatthew G. Knepley $  |     |     |
83326492d91SMatthew G. Knepley $ 16--4-17--5--18
83426492d91SMatthew G. Knepley 
83526492d91SMatthew G. Knepley   Level: beginner
83626492d91SMatthew G. Knepley 
83726492d91SMatthew G. Knepley .keywords: DM, create
83826492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
83926492d91SMatthew G. Knepley @*/
840bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
841552f7358SJed Brown {
842*3dfda0b1SToby Isaac   PetscReal      lower3[3], upper3[3];
843*3dfda0b1SToby Isaac   PetscInt       edges3[3];
844552f7358SJed Brown   PetscErrorCode ierr;
845552f7358SJed Brown 
846552f7358SJed Brown   PetscFunctionBegin;
847*3dfda0b1SToby Isaac   lower3[0] = lower[0]; lower3[1] = lower[1]; lower3[2] = 0.;
848*3dfda0b1SToby Isaac   upper3[0] = upper[0]; upper3[1] = upper[1]; upper3[2] = 0.;
849*3dfda0b1SToby Isaac   edges3[0] = edges[0]; edges3[1] = edges[1]; edges3[2] = 0;
850*3dfda0b1SToby Isaac   ierr = DMPlexCreateCubeMesh_Internal(dm, lower3, upper3, edges3, bdX, bdY, DM_BOUNDARY_NONE);CHKERRQ(ierr);
851552f7358SJed Brown   PetscFunctionReturn(0);
852552f7358SJed Brown }
853552f7358SJed Brown 
854552f7358SJed Brown #undef __FUNCT__
855552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
8566b75479cSMatthew G Knepley /*@
85726492d91SMatthew G. Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
8586b75479cSMatthew G Knepley 
8596b75479cSMatthew G Knepley   Collective on MPI_Comm
8606b75479cSMatthew G Knepley 
8616b75479cSMatthew G Knepley   Input Parameters:
8626b75479cSMatthew G Knepley + comm - The communicator for the DM object
8636b75479cSMatthew G Knepley . dim - The spatial dimension
8646b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
8656b75479cSMatthew G Knepley 
8666b75479cSMatthew G Knepley   Output Parameter:
8676b75479cSMatthew G Knepley . dm  - The DM object
8686b75479cSMatthew G Knepley 
8696b75479cSMatthew G Knepley   Level: beginner
8706b75479cSMatthew G Knepley 
8716b75479cSMatthew G Knepley .keywords: DM, create
87226492d91SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
8736b75479cSMatthew G Knepley @*/
8740adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
8750adebc6cSBarry Smith {
876552f7358SJed Brown   DM             boundary;
877552f7358SJed Brown   PetscErrorCode ierr;
878552f7358SJed Brown 
879552f7358SJed Brown   PetscFunctionBegin;
880552f7358SJed Brown   PetscValidPointer(dm, 4);
881552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
882552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
883552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
884c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr);
885552f7358SJed Brown   switch (dim) {
886552f7358SJed Brown   case 2:
887552f7358SJed Brown   {
888552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
889552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
890552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
891552f7358SJed Brown 
892552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
893552f7358SJed Brown     break;
894552f7358SJed Brown   }
895552f7358SJed Brown   case 3:
896552f7358SJed Brown   {
897552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
898552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
899552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
900552f7358SJed Brown 
901552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
902552f7358SJed Brown     break;
903552f7358SJed Brown   }
904552f7358SJed Brown   default:
905552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
906552f7358SJed Brown   }
9070298fd71SBarry Smith   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
908552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
909552f7358SJed Brown   PetscFunctionReturn(0);
910552f7358SJed Brown }
911552f7358SJed Brown 
912552f7358SJed Brown #undef __FUNCT__
913552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
91426492d91SMatthew G. Knepley /*@
91526492d91SMatthew G. Knepley   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
91626492d91SMatthew G. Knepley 
91726492d91SMatthew G. Knepley   Collective on MPI_Comm
91826492d91SMatthew G. Knepley 
91926492d91SMatthew G. Knepley   Input Parameters:
92026492d91SMatthew G. Knepley + comm  - The communicator for the DM object
92126492d91SMatthew G. Knepley . dim   - The spatial dimension
922fbf5b483SMatthew G. Knepley . periodicX - The boundary type for the X direction
923fbf5b483SMatthew G. Knepley . periodicY - The boundary type for the Y direction
924fbf5b483SMatthew G. Knepley . periodicZ - The boundary type for the Z direction
92526492d91SMatthew G. Knepley - cells - The number of cells in each direction
92626492d91SMatthew G. Knepley 
92726492d91SMatthew G. Knepley   Output Parameter:
92826492d91SMatthew G. Knepley . dm  - The DM object
92926492d91SMatthew G. Knepley 
93026492d91SMatthew G. Knepley   Level: beginner
93126492d91SMatthew G. Knepley 
93226492d91SMatthew G. Knepley .keywords: DM, create
93326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
93426492d91SMatthew G. Knepley @*/
935bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
936a6dfd86eSKarl Rupp {
937552f7358SJed Brown   PetscErrorCode ierr;
938552f7358SJed Brown 
939552f7358SJed Brown   PetscFunctionBegin;
940552f7358SJed Brown   PetscValidPointer(dm, 4);
941552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
942552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
943552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
944c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
945552f7358SJed Brown   switch (dim) {
946552f7358SJed Brown   case 2:
947552f7358SJed Brown   {
948552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
949552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
950552f7358SJed Brown 
951becd5721SMatthew G. Knepley     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
952552f7358SJed Brown     break;
953552f7358SJed Brown   }
954552f7358SJed Brown   case 3:
955552f7358SJed Brown   {
956552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
957552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
958552f7358SJed Brown 
959*3dfda0b1SToby Isaac     ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
960552f7358SJed Brown     break;
961552f7358SJed Brown   }
962552f7358SJed Brown   default:
963552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
964552f7358SJed Brown   }
965552f7358SJed Brown   PetscFunctionReturn(0);
966552f7358SJed Brown }
967552f7358SJed Brown 
968552f7358SJed Brown /* External function declarations here */
969552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
970bceba477SMatthew G. Knepley extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx);
971fd59a867SMatthew G. Knepley extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
972b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
973552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
974552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
9750a6ba040SMatthew G. Knepley extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
9760a6ba040SMatthew G. Knepley extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
9779a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
978552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
979552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
980552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
9812c40f234SMatthew G. Knepley extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
982552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
983552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
984552f7358SJed Brown 
985552f7358SJed Brown #undef __FUNCT__
9860a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexReplace_Static"
9870a6ba040SMatthew G. Knepley /* Replace dm with the contents of dmNew
9880a6ba040SMatthew G. Knepley    - Share the DM_Plex structure
9890a6ba040SMatthew G. Knepley    - Share the coordinates
990d7973521SMatthew G. Knepley    - Share the SF
9910a6ba040SMatthew G. Knepley */
9920a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
9930a6ba040SMatthew G. Knepley {
994d7973521SMatthew G. Knepley   PetscSF          sf;
99555fbe3e3SMatthew G. Knepley   DM               coordDM;
9960a6ba040SMatthew G. Knepley   Vec              coords;
99755fbe3e3SMatthew G. Knepley   const PetscReal *maxCell, *L;
9980a6ba040SMatthew G. Knepley   PetscErrorCode   ierr;
9990a6ba040SMatthew G. Knepley 
10000a6ba040SMatthew G. Knepley   PetscFunctionBegin;
1001d7973521SMatthew G. Knepley   ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr);
1002d7973521SMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
100355fbe3e3SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr);
10040a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
100555fbe3e3SMatthew G. Knepley   ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr);
10060a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
100755fbe3e3SMatthew G. Knepley   ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
100855fbe3e3SMatthew G. Knepley   if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L);CHKERRQ(ierr);}
10090a6ba040SMatthew G. Knepley   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
10100a6ba040SMatthew G. Knepley   dm->data = dmNew->data;
10110a6ba040SMatthew G. Knepley   ((DM_Plex *) dmNew->data)->refct++;
10120a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
10130a6ba040SMatthew G. Knepley }
10140a6ba040SMatthew G. Knepley 
10150a6ba040SMatthew G. Knepley #undef __FUNCT__
10160a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexSwap_Static"
10170a6ba040SMatthew G. Knepley /* Swap dm with the contents of dmNew
10180a6ba040SMatthew G. Knepley    - Swap the DM_Plex structure
10190a6ba040SMatthew G. Knepley    - Swap the coordinates
10200a6ba040SMatthew G. Knepley */
10210a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
10220a6ba040SMatthew G. Knepley {
10230a6ba040SMatthew G. Knepley   DM             coordDMA, coordDMB;
10240a6ba040SMatthew G. Knepley   Vec            coordsA,  coordsB;
10250a6ba040SMatthew G. Knepley   void          *tmp;
10260a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
10270a6ba040SMatthew G. Knepley 
10280a6ba040SMatthew G. Knepley   PetscFunctionBegin;
10290a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
10300a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
10310a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
10320a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
10330a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
10340a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
10350a6ba040SMatthew G. Knepley 
10360a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
10370a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
10380a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
10390a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
10400a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
10410a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
10420a6ba040SMatthew G. Knepley   tmp       = dmA->data;
10430a6ba040SMatthew G. Knepley   dmA->data = dmB->data;
10440a6ba040SMatthew G. Knepley   dmB->data = tmp;
10450a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
10460a6ba040SMatthew G. Knepley }
10470a6ba040SMatthew G. Knepley 
10480a6ba040SMatthew G. Knepley #undef __FUNCT__
104968d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex"
105068d4fef7SMatthew G. Knepley PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(DM dm)
10510a6ba040SMatthew G. Knepley {
10520a6ba040SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
10536c73c22cSMatthew G. Knepley   DMBoundary     b;
10540a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
10550a6ba040SMatthew G. Knepley 
10560a6ba040SMatthew G. Knepley   PetscFunctionBegin;
10576c73c22cSMatthew G. Knepley   /* Handle boundary conditions */
10586c73c22cSMatthew G. Knepley   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr);
10596c73c22cSMatthew G. Knepley   for (b = mesh->boundary; b; b = b->next) {
10606c73c22cSMatthew G. Knepley     char      optname[1024];
10616c73c22cSMatthew G. Knepley     PetscInt  ids[1024], len = 1024, i;
10626c73c22cSMatthew G. Knepley     PetscBool flg;
10636c73c22cSMatthew G. Knepley 
10646c73c22cSMatthew G. Knepley     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
10656c73c22cSMatthew G. Knepley     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
10666c73c22cSMatthew G. Knepley     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
10676c73c22cSMatthew G. Knepley     if (flg) {
10686c73c22cSMatthew G. Knepley       DMLabel label;
10696c73c22cSMatthew G. Knepley 
107063d5297fSMatthew G. Knepley       ierr = DMPlexGetLabel(dm, b->labelname, &label);CHKERRQ(ierr);
10716c73c22cSMatthew G. Knepley       for (i = 0; i < len; ++i) {
10726c73c22cSMatthew G. Knepley         PetscBool has;
10736c73c22cSMatthew G. Knepley 
10746c73c22cSMatthew G. Knepley         ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr);
10756c73c22cSMatthew G. Knepley         if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name);
10766c73c22cSMatthew G. Knepley       }
10776c73c22cSMatthew G. Knepley       b->numids = len;
10786c73c22cSMatthew G. Knepley       ierr = PetscFree(b->ids);CHKERRQ(ierr);
10796c73c22cSMatthew G. Knepley       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
10806c73c22cSMatthew G. Knepley       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
10816c73c22cSMatthew G. Knepley     }
10826c73c22cSMatthew G. Knepley   }
10836c73c22cSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
10840a6ba040SMatthew G. Knepley   /* Handle viewing */
10850a6ba040SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
10860a6ba040SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
10870a6ba040SMatthew G. Knepley   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr);
108868d4fef7SMatthew G. Knepley   PetscFunctionReturn(0);
108968d4fef7SMatthew G. Knepley }
109068d4fef7SMatthew G. Knepley 
109168d4fef7SMatthew G. Knepley #undef __FUNCT__
109268d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_Plex"
109368d4fef7SMatthew G. Knepley PetscErrorCode  DMSetFromOptions_Plex(DM dm)
109468d4fef7SMatthew G. Knepley {
109568d4fef7SMatthew G. Knepley   PetscInt       refine = 0, r;
109668d4fef7SMatthew G. Knepley   PetscBool      isHierarchy;
109768d4fef7SMatthew G. Knepley   PetscErrorCode ierr;
109868d4fef7SMatthew G. Knepley 
109968d4fef7SMatthew G. Knepley   PetscFunctionBegin;
110068d4fef7SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
110168d4fef7SMatthew G. Knepley   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
110268d4fef7SMatthew G. Knepley   /* Handle DMPlex refinement */
110368d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
110468d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
1105b6a0289aSMatthew G. Knepley   if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);}
110668d4fef7SMatthew G. Knepley   if (refine && isHierarchy) {
110768d4fef7SMatthew G. Knepley     DM *dms;
110868d4fef7SMatthew G. Knepley 
110968d4fef7SMatthew G. Knepley     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
111068d4fef7SMatthew G. Knepley     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
111168d4fef7SMatthew G. Knepley     /* Total hack since we do not pass in a pointer */
111268d4fef7SMatthew G. Knepley     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
111368d4fef7SMatthew G. Knepley     if (refine == 1) {
111468d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
111568d4fef7SMatthew G. Knepley     } else {
111668d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
111768d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
111868d4fef7SMatthew G. Knepley     }
111968d4fef7SMatthew G. Knepley     /* Free DMs */
112068d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
112168d4fef7SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
112268d4fef7SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
112368d4fef7SMatthew G. Knepley     }
112468d4fef7SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
112568d4fef7SMatthew G. Knepley   } else {
112668d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
112768d4fef7SMatthew G. Knepley       DM refinedMesh;
112868d4fef7SMatthew G. Knepley 
112968d4fef7SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
113068d4fef7SMatthew G. Knepley       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
113168d4fef7SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
113268d4fef7SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
113368d4fef7SMatthew G. Knepley       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
113468d4fef7SMatthew G. Knepley     }
113568d4fef7SMatthew G. Knepley   }
113668d4fef7SMatthew G. Knepley   ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
11370a6ba040SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
11380a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
11390a6ba040SMatthew G. Knepley }
11400a6ba040SMatthew G. Knepley 
11410a6ba040SMatthew G. Knepley #undef __FUNCT__
1142552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
1143552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
1144552f7358SJed Brown {
1145552f7358SJed Brown   PetscErrorCode ierr;
1146552f7358SJed Brown 
1147552f7358SJed Brown   PetscFunctionBegin;
1148552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
1149552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
1150552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr);
11512c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr);
1152552f7358SJed Brown   PetscFunctionReturn(0);
1153552f7358SJed Brown }
1154552f7358SJed Brown 
1155552f7358SJed Brown #undef __FUNCT__
1156552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
1157552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
1158552f7358SJed Brown {
1159552f7358SJed Brown   PetscErrorCode ierr;
1160552f7358SJed Brown 
1161552f7358SJed Brown   PetscFunctionBegin;
1162552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
1163552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
11642c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr);
1165552f7358SJed Brown   PetscFunctionReturn(0);
1166552f7358SJed Brown }
1167552f7358SJed Brown 
116838221697SMatthew G. Knepley #undef __FUNCT__
1169793f3fe5SMatthew G. Knepley #define __FUNCT__ "DMGetDimPoints_Plex"
1170793f3fe5SMatthew G. Knepley static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
1171793f3fe5SMatthew G. Knepley {
1172793f3fe5SMatthew G. Knepley   PetscInt       depth, d;
1173793f3fe5SMatthew G. Knepley   PetscErrorCode ierr;
1174793f3fe5SMatthew G. Knepley 
1175793f3fe5SMatthew G. Knepley   PetscFunctionBegin;
1176793f3fe5SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1177793f3fe5SMatthew G. Knepley   if (depth == 1) {
1178793f3fe5SMatthew G. Knepley     ierr = DMGetDimension(dm, &d);CHKERRQ(ierr);
1179793f3fe5SMatthew G. Knepley     if (dim == 0)      {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);}
1180793f3fe5SMatthew G. Knepley     else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);}
1181793f3fe5SMatthew G. Knepley     else               {*pStart = 0; *pEnd = 0;}
1182793f3fe5SMatthew G. Knepley   } else {
1183793f3fe5SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);
1184793f3fe5SMatthew G. Knepley   }
1185793f3fe5SMatthew G. Knepley   PetscFunctionReturn(0);
1186793f3fe5SMatthew G. Knepley }
1187793f3fe5SMatthew G. Knepley 
1188793f3fe5SMatthew G. Knepley #undef __FUNCT__
1189552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
1190552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
1191552f7358SJed Brown {
1192552f7358SJed Brown   PetscFunctionBegin;
1193552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
11942c40f234SMatthew G. Knepley   dm->ops->load                            = DMLoad_Plex;
1195552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
119638221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
1197552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
1198fd59a867SMatthew G. Knepley   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
1199552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
1200552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
1201184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
12020298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
1203552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
12040a6ba040SMatthew G. Knepley   dm->ops->getcoloring                     = NULL;
1205552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
1206bceba477SMatthew G. Knepley   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
1207bceba477SMatthew G. Knepley   dm->ops->getaggregates                   = NULL;
1208bceba477SMatthew G. Knepley   dm->ops->getinjection                    = DMCreateInjection_Plex;
1209552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
12100a6ba040SMatthew G. Knepley   dm->ops->coarsen                         = DMCoarsen_Plex;
12110a6ba040SMatthew G. Knepley   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
12120a6ba040SMatthew G. Knepley   dm->ops->coarsenhierarchy                = NULL;
12130298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
12140298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
12150298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
12160298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
1217552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
1218552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
1219793f3fe5SMatthew G. Knepley   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
1220552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
1221552f7358SJed Brown   PetscFunctionReturn(0);
1222552f7358SJed Brown }
1223552f7358SJed Brown 
122463a16f15SMatthew G. Knepley #undef __FUNCT__
122563a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex"
122663a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
122763a16f15SMatthew G. Knepley {
122863a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
122963a16f15SMatthew G. Knepley   PetscErrorCode ierr;
123063a16f15SMatthew G. Knepley 
123163a16f15SMatthew G. Knepley   PetscFunctionBegin;
123263a16f15SMatthew G. Knepley   mesh->refct++;
123363a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
123463a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
123563a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
123663a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
123763a16f15SMatthew G. Knepley }
123863a16f15SMatthew G. Knepley 
12398818961aSMatthew G Knepley /*MC
12408818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
12418818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
12428818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
12438818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
12448818961aSMatthew G Knepley 
12458818961aSMatthew G Knepley   Level: intermediate
12468818961aSMatthew G Knepley 
12478818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
12488818961aSMatthew G Knepley M*/
12498818961aSMatthew G Knepley 
1250552f7358SJed Brown #undef __FUNCT__
1251552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
12528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1253552f7358SJed Brown {
1254552f7358SJed Brown   DM_Plex        *mesh;
1255770b213bSMatthew G Knepley   PetscInt       unit, d;
1256552f7358SJed Brown   PetscErrorCode ierr;
1257552f7358SJed Brown 
1258552f7358SJed Brown   PetscFunctionBegin;
1259552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1260b00a9115SJed Brown   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
1261c73cfb54SMatthew G. Knepley   dm->dim  = 0;
1262552f7358SJed Brown   dm->data = mesh;
1263552f7358SJed Brown 
1264552f7358SJed Brown   mesh->refct             = 1;
126582f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
1266552f7358SJed Brown   mesh->maxConeSize       = 0;
12670298fd71SBarry Smith   mesh->cones             = NULL;
12680298fd71SBarry Smith   mesh->coneOrientations  = NULL;
126982f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
1270552f7358SJed Brown   mesh->maxSupportSize    = 0;
12710298fd71SBarry Smith   mesh->supports          = NULL;
1272552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
1273552f7358SJed Brown   mesh->refinementLimit   = -1.0;
1274552f7358SJed Brown 
12750298fd71SBarry Smith   mesh->facesTmp = NULL;
1276552f7358SJed Brown 
12770298fd71SBarry Smith   mesh->subpointMap = NULL;
1278552f7358SJed Brown 
12798865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
1280552f7358SJed Brown 
12810298fd71SBarry Smith   mesh->labels              = NULL;
12828af19771SMatthew G. Knepley   mesh->depthLabel          = NULL;
12830298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
12840298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
12858865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1286552f7358SJed Brown   mesh->vtkCellHeight       = 0;
128770034214SMatthew G. Knepley   mesh->useCone             = PETSC_FALSE;
128870034214SMatthew G. Knepley   mesh->useClosure          = PETSC_TRUE;
1289552f7358SJed Brown 
1290552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
1291552f7358SJed Brown   mesh->printFEM       = 0;
12926113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
1293552f7358SJed Brown 
1294552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
1295552f7358SJed Brown   PetscFunctionReturn(0);
1296552f7358SJed Brown }
1297552f7358SJed Brown 
1298552f7358SJed Brown #undef __FUNCT__
1299552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
1300552f7358SJed Brown /*@
1301552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1302552f7358SJed Brown 
1303552f7358SJed Brown   Collective on MPI_Comm
1304552f7358SJed Brown 
1305552f7358SJed Brown   Input Parameter:
1306552f7358SJed Brown . comm - The communicator for the DMPlex object
1307552f7358SJed Brown 
1308552f7358SJed Brown   Output Parameter:
1309552f7358SJed Brown . mesh  - The DMPlex object
1310552f7358SJed Brown 
1311552f7358SJed Brown   Level: beginner
1312552f7358SJed Brown 
1313552f7358SJed Brown .keywords: DMPlex, create
1314552f7358SJed Brown @*/
1315552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1316552f7358SJed Brown {
1317552f7358SJed Brown   PetscErrorCode ierr;
1318552f7358SJed Brown 
1319552f7358SJed Brown   PetscFunctionBegin;
1320552f7358SJed Brown   PetscValidPointer(mesh,2);
1321552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1322552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1323552f7358SJed Brown   PetscFunctionReturn(0);
1324552f7358SJed Brown }
1325552f7358SJed Brown 
1326552f7358SJed Brown #undef __FUNCT__
13279298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private"
13289298eaa6SMatthew G Knepley /*
13299298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
13309298eaa6SMatthew G Knepley */
13319298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
13329298eaa6SMatthew G Knepley {
13339298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
13349298eaa6SMatthew G Knepley   PetscErrorCode ierr;
13359298eaa6SMatthew G Knepley 
13369298eaa6SMatthew G Knepley   PetscFunctionBegin;
13379298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
13389298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
13399298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
13409298eaa6SMatthew G Knepley   }
13419298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
13429298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
13439298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
13449298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
13459298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
13469298eaa6SMatthew G Knepley     }
13479298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
13489298eaa6SMatthew G Knepley   }
13499298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
13509298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
13519298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
13529298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
13539298eaa6SMatthew G Knepley }
13549298eaa6SMatthew G Knepley 
13559298eaa6SMatthew G Knepley #undef __FUNCT__
13569298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private"
13579298eaa6SMatthew G Knepley /*
13589298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
13599298eaa6SMatthew G Knepley */
13609298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
13619298eaa6SMatthew G Knepley {
13629298eaa6SMatthew G Knepley   PetscSection   coordSection;
13639298eaa6SMatthew G Knepley   Vec            coordinates;
13649298eaa6SMatthew G Knepley   PetscScalar   *coords;
13659298eaa6SMatthew G Knepley   PetscInt       coordSize, v, d;
13669298eaa6SMatthew G Knepley   PetscErrorCode ierr;
13679298eaa6SMatthew G Knepley 
13689298eaa6SMatthew G Knepley   PetscFunctionBegin;
1369c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
13709298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
13719298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
13729298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
13739298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
13749298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
13759298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
13769298eaa6SMatthew G Knepley   }
13779298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
13789298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
13799298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
13804b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
13819298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
13829298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
13832eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
13849298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
13859298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
13869298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
13879298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
13889298eaa6SMatthew G Knepley     }
13899298eaa6SMatthew G Knepley   }
13909298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
13919298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
13929298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
13939298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
13949298eaa6SMatthew G Knepley }
13959298eaa6SMatthew G Knepley 
13969298eaa6SMatthew G Knepley #undef __FUNCT__
13979298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList"
13989298eaa6SMatthew G Knepley /*@C
13999298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
14009298eaa6SMatthew G Knepley 
14019298eaa6SMatthew G Knepley   Input Parameters:
14029298eaa6SMatthew G Knepley + comm - The communicator
14039298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
14049298eaa6SMatthew G Knepley . numCells - The number of cells
14059298eaa6SMatthew G Knepley . numVertices - The number of vertices
14069298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
14079298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
14089298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
14099298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
14109298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
14119298eaa6SMatthew G Knepley 
14129298eaa6SMatthew G Knepley   Output Parameter:
14139298eaa6SMatthew G Knepley . dm - The DM
14149298eaa6SMatthew G Knepley 
14159298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
14169298eaa6SMatthew G Knepley $
14179298eaa6SMatthew G Knepley $        2
14189298eaa6SMatthew G Knepley $      / | \
14199298eaa6SMatthew G Knepley $     /  |  \
14209298eaa6SMatthew G Knepley $    /   |   \
14219298eaa6SMatthew G Knepley $   0  0 | 1  3
14229298eaa6SMatthew G Knepley $    \   |   /
14239298eaa6SMatthew G Knepley $     \  |  /
14249298eaa6SMatthew G Knepley $      \ | /
14259298eaa6SMatthew G Knepley $        1
14269298eaa6SMatthew G Knepley would have input
14279298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
14289298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
14299298eaa6SMatthew G Knepley $
14309298eaa6SMatthew G Knepley which would result in the DMPlex
14319298eaa6SMatthew G Knepley $
14329298eaa6SMatthew G Knepley $        4
14339298eaa6SMatthew G Knepley $      / | \
14349298eaa6SMatthew G Knepley $     /  |  \
14359298eaa6SMatthew G Knepley $    /   |   \
14369298eaa6SMatthew G Knepley $   2  0 | 1  5
14379298eaa6SMatthew G Knepley $    \   |   /
14389298eaa6SMatthew G Knepley $     \  |  /
14399298eaa6SMatthew G Knepley $      \ | /
14409298eaa6SMatthew G Knepley $        3
14419298eaa6SMatthew G Knepley 
14429298eaa6SMatthew G Knepley   Level: beginner
14439298eaa6SMatthew G Knepley 
1444939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
14459298eaa6SMatthew G Knepley @*/
14469298eaa6SMatthew G Knepley 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)
14479298eaa6SMatthew G Knepley {
14489298eaa6SMatthew G Knepley   PetscErrorCode ierr;
14499298eaa6SMatthew G Knepley 
14509298eaa6SMatthew G Knepley   PetscFunctionBegin;
14519298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
14529298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1453c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
14549298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
14559298eaa6SMatthew G Knepley   if (interpolate) {
1456e56d480eSMatthew G. Knepley     DM idm = NULL;
14579298eaa6SMatthew G Knepley 
14589298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
14599298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
14609298eaa6SMatthew G Knepley     *dm  = idm;
14619298eaa6SMatthew G Knepley   }
14629298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
14639298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
14649298eaa6SMatthew G Knepley }
14659298eaa6SMatthew G Knepley 
14669298eaa6SMatthew G Knepley #undef __FUNCT__
14679298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG"
1468939f6067SMatthew G. Knepley /*@
1469939f6067SMatthew G. Knepley   DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM
1470939f6067SMatthew G. Knepley 
1471939f6067SMatthew G. Knepley   Input Parameters:
1472c73cfb54SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
1473939f6067SMatthew G. Knepley . depth - The depth of the DAG
1474939f6067SMatthew G. Knepley . numPoints - The number of points at each depth
1475939f6067SMatthew G. Knepley . coneSize - The cone size of each point
1476939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1477939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point
1478939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1479939f6067SMatthew G. Knepley 
1480939f6067SMatthew G. Knepley   Output Parameter:
1481939f6067SMatthew G. Knepley . dm - The DM
1482939f6067SMatthew G. Knepley 
1483939f6067SMatthew G. Knepley   Note: Two triangles sharing a face would have input
1484939f6067SMatthew G. Knepley $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1485939f6067SMatthew G. Knepley $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1486939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1487939f6067SMatthew G. Knepley $
1488939f6067SMatthew G. Knepley which would result in the DMPlex
1489939f6067SMatthew G. Knepley $
1490939f6067SMatthew G. Knepley $        4
1491939f6067SMatthew G. Knepley $      / | \
1492939f6067SMatthew G. Knepley $     /  |  \
1493939f6067SMatthew G. Knepley $    /   |   \
1494939f6067SMatthew G. Knepley $   2  0 | 1  5
1495939f6067SMatthew G. Knepley $    \   |   /
1496939f6067SMatthew G. Knepley $     \  |  /
1497939f6067SMatthew G. Knepley $      \ | /
1498939f6067SMatthew G. Knepley $        3
1499939f6067SMatthew G. Knepley $
1500939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1501939f6067SMatthew G. Knepley 
1502939f6067SMatthew G. Knepley   Level: advanced
1503939f6067SMatthew G. Knepley 
1504939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1505939f6067SMatthew G. Knepley @*/
15069298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
15079298eaa6SMatthew G Knepley {
15089298eaa6SMatthew G Knepley   Vec            coordinates;
15099298eaa6SMatthew G Knepley   PetscSection   coordSection;
15109298eaa6SMatthew G Knepley   PetscScalar    *coords;
15119298eaa6SMatthew G Knepley   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
15129298eaa6SMatthew G Knepley   PetscErrorCode ierr;
15139298eaa6SMatthew G Knepley 
15149298eaa6SMatthew G Knepley   PetscFunctionBegin;
1515c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
15169298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
15179298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
15189298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
15199298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
15209298eaa6SMatthew G Knepley   }
15219298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
15229298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
15239298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
15249298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
15259298eaa6SMatthew G Knepley   }
15269298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
15279298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
15289298eaa6SMatthew G Knepley   /* Build coordinates */
1529c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
15309298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
15319298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
15329298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
15339298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
15349298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
15356f8cbbeeSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
15369298eaa6SMatthew G Knepley   }
15379298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
15389298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
15399298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
15404b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
15416f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
15429298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
15432eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
15449298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
15459298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
15469298eaa6SMatthew G Knepley     PetscInt off;
15479298eaa6SMatthew G Knepley 
15489298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
15499298eaa6SMatthew G Knepley     for (d = 0; d < dim; ++d) {
15509298eaa6SMatthew G Knepley       coords[off+d] = vertexCoords[v*dim+d];
15519298eaa6SMatthew G Knepley     }
15529298eaa6SMatthew G Knepley   }
15539298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
15549298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
15559298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
15569298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
15579298eaa6SMatthew G Knepley }
15588415267dSToby Isaac 
15598415267dSToby Isaac #undef __FUNCT__
15608415267dSToby Isaac #define __FUNCT__ "DMPlexCreateReferenceCell"
15618415267dSToby Isaac /*@
15628415267dSToby Isaac   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
15638415267dSToby Isaac 
15648415267dSToby Isaac   Collective on comm
15658415267dSToby Isaac 
15668415267dSToby Isaac   Input Parameters:
15678415267dSToby Isaac + comm    - The communicator
15688415267dSToby Isaac . dim     - The spatial dimension
15698415267dSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
15708415267dSToby Isaac 
15718415267dSToby Isaac   Output Parameter:
15728415267dSToby Isaac . refdm - The reference cell
15738415267dSToby Isaac 
15748415267dSToby Isaac   Level: intermediate
15758415267dSToby Isaac 
15768415267dSToby Isaac .keywords: reference cell
15778415267dSToby Isaac .seealso:
15788415267dSToby Isaac @*/
15798415267dSToby Isaac PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
15808415267dSToby Isaac {
15818415267dSToby Isaac   DM             rdm;
15828415267dSToby Isaac   PetscErrorCode ierr;
15838415267dSToby Isaac 
15848415267dSToby Isaac   PetscFunctionBeginUser;
15858415267dSToby Isaac   ierr = DMCreate(comm, &rdm);CHKERRQ(ierr);
15868415267dSToby Isaac   ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr);
15878415267dSToby Isaac   ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr);
15888415267dSToby Isaac   switch (dim) {
15898415267dSToby Isaac   case 0:
15908415267dSToby Isaac   {
15918415267dSToby Isaac     PetscInt    numPoints[1]        = {1};
15928415267dSToby Isaac     PetscInt    coneSize[1]         = {0};
15938415267dSToby Isaac     PetscInt    cones[1]            = {0};
15948415267dSToby Isaac     PetscInt    coneOrientations[1] = {0};
15958415267dSToby Isaac     PetscScalar vertexCoords[1]     = {0.0};
15968415267dSToby Isaac 
15978415267dSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
15988415267dSToby Isaac   }
15998415267dSToby Isaac   break;
16008415267dSToby Isaac   case 1:
16018415267dSToby Isaac   {
16028415267dSToby Isaac     PetscInt    numPoints[2]        = {2, 1};
16038415267dSToby Isaac     PetscInt    coneSize[3]         = {2, 0, 0};
16048415267dSToby Isaac     PetscInt    cones[2]            = {1, 2};
16058415267dSToby Isaac     PetscInt    coneOrientations[2] = {0, 0};
16068415267dSToby Isaac     PetscScalar vertexCoords[2]     = {-1.0,  1.0};
16078415267dSToby Isaac 
16088415267dSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
16098415267dSToby Isaac   }
16108415267dSToby Isaac   break;
16118415267dSToby Isaac   case 2:
16128415267dSToby Isaac     if (simplex) {
16138415267dSToby Isaac       PetscInt    numPoints[2]        = {3, 1};
16148415267dSToby Isaac       PetscInt    coneSize[4]         = {3, 0, 0, 0};
16158415267dSToby Isaac       PetscInt    cones[3]            = {1, 2, 3};
16168415267dSToby Isaac       PetscInt    coneOrientations[3] = {0, 0, 0};
16178415267dSToby Isaac       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};
16188415267dSToby Isaac 
16198415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
16208415267dSToby Isaac     } else {
16218415267dSToby Isaac       PetscInt    numPoints[2]        = {4, 1};
16228415267dSToby Isaac       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
16238415267dSToby Isaac       PetscInt    cones[4]            = {1, 2, 3, 4};
16248415267dSToby Isaac       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
16258415267dSToby Isaac       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};
16268415267dSToby Isaac 
16278415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
16288415267dSToby Isaac     }
16298415267dSToby Isaac   break;
16308415267dSToby Isaac   case 3:
16318415267dSToby Isaac     if (simplex) {
16328415267dSToby Isaac       PetscInt    numPoints[2]        = {4, 1};
16338415267dSToby Isaac       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
16348415267dSToby Isaac       PetscInt    cones[4]            = {1, 3, 2, 4};
16358415267dSToby Isaac       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
16368415267dSToby Isaac       PetscScalar vertexCoords[12]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,  -1.0, -1.0, 1.0};
16378415267dSToby Isaac 
16388415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
16398415267dSToby Isaac     } else {
16408415267dSToby Isaac       PetscInt    numPoints[2]        = {8, 1};
16418415267dSToby Isaac       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
16428415267dSToby Isaac       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
16438415267dSToby Isaac       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
16448415267dSToby Isaac       PetscScalar vertexCoords[24]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  1.0, 1.0, -1.0,  -1.0, 1.0, -1.0,
16458415267dSToby Isaac                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};
16468415267dSToby Isaac 
16478415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
16488415267dSToby Isaac     }
16498415267dSToby Isaac   break;
16508415267dSToby Isaac   default:
16518415267dSToby Isaac     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
16528415267dSToby Isaac   }
16538415267dSToby Isaac   *refdm = NULL;
16548415267dSToby Isaac   ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr);
16558415267dSToby Isaac   ierr = DMPlexCopyCoordinates(rdm, *refdm);CHKERRQ(ierr);
16568415267dSToby Isaac   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
16578415267dSToby Isaac   PetscFunctionReturn(0);
16588415267dSToby Isaac }
1659