xref: /petsc/src/dm/impls/plex/plexcreate.c (revision 46919275265691675efdb415f86296e75faf7ea1)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <petscdmda.h>
40c312b8eSJed Brown #include <petscsf.h>
5552f7358SJed Brown 
61df5d5c5SMatthew G. Knepley /*@
71df5d5c5SMatthew G. Knepley   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
81df5d5c5SMatthew G. Knepley 
91df5d5c5SMatthew G. Knepley   Collective on MPI_Comm
101df5d5c5SMatthew G. Knepley 
111df5d5c5SMatthew G. Knepley   Input Parameters:
121df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object
131df5d5c5SMatthew G. Knepley . dim - The spatial dimension
141df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells
151df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces)
161df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement
171df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
181df5d5c5SMatthew G. Knepley 
191df5d5c5SMatthew G. Knepley   Output Parameter:
201df5d5c5SMatthew G. Knepley . dm  - The DM object
211df5d5c5SMatthew G. Knepley 
221df5d5c5SMatthew G. Knepley   Level: beginner
231df5d5c5SMatthew G. Knepley 
241df5d5c5SMatthew G. Knepley .keywords: DM, create
251df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate()
261df5d5c5SMatthew G. Knepley @*/
271df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
281df5d5c5SMatthew G. Knepley {
291df5d5c5SMatthew G. Knepley   DM             dm;
301df5d5c5SMatthew G. Knepley   PetscInt       p;
311df5d5c5SMatthew G. Knepley   PetscMPIInt    rank;
321df5d5c5SMatthew G. Knepley   PetscErrorCode ierr;
331df5d5c5SMatthew G. Knepley 
341df5d5c5SMatthew G. Knepley   PetscFunctionBegin;
351df5d5c5SMatthew G. Knepley   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
361df5d5c5SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
37c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
381df5d5c5SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
39ce78fa2fSMatthew G. Knepley   switch (dim) {
40ce78fa2fSMatthew G. Knepley   case 2:
41ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
42ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
43ce78fa2fSMatthew G. Knepley     break;
44ce78fa2fSMatthew G. Knepley   case 3:
45ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
46ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
47ce78fa2fSMatthew G. Knepley     break;
48ce78fa2fSMatthew G. Knepley   default:
49ce78fa2fSMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
50ce78fa2fSMatthew G. Knepley   }
511df5d5c5SMatthew G. Knepley   if (rank) {
521df5d5c5SMatthew G. Knepley     PetscInt numPoints[2] = {0, 0};
531df5d5c5SMatthew G. Knepley     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
541df5d5c5SMatthew G. Knepley   } else {
551df5d5c5SMatthew G. Knepley     switch (dim) {
561df5d5c5SMatthew G. Knepley     case 2:
571df5d5c5SMatthew G. Knepley       if (simplex) {
581df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {4, 2};
591df5d5c5SMatthew G. Knepley         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
601df5d5c5SMatthew G. Knepley         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
611df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
621df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
631df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
641df5d5c5SMatthew G. Knepley 
651df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
66c58f1c22SToby Isaac         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
671df5d5c5SMatthew G. Knepley       } else {
681df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {6, 2};
691df5d5c5SMatthew G. Knepley         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
701df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
711df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
721df5d5c5SMatthew 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};
731df5d5c5SMatthew G. Knepley 
741df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
751df5d5c5SMatthew G. Knepley       }
761df5d5c5SMatthew G. Knepley       break;
771df5d5c5SMatthew G. Knepley     case 3:
781df5d5c5SMatthew G. Knepley       if (simplex) {
791df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {5, 2};
801df5d5c5SMatthew G. Knepley         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
811df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
821df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
831df5d5c5SMatthew 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};
841df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
851df5d5c5SMatthew G. Knepley 
861df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
87c58f1c22SToby Isaac         for (p = 0; p < 5; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
881df5d5c5SMatthew G. Knepley       } else {
891df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]         = {12, 2};
901df5d5c5SMatthew G. Knepley         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
911df5d5c5SMatthew G. Knepley         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
921df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
931df5d5c5SMatthew 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,
941df5d5c5SMatthew 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,
951df5d5c5SMatthew 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};
961df5d5c5SMatthew G. Knepley 
971df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
981df5d5c5SMatthew G. Knepley       }
991df5d5c5SMatthew G. Knepley       break;
1001df5d5c5SMatthew G. Knepley     default:
1011df5d5c5SMatthew G. Knepley       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
1021df5d5c5SMatthew G. Knepley     }
1031df5d5c5SMatthew G. Knepley   }
1041df5d5c5SMatthew G. Knepley   *newdm = dm;
1051df5d5c5SMatthew G. Knepley   if (refinementLimit > 0.0) {
1061df5d5c5SMatthew G. Knepley     DM rdm;
1071df5d5c5SMatthew G. Knepley     const char *name;
1081df5d5c5SMatthew G. Knepley 
1091df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
1101df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
1111df5d5c5SMatthew G. Knepley     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
1121df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1131df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
1141df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1151df5d5c5SMatthew G. Knepley     *newdm = rdm;
1161df5d5c5SMatthew G. Knepley   }
1171df5d5c5SMatthew G. Knepley   if (interpolate) {
118e56d480eSMatthew G. Knepley     DM idm = NULL;
1191df5d5c5SMatthew G. Knepley     const char *name;
1201df5d5c5SMatthew G. Knepley 
1211df5d5c5SMatthew G. Knepley     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
1221df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1231df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
1241df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
125c58f1c22SToby Isaac     ierr = DMCopyLabels(*newdm, idm);CHKERRQ(ierr);
1261df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1271df5d5c5SMatthew G. Knepley     *newdm = idm;
1281df5d5c5SMatthew G. Knepley   }
1291df5d5c5SMatthew G. Knepley   {
1301df5d5c5SMatthew G. Knepley     DM refinedMesh     = NULL;
1311df5d5c5SMatthew G. Knepley     DM distributedMesh = NULL;
1321df5d5c5SMatthew G. Knepley 
1331df5d5c5SMatthew G. Knepley     /* Distribute mesh over processes */
13480cf41d5SMatthew G. Knepley     ierr = DMPlexDistribute(*newdm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
1351df5d5c5SMatthew G. Knepley     if (distributedMesh) {
1361df5d5c5SMatthew G. Knepley       ierr = DMDestroy(newdm);CHKERRQ(ierr);
1371df5d5c5SMatthew G. Knepley       *newdm = distributedMesh;
1381df5d5c5SMatthew G. Knepley     }
1391df5d5c5SMatthew G. Knepley     if (refinementUniform) {
1401df5d5c5SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
1411df5d5c5SMatthew G. Knepley       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
1421df5d5c5SMatthew G. Knepley       if (refinedMesh) {
1431df5d5c5SMatthew G. Knepley         ierr = DMDestroy(newdm);CHKERRQ(ierr);
1441df5d5c5SMatthew G. Knepley         *newdm = refinedMesh;
1451df5d5c5SMatthew G. Knepley       }
1461df5d5c5SMatthew G. Knepley     }
1471df5d5c5SMatthew G. Knepley   }
1481df5d5c5SMatthew G. Knepley   PetscFunctionReturn(0);
1491df5d5c5SMatthew G. Knepley }
1501df5d5c5SMatthew G. Knepley 
15126492d91SMatthew G. Knepley /*@
15226492d91SMatthew G. Knepley   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
153552f7358SJed Brown 
15426492d91SMatthew G. Knepley   Collective on MPI_Comm
15526492d91SMatthew G. Knepley 
15626492d91SMatthew G. Knepley   Input Parameters:
15726492d91SMatthew G. Knepley + comm  - The communicator for the DM object
15826492d91SMatthew G. Knepley . lower - The lower left corner coordinates
15926492d91SMatthew G. Knepley . upper - The upper right corner coordinates
16026492d91SMatthew G. Knepley - edges - The number of cells in each direction
16126492d91SMatthew G. Knepley 
16226492d91SMatthew G. Knepley   Output Parameter:
16326492d91SMatthew G. Knepley . dm  - The DM object
16426492d91SMatthew G. Knepley 
16526492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
16626492d91SMatthew G. Knepley $ 18--5-17--4--16
16726492d91SMatthew G. Knepley $  |     |     |
16826492d91SMatthew G. Knepley $  6    10     3
16926492d91SMatthew G. Knepley $  |     |     |
17026492d91SMatthew G. Knepley $ 19-11-20--9--15
17126492d91SMatthew G. Knepley $  |     |     |
17226492d91SMatthew G. Knepley $  7     8     2
17326492d91SMatthew G. Knepley $  |     |     |
17426492d91SMatthew G. Knepley $ 12--0-13--1--14
17526492d91SMatthew G. Knepley 
17626492d91SMatthew G. Knepley   Level: beginner
17726492d91SMatthew G. Knepley 
17826492d91SMatthew G. Knepley .keywords: DM, create
17926492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
18026492d91SMatthew G. Knepley @*/
181552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
182552f7358SJed Brown {
1831df21d24SMatthew G. Knepley   const PetscInt numVertices    = (edges[0]+1)*(edges[1]+1);
1841df21d24SMatthew G. Knepley   const PetscInt numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
185552f7358SJed Brown   PetscInt       markerTop      = 1;
186552f7358SJed Brown   PetscInt       markerBottom   = 1;
187552f7358SJed Brown   PetscInt       markerRight    = 1;
188552f7358SJed Brown   PetscInt       markerLeft     = 1;
189552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
190552f7358SJed Brown   Vec            coordinates;
191552f7358SJed Brown   PetscSection   coordSection;
192552f7358SJed Brown   PetscScalar    *coords;
193552f7358SJed Brown   PetscInt       coordSize;
194552f7358SJed Brown   PetscMPIInt    rank;
195552f7358SJed Brown   PetscInt       v, vx, vy;
196552f7358SJed Brown   PetscErrorCode ierr;
197552f7358SJed Brown 
198552f7358SJed Brown   PetscFunctionBegin;
199c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
200552f7358SJed Brown   if (markerSeparate) {
2011df21d24SMatthew G. Knepley     markerTop    = 3;
2021df21d24SMatthew G. Knepley     markerBottom = 1;
2031df21d24SMatthew G. Knepley     markerRight  = 2;
2041df21d24SMatthew G. Knepley     markerLeft   = 4;
205552f7358SJed Brown   }
20682f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
207552f7358SJed Brown   if (!rank) {
208552f7358SJed Brown     PetscInt e, ex, ey;
209552f7358SJed Brown 
210552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
211552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
212552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
213552f7358SJed Brown     }
214552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
215552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
216552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
217552f7358SJed Brown         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
218552f7358SJed Brown         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
219da80777bSKarl Rupp         PetscInt cone[2];
220552f7358SJed Brown 
221da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
222552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
223552f7358SJed Brown         if (vx == edges[0]) {
224c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
225c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
226552f7358SJed Brown           if (ey == edges[1]-1) {
227c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
228552f7358SJed Brown           }
229552f7358SJed Brown         } else if (vx == 0) {
230c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
231c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
232552f7358SJed Brown           if (ey == edges[1]-1) {
233c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
234552f7358SJed Brown           }
235552f7358SJed Brown         }
236552f7358SJed Brown       }
237552f7358SJed Brown     }
238552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
239552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
240552f7358SJed Brown         PetscInt edge   = vy*edges[0]     + ex;
241552f7358SJed Brown         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
242da80777bSKarl Rupp         PetscInt cone[2];
243552f7358SJed Brown 
244da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
245552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
246552f7358SJed Brown         if (vy == edges[1]) {
247c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
248c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
249552f7358SJed Brown           if (ex == edges[0]-1) {
250c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
251552f7358SJed Brown           }
252552f7358SJed Brown         } else if (vy == 0) {
253c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
254c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
255552f7358SJed Brown           if (ex == edges[0]-1) {
256c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
257552f7358SJed Brown           }
258552f7358SJed Brown         }
259552f7358SJed Brown       }
260552f7358SJed Brown     }
261552f7358SJed Brown   }
262552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
263552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
264552f7358SJed Brown   /* Build coordinates */
265c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
266972bc18aSToby Isaac   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
267552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
268972bc18aSToby Isaac   ierr = PetscSectionSetFieldComponents(coordSection, 0, 2);CHKERRQ(ierr);
269552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
270552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
271972bc18aSToby Isaac     ierr = PetscSectionSetFieldDof(coordSection, v, 0, 2);CHKERRQ(ierr);
272552f7358SJed Brown   }
273552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
274552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
2758b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
276da16285aSMichael Lange   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
277552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2788b9ced59SLisandro Dalcin   ierr = VecSetBlockSize(coordinates, 2);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 
29326492d91SMatthew G. Knepley /*@
29426492d91SMatthew G. Knepley   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
295552f7358SJed Brown 
29626492d91SMatthew G. Knepley   Collective on MPI_Comm
29726492d91SMatthew G. Knepley 
29826492d91SMatthew G. Knepley   Input Parameters:
29926492d91SMatthew G. Knepley + comm  - The communicator for the DM object
30026492d91SMatthew G. Knepley . lower - The lower left front corner coordinates
30126492d91SMatthew G. Knepley . upper - The upper right back corner coordinates
30226492d91SMatthew G. Knepley - edges - The number of cells in each direction
30326492d91SMatthew G. Knepley 
30426492d91SMatthew G. Knepley   Output Parameter:
30526492d91SMatthew G. Knepley . dm  - The DM object
30626492d91SMatthew G. Knepley 
30726492d91SMatthew G. Knepley   Level: beginner
30826492d91SMatthew G. Knepley 
30926492d91SMatthew G. Knepley .keywords: DM, create
31026492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
31126492d91SMatthew G. Knepley @*/
312552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
313552f7358SJed Brown {
3149e8abbc3SMichael Lange   PetscInt       vertices[3], numVertices;
3157b59f5a9SMichael Lange   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
316552f7358SJed Brown   Vec            coordinates;
317552f7358SJed Brown   PetscSection   coordSection;
318552f7358SJed Brown   PetscScalar    *coords;
319552f7358SJed Brown   PetscInt       coordSize;
320552f7358SJed Brown   PetscMPIInt    rank;
321552f7358SJed Brown   PetscInt       v, vx, vy, vz;
3227b59f5a9SMichael Lange   PetscInt       voffset, iface=0, cone[4];
323552f7358SJed Brown   PetscErrorCode ierr;
324552f7358SJed Brown 
325552f7358SJed Brown   PetscFunctionBegin;
32682f516ccSBarry 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");
32782f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
3289e8abbc3SMichael Lange   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
3299e8abbc3SMichael Lange   numVertices = vertices[0]*vertices[1]*vertices[2];
330552f7358SJed Brown   if (!rank) {
331552f7358SJed Brown     PetscInt f;
332552f7358SJed Brown 
333552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
334552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
335552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
336552f7358SJed Brown     }
337552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
338552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
339c58f1c22SToby Isaac       ierr = DMSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
340552f7358SJed Brown     }
3417b59f5a9SMichael Lange 
3427b59f5a9SMichael Lange     /* Side 0 (Top) */
3437b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3447b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3457b59f5a9SMichael Lange         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
3467b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
3477b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3487b59f5a9SMichael Lange         iface++;
349552f7358SJed Brown       }
350552f7358SJed Brown     }
3517b59f5a9SMichael Lange 
3527b59f5a9SMichael Lange     /* Side 1 (Bottom) */
3537b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3547b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3557b59f5a9SMichael Lange         voffset = numFaces + vy*(faces[0]+1) + vx;
3567b59f5a9SMichael Lange         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
3577b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3587b59f5a9SMichael Lange         iface++;
359552f7358SJed Brown       }
360552f7358SJed Brown     }
3617b59f5a9SMichael Lange 
3627b59f5a9SMichael Lange     /* Side 2 (Front) */
3637b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3647b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3657b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
3667b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
3677b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3687b59f5a9SMichael Lange         iface++;
369552f7358SJed Brown       }
3707b59f5a9SMichael Lange     }
3717b59f5a9SMichael Lange 
3727b59f5a9SMichael Lange     /* Side 3 (Back) */
3737b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3747b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3757b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
3767b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
3777b59f5a9SMichael Lange         cone[2] = voffset+1; cone[3] = voffset;
3787b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3797b59f5a9SMichael Lange         iface++;
3807b59f5a9SMichael Lange       }
3817b59f5a9SMichael Lange     }
3827b59f5a9SMichael Lange 
3837b59f5a9SMichael Lange     /* Side 4 (Left) */
3847b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3857b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3867b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
3877b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
3887b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
3897b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3907b59f5a9SMichael Lange         iface++;
3917b59f5a9SMichael Lange       }
3927b59f5a9SMichael Lange     }
3937b59f5a9SMichael Lange 
3947b59f5a9SMichael Lange     /* Side 5 (Right) */
3957b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3967b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3977b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx;
3987b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
3997b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
4007b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
4017b59f5a9SMichael Lange         iface++;
4027b59f5a9SMichael Lange       }
403552f7358SJed Brown     }
404552f7358SJed Brown   }
405552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
406552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
407552f7358SJed Brown   /* Build coordinates */
408c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
409552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
410552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
411552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
412552f7358SJed Brown   }
413552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
414552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
4158b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
416552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
417552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4188b9ced59SLisandro Dalcin   ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr);
4192eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
420552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
421552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
422552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
423552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
424552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
425552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
426552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
427552f7358SJed Brown       }
428552f7358SJed Brown     }
429552f7358SJed Brown   }
430552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
431552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
432552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
433552f7358SJed Brown   PetscFunctionReturn(0);
434552f7358SJed Brown }
435552f7358SJed Brown 
436d6218766SMatthew G. Knepley /*@
437d6218766SMatthew G. Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
438d6218766SMatthew G. Knepley 
439d6218766SMatthew G. Knepley   Collective on MPI_Comm
440d6218766SMatthew G. Knepley 
441d6218766SMatthew G. Knepley   Input Parameters:
442d6218766SMatthew G. Knepley + comm - The communicator for the DM object
443d6218766SMatthew G. Knepley . dim - The spatial dimension
444d6218766SMatthew G. Knepley . numFaces - Number of faces per dimension
445d6218766SMatthew G. Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
446d6218766SMatthew G. Knepley 
447d6218766SMatthew G. Knepley   Output Parameter:
448d6218766SMatthew G. Knepley . dm  - The DM object
449d6218766SMatthew G. Knepley 
450d6218766SMatthew G. Knepley   Level: beginner
451d6218766SMatthew G. Knepley 
452d6218766SMatthew G. Knepley .keywords: DM, create
453d6218766SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
454d6218766SMatthew G. Knepley @*/
455d6218766SMatthew G. Knepley PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm)
456d6218766SMatthew G. Knepley {
457d6218766SMatthew G. Knepley   DM             boundary;
458d6218766SMatthew G. Knepley   PetscErrorCode ierr;
459d6218766SMatthew G. Knepley 
460d6218766SMatthew G. Knepley   PetscFunctionBegin;
461d6218766SMatthew G. Knepley   PetscValidPointer(dm, 4);
462d6218766SMatthew G. Knepley   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
463d6218766SMatthew G. Knepley   PetscValidLogicalCollectiveInt(boundary,dim,2);
464d6218766SMatthew G. Knepley   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
465d6218766SMatthew G. Knepley   ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr);
466d6218766SMatthew G. Knepley   ierr = DMSetCoordinateDim(boundary, dim);CHKERRQ(ierr);
467d6218766SMatthew G. Knepley   switch (dim) {
468d6218766SMatthew G. Knepley   case 2:
469d6218766SMatthew G. Knepley   {
470d6218766SMatthew G. Knepley     PetscReal lower[2] = {0.0, 0.0};
471d6218766SMatthew G. Knepley     PetscReal upper[2] = {1.0, 1.0};
472d6218766SMatthew G. Knepley     PetscInt  edges[2];
473d6218766SMatthew G. Knepley 
474d6218766SMatthew G. Knepley     edges[0] = numFaces; edges[1] = numFaces;
475d6218766SMatthew G. Knepley     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
476d6218766SMatthew G. Knepley     break;
477d6218766SMatthew G. Knepley   }
478d6218766SMatthew G. Knepley   case 3:
479d6218766SMatthew G. Knepley   {
480d6218766SMatthew G. Knepley     PetscReal lower[3] = {0.0, 0.0, 0.0};
481d6218766SMatthew G. Knepley     PetscReal upper[3] = {1.0, 1.0, 1.0};
482d6218766SMatthew G. Knepley     PetscInt  faces[3];
483d6218766SMatthew G. Knepley 
484d6218766SMatthew G. Knepley     faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces;
485d6218766SMatthew G. Knepley     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
486d6218766SMatthew G. Knepley     break;
487d6218766SMatthew G. Knepley   }
488d6218766SMatthew G. Knepley   default:
489d6218766SMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
490d6218766SMatthew G. Knepley   }
491d6218766SMatthew G. Knepley   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
492d6218766SMatthew G. Knepley   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
493d6218766SMatthew G. Knepley   PetscFunctionReturn(0);
494d6218766SMatthew G. Knepley }
495d6218766SMatthew G. Knepley 
4963dfda0b1SToby Isaac static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
4973dfda0b1SToby Isaac {
498f4eb4c5dSMatthew G. Knepley   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
499f4eb4c5dSMatthew G. Knepley   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
500f4eb4c5dSMatthew G. Knepley   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
501f4eb4c5dSMatthew G. Knepley   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
502f4eb4c5dSMatthew G. Knepley   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
503f4eb4c5dSMatthew G. Knepley   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
5043dfda0b1SToby Isaac   PetscInt       dim;
505d8211ee3SMatthew G. Knepley   PetscBool      markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
5063dfda0b1SToby Isaac   PetscMPIInt    rank;
5073dfda0b1SToby Isaac   PetscErrorCode ierr;
5083dfda0b1SToby Isaac 
5093dfda0b1SToby Isaac   PetscFunctionBegin;
510f0226e14SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
5113dfda0b1SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
51250ae33c3SToby Isaac   ierr = DMCreateLabel(dm,"marker");CHKERRQ(ierr);
51350ae33c3SToby Isaac   ierr = DMCreateLabel(dm,"Face Sets");CHKERRQ(ierr);
514d8211ee3SMatthew G. Knepley   if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ||
515d8211ee3SMatthew G. Knepley       bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ||
516d8211ee3SMatthew G. Knepley       bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {
517*46919275SMatthew G. Knepley     ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);CHKERRQ(ierr);
518*46919275SMatthew G. Knepley     if (cutMarker) {ierr = DMCreateLabel(dm,"periodic_cut");CHKERRQ(ierr);}
519d8211ee3SMatthew G. Knepley   }
5203dfda0b1SToby Isaac   switch (dim) {
5213dfda0b1SToby Isaac   case 2:
522f4eb4c5dSMatthew G. Knepley     faceMarkerTop    = 3;
523f4eb4c5dSMatthew G. Knepley     faceMarkerBottom = 1;
524f4eb4c5dSMatthew G. Knepley     faceMarkerRight  = 2;
525f4eb4c5dSMatthew G. Knepley     faceMarkerLeft   = 4;
5263dfda0b1SToby Isaac     break;
5273dfda0b1SToby Isaac   case 3:
528f4eb4c5dSMatthew G. Knepley     faceMarkerBottom = 1;
529f4eb4c5dSMatthew G. Knepley     faceMarkerTop    = 2;
530f4eb4c5dSMatthew G. Knepley     faceMarkerFront  = 3;
531f4eb4c5dSMatthew G. Knepley     faceMarkerBack   = 4;
532f4eb4c5dSMatthew G. Knepley     faceMarkerRight  = 5;
533f4eb4c5dSMatthew G. Knepley     faceMarkerLeft   = 6;
5343dfda0b1SToby Isaac     break;
5353dfda0b1SToby Isaac   default:
5363dfda0b1SToby Isaac     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim);
5373dfda0b1SToby Isaac     break;
5383dfda0b1SToby Isaac   }
539c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
540f4eb4c5dSMatthew G. Knepley   if (markerSeparate) {
541f4eb4c5dSMatthew G. Knepley     markerBottom = faceMarkerBottom;
542f4eb4c5dSMatthew G. Knepley     markerTop    = faceMarkerTop;
543f4eb4c5dSMatthew G. Knepley     markerFront  = faceMarkerFront;
544f4eb4c5dSMatthew G. Knepley     markerBack   = faceMarkerBack;
545f4eb4c5dSMatthew G. Knepley     markerRight  = faceMarkerRight;
546f4eb4c5dSMatthew G. Knepley     markerLeft   = faceMarkerLeft;
5473dfda0b1SToby Isaac   }
5483dfda0b1SToby Isaac   {
5493dfda0b1SToby Isaac     const PetscInt numXEdges    = !rank ? edges[0] : 0;
5503dfda0b1SToby Isaac     const PetscInt numYEdges    = !rank ? edges[1] : 0;
5513dfda0b1SToby Isaac     const PetscInt numZEdges    = !rank ? edges[2] : 0;
5523dfda0b1SToby Isaac     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
5533dfda0b1SToby Isaac     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
55442206facSLisandro Dalcin     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
5553dfda0b1SToby Isaac     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
5563dfda0b1SToby Isaac     const PetscInt numXFaces    = numYEdges*numZEdges;
5573dfda0b1SToby Isaac     const PetscInt numYFaces    = numXEdges*numZEdges;
5583dfda0b1SToby Isaac     const PetscInt numZFaces    = numXEdges*numYEdges;
5593dfda0b1SToby Isaac     const PetscInt numTotXFaces = numXVertices*numXFaces;
5603dfda0b1SToby Isaac     const PetscInt numTotYFaces = numYVertices*numYFaces;
5613dfda0b1SToby Isaac     const PetscInt numTotZFaces = numZVertices*numZFaces;
5623dfda0b1SToby Isaac     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
5633dfda0b1SToby Isaac     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
5643dfda0b1SToby Isaac     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
5653dfda0b1SToby Isaac     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
5663dfda0b1SToby Isaac     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
5673dfda0b1SToby Isaac     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
5683dfda0b1SToby Isaac     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
5693dfda0b1SToby Isaac     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
5703dfda0b1SToby Isaac     const PetscInt firstYFace   = firstXFace + numTotXFaces;
5713dfda0b1SToby Isaac     const PetscInt firstZFace   = firstYFace + numTotYFaces;
5723dfda0b1SToby Isaac     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
5733dfda0b1SToby Isaac     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
5743dfda0b1SToby Isaac     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
5753dfda0b1SToby Isaac     Vec            coordinates;
5763dfda0b1SToby Isaac     PetscSection   coordSection;
5773dfda0b1SToby Isaac     PetscScalar   *coords;
5783dfda0b1SToby Isaac     PetscInt       coordSize;
5793dfda0b1SToby Isaac     PetscInt       v, vx, vy, vz;
5803dfda0b1SToby Isaac     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;
5813dfda0b1SToby Isaac 
5823dfda0b1SToby Isaac     ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr);
5833dfda0b1SToby Isaac     for (c = 0; c < numCells; c++) {
5843dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr);
5853dfda0b1SToby Isaac     }
5863dfda0b1SToby Isaac     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
5873dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
5883dfda0b1SToby Isaac     }
5893dfda0b1SToby Isaac     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
5903dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
5913dfda0b1SToby Isaac     }
5923dfda0b1SToby Isaac     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
5933dfda0b1SToby Isaac     /* Build cells */
5943dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
5953dfda0b1SToby Isaac       for (fy = 0; fy < numYEdges; ++fy) {
5963dfda0b1SToby Isaac         for (fx = 0; fx < numXEdges; ++fx) {
5973dfda0b1SToby Isaac           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
5983dfda0b1SToby Isaac           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
5993dfda0b1SToby Isaac           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
6003dfda0b1SToby Isaac           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
6013dfda0b1SToby Isaac           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
6023dfda0b1SToby Isaac           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
6033dfda0b1SToby Isaac           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
6043dfda0b1SToby Isaac                             /* B,  T,  F,  K,  R,  L */
60542206facSLisandro Dalcin           PetscInt ornt[6] = {-4,  0,  0, -1,  0, -4}; /* ??? */
60642206facSLisandro Dalcin           PetscInt cone[6];
6073dfda0b1SToby Isaac 
6083dfda0b1SToby Isaac           /* no boundary twisting in 3D */
6093dfda0b1SToby Isaac           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
6103dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr);
6113dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr);
6123dfda0b1SToby Isaac         }
6133dfda0b1SToby Isaac       }
6143dfda0b1SToby Isaac     }
6153dfda0b1SToby Isaac     /* Build x faces */
6163dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
6173dfda0b1SToby Isaac       for (fy = 0; fy < numYEdges; ++fy) {
6183dfda0b1SToby Isaac         for (fx = 0; fx < numXVertices; ++fx) {
6193dfda0b1SToby Isaac           PetscInt face    = firstXFace + (fz*numYEdges+fy)*numXVertices + fx;
6203dfda0b1SToby Isaac           PetscInt edgeL   = firstZEdge + (  fy*                 numXVertices+fx)*numZEdges + fz;
6213dfda0b1SToby Isaac           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
6223dfda0b1SToby Isaac           PetscInt edgeB   = firstYEdge + (  fz*                 numXVertices+fx)*numYEdges + fy;
6233dfda0b1SToby Isaac           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
6243dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
6253dfda0b1SToby Isaac           PetscInt cone[4];
6263dfda0b1SToby Isaac 
6273dfda0b1SToby Isaac           if (dim == 3) {
6283dfda0b1SToby Isaac             /* markers */
6293dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
6303dfda0b1SToby Isaac               if (fx == numXVertices-1) {
631c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);CHKERRQ(ierr);
632c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr);
6333dfda0b1SToby Isaac               }
6343dfda0b1SToby Isaac               else if (fx == 0) {
635c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);CHKERRQ(ierr);
636c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr);
6373dfda0b1SToby Isaac               }
6383dfda0b1SToby Isaac             }
6393dfda0b1SToby Isaac           }
6403dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
6413dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
6423dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
6433dfda0b1SToby Isaac         }
6443dfda0b1SToby Isaac       }
6453dfda0b1SToby Isaac     }
6463dfda0b1SToby Isaac     /* Build y faces */
6473dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
64842206facSLisandro Dalcin       for (fx = 0; fx < numXEdges; ++fx) {
6493dfda0b1SToby Isaac         for (fy = 0; fy < numYVertices; ++fy) {
6503dfda0b1SToby Isaac           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
6513dfda0b1SToby Isaac           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx                 )*numZEdges + fz;
6523dfda0b1SToby Isaac           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
6533dfda0b1SToby Isaac           PetscInt edgeB   = firstXEdge + (  fz                 *numYVertices+fy)*numXEdges + fx;
6543dfda0b1SToby Isaac           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
6553dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
6563dfda0b1SToby Isaac           PetscInt cone[4];
6573dfda0b1SToby Isaac 
6583dfda0b1SToby Isaac           if (dim == 3) {
6593dfda0b1SToby Isaac             /* markers */
6603dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
6613dfda0b1SToby Isaac               if (fy == numYVertices-1) {
662c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);CHKERRQ(ierr);
663c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr);
6643dfda0b1SToby Isaac               }
6653dfda0b1SToby Isaac               else if (fy == 0) {
666c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);CHKERRQ(ierr);
667c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr);
6683dfda0b1SToby Isaac               }
6693dfda0b1SToby Isaac             }
6703dfda0b1SToby Isaac           }
6713dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
6723dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
6733dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
6743dfda0b1SToby Isaac         }
6753dfda0b1SToby Isaac       }
6763dfda0b1SToby Isaac     }
6773dfda0b1SToby Isaac     /* Build z faces */
6783dfda0b1SToby Isaac     for (fy = 0; fy < numYEdges; ++fy) {
6793dfda0b1SToby Isaac       for (fx = 0; fx < numXEdges; ++fx) {
6803dfda0b1SToby Isaac         for (fz = 0; fz < numZVertices; fz++) {
6813dfda0b1SToby Isaac           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
6823dfda0b1SToby Isaac           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx                 )*numYEdges + fy;
6833dfda0b1SToby Isaac           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
6843dfda0b1SToby Isaac           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy                 )*numXEdges + fx;
6853dfda0b1SToby Isaac           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
6863dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
6873dfda0b1SToby Isaac           PetscInt cone[4];
6883dfda0b1SToby Isaac 
6893dfda0b1SToby Isaac           if (dim == 2) {
6903dfda0b1SToby Isaac             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
6913dfda0b1SToby Isaac             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
6923dfda0b1SToby Isaac           }
6933dfda0b1SToby Isaac           else {
6943dfda0b1SToby Isaac             /* markers */
6953dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
6963dfda0b1SToby Isaac               if (fz == numZVertices-1) {
697c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);CHKERRQ(ierr);
698c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr);
6993dfda0b1SToby Isaac               }
7003dfda0b1SToby Isaac               else if (fz == 0) {
701c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);CHKERRQ(ierr);
702c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr);
7033dfda0b1SToby Isaac               }
7043dfda0b1SToby Isaac             }
7053dfda0b1SToby Isaac           }
7063dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
7073dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
7083dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
7093dfda0b1SToby Isaac         }
7103dfda0b1SToby Isaac       }
7113dfda0b1SToby Isaac     }
7123dfda0b1SToby Isaac     /* Build Z edges*/
7133dfda0b1SToby Isaac     for (vy = 0; vy < numYVertices; vy++) {
7143dfda0b1SToby Isaac       for (vx = 0; vx < numXVertices; vx++) {
7153dfda0b1SToby Isaac         for (ez = 0; ez < numZEdges; ez++) {
7163dfda0b1SToby Isaac           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
7173dfda0b1SToby Isaac           const PetscInt vertexB = firstVertex + (  ez                 *numYVertices+vy)*numXVertices + vx;
7183dfda0b1SToby Isaac           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
7193dfda0b1SToby Isaac           PetscInt       cone[2];
7203dfda0b1SToby Isaac 
7213dfda0b1SToby Isaac           if (dim == 3) {
7223dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
7233dfda0b1SToby Isaac               if (vx == numXVertices-1) {
724c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
7253dfda0b1SToby Isaac               }
7263dfda0b1SToby Isaac               else if (vx == 0) {
727c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
7283dfda0b1SToby Isaac               }
7293dfda0b1SToby Isaac             }
7303dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
7313dfda0b1SToby Isaac               if (vy == numYVertices-1) {
732c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
7333dfda0b1SToby Isaac               }
7343dfda0b1SToby Isaac               else if (vy == 0) {
735c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
7363dfda0b1SToby Isaac               }
7373dfda0b1SToby Isaac             }
7383dfda0b1SToby Isaac           }
7393dfda0b1SToby Isaac           cone[0] = vertexB; cone[1] = vertexT;
7403dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
7413dfda0b1SToby Isaac         }
7423dfda0b1SToby Isaac       }
7433dfda0b1SToby Isaac     }
7443dfda0b1SToby Isaac     /* Build Y edges*/
7453dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; vz++) {
7463dfda0b1SToby Isaac       for (vx = 0; vx < numXVertices; vx++) {
7473dfda0b1SToby Isaac         for (ey = 0; ey < numYEdges; ey++) {
7483dfda0b1SToby Isaac           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
7493dfda0b1SToby Isaac           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
7503dfda0b1SToby Isaac           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
7513dfda0b1SToby Isaac           const PetscInt vertexK = firstVertex + nextv;
7523dfda0b1SToby Isaac           PetscInt       cone[2];
7533dfda0b1SToby Isaac 
7543dfda0b1SToby Isaac           cone[0] = vertexF; cone[1] = vertexK;
7553dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
7563dfda0b1SToby Isaac           if (dim == 2) {
7573dfda0b1SToby Isaac             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
7583dfda0b1SToby Isaac               if (vx == numXVertices-1) {
759c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);CHKERRQ(ierr);
760c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
761c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
7623dfda0b1SToby Isaac                 if (ey == numYEdges-1) {
763c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
7643dfda0b1SToby Isaac                 }
765d8211ee3SMatthew G. Knepley               } else if (vx == 0) {
766c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);CHKERRQ(ierr);
767c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
768c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
7693dfda0b1SToby Isaac                 if (ey == numYEdges-1) {
770c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
7713dfda0b1SToby Isaac                 }
7723dfda0b1SToby Isaac               }
773d8211ee3SMatthew G. Knepley             } else {
774d8211ee3SMatthew G. Knepley               if (vx == 0 && cutMarker) {
775*46919275SMatthew G. Knepley                 ierr = DMSetLabelValue(dm, "periodic_cut", edge,    1);CHKERRQ(ierr);
776d8211ee3SMatthew G. Knepley                 ierr = DMSetLabelValue(dm, "periodic_cut", cone[0], 1);CHKERRQ(ierr);
777d8211ee3SMatthew G. Knepley                 if (ey == numYEdges-1) {
778d8211ee3SMatthew G. Knepley                   ierr = DMSetLabelValue(dm, "periodic_cut", cone[1], 1);CHKERRQ(ierr);
7793dfda0b1SToby Isaac                 }
7803dfda0b1SToby Isaac               }
781d8211ee3SMatthew G. Knepley             }
782d8211ee3SMatthew G. Knepley           } else {
7833dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
7843dfda0b1SToby Isaac               if (vx == numXVertices-1) {
785c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
786d8211ee3SMatthew G. Knepley               } else if (vx == 0) {
787c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
7883dfda0b1SToby Isaac               }
7893dfda0b1SToby Isaac             }
7903dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
7913dfda0b1SToby Isaac               if (vz == numZVertices-1) {
792c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
793d8211ee3SMatthew G. Knepley               } else if (vz == 0) {
794c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
7953dfda0b1SToby Isaac               }
7963dfda0b1SToby Isaac             }
7973dfda0b1SToby Isaac           }
7983dfda0b1SToby Isaac         }
7993dfda0b1SToby Isaac       }
8003dfda0b1SToby Isaac     }
8013dfda0b1SToby Isaac     /* Build X edges*/
8023dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; vz++) {
8033dfda0b1SToby Isaac       for (vy = 0; vy < numYVertices; vy++) {
8043dfda0b1SToby Isaac         for (ex = 0; ex < numXEdges; ex++) {
8053dfda0b1SToby Isaac           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
8063dfda0b1SToby Isaac           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
8073dfda0b1SToby Isaac           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
8083dfda0b1SToby Isaac           const PetscInt vertexR = firstVertex + nextv;
8093dfda0b1SToby Isaac           PetscInt       cone[2];
8103dfda0b1SToby Isaac 
8113dfda0b1SToby Isaac           cone[0] = vertexL; cone[1] = vertexR;
8123dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
8133dfda0b1SToby Isaac           if (dim == 2) {
8143dfda0b1SToby Isaac             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
8153dfda0b1SToby Isaac               if (vy == numYVertices-1) {
816c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);CHKERRQ(ierr);
817c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
818c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
8193dfda0b1SToby Isaac                 if (ex == numXEdges-1) {
820c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
8213dfda0b1SToby Isaac                 }
822d8211ee3SMatthew G. Knepley               } else if (vy == 0) {
823c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);CHKERRQ(ierr);
824c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
825c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
8263dfda0b1SToby Isaac                 if (ex == numXEdges-1) {
827c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
8283dfda0b1SToby Isaac                 }
8293dfda0b1SToby Isaac               }
830d8211ee3SMatthew G. Knepley             } else {
831d8211ee3SMatthew G. Knepley               if (vy == 0 && cutMarker) {
832*46919275SMatthew G. Knepley                 ierr = DMSetLabelValue(dm, "periodic_cut", edge,    1);CHKERRQ(ierr);
833d8211ee3SMatthew G. Knepley                 ierr = DMSetLabelValue(dm, "periodic_cut", cone[0], 1);CHKERRQ(ierr);
834d8211ee3SMatthew G. Knepley                 if (ex == numXEdges-1) {
835d8211ee3SMatthew G. Knepley                   ierr = DMSetLabelValue(dm, "periodic_cut", cone[1], 1);CHKERRQ(ierr);
8363dfda0b1SToby Isaac                 }
8373dfda0b1SToby Isaac               }
838d8211ee3SMatthew G. Knepley             }
839d8211ee3SMatthew G. Knepley           } else {
8403dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
8413dfda0b1SToby Isaac               if (vy == numYVertices-1) {
842c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
8433dfda0b1SToby Isaac               }
8443dfda0b1SToby Isaac               else if (vy == 0) {
845c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
8463dfda0b1SToby Isaac               }
8473dfda0b1SToby Isaac             }
8483dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
8493dfda0b1SToby Isaac               if (vz == numZVertices-1) {
850c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
8513dfda0b1SToby Isaac               }
8523dfda0b1SToby Isaac               else if (vz == 0) {
853c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
8543dfda0b1SToby Isaac               }
8553dfda0b1SToby Isaac             }
8563dfda0b1SToby Isaac           }
8573dfda0b1SToby Isaac         }
8583dfda0b1SToby Isaac       }
8593dfda0b1SToby Isaac     }
8603dfda0b1SToby Isaac     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
8613dfda0b1SToby Isaac     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
8623dfda0b1SToby Isaac     /* Build coordinates */
8633dfda0b1SToby Isaac     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
8643dfda0b1SToby Isaac     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
8653dfda0b1SToby Isaac     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
8663dfda0b1SToby Isaac     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
8673dfda0b1SToby Isaac     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
8683dfda0b1SToby Isaac       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
8693dfda0b1SToby Isaac       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
8703dfda0b1SToby Isaac     }
8713dfda0b1SToby Isaac     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
8723dfda0b1SToby Isaac     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
8738b9ced59SLisandro Dalcin     ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
874da16285aSMichael Lange     ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
8753dfda0b1SToby Isaac     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
8768b9ced59SLisandro Dalcin     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
8773dfda0b1SToby Isaac     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
8783dfda0b1SToby Isaac     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
8793dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; ++vz) {
8803dfda0b1SToby Isaac       for (vy = 0; vy < numYVertices; ++vy) {
8813dfda0b1SToby Isaac         for (vx = 0; vx < numXVertices; ++vx) {
8823dfda0b1SToby Isaac           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
8833dfda0b1SToby Isaac           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
8843dfda0b1SToby Isaac           if (dim == 3) {
8853dfda0b1SToby Isaac             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
8863dfda0b1SToby Isaac           }
8873dfda0b1SToby Isaac         }
8883dfda0b1SToby Isaac       }
8893dfda0b1SToby Isaac     }
8903dfda0b1SToby Isaac     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
8913dfda0b1SToby Isaac     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
8923dfda0b1SToby Isaac     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
8933dfda0b1SToby Isaac   }
8943dfda0b1SToby Isaac   PetscFunctionReturn(0);
8953dfda0b1SToby Isaac }
8963dfda0b1SToby Isaac 
89726492d91SMatthew G. Knepley /*@
89826492d91SMatthew G. Knepley   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
89926492d91SMatthew G. Knepley 
90026492d91SMatthew G. Knepley   Collective on MPI_Comm
90126492d91SMatthew G. Knepley 
90226492d91SMatthew G. Knepley   Input Parameters:
90326492d91SMatthew G. Knepley + comm  - The communicator for the DM object
90426492d91SMatthew G. Knepley . dim   - The spatial dimension
905fbf5b483SMatthew G. Knepley . periodicX - The boundary type for the X direction
906fbf5b483SMatthew G. Knepley . periodicY - The boundary type for the Y direction
907fbf5b483SMatthew G. Knepley . periodicZ - The boundary type for the Z direction
90826492d91SMatthew G. Knepley - cells - The number of cells in each direction
90926492d91SMatthew G. Knepley 
91026492d91SMatthew G. Knepley   Output Parameter:
91126492d91SMatthew G. Knepley . dm  - The DM object
91226492d91SMatthew G. Knepley 
913d6218766SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
914d6218766SMatthew G. Knepley $ 22--8-23--9--24
915d6218766SMatthew G. Knepley $  |     |     |
916d6218766SMatthew G. Knepley $ 13  2 14  3  15
917d6218766SMatthew G. Knepley $  |     |     |
918d6218766SMatthew G. Knepley $ 19--6-20--7--21
919d6218766SMatthew G. Knepley $  |     |     |
920d6218766SMatthew G. Knepley $ 10  0 11  1 12
921d6218766SMatthew G. Knepley $  |     |     |
922d6218766SMatthew G. Knepley $ 16--4-17--5--18
923d6218766SMatthew G. Knepley 
92426492d91SMatthew G. Knepley   Level: beginner
92526492d91SMatthew G. Knepley 
92626492d91SMatthew G. Knepley .keywords: DM, create
92726492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
92826492d91SMatthew G. Knepley @*/
929bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
930a6dfd86eSKarl Rupp {
931c8c68bd8SToby Isaac   PetscInt       i;
932552f7358SJed Brown   PetscErrorCode ierr;
933552f7358SJed Brown 
934552f7358SJed Brown   PetscFunctionBegin;
935d67d2e29SLisandro Dalcin   PetscValidPointer(dm, 7);
936552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
937552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
938552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
939c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
940552f7358SJed Brown   switch (dim) {
941552f7358SJed Brown   case 2:
942552f7358SJed Brown   {
943d6218766SMatthew G. Knepley     PetscReal lower[3] = {0.0, 0.0, 0.0};
944d6218766SMatthew G. Knepley     PetscReal upper[3] = {1.0, 1.0, 0.0};
945d35a276dSMatthew G. Knepley     PetscInt  edges[3] = {cells[0], cells[1], 0};
946552f7358SJed Brown 
947d35a276dSMatthew G. Knepley     ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, edges, periodicX, periodicY, DM_BOUNDARY_NONE);CHKERRQ(ierr);
948c8c68bd8SToby Isaac     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
949c8c68bd8SToby Isaac         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST) {
950c8c68bd8SToby Isaac       PetscReal      L[2];
951c8c68bd8SToby Isaac       PetscReal      maxCell[2];
95270012f4bSToby Isaac       DMBoundaryType bdType[2];
953c8c68bd8SToby Isaac 
95470012f4bSToby Isaac       bdType[0] = periodicX;
95570012f4bSToby Isaac       bdType[1] = periodicY;
956c8c68bd8SToby Isaac       for (i = 0; i < dim; i++) {
957c8c68bd8SToby Isaac         L[i]       = upper[i] - lower[i];
958c8c68bd8SToby Isaac         maxCell[i] = 1.1 * (L[i] / cells[i]);
959c8c68bd8SToby Isaac       }
960c8c68bd8SToby Isaac 
961c8c68bd8SToby Isaac       ierr = DMSetPeriodicity(*dm,maxCell,L,bdType);CHKERRQ(ierr);
962c8c68bd8SToby Isaac     }
963552f7358SJed Brown     break;
964552f7358SJed Brown   }
965552f7358SJed Brown   case 3:
966552f7358SJed Brown   {
967552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
968552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
969552f7358SJed Brown 
9703dfda0b1SToby Isaac     ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
971c8c68bd8SToby Isaac     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
972c8c68bd8SToby Isaac         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST ||
973c8c68bd8SToby Isaac         periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
974c8c68bd8SToby Isaac       PetscReal      L[3];
975c8c68bd8SToby Isaac       PetscReal      maxCell[3];
97635e84bc7SToby Isaac       DMBoundaryType bdType[3];
977c8c68bd8SToby Isaac 
97835e84bc7SToby Isaac       bdType[0] = periodicX;
97935e84bc7SToby Isaac       bdType[1] = periodicY;
98035e84bc7SToby Isaac       bdType[2] = periodicZ;
981c8c68bd8SToby Isaac       for (i = 0; i < dim; i++) {
982c8c68bd8SToby Isaac         L[i]       = upper[i] - lower[i];
983c8c68bd8SToby Isaac         maxCell[i] = 1.1 * (L[i] / cells[i]);
984c8c68bd8SToby Isaac       }
985c8c68bd8SToby Isaac 
986c8c68bd8SToby Isaac       ierr = DMSetPeriodicity(*dm,maxCell,L,bdType);CHKERRQ(ierr);
987c8c68bd8SToby Isaac     }
988552f7358SJed Brown     break;
989552f7358SJed Brown   }
990552f7358SJed Brown   default:
991552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
992552f7358SJed Brown   }
993552f7358SJed Brown   PetscFunctionReturn(0);
994552f7358SJed Brown }
995552f7358SJed Brown 
996a9074c1eSMatthew G. Knepley /*@C
997a9074c1eSMatthew G. Knepley   DMPlexSetOptionsPrefix - Sets the prefix used for searching for all DM options in the database.
998a9074c1eSMatthew G. Knepley 
999a9074c1eSMatthew G. Knepley   Logically Collective on DM
1000a9074c1eSMatthew G. Knepley 
1001a9074c1eSMatthew G. Knepley   Input Parameters:
1002a9074c1eSMatthew G. Knepley + dm - the DM context
1003a9074c1eSMatthew G. Knepley - prefix - the prefix to prepend to all option names
1004a9074c1eSMatthew G. Knepley 
1005a9074c1eSMatthew G. Knepley   Notes:
1006a9074c1eSMatthew G. Knepley   A hyphen (-) must NOT be given at the beginning of the prefix name.
1007a9074c1eSMatthew G. Knepley   The first character of all runtime options is AUTOMATICALLY the hyphen.
1008a9074c1eSMatthew G. Knepley 
1009a9074c1eSMatthew G. Knepley   Level: advanced
1010a9074c1eSMatthew G. Knepley 
1011a9074c1eSMatthew G. Knepley .seealso: SNESSetFromOptions()
1012a9074c1eSMatthew G. Knepley @*/
1013a9074c1eSMatthew G. Knepley PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1014a9074c1eSMatthew G. Knepley {
1015a9074c1eSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
1016a9074c1eSMatthew G. Knepley   PetscErrorCode ierr;
1017a9074c1eSMatthew G. Knepley 
1018a9074c1eSMatthew G. Knepley   PetscFunctionBegin;
1019a9074c1eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1020a9074c1eSMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);CHKERRQ(ierr);
1021a9074c1eSMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);CHKERRQ(ierr);
1022a9074c1eSMatthew G. Knepley   PetscFunctionReturn(0);
1023a9074c1eSMatthew G. Knepley }
1024a9074c1eSMatthew G. Knepley 
10250510c589SMatthew G. Knepley /*@
10260510c589SMatthew G. Knepley   DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra.
10270510c589SMatthew G. Knepley 
10280510c589SMatthew G. Knepley   Collective on MPI_Comm
10290510c589SMatthew G. Knepley 
10300510c589SMatthew G. Knepley   Input Parameters:
10310510c589SMatthew G. Knepley + comm      - The communicator for the DM object
1032dbc1dc17SMatthew G. Knepley . numRefine - The number of regular refinements to the basic 5 cell structure
10330510c589SMatthew G. Knepley - periodicZ - The boundary type for the Z direction
10340510c589SMatthew G. Knepley 
10350510c589SMatthew G. Knepley   Output Parameter:
10360510c589SMatthew G. Knepley . dm  - The DM object
10370510c589SMatthew G. Knepley 
10380510c589SMatthew G. Knepley   Note: Here is the output numbering looking from the bottom of the cylinder:
10390510c589SMatthew G. Knepley $       17-----14
10400510c589SMatthew G. Knepley $        |     |
10410510c589SMatthew G. Knepley $        |  2  |
10420510c589SMatthew G. Knepley $        |     |
10430510c589SMatthew G. Knepley $ 17-----8-----7-----14
10440510c589SMatthew G. Knepley $  |     |     |     |
10450510c589SMatthew G. Knepley $  |  3  |  0  |  1  |
10460510c589SMatthew G. Knepley $  |     |     |     |
10470510c589SMatthew G. Knepley $ 19-----5-----6-----13
10480510c589SMatthew G. Knepley $        |     |
10490510c589SMatthew G. Knepley $        |  4  |
10500510c589SMatthew G. Knepley $        |     |
10510510c589SMatthew G. Knepley $       19-----13
10520510c589SMatthew G. Knepley $
10530510c589SMatthew G. Knepley $ and up through the top
10540510c589SMatthew G. Knepley $
10550510c589SMatthew G. Knepley $       18-----16
10560510c589SMatthew G. Knepley $        |     |
10570510c589SMatthew G. Knepley $        |  2  |
10580510c589SMatthew G. Knepley $        |     |
10590510c589SMatthew G. Knepley $ 18----10----11-----16
10600510c589SMatthew G. Knepley $  |     |     |     |
10610510c589SMatthew G. Knepley $  |  3  |  0  |  1  |
10620510c589SMatthew G. Knepley $  |     |     |     |
10630510c589SMatthew G. Knepley $ 20-----9----12-----15
10640510c589SMatthew G. Knepley $        |     |
10650510c589SMatthew G. Knepley $        |  4  |
10660510c589SMatthew G. Knepley $        |     |
10670510c589SMatthew G. Knepley $       20-----15
10680510c589SMatthew G. Knepley 
10690510c589SMatthew G. Knepley   Level: beginner
10700510c589SMatthew G. Knepley 
10710510c589SMatthew G. Knepley .keywords: DM, create
10720510c589SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
10730510c589SMatthew G. Knepley @*/
1074dbc1dc17SMatthew G. Knepley PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
10750510c589SMatthew G. Knepley {
10760510c589SMatthew G. Knepley   const PetscInt dim = 3;
1077006a8963SMatthew G. Knepley   PetscInt       numCells, numVertices, r;
10780510c589SMatthew G. Knepley   PetscErrorCode ierr;
10790510c589SMatthew G. Knepley 
10800510c589SMatthew G. Knepley   PetscFunctionBegin;
10810510c589SMatthew G. Knepley   PetscValidPointer(dm, 4);
1082dbc1dc17SMatthew G. Knepley   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
10830510c589SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
10840510c589SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
10850510c589SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
10860510c589SMatthew G. Knepley   /* Create topology */
10870510c589SMatthew G. Knepley   {
10880510c589SMatthew G. Knepley     PetscInt cone[8], c;
10890510c589SMatthew G. Knepley 
1090006a8963SMatthew G. Knepley     numCells    = 5;
1091006a8963SMatthew G. Knepley     numVertices = 16;
1092006a8963SMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1093ae8bcbbbSMatthew G. Knepley       numCells   *= 3;
1094ae8bcbbbSMatthew G. Knepley       numVertices = 24;
1095006a8963SMatthew G. Knepley     }
10960510c589SMatthew G. Knepley     ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
10970510c589SMatthew G. Knepley     for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 8);CHKERRQ(ierr);}
10980510c589SMatthew G. Knepley     ierr = DMSetUp(*dm);CHKERRQ(ierr);
1099006a8963SMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1100ae8bcbbbSMatthew G. Knepley       cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1101ae8bcbbbSMatthew G. Knepley       cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1102006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr);
1103ae8bcbbbSMatthew G. Knepley       cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1104ae8bcbbbSMatthew G. Knepley       cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1105006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr);
1106ae8bcbbbSMatthew G. Knepley       cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1107ae8bcbbbSMatthew G. Knepley       cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1108006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr);
1109ae8bcbbbSMatthew G. Knepley       cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1110ae8bcbbbSMatthew G. Knepley       cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1111006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr);
1112ae8bcbbbSMatthew G. Knepley       cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1113ae8bcbbbSMatthew G. Knepley       cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1114006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr);
1115006a8963SMatthew G. Knepley 
1116ae8bcbbbSMatthew G. Knepley       cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1117ae8bcbbbSMatthew G. Knepley       cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1118006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr);
1119ae8bcbbbSMatthew G. Knepley       cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1120ae8bcbbbSMatthew G. Knepley       cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1121006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 6, cone);CHKERRQ(ierr);
1122ae8bcbbbSMatthew G. Knepley       cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1123ae8bcbbbSMatthew G. Knepley       cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1124006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 7, cone);CHKERRQ(ierr);
1125ae8bcbbbSMatthew G. Knepley       cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1126ae8bcbbbSMatthew G. Knepley       cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1127006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 8, cone);CHKERRQ(ierr);
1128ae8bcbbbSMatthew G. Knepley       cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1129ae8bcbbbSMatthew G. Knepley       cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1130006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 9, cone);CHKERRQ(ierr);
1131ae8bcbbbSMatthew G. Knepley 
1132ae8bcbbbSMatthew G. Knepley       cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1133ae8bcbbbSMatthew G. Knepley       cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1134ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 10, cone);CHKERRQ(ierr);
1135ae8bcbbbSMatthew G. Knepley       cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1136ae8bcbbbSMatthew G. Knepley       cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1137ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 11, cone);CHKERRQ(ierr);
1138ae8bcbbbSMatthew G. Knepley       cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1139ae8bcbbbSMatthew G. Knepley       cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1140ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 12, cone);CHKERRQ(ierr);
1141ae8bcbbbSMatthew G. Knepley       cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1142ae8bcbbbSMatthew G. Knepley       cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1143ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 13, cone);CHKERRQ(ierr);
1144ae8bcbbbSMatthew G. Knepley       cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1145ae8bcbbbSMatthew G. Knepley       cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1146ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr);
1147006a8963SMatthew G. Knepley     } else {
114810c6f908SMatthew G. Knepley       cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
114910c6f908SMatthew G. Knepley       cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
11500510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr);
115110c6f908SMatthew G. Knepley       cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
115210c6f908SMatthew G. Knepley       cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
11530510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr);
115410c6f908SMatthew G. Knepley       cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
115510c6f908SMatthew G. Knepley       cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
11560510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr);
115710c6f908SMatthew G. Knepley       cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
115810c6f908SMatthew G. Knepley       cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
11590510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr);
116010c6f908SMatthew G. Knepley       cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
116110c6f908SMatthew G. Knepley       cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
11620510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr);
1163006a8963SMatthew G. Knepley     }
11640510c589SMatthew G. Knepley     ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
11650510c589SMatthew G. Knepley     ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
11660510c589SMatthew G. Knepley   }
11670510c589SMatthew G. Knepley   /* Interpolate */
11680510c589SMatthew G. Knepley   {
11690510c589SMatthew G. Knepley     DM idm = NULL;
11700510c589SMatthew G. Knepley 
11710510c589SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
11720510c589SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
11730510c589SMatthew G. Knepley     *dm  = idm;
11740510c589SMatthew G. Knepley   }
1175dbc1dc17SMatthew G. Knepley   /* Create cube geometry */
11760510c589SMatthew G. Knepley   {
11770510c589SMatthew G. Knepley     Vec             coordinates;
11780510c589SMatthew G. Knepley     PetscSection    coordSection;
11790510c589SMatthew G. Knepley     PetscScalar    *coords;
11800510c589SMatthew G. Knepley     PetscInt        coordSize, v;
11810510c589SMatthew G. Knepley     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
11820510c589SMatthew G. Knepley     const PetscReal ds2 = dis/2.0;
11830510c589SMatthew G. Knepley 
11840510c589SMatthew G. Knepley     /* Build coordinates */
11850510c589SMatthew G. Knepley     ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
11860510c589SMatthew G. Knepley     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
11870510c589SMatthew G. Knepley     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
11880510c589SMatthew G. Knepley     ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr);
11890510c589SMatthew G. Knepley     for (v = numCells; v < numCells+numVertices; ++v) {
11900510c589SMatthew G. Knepley       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
11910510c589SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
11920510c589SMatthew G. Knepley     }
11930510c589SMatthew G. Knepley     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
11940510c589SMatthew G. Knepley     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
11950510c589SMatthew G. Knepley     ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
11960510c589SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
11970510c589SMatthew G. Knepley     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
11980510c589SMatthew G. Knepley     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
11990510c589SMatthew G. Knepley     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
12000510c589SMatthew G. Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
12010510c589SMatthew G. Knepley     coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
12020510c589SMatthew G. Knepley     coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
12030510c589SMatthew G. Knepley     coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
12040510c589SMatthew G. Knepley     coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
12050510c589SMatthew G. Knepley     coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
12060510c589SMatthew G. Knepley     coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
12070510c589SMatthew G. Knepley     coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
12080510c589SMatthew G. Knepley     coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
12090510c589SMatthew G. Knepley     coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
12100510c589SMatthew G. Knepley     coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
12110510c589SMatthew G. Knepley     coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
12120510c589SMatthew G. Knepley     coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
12130510c589SMatthew G. Knepley     coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
12140510c589SMatthew G. Knepley     coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
12150510c589SMatthew G. Knepley     coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
12160510c589SMatthew G. Knepley     coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1217ae8bcbbbSMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1218ae8bcbbbSMatthew G. Knepley       /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1219ae8bcbbbSMatthew G. Knepley       /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1220ae8bcbbbSMatthew G. Knepley       /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1221ae8bcbbbSMatthew G. Knepley       /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1222ae8bcbbbSMatthew G. Knepley       /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1223ae8bcbbbSMatthew G. Knepley       /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1224ae8bcbbbSMatthew G. Knepley       /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1225ae8bcbbbSMatthew G. Knepley       /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1226ae8bcbbbSMatthew G. Knepley     }
12270510c589SMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
12280510c589SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
12290510c589SMatthew G. Knepley     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
12300510c589SMatthew G. Knepley   }
1231006a8963SMatthew G. Knepley   /* Create periodicity */
1232006a8963SMatthew G. Knepley   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1233006a8963SMatthew G. Knepley     PetscReal      L[3];
1234006a8963SMatthew G. Knepley     PetscReal      maxCell[3];
1235006a8963SMatthew G. Knepley     DMBoundaryType bdType[3];
1236006a8963SMatthew G. Knepley     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1237ae8bcbbbSMatthew G. Knepley     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1238ae8bcbbbSMatthew G. Knepley     PetscInt       i, numZCells = 3;
1239006a8963SMatthew G. Knepley 
1240006a8963SMatthew G. Knepley     bdType[0] = DM_BOUNDARY_NONE;
1241006a8963SMatthew G. Knepley     bdType[1] = DM_BOUNDARY_NONE;
1242006a8963SMatthew G. Knepley     bdType[2] = periodicZ;
1243006a8963SMatthew G. Knepley     for (i = 0; i < dim; i++) {
1244006a8963SMatthew G. Knepley       L[i]       = upper[i] - lower[i];
1245006a8963SMatthew G. Knepley       maxCell[i] = 1.1 * (L[i] / numZCells);
1246006a8963SMatthew G. Knepley     }
1247006a8963SMatthew G. Knepley     ierr = DMSetPeriodicity(*dm, maxCell, L, bdType);CHKERRQ(ierr);
1248006a8963SMatthew G. Knepley   }
1249dbc1dc17SMatthew G. Knepley   /* Refine topology */
1250dbc1dc17SMatthew G. Knepley   for (r = 0; r < numRefine; ++r) {
1251dbc1dc17SMatthew G. Knepley     DM rdm = NULL;
1252dbc1dc17SMatthew G. Knepley 
1253dbc1dc17SMatthew G. Knepley     ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr);
1254dbc1dc17SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
1255dbc1dc17SMatthew G. Knepley     *dm  = rdm;
1256dbc1dc17SMatthew G. Knepley   }
1257dbc1dc17SMatthew G. Knepley   /* Remap geometry to cylinder
1258dbc1dc17SMatthew G. Knepley        Interior square: Linear interpolation is correct
1259dbc1dc17SMatthew G. Knepley        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1260dbc1dc17SMatthew G. Knepley        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance
1261dbc1dc17SMatthew G. Knepley 
1262dbc1dc17SMatthew G. Knepley          phi     = arctan(y/x)
1263dbc1dc17SMatthew G. Knepley          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1264dbc1dc17SMatthew G. Knepley          d_far   = sqrt(1/2 + sin^2(phi))
1265dbc1dc17SMatthew G. Knepley 
1266dbc1dc17SMatthew G. Knepley        so we remap them using
1267dbc1dc17SMatthew G. Knepley 
1268dbc1dc17SMatthew G. Knepley          x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close)
1269dbc1dc17SMatthew G. Knepley          y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close)
1270dbc1dc17SMatthew G. Knepley 
1271dbc1dc17SMatthew G. Knepley        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1272dbc1dc17SMatthew G. Knepley   */
1273dbc1dc17SMatthew G. Knepley   {
1274dbc1dc17SMatthew G. Knepley     Vec           coordinates;
1275dbc1dc17SMatthew G. Knepley     PetscSection  coordSection;
1276dbc1dc17SMatthew G. Knepley     PetscScalar  *coords;
1277dbc1dc17SMatthew G. Knepley     PetscInt      vStart, vEnd, v;
1278dbc1dc17SMatthew G. Knepley     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1279dbc1dc17SMatthew G. Knepley     const PetscReal ds2 = 0.5*dis;
1280dbc1dc17SMatthew G. Knepley 
1281dbc1dc17SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1282dbc1dc17SMatthew G. Knepley     ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1283dbc1dc17SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr);
1284dbc1dc17SMatthew G. Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1285dbc1dc17SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
128671752167SMatthew G. Knepley       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1287dbc1dc17SMatthew G. Knepley       PetscInt  off;
1288dbc1dc17SMatthew G. Knepley 
1289dbc1dc17SMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
129071752167SMatthew G. Knepley       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
129171752167SMatthew G. Knepley       x    = PetscRealPart(coords[off]);
129271752167SMatthew G. Knepley       y    = PetscRealPart(coords[off+1]);
129371752167SMatthew G. Knepley       phi  = PetscAtan2Real(y, x);
1294dbc1dc17SMatthew G. Knepley       sinp = PetscSinReal(phi);
1295dbc1dc17SMatthew G. Knepley       cosp = PetscCosReal(phi);
1296dbc1dc17SMatthew G. Knepley       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
129771752167SMatthew G. Knepley         dc = PetscAbsReal(ds2/sinp);
129871752167SMatthew G. Knepley         df = PetscAbsReal(dis/sinp);
129971752167SMatthew G. Knepley         xc = ds2*x/PetscAbsScalar(y);
130071752167SMatthew G. Knepley         yc = ds2*PetscSignReal(y);
1301dbc1dc17SMatthew G. Knepley       } else {
130271752167SMatthew G. Knepley         dc = PetscAbsReal(ds2/cosp);
130371752167SMatthew G. Knepley         df = PetscAbsReal(dis/cosp);
130471752167SMatthew G. Knepley         xc = ds2*PetscSignReal(x);
130571752167SMatthew G. Knepley         yc = ds2*y/PetscAbsScalar(x);
1306dbc1dc17SMatthew G. Knepley       }
1307dbc1dc17SMatthew G. Knepley       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1308dbc1dc17SMatthew G. Knepley       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1309dbc1dc17SMatthew G. Knepley     }
1310dbc1dc17SMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
13110510c589SMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1312006a8963SMatthew G. Knepley       ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
13130510c589SMatthew G. Knepley     }
13140510c589SMatthew G. Knepley   }
13150510c589SMatthew G. Knepley   PetscFunctionReturn(0);
13160510c589SMatthew G. Knepley }
13170510c589SMatthew G. Knepley 
1318552f7358SJed Brown /* External function declarations here */
1319552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
13206dbf9973SLawrence Mitchell extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
1321fd59a867SMatthew G. Knepley extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
132266ad2231SToby Isaac extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
1323b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
1324552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
1325d9879d6bSMatthew G. Knepley PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
1326552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
1327552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
1328552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
13292c40f234SMatthew G. Knepley extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
1330552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
1331552f7358SJed Brown 
13320a6ba040SMatthew G. Knepley /* Replace dm with the contents of dmNew
13330a6ba040SMatthew G. Knepley    - Share the DM_Plex structure
13340a6ba040SMatthew G. Knepley    - Share the coordinates
1335d7973521SMatthew G. Knepley    - Share the SF
13360a6ba040SMatthew G. Knepley */
13370a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
13380a6ba040SMatthew G. Knepley {
1339d7973521SMatthew G. Knepley   PetscSF          sf;
1340acdc6f61SToby Isaac   DM               coordDM, coarseDM;
13410a6ba040SMatthew G. Knepley   Vec              coords;
134255fbe3e3SMatthew G. Knepley   const PetscReal *maxCell, *L;
13435dc8c3f7SMatthew G. Knepley   const DMBoundaryType *bd;
13440a6ba040SMatthew G. Knepley   PetscErrorCode   ierr;
13450a6ba040SMatthew G. Knepley 
13460a6ba040SMatthew G. Knepley   PetscFunctionBegin;
1347d7973521SMatthew G. Knepley   ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr);
1348d7973521SMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
134955fbe3e3SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr);
13500a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
135155fbe3e3SMatthew G. Knepley   ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr);
13520a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
13535dc8c3f7SMatthew G. Knepley   ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr);
13545dc8c3f7SMatthew G. Knepley   if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L, bd);CHKERRQ(ierr);}
13550a6ba040SMatthew G. Knepley   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
13560a6ba040SMatthew G. Knepley   dm->data = dmNew->data;
13570a6ba040SMatthew G. Knepley   ((DM_Plex *) dmNew->data)->refct++;
1358c58f1c22SToby Isaac   dmNew->labels->refct++;
1359c58f1c22SToby Isaac   if (!--(dm->labels->refct)) {
1360c58f1c22SToby Isaac     DMLabelLink next = dm->labels->next;
1361c58f1c22SToby Isaac 
1362c58f1c22SToby Isaac     /* destroy the labels */
1363c58f1c22SToby Isaac     while (next) {
1364c58f1c22SToby Isaac       DMLabelLink tmp = next->next;
1365c58f1c22SToby Isaac 
1366c58f1c22SToby Isaac       ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr);
1367c58f1c22SToby Isaac       ierr = PetscFree(next);CHKERRQ(ierr);
1368c58f1c22SToby Isaac       next = tmp;
1369c58f1c22SToby Isaac     }
1370c58f1c22SToby Isaac     ierr = PetscFree(dm->labels);CHKERRQ(ierr);
1371c58f1c22SToby Isaac   }
1372c58f1c22SToby Isaac   dm->labels = dmNew->labels;
1373c58f1c22SToby Isaac   dm->depthLabel = dmNew->depthLabel;
1374acdc6f61SToby Isaac   ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr);
1375acdc6f61SToby Isaac   ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr);
13760a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
13770a6ba040SMatthew G. Knepley }
13780a6ba040SMatthew G. Knepley 
13790a6ba040SMatthew G. Knepley /* Swap dm with the contents of dmNew
13800a6ba040SMatthew G. Knepley    - Swap the DM_Plex structure
13810a6ba040SMatthew G. Knepley    - Swap the coordinates
138279a015ccSMatthew G. Knepley    - Swap the point PetscSF
13830a6ba040SMatthew G. Knepley */
13840a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
13850a6ba040SMatthew G. Knepley {
13860a6ba040SMatthew G. Knepley   DM              coordDMA, coordDMB;
13870a6ba040SMatthew G. Knepley   Vec             coordsA,  coordsB;
138879a015ccSMatthew G. Knepley   PetscSF         sfA,      sfB;
13890a6ba040SMatthew G. Knepley   void            *tmp;
1390c58f1c22SToby Isaac   DMLabelLinkList listTmp;
1391c58f1c22SToby Isaac   DMLabel         depthTmp;
13920a6ba040SMatthew G. Knepley   PetscErrorCode  ierr;
13930a6ba040SMatthew G. Knepley 
13940a6ba040SMatthew G. Knepley   PetscFunctionBegin;
139579a015ccSMatthew G. Knepley   ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr);
139679a015ccSMatthew G. Knepley   ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr);
139779a015ccSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr);
139879a015ccSMatthew G. Knepley   ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr);
139979a015ccSMatthew G. Knepley   ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr);
140079a015ccSMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr);
140179a015ccSMatthew G. Knepley 
14020a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
14030a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
14040a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
14050a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
14060a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
14070a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
14080a6ba040SMatthew G. Knepley 
14090a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
14100a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
14110a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
14120a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
14130a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
14140a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
1415acdc6f61SToby Isaac 
14160a6ba040SMatthew G. Knepley   tmp       = dmA->data;
14170a6ba040SMatthew G. Knepley   dmA->data = dmB->data;
14180a6ba040SMatthew G. Knepley   dmB->data = tmp;
1419c58f1c22SToby Isaac   listTmp   = dmA->labels;
1420c58f1c22SToby Isaac   dmA->labels = dmB->labels;
1421c58f1c22SToby Isaac   dmB->labels = listTmp;
1422c58f1c22SToby Isaac   depthTmp  = dmA->depthLabel;
1423c58f1c22SToby Isaac   dmA->depthLabel = dmB->depthLabel;
1424c58f1c22SToby Isaac   dmB->depthLabel = depthTmp;
14250a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
14260a6ba040SMatthew G. Knepley }
14270a6ba040SMatthew G. Knepley 
14284416b707SBarry Smith PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
14290a6ba040SMatthew G. Knepley {
14300a6ba040SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
14310a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
14320a6ba040SMatthew G. Knepley 
14330a6ba040SMatthew G. Knepley   PetscFunctionBegin;
14340a6ba040SMatthew G. Knepley   /* Handle viewing */
14350a6ba040SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
14360a6ba040SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
143794ae4db5SBarry Smith   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr);
1438953fc75cSMatthew G. Knepley   /* Point Location */
1439953fc75cSMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr);
14402e62ab5aSMatthew G. Knepley   /* Generation and remeshing */
14412e62ab5aSMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr);
1442b29cfa1cSToby Isaac   /* Projection behavior */
1443b29cfa1cSToby Isaac   ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr);
144464141f95SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr);
14454f3833eaSMatthew G. Knepley 
14464f3833eaSMatthew G. Knepley   ierr = PetscPartitionerSetFromOptions(mesh->partitioner);CHKERRQ(ierr);
144768d4fef7SMatthew G. Knepley   PetscFunctionReturn(0);
144868d4fef7SMatthew G. Knepley }
144968d4fef7SMatthew G. Knepley 
145046fa42a0SMatthew G. Knepley static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
145168d4fef7SMatthew G. Knepley {
1452b653a561SMatthew G. Knepley   PetscInt       refine = 0, coarsen = 0, r;
1453b653a561SMatthew G. Knepley   PetscBool      isHierarchy;
145468d4fef7SMatthew G. Knepley   PetscErrorCode ierr;
145568d4fef7SMatthew G. Knepley 
145668d4fef7SMatthew G. Knepley   PetscFunctionBegin;
145768d4fef7SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14581a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr);
145968d4fef7SMatthew G. Knepley   /* Handle DMPlex refinement */
146068d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
146168d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
1462b6a0289aSMatthew G. Knepley   if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);}
146368d4fef7SMatthew G. Knepley   if (refine && isHierarchy) {
1464acdc6f61SToby Isaac     DM *dms, coarseDM;
146568d4fef7SMatthew G. Knepley 
1466acdc6f61SToby Isaac     ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr);
1467acdc6f61SToby Isaac     ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
146868d4fef7SMatthew G. Knepley     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
146968d4fef7SMatthew G. Knepley     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
147068d4fef7SMatthew G. Knepley     /* Total hack since we do not pass in a pointer */
147168d4fef7SMatthew G. Knepley     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
147268d4fef7SMatthew G. Knepley     if (refine == 1) {
1473a8fb8f29SToby Isaac       ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
14740aef6b92SMatthew G. Knepley       ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr);
147568d4fef7SMatthew G. Knepley     } else {
1476a8fb8f29SToby Isaac       ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
14770aef6b92SMatthew G. Knepley       ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr);
1478a8fb8f29SToby Isaac       ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
14790aef6b92SMatthew G. Knepley       ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr);
148068d4fef7SMatthew G. Knepley     }
1481acdc6f61SToby Isaac     ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr);
1482acdc6f61SToby Isaac     ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr);
148368d4fef7SMatthew G. Knepley     /* Free DMs */
148468d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
1485547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr);
148668d4fef7SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
148768d4fef7SMatthew G. Knepley     }
148868d4fef7SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
148968d4fef7SMatthew G. Knepley   } else {
149068d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
149168d4fef7SMatthew G. Knepley       DM refinedMesh;
149268d4fef7SMatthew G. Knepley 
14931a1499c8SBarry Smith       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
149468d4fef7SMatthew G. Knepley       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
149568d4fef7SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
149668d4fef7SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
1497547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
149868d4fef7SMatthew G. Knepley       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
149968d4fef7SMatthew G. Knepley     }
150068d4fef7SMatthew G. Knepley   }
15013cf6fe12SMatthew G. Knepley   /* Handle DMPlex coarsening */
1502b653a561SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr);
1503b653a561SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr);
1504b653a561SMatthew G. Knepley   if (coarsen && isHierarchy) {
1505b653a561SMatthew G. Knepley     DM *dms;
1506b653a561SMatthew G. Knepley 
1507b653a561SMatthew G. Knepley     ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr);
1508b653a561SMatthew G. Knepley     ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr);
1509b653a561SMatthew G. Knepley     /* Free DMs */
1510b653a561SMatthew G. Knepley     for (r = 0; r < coarsen; ++r) {
1511547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr);
1512b653a561SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
1513b653a561SMatthew G. Knepley     }
1514b653a561SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
1515b653a561SMatthew G. Knepley   } else {
1516b653a561SMatthew G. Knepley     for (r = 0; r < coarsen; ++r) {
15173cf6fe12SMatthew G. Knepley       DM coarseMesh;
15183cf6fe12SMatthew G. Knepley 
15193cf6fe12SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
15203cf6fe12SMatthew G. Knepley       ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr);
15213cf6fe12SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
15223cf6fe12SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr);
1523547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
15243cf6fe12SMatthew G. Knepley       ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr);
15253cf6fe12SMatthew G. Knepley     }
1526b653a561SMatthew G. Knepley   }
15273cf6fe12SMatthew G. Knepley   /* Handle */
15281a1499c8SBarry Smith   ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
15290a6ba040SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
15300a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
15310a6ba040SMatthew G. Knepley }
15320a6ba040SMatthew G. Knepley 
1533552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
1534552f7358SJed Brown {
1535552f7358SJed Brown   PetscErrorCode ierr;
1536552f7358SJed Brown 
1537552f7358SJed Brown   PetscFunctionBegin;
1538552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
1539552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
1540552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr);
1541d930f514SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr);
15422c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr);
1543d930f514SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr);
1544552f7358SJed Brown   PetscFunctionReturn(0);
1545552f7358SJed Brown }
1546552f7358SJed Brown 
1547552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
1548552f7358SJed Brown {
1549552f7358SJed Brown   PetscErrorCode ierr;
1550552f7358SJed Brown 
1551552f7358SJed Brown   PetscFunctionBegin;
1552552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
1553552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
15542c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr);
1555552f7358SJed Brown   PetscFunctionReturn(0);
1556552f7358SJed Brown }
1557552f7358SJed Brown 
1558793f3fe5SMatthew G. Knepley static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
1559793f3fe5SMatthew G. Knepley {
1560793f3fe5SMatthew G. Knepley   PetscInt       depth, d;
1561793f3fe5SMatthew G. Knepley   PetscErrorCode ierr;
1562793f3fe5SMatthew G. Knepley 
1563793f3fe5SMatthew G. Knepley   PetscFunctionBegin;
1564793f3fe5SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1565793f3fe5SMatthew G. Knepley   if (depth == 1) {
1566793f3fe5SMatthew G. Knepley     ierr = DMGetDimension(dm, &d);CHKERRQ(ierr);
1567793f3fe5SMatthew G. Knepley     if (dim == 0)      {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);}
1568793f3fe5SMatthew G. Knepley     else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);}
1569793f3fe5SMatthew G. Knepley     else               {*pStart = 0; *pEnd = 0;}
1570793f3fe5SMatthew G. Knepley   } else {
1571793f3fe5SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);
1572793f3fe5SMatthew G. Knepley   }
1573793f3fe5SMatthew G. Knepley   PetscFunctionReturn(0);
1574793f3fe5SMatthew G. Knepley }
1575793f3fe5SMatthew G. Knepley 
15762f356facSMatthew G. Knepley static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
1577502a2867SDave May {
1578502a2867SDave May   PetscSF        sf;
15792f356facSMatthew G. Knepley   PetscErrorCode ierr;
1580502a2867SDave May 
15812f356facSMatthew G. Knepley   PetscFunctionBegin;
1582502a2867SDave May   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
1583502a2867SDave May   ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr);
1584502a2867SDave May   PetscFunctionReturn(0);
1585502a2867SDave May }
1586502a2867SDave May 
158746fa42a0SMatthew G. Knepley static PetscErrorCode DMInitialize_Plex(DM dm)
1588552f7358SJed Brown {
1589713918a9SToby Isaac   PetscErrorCode ierr;
1590713918a9SToby Isaac 
1591552f7358SJed Brown   PetscFunctionBegin;
1592552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
15932c40f234SMatthew G. Knepley   dm->ops->load                            = DMLoad_Plex;
1594552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
159538221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
1596552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
1597fd59a867SMatthew G. Knepley   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
159866ad2231SToby Isaac   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
1599552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
1600552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
1601184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
16020298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
1603552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
16040a6ba040SMatthew G. Knepley   dm->ops->getcoloring                     = NULL;
1605552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
1606bceba477SMatthew G. Knepley   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
1607bceba477SMatthew G. Knepley   dm->ops->getaggregates                   = NULL;
1608bceba477SMatthew G. Knepley   dm->ops->getinjection                    = DMCreateInjection_Plex;
1609552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
16100a6ba040SMatthew G. Knepley   dm->ops->coarsen                         = DMCoarsen_Plex;
16110a6ba040SMatthew G. Knepley   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
1612b653a561SMatthew G. Knepley   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
16130298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
16140298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
16150298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
16160298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
1617552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
1618552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
1619793f3fe5SMatthew G. Knepley   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
1620552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
16210709b2feSToby Isaac   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
16220709b2feSToby Isaac   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
1623bfc4295aSToby Isaac   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
16248c6c5593SMatthew G. Knepley   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
16250709b2feSToby Isaac   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
1626b698f381SToby Isaac   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
16272a16baeaSToby Isaac   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
1628502a2867SDave May   dm->ops->getneighbors                    = DMGetNeighors_Plex;
1629713918a9SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Plex);CHKERRQ(ierr);
1630552f7358SJed Brown   PetscFunctionReturn(0);
1631552f7358SJed Brown }
1632552f7358SJed Brown 
163346fa42a0SMatthew G. Knepley PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
163463a16f15SMatthew G. Knepley {
163563a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
163663a16f15SMatthew G. Knepley   PetscErrorCode ierr;
163763a16f15SMatthew G. Knepley 
163863a16f15SMatthew G. Knepley   PetscFunctionBegin;
163963a16f15SMatthew G. Knepley   mesh->refct++;
164063a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
164163a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
164263a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
164363a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
164463a16f15SMatthew G. Knepley }
164563a16f15SMatthew G. Knepley 
16468818961aSMatthew G Knepley /*MC
16478818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
16488818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
16498818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
16508818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
16518818961aSMatthew G Knepley 
16528818961aSMatthew G Knepley   Level: intermediate
16538818961aSMatthew G Knepley 
16548818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
16558818961aSMatthew G Knepley M*/
16568818961aSMatthew G Knepley 
16578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1658552f7358SJed Brown {
1659552f7358SJed Brown   DM_Plex        *mesh;
1660770b213bSMatthew G Knepley   PetscInt       unit, d;
1661552f7358SJed Brown   PetscErrorCode ierr;
1662552f7358SJed Brown 
1663552f7358SJed Brown   PetscFunctionBegin;
1664552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1665b00a9115SJed Brown   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
1666c73cfb54SMatthew G. Knepley   dm->dim  = 0;
1667552f7358SJed Brown   dm->data = mesh;
1668552f7358SJed Brown 
1669552f7358SJed Brown   mesh->refct             = 1;
167082f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
1671552f7358SJed Brown   mesh->maxConeSize       = 0;
16720298fd71SBarry Smith   mesh->cones             = NULL;
16730298fd71SBarry Smith   mesh->coneOrientations  = NULL;
167482f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
1675552f7358SJed Brown   mesh->maxSupportSize    = 0;
16760298fd71SBarry Smith   mesh->supports          = NULL;
1677552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
1678552f7358SJed Brown   mesh->refinementLimit   = -1.0;
1679552f7358SJed Brown 
16800298fd71SBarry Smith   mesh->facesTmp = NULL;
1681552f7358SJed Brown 
1682d9deefdfSMatthew G. Knepley   mesh->tetgenOpts   = NULL;
1683d9deefdfSMatthew G. Knepley   mesh->triangleOpts = NULL;
168477623264SMatthew G. Knepley   ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr);
16852e62ab5aSMatthew G. Knepley   mesh->remeshBd     = PETSC_FALSE;
1686d9deefdfSMatthew G. Knepley 
16870298fd71SBarry Smith   mesh->subpointMap = NULL;
1688552f7358SJed Brown 
16898865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
1690552f7358SJed Brown 
16910aef6b92SMatthew G. Knepley   mesh->regularRefinement   = PETSC_FALSE;
1692df0420ecSMatthew G. Knepley   mesh->depthState          = -1;
16930298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
16940298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
1695a68b90caSToby Isaac   mesh->anchorSection       = NULL;
1696a68b90caSToby Isaac   mesh->anchorIS            = NULL;
169741e6d900SToby Isaac   mesh->createanchors       = NULL;
1698fa73a4e1SToby Isaac   mesh->computeanchormatrix = NULL;
1699d961a43aSToby Isaac   mesh->parentSection       = NULL;
1700d961a43aSToby Isaac   mesh->parents             = NULL;
1701d961a43aSToby Isaac   mesh->childIDs            = NULL;
1702d961a43aSToby Isaac   mesh->childSection        = NULL;
1703d961a43aSToby Isaac   mesh->children            = NULL;
1704d6a7ad0dSToby Isaac   mesh->referenceTree       = NULL;
1705dcbd3bf7SToby Isaac   mesh->getchildsymmetry    = NULL;
17068865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1707552f7358SJed Brown   mesh->vtkCellHeight       = 0;
170870034214SMatthew G. Knepley   mesh->useCone             = PETSC_FALSE;
170970034214SMatthew G. Knepley   mesh->useClosure          = PETSC_TRUE;
1710e228b242SToby Isaac   mesh->useAnchors          = PETSC_FALSE;
1711552f7358SJed Brown 
1712b29cfa1cSToby Isaac   mesh->maxProjectionHeight = 0;
1713b29cfa1cSToby Isaac 
1714552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
1715552f7358SJed Brown   mesh->printFEM       = 0;
17166113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
1717552f7358SJed Brown 
1718552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
1719552f7358SJed Brown   PetscFunctionReturn(0);
1720552f7358SJed Brown }
1721552f7358SJed Brown 
1722552f7358SJed Brown /*@
1723552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1724552f7358SJed Brown 
1725552f7358SJed Brown   Collective on MPI_Comm
1726552f7358SJed Brown 
1727552f7358SJed Brown   Input Parameter:
1728552f7358SJed Brown . comm - The communicator for the DMPlex object
1729552f7358SJed Brown 
1730552f7358SJed Brown   Output Parameter:
1731552f7358SJed Brown . mesh  - The DMPlex object
1732552f7358SJed Brown 
1733552f7358SJed Brown   Level: beginner
1734552f7358SJed Brown 
1735552f7358SJed Brown .keywords: DMPlex, create
1736552f7358SJed Brown @*/
1737552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1738552f7358SJed Brown {
1739552f7358SJed Brown   PetscErrorCode ierr;
1740552f7358SJed Brown 
1741552f7358SJed Brown   PetscFunctionBegin;
1742552f7358SJed Brown   PetscValidPointer(mesh,2);
1743552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1744552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1745552f7358SJed Brown   PetscFunctionReturn(0);
1746552f7358SJed Brown }
1747552f7358SJed Brown 
1748a47d0d45SMatthew G. Knepley /*
1749a47d0d45SMatthew G. Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell, but vertex numbers are global and an SF is built for them
1750a47d0d45SMatthew G. Knepley */
1751a47d0d45SMatthew G. Knepley static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert)
1752a47d0d45SMatthew G. Knepley {
1753a47d0d45SMatthew G. Knepley   PetscSF         sfPoint;
1754a47d0d45SMatthew G. Knepley   PetscLayout     vLayout;
1755a47d0d45SMatthew G. Knepley   PetscHashI      vhash;
1756a47d0d45SMatthew G. Knepley   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
1757a47d0d45SMatthew G. Knepley   const PetscInt *vrange;
1758a47d0d45SMatthew G. Knepley   PetscInt        numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
1759a47d0d45SMatthew G. Knepley   PETSC_UNUSED PetscHashIIter ret, iter;
17609852e123SBarry Smith   PetscMPIInt     rank, size;
1761a47d0d45SMatthew G. Knepley   PetscErrorCode  ierr;
1762a47d0d45SMatthew G. Knepley 
1763a47d0d45SMatthew G. Knepley   PetscFunctionBegin;
1764a47d0d45SMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
17659852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
1766a47d0d45SMatthew G. Knepley   /* Partition vertices */
1767a47d0d45SMatthew G. Knepley   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr);
1768a47d0d45SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr);
1769a47d0d45SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
1770a47d0d45SMatthew G. Knepley   ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
1771a47d0d45SMatthew G. Knepley   ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr);
1772a47d0d45SMatthew G. Knepley   /* Count vertices and map them to procs */
1773a47d0d45SMatthew G. Knepley   PetscHashICreate(vhash);
1774a47d0d45SMatthew G. Knepley   for (c = 0; c < numCells; ++c) {
1775a47d0d45SMatthew G. Knepley     for (p = 0; p < numCorners; ++p) {
1776a47d0d45SMatthew G. Knepley       PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter);
1777a47d0d45SMatthew G. Knepley     }
1778a47d0d45SMatthew G. Knepley   }
1779a47d0d45SMatthew G. Knepley   PetscHashISize(vhash, numVerticesAdj);
1780a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr);
1781324a3438SMichael Lange   off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr);
1782a47d0d45SMatthew G. Knepley   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
17838cd1fd5cSMatthew G. Knepley   ierr = PetscSortInt(numVerticesAdj, verticesAdj);CHKERRQ(ierr);
1784a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr);
1785a47d0d45SMatthew G. Knepley   for (v = 0; v < numVerticesAdj; ++v) {
1786a47d0d45SMatthew G. Knepley     const PetscInt gv = verticesAdj[v];
1787a47d0d45SMatthew G. Knepley     PetscInt       vrank;
1788a47d0d45SMatthew G. Knepley 
17899852e123SBarry Smith     ierr = PetscFindInt(gv, size+1, vrange, &vrank);CHKERRQ(ierr);
1790a47d0d45SMatthew G. Knepley     vrank = vrank < 0 ? -(vrank+2) : vrank;
1791a47d0d45SMatthew G. Knepley     remoteVerticesAdj[v].index = gv - vrange[vrank];
1792a47d0d45SMatthew G. Knepley     remoteVerticesAdj[v].rank  = vrank;
1793a47d0d45SMatthew G. Knepley   }
1794a47d0d45SMatthew G. Knepley   /* Create cones */
1795a47d0d45SMatthew G. Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr);
1796a47d0d45SMatthew G. Knepley   for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);}
1797a47d0d45SMatthew G. Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
1798a47d0d45SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1799a47d0d45SMatthew G. Knepley   for (c = 0; c < numCells; ++c) {
1800a47d0d45SMatthew G. Knepley     for (p = 0; p < numCorners; ++p) {
1801a47d0d45SMatthew G. Knepley       const PetscInt gv = cells[c*numCorners+p];
1802a47d0d45SMatthew G. Knepley       PetscInt       lv;
1803a47d0d45SMatthew G. Knepley 
1804a47d0d45SMatthew G. Knepley       ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr);
1805a47d0d45SMatthew G. Knepley       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
1806a47d0d45SMatthew G. Knepley       cone[p] = lv+numCells;
1807a47d0d45SMatthew G. Knepley     }
1808a47d0d45SMatthew G. Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
1809a47d0d45SMatthew G. Knepley   }
1810a47d0d45SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
1811a47d0d45SMatthew G. Knepley   /* Create SF for vertices */
1812a47d0d45SMatthew G. Knepley   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr);
1813a47d0d45SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr);
1814a47d0d45SMatthew G. Knepley   ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr);
1815a47d0d45SMatthew G. Knepley   ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr);
1816a47d0d45SMatthew G. Knepley   ierr = PetscFree(verticesAdj);CHKERRQ(ierr);
1817a47d0d45SMatthew G. Knepley   /* Build pointSF */
1818a47d0d45SMatthew G. Knepley   ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr);
1819a47d0d45SMatthew G. Knepley   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
1820a47d0d45SMatthew G. Knepley   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
1821a47d0d45SMatthew G. Knepley   ierr = PetscSFReduceBegin(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr);
1822a47d0d45SMatthew G. Knepley   ierr = PetscSFReduceEnd(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr);
1823a47d0d45SMatthew G. Knepley   for (v = 0; v < numVertices;    ++v) if (vertexOwner[v].rank < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global vertex %d on rank %d was unclaimed", v + vrange[rank], rank);
1824a47d0d45SMatthew G. Knepley   ierr = PetscSFBcastBegin(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr);
1825a47d0d45SMatthew G. Knepley   ierr = PetscSFBcastEnd(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr);
1826a47d0d45SMatthew G. Knepley   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
1827a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr);
1828a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr);
1829a47d0d45SMatthew G. Knepley   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
1830a47d0d45SMatthew G. Knepley     if (vertexLocal[v].rank != rank) {
1831a47d0d45SMatthew G. Knepley       localVertex[g]        = v+numCells;
1832a47d0d45SMatthew G. Knepley       remoteVertex[g].index = vertexLocal[v].index;
1833a47d0d45SMatthew G. Knepley       remoteVertex[g].rank  = vertexLocal[v].rank;
1834a47d0d45SMatthew G. Knepley       ++g;
1835a47d0d45SMatthew G. Knepley     }
1836a47d0d45SMatthew G. Knepley   }
1837a47d0d45SMatthew G. Knepley   ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr);
1838a47d0d45SMatthew G. Knepley   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
1839a47d0d45SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1840a47d0d45SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr);
1841a47d0d45SMatthew G. Knepley   ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr);
1842a47d0d45SMatthew G. Knepley   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
1843a47d0d45SMatthew G. Knepley   PetscHashIDestroy(vhash);
1844a47d0d45SMatthew G. Knepley   /* Fill in the rest of the topology structure */
1845a47d0d45SMatthew G. Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
1846a47d0d45SMatthew G. Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
1847a47d0d45SMatthew G. Knepley   PetscFunctionReturn(0);
1848a47d0d45SMatthew G. Knepley }
1849a47d0d45SMatthew G. Knepley 
1850a47d0d45SMatthew G. Knepley /*
1851a47d0d45SMatthew G. Knepley   This takes as input the coordinates for each owned vertex
1852a47d0d45SMatthew G. Knepley */
185321016a8bSBarry Smith static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[])
1854a47d0d45SMatthew G. Knepley {
1855a47d0d45SMatthew G. Knepley   PetscSection   coordSection;
1856a47d0d45SMatthew G. Knepley   Vec            coordinates;
1857a47d0d45SMatthew G. Knepley   PetscScalar   *coords;
1858c5e08f72SMatthew G. Knepley   PetscInt       numVertices, numVerticesAdj, coordSize, v;
1859a47d0d45SMatthew G. Knepley   PetscErrorCode ierr;
1860a47d0d45SMatthew G. Knepley 
1861a47d0d45SMatthew G. Knepley   PetscFunctionBegin;
1862a47d0d45SMatthew G. Knepley   ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr);
1863a47d0d45SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1864a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1865a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
1866a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr);
1867a47d0d45SMatthew G. Knepley   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
1868a47d0d45SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
1869a47d0d45SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
1870a47d0d45SMatthew G. Knepley   }
1871a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1872a47d0d45SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
1873a47d0d45SMatthew G. Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
1874a47d0d45SMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
1875a47d0d45SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1876a47d0d45SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1877a47d0d45SMatthew G. Knepley   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1878a47d0d45SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1879a47d0d45SMatthew G. Knepley   {
1880a47d0d45SMatthew G. Knepley     MPI_Datatype coordtype;
1881a47d0d45SMatthew G. Knepley 
1882a47d0d45SMatthew G. Knepley     /* Need a temp buffer for coords if we have complex/single */
188321016a8bSBarry Smith     ierr = MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);CHKERRQ(ierr);
1884a47d0d45SMatthew G. Knepley     ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr);
188521016a8bSBarry Smith #if defined(PETSC_USE_COMPLEX)
188621016a8bSBarry Smith     {
188721016a8bSBarry Smith     PetscScalar *svertexCoords;
188821016a8bSBarry Smith     PetscInt    i;
188921016a8bSBarry Smith     ierr = PetscMalloc1(numV*spaceDim,&svertexCoords);CHKERRQ(ierr);
189021016a8bSBarry Smith     for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
189121016a8bSBarry Smith     ierr = PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr);
189221016a8bSBarry Smith     ierr = PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr);
189321016a8bSBarry Smith     ierr = PetscFree(svertexCoords);CHKERRQ(ierr);
189421016a8bSBarry Smith     }
189521016a8bSBarry Smith #else
1896a47d0d45SMatthew G. Knepley     ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr);
1897a47d0d45SMatthew G. Knepley     ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr);
189821016a8bSBarry Smith #endif
1899a47d0d45SMatthew G. Knepley     ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr);
1900a47d0d45SMatthew G. Knepley   }
1901a47d0d45SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1902a47d0d45SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1903a47d0d45SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1904a47d0d45SMatthew G. Knepley   PetscFunctionReturn(0);
1905a47d0d45SMatthew G. Knepley }
1906a47d0d45SMatthew G. Knepley 
1907a47d0d45SMatthew G. Knepley /*@C
1908a47d0d45SMatthew G. Knepley   DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
1909a47d0d45SMatthew G. Knepley 
1910a47d0d45SMatthew G. Knepley   Input Parameters:
1911a47d0d45SMatthew G. Knepley + comm - The communicator
1912a47d0d45SMatthew G. Knepley . dim - The topological dimension of the mesh
1913a47d0d45SMatthew G. Knepley . numCells - The number of cells owned by this process
1914a47d0d45SMatthew G. Knepley . numVertices - The number of vertices owned by this process
1915a47d0d45SMatthew G. Knepley . numCorners - The number of vertices for each cell
1916a47d0d45SMatthew G. Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1917a47d0d45SMatthew G. Knepley . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
1918a47d0d45SMatthew G. Knepley . spaceDim - The spatial dimension used for coordinates
1919a47d0d45SMatthew G. Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
1920a47d0d45SMatthew G. Knepley 
1921a47d0d45SMatthew G. Knepley   Output Parameter:
192218d54ad4SMichael Lange + dm - The DM
192318d54ad4SMichael Lange - vertexSF - Optional, SF describing complete vertex ownership
1924a47d0d45SMatthew G. Knepley 
1925a47d0d45SMatthew G. Knepley   Note: Two triangles sharing a face
1926a47d0d45SMatthew G. Knepley $
1927a47d0d45SMatthew G. Knepley $        2
1928a47d0d45SMatthew G. Knepley $      / | \
1929a47d0d45SMatthew G. Knepley $     /  |  \
1930a47d0d45SMatthew G. Knepley $    /   |   \
1931a47d0d45SMatthew G. Knepley $   0  0 | 1  3
1932a47d0d45SMatthew G. Knepley $    \   |   /
1933a47d0d45SMatthew G. Knepley $     \  |  /
1934a47d0d45SMatthew G. Knepley $      \ | /
1935a47d0d45SMatthew G. Knepley $        1
1936a47d0d45SMatthew G. Knepley would have input
1937a47d0d45SMatthew G. Knepley $  numCells = 2, numVertices = 4
1938a47d0d45SMatthew G. Knepley $  cells = [0 1 2  1 3 2]
1939a47d0d45SMatthew G. Knepley $
1940a47d0d45SMatthew G. Knepley which would result in the DMPlex
1941a47d0d45SMatthew G. Knepley $
1942a47d0d45SMatthew G. Knepley $        4
1943a47d0d45SMatthew G. Knepley $      / | \
1944a47d0d45SMatthew G. Knepley $     /  |  \
1945a47d0d45SMatthew G. Knepley $    /   |   \
1946a47d0d45SMatthew G. Knepley $   2  0 | 1  5
1947a47d0d45SMatthew G. Knepley $    \   |   /
1948a47d0d45SMatthew G. Knepley $     \  |  /
1949a47d0d45SMatthew G. Knepley $      \ | /
1950a47d0d45SMatthew G. Knepley $        3
1951a47d0d45SMatthew G. Knepley 
1952a47d0d45SMatthew G. Knepley   Level: beginner
1953a47d0d45SMatthew G. Knepley 
1954a47d0d45SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
1955a47d0d45SMatthew G. Knepley @*/
1956fba955ccSBarry Smith PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, DM *dm)
1957a47d0d45SMatthew G. Knepley {
1958a47d0d45SMatthew G. Knepley   PetscSF        sfVert;
1959a47d0d45SMatthew G. Knepley   PetscErrorCode ierr;
1960a47d0d45SMatthew G. Knepley 
1961a47d0d45SMatthew G. Knepley   PetscFunctionBegin;
1962a47d0d45SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1963a47d0d45SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1964a47d0d45SMatthew G. Knepley   PetscValidLogicalCollectiveInt(*dm, dim, 2);
1965a47d0d45SMatthew G. Knepley   PetscValidLogicalCollectiveInt(*dm, spaceDim, 8);
1966a47d0d45SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
1967a47d0d45SMatthew G. Knepley   ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr);
1968a47d0d45SMatthew G. Knepley   if (interpolate) {
1969a47d0d45SMatthew G. Knepley     DM idm = NULL;
1970a47d0d45SMatthew G. Knepley 
1971a47d0d45SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1972a47d0d45SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
1973a47d0d45SMatthew G. Knepley     *dm  = idm;
1974a47d0d45SMatthew G. Knepley   }
197521016a8bSBarry Smith   ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);CHKERRQ(ierr);
197618d54ad4SMichael Lange   if (vertexSF) *vertexSF = sfVert;
1977fba955ccSBarry Smith   else {ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);}
1978a47d0d45SMatthew G. Knepley   PetscFunctionReturn(0);
1979a47d0d45SMatthew G. Knepley }
1980a47d0d45SMatthew G. Knepley 
19819298eaa6SMatthew G Knepley /*
19829298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
19839298eaa6SMatthew G Knepley */
19848f767406SMatthew G. Knepley static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
19859298eaa6SMatthew G Knepley {
19869298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
19879298eaa6SMatthew G Knepley   PetscErrorCode ierr;
19889298eaa6SMatthew G Knepley 
19899298eaa6SMatthew G Knepley   PetscFunctionBegin;
19909298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
19919298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
19929298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
19939298eaa6SMatthew G Knepley   }
19949298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
19959298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
19969298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
19979298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
19989298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
19999298eaa6SMatthew G Knepley     }
20009298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
20019298eaa6SMatthew G Knepley   }
20029298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
20039298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
20049298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
20059298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
20069298eaa6SMatthew G Knepley }
20079298eaa6SMatthew G Knepley 
20089298eaa6SMatthew G Knepley /*
20099298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
20109298eaa6SMatthew G Knepley */
20118f767406SMatthew G. Knepley static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
20129298eaa6SMatthew G Knepley {
20139298eaa6SMatthew G Knepley   PetscSection   coordSection;
20149298eaa6SMatthew G Knepley   Vec            coordinates;
20158b9ced59SLisandro Dalcin   DM             cdm;
20169298eaa6SMatthew G Knepley   PetscScalar   *coords;
20178b9ced59SLisandro Dalcin   PetscInt       v, d;
20189298eaa6SMatthew G Knepley   PetscErrorCode ierr;
20199298eaa6SMatthew G Knepley 
20209298eaa6SMatthew G Knepley   PetscFunctionBegin;
2021c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
20229298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
20239298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
20249298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
20259298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
20269298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
20279298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
20289298eaa6SMatthew G Knepley   }
20299298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
20308b9ced59SLisandro Dalcin 
20318b9ced59SLisandro Dalcin   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
20328b9ced59SLisandro Dalcin   ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr);
20334b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
20349298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
20359298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
20369298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
20379298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
20389298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
20399298eaa6SMatthew G Knepley     }
20409298eaa6SMatthew G Knepley   }
20419298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
20429298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
20439298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
20449298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
20459298eaa6SMatthew G Knepley }
20469298eaa6SMatthew G Knepley 
20479298eaa6SMatthew G Knepley /*@C
20489298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
20499298eaa6SMatthew G Knepley 
20509298eaa6SMatthew G Knepley   Input Parameters:
20519298eaa6SMatthew G Knepley + comm - The communicator
20529298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
20539298eaa6SMatthew G Knepley . numCells - The number of cells
20549298eaa6SMatthew G Knepley . numVertices - The number of vertices
20559298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
20569298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
20579298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
20589298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
20599298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
20609298eaa6SMatthew G Knepley 
20619298eaa6SMatthew G Knepley   Output Parameter:
20629298eaa6SMatthew G Knepley . dm - The DM
20639298eaa6SMatthew G Knepley 
20649298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
20659298eaa6SMatthew G Knepley $
20669298eaa6SMatthew G Knepley $        2
20679298eaa6SMatthew G Knepley $      / | \
20689298eaa6SMatthew G Knepley $     /  |  \
20699298eaa6SMatthew G Knepley $    /   |   \
20709298eaa6SMatthew G Knepley $   0  0 | 1  3
20719298eaa6SMatthew G Knepley $    \   |   /
20729298eaa6SMatthew G Knepley $     \  |  /
20739298eaa6SMatthew G Knepley $      \ | /
20749298eaa6SMatthew G Knepley $        1
20759298eaa6SMatthew G Knepley would have input
20769298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
20779298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
20789298eaa6SMatthew G Knepley $
20799298eaa6SMatthew G Knepley which would result in the DMPlex
20809298eaa6SMatthew G Knepley $
20819298eaa6SMatthew G Knepley $        4
20829298eaa6SMatthew G Knepley $      / | \
20839298eaa6SMatthew G Knepley $     /  |  \
20849298eaa6SMatthew G Knepley $    /   |   \
20859298eaa6SMatthew G Knepley $   2  0 | 1  5
20869298eaa6SMatthew G Knepley $    \   |   /
20879298eaa6SMatthew G Knepley $     \  |  /
20889298eaa6SMatthew G Knepley $      \ | /
20899298eaa6SMatthew G Knepley $        3
20909298eaa6SMatthew G Knepley 
20919298eaa6SMatthew G Knepley   Level: beginner
20929298eaa6SMatthew G Knepley 
2093939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
20949298eaa6SMatthew G Knepley @*/
20959298eaa6SMatthew 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)
20969298eaa6SMatthew G Knepley {
20979298eaa6SMatthew G Knepley   PetscErrorCode ierr;
20989298eaa6SMatthew G Knepley 
20999298eaa6SMatthew G Knepley   PetscFunctionBegin;
21009298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
21019298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2102c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
21039298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
21049298eaa6SMatthew G Knepley   if (interpolate) {
2105e56d480eSMatthew G. Knepley     DM idm = NULL;
21069298eaa6SMatthew G Knepley 
21079298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
21089298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
21099298eaa6SMatthew G Knepley     *dm  = idm;
21109298eaa6SMatthew G Knepley   }
21119298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
21129298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
21139298eaa6SMatthew G Knepley }
21149298eaa6SMatthew G Knepley 
2115939f6067SMatthew G. Knepley /*@
2116939f6067SMatthew 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
2117939f6067SMatthew G. Knepley 
2118939f6067SMatthew G. Knepley   Input Parameters:
2119c73cfb54SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
2120939f6067SMatthew G. Knepley . depth - The depth of the DAG
2121939f6067SMatthew G. Knepley . numPoints - The number of points at each depth
2122939f6067SMatthew G. Knepley . coneSize - The cone size of each point
2123939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
2124939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point
2125939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
2126939f6067SMatthew G. Knepley 
2127939f6067SMatthew G. Knepley   Output Parameter:
2128939f6067SMatthew G. Knepley . dm - The DM
2129939f6067SMatthew G. Knepley 
2130939f6067SMatthew G. Knepley   Note: Two triangles sharing a face would have input
2131939f6067SMatthew G. Knepley $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
2132939f6067SMatthew G. Knepley $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
2133939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
2134939f6067SMatthew G. Knepley $
2135939f6067SMatthew G. Knepley which would result in the DMPlex
2136939f6067SMatthew G. Knepley $
2137939f6067SMatthew G. Knepley $        4
2138939f6067SMatthew G. Knepley $      / | \
2139939f6067SMatthew G. Knepley $     /  |  \
2140939f6067SMatthew G. Knepley $    /   |   \
2141939f6067SMatthew G. Knepley $   2  0 | 1  5
2142939f6067SMatthew G. Knepley $    \   |   /
2143939f6067SMatthew G. Knepley $     \  |  /
2144939f6067SMatthew G. Knepley $      \ | /
2145939f6067SMatthew G. Knepley $        3
2146939f6067SMatthew G. Knepley $
2147939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
2148939f6067SMatthew G. Knepley 
2149939f6067SMatthew G. Knepley   Level: advanced
2150939f6067SMatthew G. Knepley 
2151939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
2152939f6067SMatthew G. Knepley @*/
21539298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
21549298eaa6SMatthew G Knepley {
21559298eaa6SMatthew G Knepley   Vec            coordinates;
21569298eaa6SMatthew G Knepley   PetscSection   coordSection;
21579298eaa6SMatthew G Knepley   PetscScalar    *coords;
2158811e8653SToby Isaac   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;
21599298eaa6SMatthew G Knepley   PetscErrorCode ierr;
21609298eaa6SMatthew G Knepley 
21619298eaa6SMatthew G Knepley   PetscFunctionBegin;
2162c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2163811e8653SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
2164811e8653SToby Isaac   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
21659298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
21669298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
21679298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
21689298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
216997e052ccSToby Isaac     if (firstVertex < 0 && !coneSize[p - pStart]) {
217097e052ccSToby Isaac       firstVertex = p - pStart;
21719298eaa6SMatthew G Knepley     }
217297e052ccSToby Isaac   }
217397e052ccSToby Isaac   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
21749298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
21759298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
21769298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
21779298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
21789298eaa6SMatthew G Knepley   }
21799298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
21809298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
21819298eaa6SMatthew G Knepley   /* Build coordinates */
2182c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
21839298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
2184811e8653SToby Isaac   ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr);
21859298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
21869298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
2187811e8653SToby Isaac     ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr);
2188811e8653SToby Isaac     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr);
21899298eaa6SMatthew G Knepley   }
21909298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
21919298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
21928b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
21936f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
21949298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
21958b9ced59SLisandro Dalcin   ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr);
21962eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
21979298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
21989298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
21999298eaa6SMatthew G Knepley     PetscInt off;
22009298eaa6SMatthew G Knepley 
22019298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
2202811e8653SToby Isaac     for (d = 0; d < dimEmbed; ++d) {
2203811e8653SToby Isaac       coords[off+d] = vertexCoords[v*dimEmbed+d];
22049298eaa6SMatthew G Knepley     }
22059298eaa6SMatthew G Knepley   }
22069298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
22079298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
22089298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
22099298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
22109298eaa6SMatthew G Knepley }
22118415267dSToby Isaac 
2212ca522641SMatthew G. Knepley /*@C
2213ca522641SMatthew G. Knepley   DMPlexCreateFromFile - This takes a filename and produces a DM
2214ca522641SMatthew G. Knepley 
2215ca522641SMatthew G. Knepley   Input Parameters:
2216ca522641SMatthew G. Knepley + comm - The communicator
2217ca522641SMatthew G. Knepley . filename - A file name
2218ca522641SMatthew G. Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
2219ca522641SMatthew G. Knepley 
2220ca522641SMatthew G. Knepley   Output Parameter:
2221ca522641SMatthew G. Knepley . dm - The DM
2222ca522641SMatthew G. Knepley 
2223ca522641SMatthew G. Knepley   Level: beginner
2224ca522641SMatthew G. Knepley 
2225ca522641SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
2226ca522641SMatthew G. Knepley @*/
2227ca522641SMatthew G. Knepley PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
2228ca522641SMatthew G. Knepley {
222996a26029SBarry Smith   const char    *extGmsh   = ".msh";
2230ca522641SMatthew G. Knepley   const char    *extCGNS   = ".cgns";
2231ca522641SMatthew G. Knepley   const char    *extExodus = ".exo";
22322f0bd6dcSMichael Lange   const char    *extFluent = ".cas";
2233cc2f8f65SMatthew G. Knepley   const char    *extHDF5   = ".h5";
2234707dd687SMichael Lange   const char    *extMed    = ".med";
2235ca522641SMatthew G. Knepley   size_t         len;
2236707dd687SMichael Lange   PetscBool      isGmsh, isCGNS, isExodus, isFluent, isHDF5, isMed;
2237ca522641SMatthew G. Knepley   PetscMPIInt    rank;
2238ca522641SMatthew G. Knepley   PetscErrorCode ierr;
2239ca522641SMatthew G. Knepley 
2240ca522641SMatthew G. Knepley   PetscFunctionBegin;
2241ca522641SMatthew G. Knepley   PetscValidPointer(filename, 2);
2242ca522641SMatthew G. Knepley   PetscValidPointer(dm, 4);
2243ca522641SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2244ca522641SMatthew G. Knepley   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
2245ca522641SMatthew G. Knepley   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
2246ca522641SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,   4, &isGmsh);CHKERRQ(ierr);
2247ca522641SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,   5, &isCGNS);CHKERRQ(ierr);
2248ca522641SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr);
22492f0bd6dcSMichael Lange   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr);
2250cc2f8f65SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,   3, &isHDF5);CHKERRQ(ierr);
2251707dd687SMichael Lange   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,    4, &isMed);CHKERRQ(ierr);
2252ca522641SMatthew G. Knepley   if (isGmsh) {
22537d282ae0SMichael Lange     ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2254ca522641SMatthew G. Knepley   } else if (isCGNS) {
2255ca522641SMatthew G. Knepley     ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2256ca522641SMatthew G. Knepley   } else if (isExodus) {
2257ca522641SMatthew G. Knepley     ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
22582f0bd6dcSMichael Lange   } else if (isFluent) {
22592f0bd6dcSMichael Lange     ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2260cc2f8f65SMatthew G. Knepley   } else if (isHDF5) {
2261cc2f8f65SMatthew G. Knepley     PetscViewer viewer;
2262cc2f8f65SMatthew G. Knepley 
2263cc2f8f65SMatthew G. Knepley     ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
2264cc2f8f65SMatthew G. Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr);
2265cc2f8f65SMatthew G. Knepley     ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
2266cc2f8f65SMatthew G. Knepley     ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
2267cc2f8f65SMatthew G. Knepley     ierr = DMCreate(comm, dm);CHKERRQ(ierr);
2268cc2f8f65SMatthew G. Knepley     ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2269cc2f8f65SMatthew G. Knepley     ierr = DMLoad(*dm, viewer);CHKERRQ(ierr);
2270cc2f8f65SMatthew G. Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2271707dd687SMichael Lange   } else if (isMed) {
2272707dd687SMichael Lange     ierr = DMPlexCreateMedFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2273ca522641SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
2274ca522641SMatthew G. Knepley   PetscFunctionReturn(0);
2275ca522641SMatthew G. Knepley }
2276ca522641SMatthew G. Knepley 
22778415267dSToby Isaac /*@
22788415267dSToby Isaac   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
22798415267dSToby Isaac 
22808415267dSToby Isaac   Collective on comm
22818415267dSToby Isaac 
22828415267dSToby Isaac   Input Parameters:
22838415267dSToby Isaac + comm    - The communicator
22848415267dSToby Isaac . dim     - The spatial dimension
22858415267dSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
22868415267dSToby Isaac 
22878415267dSToby Isaac   Output Parameter:
22888415267dSToby Isaac . refdm - The reference cell
22898415267dSToby Isaac 
22908415267dSToby Isaac   Level: intermediate
22918415267dSToby Isaac 
22928415267dSToby Isaac .keywords: reference cell
22938415267dSToby Isaac .seealso:
22948415267dSToby Isaac @*/
22958415267dSToby Isaac PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
22968415267dSToby Isaac {
22978415267dSToby Isaac   DM             rdm;
2298cdf740e0SMatthew G. Knepley   Vec            coords;
22998415267dSToby Isaac   PetscErrorCode ierr;
23008415267dSToby Isaac 
23018415267dSToby Isaac   PetscFunctionBeginUser;
23028415267dSToby Isaac   ierr = DMCreate(comm, &rdm);CHKERRQ(ierr);
23038415267dSToby Isaac   ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr);
23048415267dSToby Isaac   ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr);
23058415267dSToby Isaac   switch (dim) {
23068415267dSToby Isaac   case 0:
23078415267dSToby Isaac   {
23088415267dSToby Isaac     PetscInt    numPoints[1]        = {1};
23098415267dSToby Isaac     PetscInt    coneSize[1]         = {0};
23108415267dSToby Isaac     PetscInt    cones[1]            = {0};
23118415267dSToby Isaac     PetscInt    coneOrientations[1] = {0};
23128415267dSToby Isaac     PetscScalar vertexCoords[1]     = {0.0};
23138415267dSToby Isaac 
23148415267dSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
23158415267dSToby Isaac   }
23168415267dSToby Isaac   break;
23178415267dSToby Isaac   case 1:
23188415267dSToby Isaac   {
23198415267dSToby Isaac     PetscInt    numPoints[2]        = {2, 1};
23208415267dSToby Isaac     PetscInt    coneSize[3]         = {2, 0, 0};
23218415267dSToby Isaac     PetscInt    cones[2]            = {1, 2};
23228415267dSToby Isaac     PetscInt    coneOrientations[2] = {0, 0};
23238415267dSToby Isaac     PetscScalar vertexCoords[2]     = {-1.0,  1.0};
23248415267dSToby Isaac 
23258415267dSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
23268415267dSToby Isaac   }
23278415267dSToby Isaac   break;
23288415267dSToby Isaac   case 2:
23298415267dSToby Isaac     if (simplex) {
23308415267dSToby Isaac       PetscInt    numPoints[2]        = {3, 1};
23318415267dSToby Isaac       PetscInt    coneSize[4]         = {3, 0, 0, 0};
23328415267dSToby Isaac       PetscInt    cones[3]            = {1, 2, 3};
23338415267dSToby Isaac       PetscInt    coneOrientations[3] = {0, 0, 0};
23348415267dSToby Isaac       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};
23358415267dSToby Isaac 
23368415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
23378415267dSToby Isaac     } else {
23388415267dSToby Isaac       PetscInt    numPoints[2]        = {4, 1};
23398415267dSToby Isaac       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
23408415267dSToby Isaac       PetscInt    cones[4]            = {1, 2, 3, 4};
23418415267dSToby Isaac       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
23428415267dSToby Isaac       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};
23438415267dSToby Isaac 
23448415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
23458415267dSToby Isaac     }
23468415267dSToby Isaac   break;
23478415267dSToby Isaac   case 3:
23488415267dSToby Isaac     if (simplex) {
23498415267dSToby Isaac       PetscInt    numPoints[2]        = {4, 1};
23508415267dSToby Isaac       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
23518415267dSToby Isaac       PetscInt    cones[4]            = {1, 3, 2, 4};
23528415267dSToby Isaac       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
23538415267dSToby 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};
23548415267dSToby Isaac 
23558415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
23568415267dSToby Isaac     } else {
23578415267dSToby Isaac       PetscInt    numPoints[2]        = {8, 1};
23588415267dSToby Isaac       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
23598415267dSToby Isaac       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
23608415267dSToby Isaac       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
23618415267dSToby 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,
23628415267dSToby 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};
23638415267dSToby Isaac 
23648415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
23658415267dSToby Isaac     }
23668415267dSToby Isaac   break;
23678415267dSToby Isaac   default:
23688415267dSToby Isaac     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
23698415267dSToby Isaac   }
23708415267dSToby Isaac   *refdm = NULL;
23718415267dSToby Isaac   ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr);
2372cdf740e0SMatthew G. Knepley   if (rdm->coordinateDM) {
2373cdf740e0SMatthew G. Knepley     DM           ncdm;
2374cdf740e0SMatthew G. Knepley     PetscSection cs;
2375cdf740e0SMatthew G. Knepley     PetscInt     pEnd = -1;
2376cdf740e0SMatthew G. Knepley 
2377cdf740e0SMatthew G. Knepley     ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr);
2378cdf740e0SMatthew G. Knepley     if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);}
2379cdf740e0SMatthew G. Knepley     if (pEnd >= 0) {
2380cdf740e0SMatthew G. Knepley       ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr);
2381cdf740e0SMatthew G. Knepley       ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr);
2382a61e840bSMatthew G. Knepley       ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr);
2383cdf740e0SMatthew G. Knepley       ierr = DMDestroy(&ncdm);CHKERRQ(ierr);
2384cdf740e0SMatthew G. Knepley     }
2385cdf740e0SMatthew G. Knepley   }
2386cdf740e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr);
2387cdf740e0SMatthew G. Knepley   if (coords) {
2388cdf740e0SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr);
2389cdf740e0SMatthew G. Knepley   } else {
2390cdf740e0SMatthew G. Knepley     ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr);
2391cdf740e0SMatthew G. Knepley     if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);}
2392cdf740e0SMatthew G. Knepley   }
23938415267dSToby Isaac   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
23948415267dSToby Isaac   PetscFunctionReturn(0);
23958415267dSToby Isaac }
2396