xref: /petsc/src/dm/impls/plex/plexcreate.c (revision b29cfa1c0f7944571ef0a263c2bd25d4db888f1b)
1552f7358SJed Brown #define PETSCDM_DLL
234541f0dSBarry Smith #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <petscdmda.h>
40c312b8eSJed Brown #include <petscsf.h>
5552f7358SJed Brown 
6552f7358SJed Brown #undef __FUNCT__
71df5d5c5SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateDoublet"
81df5d5c5SMatthew G. Knepley /*@
91df5d5c5SMatthew G. Knepley   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
101df5d5c5SMatthew G. Knepley 
111df5d5c5SMatthew G. Knepley   Collective on MPI_Comm
121df5d5c5SMatthew G. Knepley 
131df5d5c5SMatthew G. Knepley   Input Parameters:
141df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object
151df5d5c5SMatthew G. Knepley . dim - The spatial dimension
161df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells
171df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces)
181df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement
191df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
201df5d5c5SMatthew G. Knepley 
211df5d5c5SMatthew G. Knepley   Output Parameter:
221df5d5c5SMatthew G. Knepley . dm  - The DM object
231df5d5c5SMatthew G. Knepley 
241df5d5c5SMatthew G. Knepley   Level: beginner
251df5d5c5SMatthew G. Knepley 
261df5d5c5SMatthew G. Knepley .keywords: DM, create
271df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate()
281df5d5c5SMatthew G. Knepley @*/
291df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
301df5d5c5SMatthew G. Knepley {
311df5d5c5SMatthew G. Knepley   DM             dm;
321df5d5c5SMatthew G. Knepley   PetscInt       p;
331df5d5c5SMatthew G. Knepley   PetscMPIInt    rank;
341df5d5c5SMatthew G. Knepley   PetscErrorCode ierr;
351df5d5c5SMatthew G. Knepley 
361df5d5c5SMatthew G. Knepley   PetscFunctionBegin;
371df5d5c5SMatthew G. Knepley   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
381df5d5c5SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
39c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
401df5d5c5SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
41ce78fa2fSMatthew G. Knepley   switch (dim) {
42ce78fa2fSMatthew G. Knepley   case 2:
43ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
44ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
45ce78fa2fSMatthew G. Knepley     break;
46ce78fa2fSMatthew G. Knepley   case 3:
47ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
48ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
49ce78fa2fSMatthew G. Knepley     break;
50ce78fa2fSMatthew G. Knepley   default:
51ce78fa2fSMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
52ce78fa2fSMatthew G. Knepley   }
531df5d5c5SMatthew G. Knepley   if (rank) {
541df5d5c5SMatthew G. Knepley     PetscInt numPoints[2] = {0, 0};
551df5d5c5SMatthew G. Knepley     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
561df5d5c5SMatthew G. Knepley   } else {
571df5d5c5SMatthew G. Knepley     switch (dim) {
581df5d5c5SMatthew G. Knepley     case 2:
591df5d5c5SMatthew G. Knepley       if (simplex) {
601df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {4, 2};
611df5d5c5SMatthew G. Knepley         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
621df5d5c5SMatthew G. Knepley         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
631df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
641df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
651df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
661df5d5c5SMatthew G. Knepley 
671df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
681df5d5c5SMatthew G. Knepley         for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
691df5d5c5SMatthew G. Knepley       } else {
701df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {6, 2};
711df5d5c5SMatthew G. Knepley         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
721df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
731df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
741df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[12]    = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  -1.0, 0.5,  1.0, -0.5,  1.0, 0.5};
751df5d5c5SMatthew G. Knepley 
761df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
771df5d5c5SMatthew G. Knepley       }
781df5d5c5SMatthew G. Knepley       break;
791df5d5c5SMatthew G. Knepley     case 3:
801df5d5c5SMatthew G. Knepley       if (simplex) {
811df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {5, 2};
821df5d5c5SMatthew G. Knepley         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
831df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
841df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
851df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  1.0, 0.0, 0.0};
861df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
871df5d5c5SMatthew G. Knepley 
881df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
891df5d5c5SMatthew G. Knepley         for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
901df5d5c5SMatthew G. Knepley       } else {
911df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]         = {12, 2};
921df5d5c5SMatthew G. Knepley         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
931df5d5c5SMatthew G. Knepley         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
941df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
951df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5,  -1.0,  0.5, -0.5,  0.0,  0.5, -0.5,   0.0, -0.5, -0.5,
961df5d5c5SMatthew G. Knepley                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
971df5d5c5SMatthew G. Knepley                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
981df5d5c5SMatthew G. Knepley 
991df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
1001df5d5c5SMatthew G. Knepley       }
1011df5d5c5SMatthew G. Knepley       break;
1021df5d5c5SMatthew G. Knepley     default:
1031df5d5c5SMatthew G. Knepley       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
1041df5d5c5SMatthew G. Knepley     }
1051df5d5c5SMatthew G. Knepley   }
1061df5d5c5SMatthew G. Knepley   *newdm = dm;
1071df5d5c5SMatthew G. Knepley   if (refinementLimit > 0.0) {
1081df5d5c5SMatthew G. Knepley     DM rdm;
1091df5d5c5SMatthew G. Knepley     const char *name;
1101df5d5c5SMatthew G. Knepley 
1111df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
1121df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
1131df5d5c5SMatthew G. Knepley     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
1141df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1151df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
1161df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1171df5d5c5SMatthew G. Knepley     *newdm = rdm;
1181df5d5c5SMatthew G. Knepley   }
1191df5d5c5SMatthew G. Knepley   if (interpolate) {
120e56d480eSMatthew G. Knepley     DM idm = NULL;
1211df5d5c5SMatthew G. Knepley     const char *name;
1221df5d5c5SMatthew G. Knepley 
1231df5d5c5SMatthew G. Knepley     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
1241df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1251df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
1261df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
1271df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr);
1281df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1291df5d5c5SMatthew G. Knepley     *newdm = idm;
1301df5d5c5SMatthew G. Knepley   }
1311df5d5c5SMatthew G. Knepley   {
1321df5d5c5SMatthew G. Knepley     DM refinedMesh     = NULL;
1331df5d5c5SMatthew G. Knepley     DM distributedMesh = NULL;
1341df5d5c5SMatthew G. Knepley 
1351df5d5c5SMatthew G. Knepley     /* Distribute mesh over processes */
13607104863SJed Brown     ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr);
1371df5d5c5SMatthew G. Knepley     if (distributedMesh) {
1381df5d5c5SMatthew G. Knepley       ierr = DMDestroy(newdm);CHKERRQ(ierr);
1391df5d5c5SMatthew G. Knepley       *newdm = distributedMesh;
1401df5d5c5SMatthew G. Knepley     }
1411df5d5c5SMatthew G. Knepley     if (refinementUniform) {
1421df5d5c5SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
1431df5d5c5SMatthew G. Knepley       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
1441df5d5c5SMatthew G. Knepley       if (refinedMesh) {
1451df5d5c5SMatthew G. Knepley         ierr = DMDestroy(newdm);CHKERRQ(ierr);
1461df5d5c5SMatthew G. Knepley         *newdm = refinedMesh;
1471df5d5c5SMatthew G. Knepley       }
1481df5d5c5SMatthew G. Knepley     }
1491df5d5c5SMatthew G. Knepley   }
1501df5d5c5SMatthew G. Knepley   PetscFunctionReturn(0);
1511df5d5c5SMatthew G. Knepley }
1521df5d5c5SMatthew G. Knepley 
1531df5d5c5SMatthew G. Knepley #undef __FUNCT__
154552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary"
15526492d91SMatthew G. Knepley /*@
15626492d91SMatthew G. Knepley   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
157552f7358SJed Brown 
15826492d91SMatthew G. Knepley   Collective on MPI_Comm
15926492d91SMatthew G. Knepley 
16026492d91SMatthew G. Knepley   Input Parameters:
16126492d91SMatthew G. Knepley + comm  - The communicator for the DM object
16226492d91SMatthew G. Knepley . lower - The lower left corner coordinates
16326492d91SMatthew G. Knepley . upper - The upper right corner coordinates
16426492d91SMatthew G. Knepley - edges - The number of cells in each direction
16526492d91SMatthew G. Knepley 
16626492d91SMatthew G. Knepley   Output Parameter:
16726492d91SMatthew G. Knepley . dm  - The DM object
16826492d91SMatthew G. Knepley 
16926492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
17026492d91SMatthew G. Knepley $ 18--5-17--4--16
17126492d91SMatthew G. Knepley $  |     |     |
17226492d91SMatthew G. Knepley $  6    10     3
17326492d91SMatthew G. Knepley $  |     |     |
17426492d91SMatthew G. Knepley $ 19-11-20--9--15
17526492d91SMatthew G. Knepley $  |     |     |
17626492d91SMatthew G. Knepley $  7     8     2
17726492d91SMatthew G. Knepley $  |     |     |
17826492d91SMatthew G. Knepley $ 12--0-13--1--14
17926492d91SMatthew G. Knepley 
18026492d91SMatthew G. Knepley   Level: beginner
18126492d91SMatthew G. Knepley 
18226492d91SMatthew G. Knepley .keywords: DM, create
18326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
18426492d91SMatthew G. Knepley @*/
185552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
186552f7358SJed Brown {
187552f7358SJed Brown   PetscInt       numVertices    = (edges[0]+1)*(edges[1]+1);
188552f7358SJed Brown   PetscInt       numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
189552f7358SJed Brown   PetscInt       markerTop      = 1;
190552f7358SJed Brown   PetscInt       markerBottom   = 1;
191552f7358SJed Brown   PetscInt       markerRight    = 1;
192552f7358SJed Brown   PetscInt       markerLeft     = 1;
193552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
194552f7358SJed Brown   Vec            coordinates;
195552f7358SJed Brown   PetscSection   coordSection;
196552f7358SJed Brown   PetscScalar    *coords;
197552f7358SJed Brown   PetscInt       coordSize;
198552f7358SJed Brown   PetscMPIInt    rank;
199552f7358SJed Brown   PetscInt       v, vx, vy;
200552f7358SJed Brown   PetscErrorCode ierr;
201552f7358SJed Brown 
202552f7358SJed Brown   PetscFunctionBegin;
2030298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
204552f7358SJed Brown   if (markerSeparate) {
205552f7358SJed Brown     markerTop    = 1;
206552f7358SJed Brown     markerBottom = 0;
207552f7358SJed Brown     markerRight  = 0;
208552f7358SJed Brown     markerLeft   = 0;
209552f7358SJed Brown   }
21082f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
211552f7358SJed Brown   if (!rank) {
212552f7358SJed Brown     PetscInt e, ex, ey;
213552f7358SJed Brown 
214552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
215552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
216552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
217552f7358SJed Brown     }
218552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
219552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
220552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
221552f7358SJed Brown         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
222552f7358SJed Brown         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
223da80777bSKarl Rupp         PetscInt cone[2];
224552f7358SJed Brown 
225da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
226552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
227552f7358SJed Brown         if (vx == edges[0]) {
228552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
229552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
230552f7358SJed Brown           if (ey == edges[1]-1) {
231552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
232552f7358SJed Brown           }
233552f7358SJed Brown         } else if (vx == 0) {
234552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
235552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
236552f7358SJed Brown           if (ey == edges[1]-1) {
237552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
238552f7358SJed Brown           }
239552f7358SJed Brown         }
240552f7358SJed Brown       }
241552f7358SJed Brown     }
242552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
243552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
244552f7358SJed Brown         PetscInt edge   = vy*edges[0]     + ex;
245552f7358SJed Brown         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
246da80777bSKarl Rupp         PetscInt cone[2];
247552f7358SJed Brown 
248da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
249552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
250552f7358SJed Brown         if (vy == edges[1]) {
251552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
252552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
253552f7358SJed Brown           if (ex == edges[0]-1) {
254552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
255552f7358SJed Brown           }
256552f7358SJed Brown         } else if (vy == 0) {
257552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
258552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
259552f7358SJed Brown           if (ex == edges[0]-1) {
260552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
261552f7358SJed Brown           }
262552f7358SJed Brown         }
263552f7358SJed Brown       }
264552f7358SJed Brown     }
265552f7358SJed Brown   }
266552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
267552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
268552f7358SJed Brown   /* Build coordinates */
269c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
270552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
271552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
272552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
273552f7358SJed Brown   }
274552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
275552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
27682f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
2774b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
278552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2792eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
280552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
281552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
282552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
283552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
284552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
285552f7358SJed Brown     }
286552f7358SJed Brown   }
287552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
288552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
289552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
290552f7358SJed Brown   PetscFunctionReturn(0);
291552f7358SJed Brown }
292552f7358SJed Brown 
293552f7358SJed Brown #undef __FUNCT__
294552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
29526492d91SMatthew G. Knepley /*@
29626492d91SMatthew G. Knepley   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
297552f7358SJed Brown 
29826492d91SMatthew G. Knepley   Collective on MPI_Comm
29926492d91SMatthew G. Knepley 
30026492d91SMatthew G. Knepley   Input Parameters:
30126492d91SMatthew G. Knepley + comm  - The communicator for the DM object
30226492d91SMatthew G. Knepley . lower - The lower left front corner coordinates
30326492d91SMatthew G. Knepley . upper - The upper right back corner coordinates
30426492d91SMatthew G. Knepley - edges - The number of cells in each direction
30526492d91SMatthew G. Knepley 
30626492d91SMatthew G. Knepley   Output Parameter:
30726492d91SMatthew G. Knepley . dm  - The DM object
30826492d91SMatthew G. Knepley 
30926492d91SMatthew G. Knepley   Level: beginner
31026492d91SMatthew G. Knepley 
31126492d91SMatthew G. Knepley .keywords: DM, create
31226492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
31326492d91SMatthew G. Knepley @*/
314552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
315552f7358SJed Brown {
3169e8abbc3SMichael Lange   PetscInt       vertices[3], numVertices;
3177b59f5a9SMichael Lange   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
318552f7358SJed Brown   Vec            coordinates;
319552f7358SJed Brown   PetscSection   coordSection;
320552f7358SJed Brown   PetscScalar    *coords;
321552f7358SJed Brown   PetscInt       coordSize;
322552f7358SJed Brown   PetscMPIInt    rank;
323552f7358SJed Brown   PetscInt       v, vx, vy, vz;
3247b59f5a9SMichael Lange   PetscInt       voffset, iface=0, cone[4];
325552f7358SJed Brown   PetscErrorCode ierr;
326552f7358SJed Brown 
327552f7358SJed Brown   PetscFunctionBegin;
32882f516ccSBarry Smith   if ((faces[0] < 1) || (faces[1] < 1) || (faces[2] < 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side");
32982f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
3309e8abbc3SMichael Lange   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
3319e8abbc3SMichael Lange   numVertices = vertices[0]*vertices[1]*vertices[2];
332552f7358SJed Brown   if (!rank) {
333552f7358SJed Brown     PetscInt f;
334552f7358SJed Brown 
335552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
336552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
337552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
338552f7358SJed Brown     }
339552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
340552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
341552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
342552f7358SJed Brown     }
3437b59f5a9SMichael Lange 
3447b59f5a9SMichael Lange     /* Side 0 (Top) */
3457b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3467b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3477b59f5a9SMichael Lange         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
3487b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
3497b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3507b59f5a9SMichael Lange         iface++;
351552f7358SJed Brown       }
352552f7358SJed Brown     }
3537b59f5a9SMichael Lange 
3547b59f5a9SMichael Lange     /* Side 1 (Bottom) */
3557b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3567b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3577b59f5a9SMichael Lange         voffset = numFaces + vy*(faces[0]+1) + vx;
3587b59f5a9SMichael Lange         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
3597b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3607b59f5a9SMichael Lange         iface++;
361552f7358SJed Brown       }
362552f7358SJed Brown     }
3637b59f5a9SMichael Lange 
3647b59f5a9SMichael Lange     /* Side 2 (Front) */
3657b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3667b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3677b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
3687b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
3697b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3707b59f5a9SMichael Lange         iface++;
371552f7358SJed Brown       }
3727b59f5a9SMichael Lange     }
3737b59f5a9SMichael Lange 
3747b59f5a9SMichael Lange     /* Side 3 (Back) */
3757b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3767b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3777b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
3787b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
3797b59f5a9SMichael Lange         cone[2] = voffset+1; cone[3] = voffset;
3807b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3817b59f5a9SMichael Lange         iface++;
3827b59f5a9SMichael Lange       }
3837b59f5a9SMichael Lange     }
3847b59f5a9SMichael Lange 
3857b59f5a9SMichael Lange     /* Side 4 (Left) */
3867b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3877b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3887b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
3897b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
3907b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
3917b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3927b59f5a9SMichael Lange         iface++;
3937b59f5a9SMichael Lange       }
3947b59f5a9SMichael Lange     }
3957b59f5a9SMichael Lange 
3967b59f5a9SMichael Lange     /* Side 5 (Right) */
3977b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3987b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3997b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx;
4007b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
4017b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
4027b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
4037b59f5a9SMichael Lange         iface++;
4047b59f5a9SMichael Lange       }
405552f7358SJed Brown     }
406552f7358SJed Brown   }
407552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
408552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
409552f7358SJed Brown   /* Build coordinates */
410c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
411552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
412552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
413552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
414552f7358SJed Brown   }
415552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
416552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
41782f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
4184b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr);
419552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
420552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4212eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
422552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
423552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
424552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
425552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
426552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
427552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
428552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
429552f7358SJed Brown       }
430552f7358SJed Brown     }
431552f7358SJed Brown   }
432552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
433552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
434552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
435552f7358SJed Brown   PetscFunctionReturn(0);
436552f7358SJed Brown }
437552f7358SJed Brown 
438552f7358SJed Brown #undef __FUNCT__
439552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
44026492d91SMatthew G. Knepley /*@
44126492d91SMatthew G. Knepley   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
442552f7358SJed Brown 
44326492d91SMatthew G. Knepley   Collective on MPI_Comm
44426492d91SMatthew G. Knepley 
44526492d91SMatthew G. Knepley   Input Parameters:
44626492d91SMatthew G. Knepley + comm  - The communicator for the DM object
44726492d91SMatthew G. Knepley . lower - The lower left corner coordinates
44826492d91SMatthew G. Knepley . upper - The upper right corner coordinates
449fbf5b483SMatthew G. Knepley . edges - The number of cells in each direction
450fbf5b483SMatthew G. Knepley . bdX   - The boundary type for the X direction
451fbf5b483SMatthew G. Knepley - bdY   - The boundary type for the Y direction
45226492d91SMatthew G. Knepley 
45326492d91SMatthew G. Knepley   Output Parameter:
45426492d91SMatthew G. Knepley . dm  - The DM object
45526492d91SMatthew G. Knepley 
45626492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
45726492d91SMatthew G. Knepley $ 22--8-23--9--24
45826492d91SMatthew G. Knepley $  |     |     |
45926492d91SMatthew G. Knepley $ 13  2 14  3  15
46026492d91SMatthew G. Knepley $  |     |     |
46126492d91SMatthew G. Knepley $ 19--6-20--7--21
46226492d91SMatthew G. Knepley $  |     |     |
46326492d91SMatthew G. Knepley $ 10  0 11  1 12
46426492d91SMatthew G. Knepley $  |     |     |
46526492d91SMatthew G. Knepley $ 16--4-17--5--18
46626492d91SMatthew G. Knepley 
46726492d91SMatthew G. Knepley   Level: beginner
46826492d91SMatthew G. Knepley 
46926492d91SMatthew G. Knepley .keywords: DM, create
47026492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
47126492d91SMatthew G. Knepley @*/
472bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
473552f7358SJed Brown {
474552f7358SJed Brown   PetscInt       markerTop      = 1;
475552f7358SJed Brown   PetscInt       markerBottom   = 1;
476552f7358SJed Brown   PetscInt       markerRight    = 1;
477552f7358SJed Brown   PetscInt       markerLeft     = 1;
478552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
479552f7358SJed Brown   PetscMPIInt    rank;
480552f7358SJed Brown   PetscErrorCode ierr;
481552f7358SJed Brown 
482552f7358SJed Brown   PetscFunctionBegin;
48382f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
4840298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
485552f7358SJed Brown   if (markerSeparate) {
486552f7358SJed Brown     markerTop    = 3;
487552f7358SJed Brown     markerBottom = 1;
488552f7358SJed Brown     markerRight  = 2;
489552f7358SJed Brown     markerLeft   = 4;
490552f7358SJed Brown   }
491552f7358SJed Brown   {
492552f7358SJed Brown     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
493552f7358SJed Brown     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
49446675b8eSMatthew G. Knepley     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
49546675b8eSMatthew G. Knepley     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
496552f7358SJed Brown     const PetscInt numTotXEdges = numXEdges*numYVertices;
497552f7358SJed Brown     const PetscInt numTotYEdges = numYEdges*numXVertices;
498552f7358SJed Brown     const PetscInt numVertices  = numXVertices*numYVertices;
499552f7358SJed Brown     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
500552f7358SJed Brown     const PetscInt numFaces     = numXEdges*numYEdges;
501552f7358SJed Brown     const PetscInt firstVertex  = numFaces;
502552f7358SJed Brown     const PetscInt firstXEdge   = numFaces + numVertices;
503552f7358SJed Brown     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
504552f7358SJed Brown     Vec            coordinates;
505552f7358SJed Brown     PetscSection   coordSection;
506552f7358SJed Brown     PetscScalar   *coords;
507552f7358SJed Brown     PetscInt       coordSize;
508552f7358SJed Brown     PetscInt       v, vx, vy;
509552f7358SJed Brown     PetscInt       f, fx, fy, e, ex, ey;
510552f7358SJed Brown 
511552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
512552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
513552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
514552f7358SJed Brown     }
515552f7358SJed Brown     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
516552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
517552f7358SJed Brown     }
518552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
519552f7358SJed Brown     /* Build faces */
5208a6604f2SMatthew G. Knepley     for (fy = 0; fy < numYEdges; ++fy) {
5218a6604f2SMatthew G. Knepley       for (fx = 0; fx < numXEdges; ++fx) {
5228a6604f2SMatthew G. Knepley         PetscInt face    = fy*numXEdges + fx;
5238a6604f2SMatthew G. Knepley         PetscInt edgeL   = firstYEdge +   fx                 *numYEdges + fy;
5248a6604f2SMatthew G. Knepley         PetscInt edgeR   = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy;
5258a6604f2SMatthew G. Knepley         PetscInt edgeB   = firstXEdge +   fy                 *numXEdges + fx;
5268a6604f2SMatthew G. Knepley         PetscInt edgeT   = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx;
5278a6604f2SMatthew G. Knepley         PetscInt ornt[4] = {0, 0, -2, -2};
528da80777bSKarl Rupp         PetscInt cone[4];
529552f7358SJed Brown 
5308a6604f2SMatthew G. Knepley         if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
5318a6604f2SMatthew G. Knepley         if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
532becd5721SMatthew G. Knepley         cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
533552f7358SJed Brown         ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
534552f7358SJed Brown         ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
535552f7358SJed Brown       }
536552f7358SJed Brown     }
537552f7358SJed Brown     /* Build Y edges*/
538552f7358SJed Brown     for (vx = 0; vx < numXVertices; vx++) {
539552f7358SJed Brown       for (ey = 0; ey < numYEdges; ey++) {
54046675b8eSMatthew G. Knepley         const PetscInt nextv   = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx;
541552f7358SJed Brown         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
542becd5721SMatthew G. Knepley         const PetscInt vertexB = firstVertex + ey*numXVertices + vx;
54346675b8eSMatthew G. Knepley         const PetscInt vertexT = firstVertex + nextv;
544da80777bSKarl Rupp         PetscInt       cone[2];
545552f7358SJed Brown 
546becd5721SMatthew G. Knepley         cone[0] = vertexB; cone[1] = vertexT;
547552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
548c4b39bc1SMatthew G. Knepley         if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
549552f7358SJed Brown           if (vx == numXVertices-1) {
550552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
551552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
552552f7358SJed Brown             if (ey == numYEdges-1) {
553552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
554552f7358SJed Brown             }
555552f7358SJed Brown           } else if (vx == 0) {
556552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
557552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
558552f7358SJed Brown             if (ey == numYEdges-1) {
559552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
560552f7358SJed Brown             }
561552f7358SJed Brown           }
562552f7358SJed Brown         }
563552f7358SJed Brown       }
564c4b39bc1SMatthew G. Knepley     }
565552f7358SJed Brown     /* Build X edges*/
566552f7358SJed Brown     for (vy = 0; vy < numYVertices; vy++) {
567552f7358SJed Brown       for (ex = 0; ex < numXEdges; ex++) {
56846675b8eSMatthew G. Knepley         const PetscInt nextv   = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices;
569552f7358SJed Brown         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
570becd5721SMatthew G. Knepley         const PetscInt vertexL = firstVertex + vy*numXVertices + ex;
57146675b8eSMatthew G. Knepley         const PetscInt vertexR = firstVertex + nextv;
572da80777bSKarl Rupp         PetscInt       cone[2];
573552f7358SJed Brown 
574becd5721SMatthew G. Knepley         cone[0] = vertexL; cone[1] = vertexR;
575552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
576c4b39bc1SMatthew G. Knepley         if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
577552f7358SJed Brown           if (vy == numYVertices-1) {
578552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
579552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
580552f7358SJed Brown             if (ex == numXEdges-1) {
581552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
582552f7358SJed Brown             }
583552f7358SJed Brown           } else if (vy == 0) {
584552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
585552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
586552f7358SJed Brown             if (ex == numXEdges-1) {
587552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
588552f7358SJed Brown             }
589552f7358SJed Brown           }
590552f7358SJed Brown         }
591552f7358SJed Brown       }
592c4b39bc1SMatthew G. Knepley     }
593552f7358SJed Brown     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
594552f7358SJed Brown     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
595552f7358SJed Brown     /* Build coordinates */
596c2166f76SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
597130924b7SMatthew G. Knepley     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
598130924b7SMatthew G. Knepley     ierr = PetscSectionSetFieldComponents(coordSection, 0, 2);CHKERRQ(ierr);
599552f7358SJed Brown     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
600552f7358SJed Brown     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
601552f7358SJed Brown       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
602130924b7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(coordSection, v, 0, 2);CHKERRQ(ierr);
603552f7358SJed Brown     }
604552f7358SJed Brown     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
605552f7358SJed Brown     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
60682f516ccSBarry Smith     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
6074b66c01dSMatthew G. Knepley     ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
608552f7358SJed Brown     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
6092eb5907fSJed Brown     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
610552f7358SJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
611552f7358SJed Brown     for (vy = 0; vy < numYVertices; ++vy) {
612552f7358SJed Brown       for (vx = 0; vx < numXVertices; ++vx) {
613552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
614552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
615552f7358SJed Brown       }
616552f7358SJed Brown     }
617552f7358SJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
618552f7358SJed Brown     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
619552f7358SJed Brown     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
620552f7358SJed Brown   }
621552f7358SJed Brown   PetscFunctionReturn(0);
622552f7358SJed Brown }
623552f7358SJed Brown 
624552f7358SJed Brown #undef __FUNCT__
625552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
6266b75479cSMatthew G Knepley /*@
62726492d91SMatthew G. Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
6286b75479cSMatthew G Knepley 
6296b75479cSMatthew G Knepley   Collective on MPI_Comm
6306b75479cSMatthew G Knepley 
6316b75479cSMatthew G Knepley   Input Parameters:
6326b75479cSMatthew G Knepley + comm - The communicator for the DM object
6336b75479cSMatthew G Knepley . dim - The spatial dimension
6346b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
6356b75479cSMatthew G Knepley 
6366b75479cSMatthew G Knepley   Output Parameter:
6376b75479cSMatthew G Knepley . dm  - The DM object
6386b75479cSMatthew G Knepley 
6396b75479cSMatthew G Knepley   Level: beginner
6406b75479cSMatthew G Knepley 
6416b75479cSMatthew G Knepley .keywords: DM, create
64226492d91SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
6436b75479cSMatthew G Knepley @*/
6440adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
6450adebc6cSBarry Smith {
646552f7358SJed Brown   DM             boundary;
647552f7358SJed Brown   PetscErrorCode ierr;
648552f7358SJed Brown 
649552f7358SJed Brown   PetscFunctionBegin;
650552f7358SJed Brown   PetscValidPointer(dm, 4);
651552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
652552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
653552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
654c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr);
655552f7358SJed Brown   switch (dim) {
656552f7358SJed Brown   case 2:
657552f7358SJed Brown   {
658552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
659552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
660552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
661552f7358SJed Brown 
662552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
663552f7358SJed Brown     break;
664552f7358SJed Brown   }
665552f7358SJed Brown   case 3:
666552f7358SJed Brown   {
667552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
668552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
669552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
670552f7358SJed Brown 
671552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
672552f7358SJed Brown     break;
673552f7358SJed Brown   }
674552f7358SJed Brown   default:
675552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
676552f7358SJed Brown   }
6770298fd71SBarry Smith   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
678552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
679552f7358SJed Brown   PetscFunctionReturn(0);
680552f7358SJed Brown }
681552f7358SJed Brown 
682552f7358SJed Brown #undef __FUNCT__
683552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
68426492d91SMatthew G. Knepley /*@
68526492d91SMatthew G. Knepley   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
68626492d91SMatthew G. Knepley 
68726492d91SMatthew G. Knepley   Collective on MPI_Comm
68826492d91SMatthew G. Knepley 
68926492d91SMatthew G. Knepley   Input Parameters:
69026492d91SMatthew G. Knepley + comm  - The communicator for the DM object
69126492d91SMatthew G. Knepley . dim   - The spatial dimension
692fbf5b483SMatthew G. Knepley . periodicX - The boundary type for the X direction
693fbf5b483SMatthew G. Knepley . periodicY - The boundary type for the Y direction
694fbf5b483SMatthew G. Knepley . periodicZ - The boundary type for the Z direction
69526492d91SMatthew G. Knepley - cells - The number of cells in each direction
69626492d91SMatthew G. Knepley 
69726492d91SMatthew G. Knepley   Output Parameter:
69826492d91SMatthew G. Knepley . dm  - The DM object
69926492d91SMatthew G. Knepley 
70026492d91SMatthew G. Knepley   Level: beginner
70126492d91SMatthew G. Knepley 
70226492d91SMatthew G. Knepley .keywords: DM, create
70326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
70426492d91SMatthew G. Knepley @*/
705bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
706a6dfd86eSKarl Rupp {
707552f7358SJed Brown   PetscErrorCode ierr;
708552f7358SJed Brown 
709552f7358SJed Brown   PetscFunctionBegin;
710552f7358SJed Brown   PetscValidPointer(dm, 4);
711552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
712552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
713552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
714c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
715552f7358SJed Brown   switch (dim) {
716552f7358SJed Brown   case 2:
717552f7358SJed Brown   {
718552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
719552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
720552f7358SJed Brown 
721becd5721SMatthew G. Knepley     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
722552f7358SJed Brown     break;
723552f7358SJed Brown   }
724552f7358SJed Brown #if 0
725552f7358SJed Brown   case 3:
726552f7358SJed Brown   {
727552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
728552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
729552f7358SJed Brown 
730becd5721SMatthew G. Knepley     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
731552f7358SJed Brown     break;
732552f7358SJed Brown   }
733552f7358SJed Brown #endif
734552f7358SJed Brown   default:
735552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
736552f7358SJed Brown   }
737552f7358SJed Brown   PetscFunctionReturn(0);
738552f7358SJed Brown }
739552f7358SJed Brown 
740552f7358SJed Brown /* External function declarations here */
741552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
742bceba477SMatthew G. Knepley extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx);
743fd59a867SMatthew G. Knepley extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
744b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
745552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
746552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
7470a6ba040SMatthew G. Knepley extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
7480a6ba040SMatthew G. Knepley extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
7499a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
750552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
751552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
752552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
7532c40f234SMatthew G. Knepley extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
754552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
755552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
756552f7358SJed Brown 
757552f7358SJed Brown #undef __FUNCT__
7580a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexReplace_Static"
7590a6ba040SMatthew G. Knepley /* Replace dm with the contents of dmNew
7600a6ba040SMatthew G. Knepley    - Share the DM_Plex structure
7610a6ba040SMatthew G. Knepley    - Share the coordinates
762d7973521SMatthew G. Knepley    - Share the SF
7630a6ba040SMatthew G. Knepley */
7640a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
7650a6ba040SMatthew G. Knepley {
766d7973521SMatthew G. Knepley   PetscSF          sf;
76755fbe3e3SMatthew G. Knepley   DM               coordDM;
7680a6ba040SMatthew G. Knepley   Vec              coords;
76955fbe3e3SMatthew G. Knepley   const PetscReal *maxCell, *L;
7700a6ba040SMatthew G. Knepley   PetscErrorCode   ierr;
7710a6ba040SMatthew G. Knepley 
7720a6ba040SMatthew G. Knepley   PetscFunctionBegin;
773d7973521SMatthew G. Knepley   ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr);
774d7973521SMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
77555fbe3e3SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr);
7760a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
77755fbe3e3SMatthew G. Knepley   ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr);
7780a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
77955fbe3e3SMatthew G. Knepley   ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
78055fbe3e3SMatthew G. Knepley   if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L);CHKERRQ(ierr);}
7810a6ba040SMatthew G. Knepley   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
7820a6ba040SMatthew G. Knepley   dm->data = dmNew->data;
7830a6ba040SMatthew G. Knepley   ((DM_Plex *) dmNew->data)->refct++;
7840a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
7850a6ba040SMatthew G. Knepley }
7860a6ba040SMatthew G. Knepley 
7870a6ba040SMatthew G. Knepley #undef __FUNCT__
7880a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexSwap_Static"
7890a6ba040SMatthew G. Knepley /* Swap dm with the contents of dmNew
7900a6ba040SMatthew G. Knepley    - Swap the DM_Plex structure
7910a6ba040SMatthew G. Knepley    - Swap the coordinates
7920a6ba040SMatthew G. Knepley */
7930a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
7940a6ba040SMatthew G. Knepley {
7950a6ba040SMatthew G. Knepley   DM             coordDMA, coordDMB;
7960a6ba040SMatthew G. Knepley   Vec            coordsA,  coordsB;
7970a6ba040SMatthew G. Knepley   void          *tmp;
7980a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
7990a6ba040SMatthew G. Knepley 
8000a6ba040SMatthew G. Knepley   PetscFunctionBegin;
8010a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
8020a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
8030a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
8040a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
8050a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
8060a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
8070a6ba040SMatthew G. Knepley 
8080a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
8090a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
8100a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
8110a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
8120a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
8130a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
8140a6ba040SMatthew G. Knepley   tmp       = dmA->data;
8150a6ba040SMatthew G. Knepley   dmA->data = dmB->data;
8160a6ba040SMatthew G. Knepley   dmB->data = tmp;
8170a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
8180a6ba040SMatthew G. Knepley }
8190a6ba040SMatthew G. Knepley 
8200a6ba040SMatthew G. Knepley #undef __FUNCT__
82168d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex"
82268d4fef7SMatthew G. Knepley PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(DM dm)
8230a6ba040SMatthew G. Knepley {
8240a6ba040SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
8256c73c22cSMatthew G. Knepley   DMBoundary     b;
8260a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
8270a6ba040SMatthew G. Knepley 
8280a6ba040SMatthew G. Knepley   PetscFunctionBegin;
8296c73c22cSMatthew G. Knepley   /* Handle boundary conditions */
8306c73c22cSMatthew G. Knepley   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr);
8316c73c22cSMatthew G. Knepley   for (b = mesh->boundary; b; b = b->next) {
8326c73c22cSMatthew G. Knepley     char      optname[1024];
8336c73c22cSMatthew G. Knepley     PetscInt  ids[1024], len = 1024, i;
8346c73c22cSMatthew G. Knepley     PetscBool flg;
8356c73c22cSMatthew G. Knepley 
8366c73c22cSMatthew G. Knepley     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
8376c73c22cSMatthew G. Knepley     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
8386c73c22cSMatthew G. Knepley     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
8396c73c22cSMatthew G. Knepley     if (flg) {
8406c73c22cSMatthew G. Knepley       DMLabel label;
8416c73c22cSMatthew G. Knepley 
84263d5297fSMatthew G. Knepley       ierr = DMPlexGetLabel(dm, b->labelname, &label);CHKERRQ(ierr);
8436c73c22cSMatthew G. Knepley       for (i = 0; i < len; ++i) {
8446c73c22cSMatthew G. Knepley         PetscBool has;
8456c73c22cSMatthew G. Knepley 
8466c73c22cSMatthew G. Knepley         ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr);
8476c73c22cSMatthew G. Knepley         if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name);
8486c73c22cSMatthew G. Knepley       }
8496c73c22cSMatthew G. Knepley       b->numids = len;
8506c73c22cSMatthew G. Knepley       ierr = PetscFree(b->ids);CHKERRQ(ierr);
8516c73c22cSMatthew G. Knepley       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
8526c73c22cSMatthew G. Knepley       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
8536c73c22cSMatthew G. Knepley     }
8546c73c22cSMatthew G. Knepley   }
8556c73c22cSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8560a6ba040SMatthew G. Knepley   /* Handle viewing */
8570a6ba040SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
8580a6ba040SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
8590a6ba040SMatthew G. Knepley   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr);
860*b29cfa1cSToby Isaac   /* Projection behavior */
861*b29cfa1cSToby Isaac   ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr);
86268d4fef7SMatthew G. Knepley   PetscFunctionReturn(0);
86368d4fef7SMatthew G. Knepley }
86468d4fef7SMatthew G. Knepley 
86568d4fef7SMatthew G. Knepley #undef __FUNCT__
86668d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_Plex"
86768d4fef7SMatthew G. Knepley PetscErrorCode  DMSetFromOptions_Plex(DM dm)
86868d4fef7SMatthew G. Knepley {
86968d4fef7SMatthew G. Knepley   PetscInt       refine = 0, r;
87068d4fef7SMatthew G. Knepley   PetscBool      isHierarchy;
87168d4fef7SMatthew G. Knepley   PetscErrorCode ierr;
87268d4fef7SMatthew G. Knepley 
87368d4fef7SMatthew G. Knepley   PetscFunctionBegin;
87468d4fef7SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
87568d4fef7SMatthew G. Knepley   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
87668d4fef7SMatthew G. Knepley   /* Handle DMPlex refinement */
87768d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
87868d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
879b6a0289aSMatthew G. Knepley   if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);}
88068d4fef7SMatthew G. Knepley   if (refine && isHierarchy) {
88168d4fef7SMatthew G. Knepley     DM *dms;
88268d4fef7SMatthew G. Knepley 
88368d4fef7SMatthew G. Knepley     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
88468d4fef7SMatthew G. Knepley     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
88568d4fef7SMatthew G. Knepley     /* Total hack since we do not pass in a pointer */
88668d4fef7SMatthew G. Knepley     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
88768d4fef7SMatthew G. Knepley     if (refine == 1) {
88868d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
88968d4fef7SMatthew G. Knepley     } else {
89068d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
89168d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
89268d4fef7SMatthew G. Knepley     }
89368d4fef7SMatthew G. Knepley     /* Free DMs */
89468d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
89568d4fef7SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
89668d4fef7SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
89768d4fef7SMatthew G. Knepley     }
89868d4fef7SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
89968d4fef7SMatthew G. Knepley   } else {
90068d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
90168d4fef7SMatthew G. Knepley       DM refinedMesh;
90268d4fef7SMatthew G. Knepley 
90368d4fef7SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
90468d4fef7SMatthew G. Knepley       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
90568d4fef7SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
90668d4fef7SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
90768d4fef7SMatthew G. Knepley       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
90868d4fef7SMatthew G. Knepley     }
90968d4fef7SMatthew G. Knepley   }
91068d4fef7SMatthew G. Knepley   ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
9110a6ba040SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
9120a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
9130a6ba040SMatthew G. Knepley }
9140a6ba040SMatthew G. Knepley 
9150a6ba040SMatthew G. Knepley #undef __FUNCT__
916552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
917552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
918552f7358SJed Brown {
919552f7358SJed Brown   PetscErrorCode ierr;
920552f7358SJed Brown 
921552f7358SJed Brown   PetscFunctionBegin;
922552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
923552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
924552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr);
9252c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr);
926552f7358SJed Brown   PetscFunctionReturn(0);
927552f7358SJed Brown }
928552f7358SJed Brown 
929552f7358SJed Brown #undef __FUNCT__
930552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
931552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
932552f7358SJed Brown {
933552f7358SJed Brown   PetscErrorCode ierr;
934552f7358SJed Brown 
935552f7358SJed Brown   PetscFunctionBegin;
936552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
937552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
9382c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr);
939552f7358SJed Brown   PetscFunctionReturn(0);
940552f7358SJed Brown }
941552f7358SJed Brown 
94238221697SMatthew G. Knepley #undef __FUNCT__
943793f3fe5SMatthew G. Knepley #define __FUNCT__ "DMGetDimPoints_Plex"
944793f3fe5SMatthew G. Knepley static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
945793f3fe5SMatthew G. Knepley {
946793f3fe5SMatthew G. Knepley   PetscInt       depth, d;
947793f3fe5SMatthew G. Knepley   PetscErrorCode ierr;
948793f3fe5SMatthew G. Knepley 
949793f3fe5SMatthew G. Knepley   PetscFunctionBegin;
950793f3fe5SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
951793f3fe5SMatthew G. Knepley   if (depth == 1) {
952793f3fe5SMatthew G. Knepley     ierr = DMGetDimension(dm, &d);CHKERRQ(ierr);
953793f3fe5SMatthew G. Knepley     if (dim == 0)      {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);}
954793f3fe5SMatthew G. Knepley     else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);}
955793f3fe5SMatthew G. Knepley     else               {*pStart = 0; *pEnd = 0;}
956793f3fe5SMatthew G. Knepley   } else {
957793f3fe5SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);
958793f3fe5SMatthew G. Knepley   }
959793f3fe5SMatthew G. Knepley   PetscFunctionReturn(0);
960793f3fe5SMatthew G. Knepley }
961793f3fe5SMatthew G. Knepley 
962793f3fe5SMatthew G. Knepley #undef __FUNCT__
963552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
964552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
965552f7358SJed Brown {
966552f7358SJed Brown   PetscFunctionBegin;
967552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
9682c40f234SMatthew G. Knepley   dm->ops->load                            = DMLoad_Plex;
969552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
97038221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
971552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
972fd59a867SMatthew G. Knepley   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
973552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
974552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
975184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
9760298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
977552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
9780a6ba040SMatthew G. Knepley   dm->ops->getcoloring                     = NULL;
979552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
980bceba477SMatthew G. Knepley   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
981bceba477SMatthew G. Knepley   dm->ops->getaggregates                   = NULL;
982bceba477SMatthew G. Knepley   dm->ops->getinjection                    = DMCreateInjection_Plex;
983552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
9840a6ba040SMatthew G. Knepley   dm->ops->coarsen                         = DMCoarsen_Plex;
9850a6ba040SMatthew G. Knepley   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
9860a6ba040SMatthew G. Knepley   dm->ops->coarsenhierarchy                = NULL;
9870298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
9880298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
9890298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
9900298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
991552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
992552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
993793f3fe5SMatthew G. Knepley   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
994552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
995552f7358SJed Brown   PetscFunctionReturn(0);
996552f7358SJed Brown }
997552f7358SJed Brown 
99863a16f15SMatthew G. Knepley #undef __FUNCT__
99963a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex"
100063a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
100163a16f15SMatthew G. Knepley {
100263a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
100363a16f15SMatthew G. Knepley   PetscErrorCode ierr;
100463a16f15SMatthew G. Knepley 
100563a16f15SMatthew G. Knepley   PetscFunctionBegin;
100663a16f15SMatthew G. Knepley   mesh->refct++;
100763a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
100863a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
100963a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
101063a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
101163a16f15SMatthew G. Knepley }
101263a16f15SMatthew G. Knepley 
10138818961aSMatthew G Knepley /*MC
10148818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
10158818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
10168818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
10178818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
10188818961aSMatthew G Knepley 
10198818961aSMatthew G Knepley   Level: intermediate
10208818961aSMatthew G Knepley 
10218818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
10228818961aSMatthew G Knepley M*/
10238818961aSMatthew G Knepley 
1024552f7358SJed Brown #undef __FUNCT__
1025552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
10268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1027552f7358SJed Brown {
1028552f7358SJed Brown   DM_Plex        *mesh;
1029770b213bSMatthew G Knepley   PetscInt       unit, d;
1030552f7358SJed Brown   PetscErrorCode ierr;
1031552f7358SJed Brown 
1032552f7358SJed Brown   PetscFunctionBegin;
1033552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1034b00a9115SJed Brown   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
1035c73cfb54SMatthew G. Knepley   dm->dim  = 0;
1036552f7358SJed Brown   dm->data = mesh;
1037552f7358SJed Brown 
1038552f7358SJed Brown   mesh->refct             = 1;
103982f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
1040552f7358SJed Brown   mesh->maxConeSize       = 0;
10410298fd71SBarry Smith   mesh->cones             = NULL;
10420298fd71SBarry Smith   mesh->coneOrientations  = NULL;
104382f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
1044552f7358SJed Brown   mesh->maxSupportSize    = 0;
10450298fd71SBarry Smith   mesh->supports          = NULL;
1046552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
1047552f7358SJed Brown   mesh->refinementLimit   = -1.0;
1048552f7358SJed Brown 
10490298fd71SBarry Smith   mesh->facesTmp = NULL;
1050552f7358SJed Brown 
10510298fd71SBarry Smith   mesh->subpointMap = NULL;
1052552f7358SJed Brown 
10538865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
1054552f7358SJed Brown 
10550298fd71SBarry Smith   mesh->labels              = NULL;
10568af19771SMatthew G. Knepley   mesh->depthLabel          = NULL;
10570298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
10580298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
10598865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1060552f7358SJed Brown   mesh->vtkCellHeight       = 0;
106170034214SMatthew G. Knepley   mesh->useCone             = PETSC_FALSE;
106270034214SMatthew G. Knepley   mesh->useClosure          = PETSC_TRUE;
1063552f7358SJed Brown 
1064*b29cfa1cSToby Isaac   mesh->maxProjectionHeight = 0;
1065*b29cfa1cSToby Isaac 
1066552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
1067552f7358SJed Brown   mesh->printFEM       = 0;
10686113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
1069552f7358SJed Brown 
1070552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
1071552f7358SJed Brown   PetscFunctionReturn(0);
1072552f7358SJed Brown }
1073552f7358SJed Brown 
1074552f7358SJed Brown #undef __FUNCT__
1075552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
1076552f7358SJed Brown /*@
1077552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1078552f7358SJed Brown 
1079552f7358SJed Brown   Collective on MPI_Comm
1080552f7358SJed Brown 
1081552f7358SJed Brown   Input Parameter:
1082552f7358SJed Brown . comm - The communicator for the DMPlex object
1083552f7358SJed Brown 
1084552f7358SJed Brown   Output Parameter:
1085552f7358SJed Brown . mesh  - The DMPlex object
1086552f7358SJed Brown 
1087552f7358SJed Brown   Level: beginner
1088552f7358SJed Brown 
1089552f7358SJed Brown .keywords: DMPlex, create
1090552f7358SJed Brown @*/
1091552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1092552f7358SJed Brown {
1093552f7358SJed Brown   PetscErrorCode ierr;
1094552f7358SJed Brown 
1095552f7358SJed Brown   PetscFunctionBegin;
1096552f7358SJed Brown   PetscValidPointer(mesh,2);
1097552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1098552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1099552f7358SJed Brown   PetscFunctionReturn(0);
1100552f7358SJed Brown }
1101552f7358SJed Brown 
1102552f7358SJed Brown #undef __FUNCT__
11039298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private"
11049298eaa6SMatthew G Knepley /*
11059298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
11069298eaa6SMatthew G Knepley */
11079298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
11089298eaa6SMatthew G Knepley {
11099298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
11109298eaa6SMatthew G Knepley   PetscErrorCode ierr;
11119298eaa6SMatthew G Knepley 
11129298eaa6SMatthew G Knepley   PetscFunctionBegin;
11139298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
11149298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
11159298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
11169298eaa6SMatthew G Knepley   }
11179298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
11189298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
11199298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
11209298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
11219298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
11229298eaa6SMatthew G Knepley     }
11239298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
11249298eaa6SMatthew G Knepley   }
11259298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
11269298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
11279298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
11289298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
11299298eaa6SMatthew G Knepley }
11309298eaa6SMatthew G Knepley 
11319298eaa6SMatthew G Knepley #undef __FUNCT__
11329298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private"
11339298eaa6SMatthew G Knepley /*
11349298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
11359298eaa6SMatthew G Knepley */
11369298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
11379298eaa6SMatthew G Knepley {
11389298eaa6SMatthew G Knepley   PetscSection   coordSection;
11399298eaa6SMatthew G Knepley   Vec            coordinates;
11409298eaa6SMatthew G Knepley   PetscScalar   *coords;
11419298eaa6SMatthew G Knepley   PetscInt       coordSize, v, d;
11429298eaa6SMatthew G Knepley   PetscErrorCode ierr;
11439298eaa6SMatthew G Knepley 
11449298eaa6SMatthew G Knepley   PetscFunctionBegin;
1145c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
11469298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
11479298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
11489298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
11499298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
11509298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
11519298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
11529298eaa6SMatthew G Knepley   }
11539298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
11549298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
11559298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
11564b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
11579298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
11589298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
11592eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
11609298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
11619298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
11629298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
11639298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
11649298eaa6SMatthew G Knepley     }
11659298eaa6SMatthew G Knepley   }
11669298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
11679298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
11689298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
11699298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
11709298eaa6SMatthew G Knepley }
11719298eaa6SMatthew G Knepley 
11729298eaa6SMatthew G Knepley #undef __FUNCT__
11739298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList"
11749298eaa6SMatthew G Knepley /*@C
11759298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
11769298eaa6SMatthew G Knepley 
11779298eaa6SMatthew G Knepley   Input Parameters:
11789298eaa6SMatthew G Knepley + comm - The communicator
11799298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
11809298eaa6SMatthew G Knepley . numCells - The number of cells
11819298eaa6SMatthew G Knepley . numVertices - The number of vertices
11829298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
11839298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
11849298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
11859298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
11869298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
11879298eaa6SMatthew G Knepley 
11889298eaa6SMatthew G Knepley   Output Parameter:
11899298eaa6SMatthew G Knepley . dm - The DM
11909298eaa6SMatthew G Knepley 
11919298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
11929298eaa6SMatthew G Knepley $
11939298eaa6SMatthew G Knepley $        2
11949298eaa6SMatthew G Knepley $      / | \
11959298eaa6SMatthew G Knepley $     /  |  \
11969298eaa6SMatthew G Knepley $    /   |   \
11979298eaa6SMatthew G Knepley $   0  0 | 1  3
11989298eaa6SMatthew G Knepley $    \   |   /
11999298eaa6SMatthew G Knepley $     \  |  /
12009298eaa6SMatthew G Knepley $      \ | /
12019298eaa6SMatthew G Knepley $        1
12029298eaa6SMatthew G Knepley would have input
12039298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
12049298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
12059298eaa6SMatthew G Knepley $
12069298eaa6SMatthew G Knepley which would result in the DMPlex
12079298eaa6SMatthew G Knepley $
12089298eaa6SMatthew G Knepley $        4
12099298eaa6SMatthew G Knepley $      / | \
12109298eaa6SMatthew G Knepley $     /  |  \
12119298eaa6SMatthew G Knepley $    /   |   \
12129298eaa6SMatthew G Knepley $   2  0 | 1  5
12139298eaa6SMatthew G Knepley $    \   |   /
12149298eaa6SMatthew G Knepley $     \  |  /
12159298eaa6SMatthew G Knepley $      \ | /
12169298eaa6SMatthew G Knepley $        3
12179298eaa6SMatthew G Knepley 
12189298eaa6SMatthew G Knepley   Level: beginner
12199298eaa6SMatthew G Knepley 
1220939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
12219298eaa6SMatthew G Knepley @*/
12229298eaa6SMatthew 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)
12239298eaa6SMatthew G Knepley {
12249298eaa6SMatthew G Knepley   PetscErrorCode ierr;
12259298eaa6SMatthew G Knepley 
12269298eaa6SMatthew G Knepley   PetscFunctionBegin;
12279298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
12289298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1229c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
12309298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
12319298eaa6SMatthew G Knepley   if (interpolate) {
1232e56d480eSMatthew G. Knepley     DM idm = NULL;
12339298eaa6SMatthew G Knepley 
12349298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
12359298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
12369298eaa6SMatthew G Knepley     *dm  = idm;
12379298eaa6SMatthew G Knepley   }
12389298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
12399298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
12409298eaa6SMatthew G Knepley }
12419298eaa6SMatthew G Knepley 
12429298eaa6SMatthew G Knepley #undef __FUNCT__
12439298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG"
1244939f6067SMatthew G. Knepley /*@
1245939f6067SMatthew 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
1246939f6067SMatthew G. Knepley 
1247939f6067SMatthew G. Knepley   Input Parameters:
1248c73cfb54SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
1249939f6067SMatthew G. Knepley . depth - The depth of the DAG
1250939f6067SMatthew G. Knepley . numPoints - The number of points at each depth
1251939f6067SMatthew G. Knepley . coneSize - The cone size of each point
1252939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1253939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point
1254939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1255939f6067SMatthew G. Knepley 
1256939f6067SMatthew G. Knepley   Output Parameter:
1257939f6067SMatthew G. Knepley . dm - The DM
1258939f6067SMatthew G. Knepley 
1259939f6067SMatthew G. Knepley   Note: Two triangles sharing a face would have input
1260939f6067SMatthew G. Knepley $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1261939f6067SMatthew G. Knepley $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1262939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1263939f6067SMatthew G. Knepley $
1264939f6067SMatthew G. Knepley which would result in the DMPlex
1265939f6067SMatthew G. Knepley $
1266939f6067SMatthew G. Knepley $        4
1267939f6067SMatthew G. Knepley $      / | \
1268939f6067SMatthew G. Knepley $     /  |  \
1269939f6067SMatthew G. Knepley $    /   |   \
1270939f6067SMatthew G. Knepley $   2  0 | 1  5
1271939f6067SMatthew G. Knepley $    \   |   /
1272939f6067SMatthew G. Knepley $     \  |  /
1273939f6067SMatthew G. Knepley $      \ | /
1274939f6067SMatthew G. Knepley $        3
1275939f6067SMatthew G. Knepley $
1276939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1277939f6067SMatthew G. Knepley 
1278939f6067SMatthew G. Knepley   Level: advanced
1279939f6067SMatthew G. Knepley 
1280939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1281939f6067SMatthew G. Knepley @*/
12829298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
12839298eaa6SMatthew G Knepley {
12849298eaa6SMatthew G Knepley   Vec            coordinates;
12859298eaa6SMatthew G Knepley   PetscSection   coordSection;
12869298eaa6SMatthew G Knepley   PetscScalar    *coords;
12879298eaa6SMatthew G Knepley   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
12889298eaa6SMatthew G Knepley   PetscErrorCode ierr;
12899298eaa6SMatthew G Knepley 
12909298eaa6SMatthew G Knepley   PetscFunctionBegin;
1291c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
12929298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
12939298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
12949298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
12959298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
12969298eaa6SMatthew G Knepley   }
12979298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
12989298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
12999298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
13009298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
13019298eaa6SMatthew G Knepley   }
13029298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
13039298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
13049298eaa6SMatthew G Knepley   /* Build coordinates */
1305c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
13069298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
13079298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
13089298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
13099298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
13109298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
13116f8cbbeeSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
13129298eaa6SMatthew G Knepley   }
13139298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
13149298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
13159298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
13164b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
13176f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
13189298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
13192eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
13209298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
13219298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
13229298eaa6SMatthew G Knepley     PetscInt off;
13239298eaa6SMatthew G Knepley 
13249298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
13259298eaa6SMatthew G Knepley     for (d = 0; d < dim; ++d) {
13269298eaa6SMatthew G Knepley       coords[off+d] = vertexCoords[v*dim+d];
13279298eaa6SMatthew G Knepley     }
13289298eaa6SMatthew G Knepley   }
13299298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
13309298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
13319298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
13329298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
13339298eaa6SMatthew G Knepley }
1334