xref: /petsc/src/dm/impls/plex/plexcreate.c (revision e56d480e16d5f8d007751dcb65a0c5438fe5de74)
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);
391df5d5c5SMatthew G. Knepley   ierr = DMPlexSetDimension(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) {
120*e56d480eSMatthew 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);
277552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2782eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
279552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
280552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
281552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
282552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
283552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
284552f7358SJed Brown     }
285552f7358SJed Brown   }
286552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
287552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
288552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
289552f7358SJed Brown   PetscFunctionReturn(0);
290552f7358SJed Brown }
291552f7358SJed Brown 
292552f7358SJed Brown #undef __FUNCT__
293552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
29426492d91SMatthew G. Knepley /*@
29526492d91SMatthew G. Knepley   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
296552f7358SJed Brown 
29726492d91SMatthew G. Knepley   Collective on MPI_Comm
29826492d91SMatthew G. Knepley 
29926492d91SMatthew G. Knepley   Input Parameters:
30026492d91SMatthew G. Knepley + comm  - The communicator for the DM object
30126492d91SMatthew G. Knepley . lower - The lower left front corner coordinates
30226492d91SMatthew G. Knepley . upper - The upper right back corner coordinates
30326492d91SMatthew G. Knepley - edges - The number of cells in each direction
30426492d91SMatthew G. Knepley 
30526492d91SMatthew G. Knepley   Output Parameter:
30626492d91SMatthew G. Knepley . dm  - The DM object
30726492d91SMatthew G. Knepley 
30826492d91SMatthew G. Knepley   Level: beginner
30926492d91SMatthew G. Knepley 
31026492d91SMatthew G. Knepley .keywords: DM, create
31126492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
31226492d91SMatthew G. Knepley @*/
313552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
314552f7358SJed Brown {
315552f7358SJed Brown   PetscInt       numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
316552f7358SJed Brown   PetscInt       numFaces    = 6;
317552f7358SJed Brown   Vec            coordinates;
318552f7358SJed Brown   PetscSection   coordSection;
319552f7358SJed Brown   PetscScalar    *coords;
320552f7358SJed Brown   PetscInt       coordSize;
321552f7358SJed Brown   PetscMPIInt    rank;
322552f7358SJed Brown   PetscInt       v, vx, vy, vz;
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   if ((faces[0] > 1) || (faces[1] > 1) || (faces[2] > 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Currently can't handle more than 1 face per side");
32882f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
329552f7358SJed Brown   if (!rank) {
330552f7358SJed Brown     PetscInt f;
331552f7358SJed Brown 
332552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
333552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
334552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
335552f7358SJed Brown     }
336552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
337552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
338552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
339552f7358SJed Brown     }
340552f7358SJed Brown     { /* Side 0 (Front) */
341da80777bSKarl Rupp       PetscInt cone[4];
342da80777bSKarl Rupp       cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6;
343552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr);
344552f7358SJed Brown     }
345552f7358SJed Brown     { /* Side 1 (Back) */
346da80777bSKarl Rupp       PetscInt cone[4];
347da80777bSKarl Rupp       cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3;
348552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr);
349552f7358SJed Brown     }
350552f7358SJed Brown     { /* Side 2 (Bottom) */
351da80777bSKarl Rupp       PetscInt cone[4];
352da80777bSKarl Rupp       cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4;
353552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr);
354552f7358SJed Brown     }
355552f7358SJed Brown     { /* Side 3 (Top) */
356da80777bSKarl Rupp       PetscInt cone[4];
357da80777bSKarl Rupp       cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2;
358552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr);
359552f7358SJed Brown     }
360552f7358SJed Brown     { /* Side 4 (Left) */
361da80777bSKarl Rupp       PetscInt cone[4];
362da80777bSKarl Rupp       cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2;
363552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr);
364552f7358SJed Brown     }
365552f7358SJed Brown     { /* Side 5 (Right) */
366da80777bSKarl Rupp       PetscInt cone[4];
367da80777bSKarl Rupp       cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7;
368552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr);
369552f7358SJed Brown     }
370552f7358SJed Brown   }
371552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
372552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
373552f7358SJed Brown   /* Build coordinates */
374c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
375552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
376552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
377552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
378552f7358SJed Brown   }
379552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
380552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
38182f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
382552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
383552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
3842eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
385552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
386552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
387552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
388552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
389552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
390552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
391552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
392552f7358SJed Brown       }
393552f7358SJed Brown     }
394552f7358SJed Brown   }
395552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
396552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
397552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
398552f7358SJed Brown   PetscFunctionReturn(0);
399552f7358SJed Brown }
400552f7358SJed Brown 
401552f7358SJed Brown #undef __FUNCT__
402552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
40326492d91SMatthew G. Knepley /*@
40426492d91SMatthew G. Knepley   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
405552f7358SJed Brown 
40626492d91SMatthew G. Knepley   Collective on MPI_Comm
40726492d91SMatthew G. Knepley 
40826492d91SMatthew G. Knepley   Input Parameters:
40926492d91SMatthew G. Knepley + comm  - The communicator for the DM object
41026492d91SMatthew G. Knepley . lower - The lower left corner coordinates
41126492d91SMatthew G. Knepley . upper - The upper right corner coordinates
412fbf5b483SMatthew G. Knepley . edges - The number of cells in each direction
413fbf5b483SMatthew G. Knepley . bdX   - The boundary type for the X direction
414fbf5b483SMatthew G. Knepley - bdY   - The boundary type for the Y direction
41526492d91SMatthew G. Knepley 
41626492d91SMatthew G. Knepley   Output Parameter:
41726492d91SMatthew G. Knepley . dm  - The DM object
41826492d91SMatthew G. Knepley 
41926492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
42026492d91SMatthew G. Knepley $ 22--8-23--9--24
42126492d91SMatthew G. Knepley $  |     |     |
42226492d91SMatthew G. Knepley $ 13  2 14  3  15
42326492d91SMatthew G. Knepley $  |     |     |
42426492d91SMatthew G. Knepley $ 19--6-20--7--21
42526492d91SMatthew G. Knepley $  |     |     |
42626492d91SMatthew G. Knepley $ 10  0 11  1 12
42726492d91SMatthew G. Knepley $  |     |     |
42826492d91SMatthew G. Knepley $ 16--4-17--5--18
42926492d91SMatthew G. Knepley 
43026492d91SMatthew G. Knepley   Level: beginner
43126492d91SMatthew G. Knepley 
43226492d91SMatthew G. Knepley .keywords: DM, create
43326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
43426492d91SMatthew G. Knepley @*/
435bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
436552f7358SJed Brown {
437552f7358SJed Brown   PetscInt       markerTop      = 1;
438552f7358SJed Brown   PetscInt       markerBottom   = 1;
439552f7358SJed Brown   PetscInt       markerRight    = 1;
440552f7358SJed Brown   PetscInt       markerLeft     = 1;
441552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
442552f7358SJed Brown   PetscMPIInt    rank;
443552f7358SJed Brown   PetscErrorCode ierr;
444552f7358SJed Brown 
445552f7358SJed Brown   PetscFunctionBegin;
44682f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
4470298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
448552f7358SJed Brown   if (markerSeparate) {
449552f7358SJed Brown     markerTop    = 3;
450552f7358SJed Brown     markerBottom = 1;
451552f7358SJed Brown     markerRight  = 2;
452552f7358SJed Brown     markerLeft   = 4;
453552f7358SJed Brown   }
454552f7358SJed Brown   {
455552f7358SJed Brown     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
456552f7358SJed Brown     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
45746675b8eSMatthew G. Knepley     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
45846675b8eSMatthew G. Knepley     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
459552f7358SJed Brown     const PetscInt numTotXEdges = numXEdges*numYVertices;
460552f7358SJed Brown     const PetscInt numTotYEdges = numYEdges*numXVertices;
461552f7358SJed Brown     const PetscInt numVertices  = numXVertices*numYVertices;
462552f7358SJed Brown     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
463552f7358SJed Brown     const PetscInt numFaces     = numXEdges*numYEdges;
464552f7358SJed Brown     const PetscInt firstVertex  = numFaces;
465552f7358SJed Brown     const PetscInt firstXEdge   = numFaces + numVertices;
466552f7358SJed Brown     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
467552f7358SJed Brown     Vec            coordinates;
468552f7358SJed Brown     PetscSection   coordSection;
469552f7358SJed Brown     PetscScalar   *coords;
470552f7358SJed Brown     PetscInt       coordSize;
471552f7358SJed Brown     PetscInt       v, vx, vy;
472552f7358SJed Brown     PetscInt       f, fx, fy, e, ex, ey;
473552f7358SJed Brown 
474552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
475552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
476552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
477552f7358SJed Brown     }
478552f7358SJed Brown     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
479552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
480552f7358SJed Brown     }
481552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
482552f7358SJed Brown     /* Build faces */
483552f7358SJed Brown     for (fy = 0; fy < numYEdges; fy++) {
484552f7358SJed Brown       for (fx = 0; fx < numXEdges; fx++) {
485552f7358SJed Brown         const PetscInt face    = fy*numXEdges + fx;
486552f7358SJed Brown         const PetscInt edgeL   = firstYEdge +   fx                 *numYEdges + fy;
487becd5721SMatthew G. Knepley         const PetscInt edgeR   = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy;
488552f7358SJed Brown         const PetscInt edgeB   = firstXEdge +   fy                 *numXEdges + fx;
489becd5721SMatthew G. Knepley         const PetscInt edgeT   = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx;
490552f7358SJed Brown         const PetscInt ornt[4] = {0, 0, -2, -2};
491da80777bSKarl Rupp         PetscInt       cone[4];
492552f7358SJed Brown 
493becd5721SMatthew G. Knepley         cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
494552f7358SJed Brown         ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
495552f7358SJed Brown         ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
496552f7358SJed Brown       }
497552f7358SJed Brown     }
498552f7358SJed Brown     /* Build Y edges*/
499552f7358SJed Brown     for (vx = 0; vx < numXVertices; vx++) {
500552f7358SJed Brown       for (ey = 0; ey < numYEdges; ey++) {
50146675b8eSMatthew G. Knepley         const PetscInt nextv   = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx;
502552f7358SJed Brown         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
503becd5721SMatthew G. Knepley         const PetscInt vertexB = firstVertex + ey*numXVertices + vx;
50446675b8eSMatthew G. Knepley         const PetscInt vertexT = firstVertex + nextv;
505da80777bSKarl Rupp         PetscInt       cone[2];
506552f7358SJed Brown 
507becd5721SMatthew G. Knepley         cone[0] = vertexB; cone[1] = vertexT;
508552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
509c4b39bc1SMatthew G. Knepley         if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
510552f7358SJed Brown           if (vx == numXVertices-1) {
511552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
512552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
513552f7358SJed Brown             if (ey == numYEdges-1) {
514552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
515552f7358SJed Brown             }
516552f7358SJed Brown           } else if (vx == 0) {
517552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
518552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
519552f7358SJed Brown             if (ey == numYEdges-1) {
520552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
521552f7358SJed Brown             }
522552f7358SJed Brown           }
523552f7358SJed Brown         }
524552f7358SJed Brown       }
525c4b39bc1SMatthew G. Knepley     }
526552f7358SJed Brown     /* Build X edges*/
527552f7358SJed Brown     for (vy = 0; vy < numYVertices; vy++) {
528552f7358SJed Brown       for (ex = 0; ex < numXEdges; ex++) {
52946675b8eSMatthew G. Knepley         const PetscInt nextv   = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices;
530552f7358SJed Brown         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
531becd5721SMatthew G. Knepley         const PetscInt vertexL = firstVertex + vy*numXVertices + ex;
53246675b8eSMatthew G. Knepley         const PetscInt vertexR = firstVertex + nextv;
533da80777bSKarl Rupp         PetscInt       cone[2];
534552f7358SJed Brown 
535becd5721SMatthew G. Knepley         cone[0] = vertexL; cone[1] = vertexR;
536552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
537c4b39bc1SMatthew G. Knepley         if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
538552f7358SJed Brown           if (vy == numYVertices-1) {
539552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
540552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
541552f7358SJed Brown             if (ex == numXEdges-1) {
542552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
543552f7358SJed Brown             }
544552f7358SJed Brown           } else if (vy == 0) {
545552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
546552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
547552f7358SJed Brown             if (ex == numXEdges-1) {
548552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
549552f7358SJed Brown             }
550552f7358SJed Brown           }
551552f7358SJed Brown         }
552552f7358SJed Brown       }
553c4b39bc1SMatthew G. Knepley     }
554552f7358SJed Brown     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
555552f7358SJed Brown     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
556552f7358SJed Brown     /* Build coordinates */
557c2166f76SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
558552f7358SJed Brown     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
559552f7358SJed Brown     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
560552f7358SJed Brown       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
561552f7358SJed Brown     }
562552f7358SJed Brown     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
563552f7358SJed Brown     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
56482f516ccSBarry Smith     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
565552f7358SJed Brown     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
5662eb5907fSJed Brown     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
567552f7358SJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
568552f7358SJed Brown     for (vy = 0; vy < numYVertices; ++vy) {
569552f7358SJed Brown       for (vx = 0; vx < numXVertices; ++vx) {
570552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
571552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
572552f7358SJed Brown       }
573552f7358SJed Brown     }
574552f7358SJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
575552f7358SJed Brown     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
576552f7358SJed Brown     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
577552f7358SJed Brown   }
578552f7358SJed Brown   PetscFunctionReturn(0);
579552f7358SJed Brown }
580552f7358SJed Brown 
581552f7358SJed Brown #undef __FUNCT__
582552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
5836b75479cSMatthew G Knepley /*@
58426492d91SMatthew G. Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
5856b75479cSMatthew G Knepley 
5866b75479cSMatthew G Knepley   Collective on MPI_Comm
5876b75479cSMatthew G Knepley 
5886b75479cSMatthew G Knepley   Input Parameters:
5896b75479cSMatthew G Knepley + comm - The communicator for the DM object
5906b75479cSMatthew G Knepley . dim - The spatial dimension
5916b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
5926b75479cSMatthew G Knepley 
5936b75479cSMatthew G Knepley   Output Parameter:
5946b75479cSMatthew G Knepley . dm  - The DM object
5956b75479cSMatthew G Knepley 
5966b75479cSMatthew G Knepley   Level: beginner
5976b75479cSMatthew G Knepley 
5986b75479cSMatthew G Knepley .keywords: DM, create
59926492d91SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
6006b75479cSMatthew G Knepley @*/
6010adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
6020adebc6cSBarry Smith {
603552f7358SJed Brown   DM             boundary;
604552f7358SJed Brown   PetscErrorCode ierr;
605552f7358SJed Brown 
606552f7358SJed Brown   PetscFunctionBegin;
607552f7358SJed Brown   PetscValidPointer(dm, 4);
608552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
609552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
610552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
611552f7358SJed Brown   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
612552f7358SJed Brown   switch (dim) {
613552f7358SJed Brown   case 2:
614552f7358SJed Brown   {
615552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
616552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
617552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
618552f7358SJed Brown 
619552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
620552f7358SJed Brown     break;
621552f7358SJed Brown   }
622552f7358SJed Brown   case 3:
623552f7358SJed Brown   {
624552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
625552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
626552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
627552f7358SJed Brown 
628552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
629552f7358SJed Brown     break;
630552f7358SJed Brown   }
631552f7358SJed Brown   default:
632552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
633552f7358SJed Brown   }
6340298fd71SBarry Smith   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
635552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
636552f7358SJed Brown   PetscFunctionReturn(0);
637552f7358SJed Brown }
638552f7358SJed Brown 
639552f7358SJed Brown #undef __FUNCT__
640552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
64126492d91SMatthew G. Knepley /*@
64226492d91SMatthew G. Knepley   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
64326492d91SMatthew G. Knepley 
64426492d91SMatthew G. Knepley   Collective on MPI_Comm
64526492d91SMatthew G. Knepley 
64626492d91SMatthew G. Knepley   Input Parameters:
64726492d91SMatthew G. Knepley + comm  - The communicator for the DM object
64826492d91SMatthew G. Knepley . dim   - The spatial dimension
649fbf5b483SMatthew G. Knepley . periodicX - The boundary type for the X direction
650fbf5b483SMatthew G. Knepley . periodicY - The boundary type for the Y direction
651fbf5b483SMatthew G. Knepley . periodicZ - The boundary type for the Z direction
65226492d91SMatthew G. Knepley - cells - The number of cells in each direction
65326492d91SMatthew G. Knepley 
65426492d91SMatthew G. Knepley   Output Parameter:
65526492d91SMatthew G. Knepley . dm  - The DM object
65626492d91SMatthew G. Knepley 
65726492d91SMatthew G. Knepley   Level: beginner
65826492d91SMatthew G. Knepley 
65926492d91SMatthew G. Knepley .keywords: DM, create
66026492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
66126492d91SMatthew G. Knepley @*/
662bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
663a6dfd86eSKarl Rupp {
664552f7358SJed Brown   PetscErrorCode ierr;
665552f7358SJed Brown 
666552f7358SJed Brown   PetscFunctionBegin;
667552f7358SJed Brown   PetscValidPointer(dm, 4);
668552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
669552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
670552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
671552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
672552f7358SJed Brown   switch (dim) {
673552f7358SJed Brown   case 2:
674552f7358SJed Brown   {
675552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
676552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
677552f7358SJed Brown 
678becd5721SMatthew G. Knepley     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
679552f7358SJed Brown     break;
680552f7358SJed Brown   }
681552f7358SJed Brown #if 0
682552f7358SJed Brown   case 3:
683552f7358SJed Brown   {
684552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
685552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
686552f7358SJed Brown 
687becd5721SMatthew G. Knepley     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
688552f7358SJed Brown     break;
689552f7358SJed Brown   }
690552f7358SJed Brown #endif
691552f7358SJed Brown   default:
692552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
693552f7358SJed Brown   }
694552f7358SJed Brown   PetscFunctionReturn(0);
695552f7358SJed Brown }
696552f7358SJed Brown 
697552f7358SJed Brown /* External function declarations here */
698552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
699bceba477SMatthew G. Knepley extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx);
700fd59a867SMatthew G. Knepley extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
701b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
702552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
703552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
7040a6ba040SMatthew G. Knepley extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
7050a6ba040SMatthew G. Knepley extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
7069a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
707552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
708552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
709552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
710552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
711552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
712552f7358SJed Brown 
713552f7358SJed Brown #undef __FUNCT__
7140a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexReplace_Static"
7150a6ba040SMatthew G. Knepley /* Replace dm with the contents of dmNew
7160a6ba040SMatthew G. Knepley    - Share the DM_Plex structure
7170a6ba040SMatthew G. Knepley    - Share the coordinates
7180a6ba040SMatthew G. Knepley */
7190a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
7200a6ba040SMatthew G. Knepley {
7210a6ba040SMatthew G. Knepley   PetscSection   coordSection;
7220a6ba040SMatthew G. Knepley   Vec            coords;
7230a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
7240a6ba040SMatthew G. Knepley 
7250a6ba040SMatthew G. Knepley   PetscFunctionBegin;
7260a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateSection(dmNew, &coordSection);CHKERRQ(ierr);
7270a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
7280a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, coordSection);CHKERRQ(ierr);
7290a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
7300a6ba040SMatthew G. Knepley   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
7310a6ba040SMatthew G. Knepley   dm->data = dmNew->data;
7320a6ba040SMatthew G. Knepley   ((DM_Plex *) dmNew->data)->refct++;
7330a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
7340a6ba040SMatthew G. Knepley }
7350a6ba040SMatthew G. Knepley 
7360a6ba040SMatthew G. Knepley #undef __FUNCT__
7370a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexSwap_Static"
7380a6ba040SMatthew G. Knepley /* Swap dm with the contents of dmNew
7390a6ba040SMatthew G. Knepley    - Swap the DM_Plex structure
7400a6ba040SMatthew G. Knepley    - Swap the coordinates
7410a6ba040SMatthew G. Knepley */
7420a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
7430a6ba040SMatthew G. Knepley {
7440a6ba040SMatthew G. Knepley   DM             coordDMA, coordDMB;
7450a6ba040SMatthew G. Knepley   Vec            coordsA,  coordsB;
7460a6ba040SMatthew G. Knepley   void          *tmp;
7470a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
7480a6ba040SMatthew G. Knepley 
7490a6ba040SMatthew G. Knepley   PetscFunctionBegin;
7500a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
7510a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
7520a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
7530a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
7540a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
7550a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
7560a6ba040SMatthew G. Knepley 
7570a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
7580a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
7590a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
7600a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
7610a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
7620a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
7630a6ba040SMatthew G. Knepley   tmp       = dmA->data;
7640a6ba040SMatthew G. Knepley   dmA->data = dmB->data;
7650a6ba040SMatthew G. Knepley   dmB->data = tmp;
7660a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
7670a6ba040SMatthew G. Knepley }
7680a6ba040SMatthew G. Knepley 
7690a6ba040SMatthew G. Knepley #undef __FUNCT__
77068d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex"
77168d4fef7SMatthew G. Knepley PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(DM dm)
7720a6ba040SMatthew G. Knepley {
7730a6ba040SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
7746c73c22cSMatthew G. Knepley   DMBoundary     b;
7750a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
7760a6ba040SMatthew G. Knepley 
7770a6ba040SMatthew G. Knepley   PetscFunctionBegin;
7786c73c22cSMatthew G. Knepley   /* Handle boundary conditions */
7796c73c22cSMatthew G. Knepley   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr);
7806c73c22cSMatthew G. Knepley   for (b = mesh->boundary; b; b = b->next) {
7816c73c22cSMatthew G. Knepley     char      optname[1024];
7826c73c22cSMatthew G. Knepley     PetscInt  ids[1024], len = 1024, i;
7836c73c22cSMatthew G. Knepley     PetscBool flg;
7846c73c22cSMatthew G. Knepley 
7856c73c22cSMatthew G. Knepley     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
7866c73c22cSMatthew G. Knepley     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
7876c73c22cSMatthew G. Knepley     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
7886c73c22cSMatthew G. Knepley     if (flg) {
7896c73c22cSMatthew G. Knepley       DMLabel label;
7906c73c22cSMatthew G. Knepley 
7916c73c22cSMatthew G. Knepley       ierr = DMPlexGetLabel(dm, b->name, &label);CHKERRQ(ierr);
7926c73c22cSMatthew G. Knepley       for (i = 0; i < len; ++i) {
7936c73c22cSMatthew G. Knepley         PetscBool has;
7946c73c22cSMatthew G. Knepley 
7956c73c22cSMatthew G. Knepley         ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr);
7966c73c22cSMatthew 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);
7976c73c22cSMatthew G. Knepley       }
7986c73c22cSMatthew G. Knepley       b->numids = len;
7996c73c22cSMatthew G. Knepley       ierr = PetscFree(b->ids);CHKERRQ(ierr);
8006c73c22cSMatthew G. Knepley       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
8016c73c22cSMatthew G. Knepley       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
8026c73c22cSMatthew G. Knepley     }
8036c73c22cSMatthew G. Knepley   }
8046c73c22cSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8050a6ba040SMatthew G. Knepley   /* Handle viewing */
8060a6ba040SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
8070a6ba040SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
8080a6ba040SMatthew G. Knepley   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr);
80968d4fef7SMatthew G. Knepley   PetscFunctionReturn(0);
81068d4fef7SMatthew G. Knepley }
81168d4fef7SMatthew G. Knepley 
81268d4fef7SMatthew G. Knepley #undef __FUNCT__
81368d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_Plex"
81468d4fef7SMatthew G. Knepley PetscErrorCode  DMSetFromOptions_Plex(DM dm)
81568d4fef7SMatthew G. Knepley {
81668d4fef7SMatthew G. Knepley   PetscInt       refine = 0, r;
81768d4fef7SMatthew G. Knepley   PetscBool      isHierarchy;
81868d4fef7SMatthew G. Knepley   PetscErrorCode ierr;
81968d4fef7SMatthew G. Knepley 
82068d4fef7SMatthew G. Knepley   PetscFunctionBegin;
82168d4fef7SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
82268d4fef7SMatthew G. Knepley   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
82368d4fef7SMatthew G. Knepley   /* Handle DMPlex refinement */
82468d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
82568d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
82668d4fef7SMatthew G. Knepley   ierr = DMPlexSetRefinementUniform(dm, refine ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
82768d4fef7SMatthew G. Knepley   if (refine && isHierarchy) {
82868d4fef7SMatthew G. Knepley     DM *dms;
82968d4fef7SMatthew G. Knepley 
83068d4fef7SMatthew G. Knepley     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
83168d4fef7SMatthew G. Knepley     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
83268d4fef7SMatthew G. Knepley     /* Total hack since we do not pass in a pointer */
83368d4fef7SMatthew G. Knepley     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
83468d4fef7SMatthew G. Knepley     if (refine == 1) {
83568d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
83668d4fef7SMatthew G. Knepley     } else {
83768d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
83868d4fef7SMatthew G. Knepley       ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
83968d4fef7SMatthew G. Knepley     }
84068d4fef7SMatthew G. Knepley     /* Free DMs */
84168d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
84268d4fef7SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
84368d4fef7SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
84468d4fef7SMatthew G. Knepley     }
84568d4fef7SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
84668d4fef7SMatthew G. Knepley   } else {
84768d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
84868d4fef7SMatthew G. Knepley       DM refinedMesh;
84968d4fef7SMatthew G. Knepley 
85068d4fef7SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
85168d4fef7SMatthew G. Knepley       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
85268d4fef7SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
85368d4fef7SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
85468d4fef7SMatthew G. Knepley       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
85568d4fef7SMatthew G. Knepley     }
85668d4fef7SMatthew G. Knepley   }
85768d4fef7SMatthew G. Knepley   ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr);
8580a6ba040SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
8590a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
8600a6ba040SMatthew G. Knepley }
8610a6ba040SMatthew G. Knepley 
8620a6ba040SMatthew G. Knepley #undef __FUNCT__
863552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
864552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
865552f7358SJed Brown {
866552f7358SJed Brown   PetscErrorCode ierr;
867552f7358SJed Brown 
868552f7358SJed Brown   PetscFunctionBegin;
869552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
870552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
871552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr);
872552f7358SJed Brown   PetscFunctionReturn(0);
873552f7358SJed Brown }
874552f7358SJed Brown 
875552f7358SJed Brown #undef __FUNCT__
876552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
877552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
878552f7358SJed Brown {
879552f7358SJed Brown   PetscErrorCode ierr;
880552f7358SJed Brown 
881552f7358SJed Brown   PetscFunctionBegin;
882552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
883552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
884552f7358SJed Brown   PetscFunctionReturn(0);
885552f7358SJed Brown }
886552f7358SJed Brown 
88738221697SMatthew G. Knepley #undef __FUNCT__
888552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
889552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
890552f7358SJed Brown {
891552f7358SJed Brown   PetscFunctionBegin;
892552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
893552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
89438221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
895552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
896fd59a867SMatthew G. Knepley   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
897552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
898552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
899184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
900184d77edSJed Brown   dm->ops->getlocaltoglobalmappingblock    = NULL;
9010298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
902552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
9030a6ba040SMatthew G. Knepley   dm->ops->getcoloring                     = NULL;
904552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
905bceba477SMatthew G. Knepley   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
906bceba477SMatthew G. Knepley   dm->ops->getaggregates                   = NULL;
907bceba477SMatthew G. Knepley   dm->ops->getinjection                    = DMCreateInjection_Plex;
908552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
9090a6ba040SMatthew G. Knepley   dm->ops->coarsen                         = DMCoarsen_Plex;
9100a6ba040SMatthew G. Knepley   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
9110a6ba040SMatthew G. Knepley   dm->ops->coarsenhierarchy                = NULL;
9120298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
9130298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
9140298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
9150298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
916552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
917552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
918552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
919552f7358SJed Brown   PetscFunctionReturn(0);
920552f7358SJed Brown }
921552f7358SJed Brown 
92263a16f15SMatthew G. Knepley #undef __FUNCT__
92363a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex"
92463a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
92563a16f15SMatthew G. Knepley {
92663a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
92763a16f15SMatthew G. Knepley   PetscErrorCode ierr;
92863a16f15SMatthew G. Knepley 
92963a16f15SMatthew G. Knepley   PetscFunctionBegin;
93063a16f15SMatthew G. Knepley   mesh->refct++;
93163a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
93263a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
93363a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
93463a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
93563a16f15SMatthew G. Knepley }
93663a16f15SMatthew G. Knepley 
9378818961aSMatthew G Knepley /*MC
9388818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
9398818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
9408818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
9418818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
9428818961aSMatthew G Knepley 
9438818961aSMatthew G Knepley   Level: intermediate
9448818961aSMatthew G Knepley 
9458818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
9468818961aSMatthew G Knepley M*/
9478818961aSMatthew G Knepley 
948552f7358SJed Brown #undef __FUNCT__
949552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
9508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
951552f7358SJed Brown {
952552f7358SJed Brown   DM_Plex        *mesh;
953770b213bSMatthew G Knepley   PetscInt       unit, d;
954552f7358SJed Brown   PetscErrorCode ierr;
955552f7358SJed Brown 
956552f7358SJed Brown   PetscFunctionBegin;
957552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
958b00a9115SJed Brown   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
959552f7358SJed Brown   dm->data = mesh;
960552f7358SJed Brown 
961552f7358SJed Brown   mesh->refct             = 1;
962552f7358SJed Brown   mesh->dim               = 0;
96382f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
964552f7358SJed Brown   mesh->maxConeSize       = 0;
9650298fd71SBarry Smith   mesh->cones             = NULL;
9660298fd71SBarry Smith   mesh->coneOrientations  = NULL;
96782f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
968552f7358SJed Brown   mesh->maxSupportSize    = 0;
9690298fd71SBarry Smith   mesh->supports          = NULL;
970552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
971552f7358SJed Brown   mesh->refinementLimit   = -1.0;
972552f7358SJed Brown 
9730298fd71SBarry Smith   mesh->facesTmp = NULL;
974552f7358SJed Brown 
9750298fd71SBarry Smith   mesh->subpointMap = NULL;
976552f7358SJed Brown 
9778865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
978552f7358SJed Brown 
9790298fd71SBarry Smith   mesh->labels              = NULL;
9808af19771SMatthew G. Knepley   mesh->depthLabel          = NULL;
9810298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
9820298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
9838865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
984552f7358SJed Brown   mesh->vtkCellHeight       = 0;
98570034214SMatthew G. Knepley   mesh->useCone             = PETSC_FALSE;
98670034214SMatthew G. Knepley   mesh->useClosure          = PETSC_TRUE;
987552f7358SJed Brown 
988552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
989552f7358SJed Brown   mesh->printFEM       = 0;
9906113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
991552f7358SJed Brown 
992552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
993552f7358SJed Brown   PetscFunctionReturn(0);
994552f7358SJed Brown }
995552f7358SJed Brown 
996552f7358SJed Brown #undef __FUNCT__
997552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
998552f7358SJed Brown /*@
999552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
1000552f7358SJed Brown 
1001552f7358SJed Brown   Collective on MPI_Comm
1002552f7358SJed Brown 
1003552f7358SJed Brown   Input Parameter:
1004552f7358SJed Brown . comm - The communicator for the DMPlex object
1005552f7358SJed Brown 
1006552f7358SJed Brown   Output Parameter:
1007552f7358SJed Brown . mesh  - The DMPlex object
1008552f7358SJed Brown 
1009552f7358SJed Brown   Level: beginner
1010552f7358SJed Brown 
1011552f7358SJed Brown .keywords: DMPlex, create
1012552f7358SJed Brown @*/
1013552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1014552f7358SJed Brown {
1015552f7358SJed Brown   PetscErrorCode ierr;
1016552f7358SJed Brown 
1017552f7358SJed Brown   PetscFunctionBegin;
1018552f7358SJed Brown   PetscValidPointer(mesh,2);
1019552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
1020552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
1021552f7358SJed Brown   PetscFunctionReturn(0);
1022552f7358SJed Brown }
1023552f7358SJed Brown 
1024552f7358SJed Brown #undef __FUNCT__
10259298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private"
10269298eaa6SMatthew G Knepley /*
10279298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
10289298eaa6SMatthew G Knepley */
10299298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
10309298eaa6SMatthew G Knepley {
10319298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
10329298eaa6SMatthew G Knepley   PetscErrorCode ierr;
10339298eaa6SMatthew G Knepley 
10349298eaa6SMatthew G Knepley   PetscFunctionBegin;
10359298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
10369298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
10379298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
10389298eaa6SMatthew G Knepley   }
10399298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
10409298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
10419298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
10429298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
10439298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
10449298eaa6SMatthew G Knepley     }
10459298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
10469298eaa6SMatthew G Knepley   }
10479298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
10489298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
10499298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
10509298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
10519298eaa6SMatthew G Knepley }
10529298eaa6SMatthew G Knepley 
10539298eaa6SMatthew G Knepley #undef __FUNCT__
10549298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private"
10559298eaa6SMatthew G Knepley /*
10569298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
10579298eaa6SMatthew G Knepley */
10589298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
10599298eaa6SMatthew G Knepley {
10609298eaa6SMatthew G Knepley   PetscSection   coordSection;
10619298eaa6SMatthew G Knepley   Vec            coordinates;
10629298eaa6SMatthew G Knepley   PetscScalar   *coords;
10639298eaa6SMatthew G Knepley   PetscInt       coordSize, v, d;
10649298eaa6SMatthew G Knepley   PetscErrorCode ierr;
10659298eaa6SMatthew G Knepley 
10669298eaa6SMatthew G Knepley   PetscFunctionBegin;
1067c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
10689298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
10699298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
10709298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
10719298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
10729298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
10739298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
10749298eaa6SMatthew G Knepley   }
10759298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
10769298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
10779298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
10789298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
10799298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
10802eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
10819298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
10829298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
10839298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
10849298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
10859298eaa6SMatthew G Knepley     }
10869298eaa6SMatthew G Knepley   }
10879298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
10889298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
10899298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
10909298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
10919298eaa6SMatthew G Knepley }
10929298eaa6SMatthew G Knepley 
10939298eaa6SMatthew G Knepley #undef __FUNCT__
10949298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList"
10959298eaa6SMatthew G Knepley /*@C
10969298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
10979298eaa6SMatthew G Knepley 
10989298eaa6SMatthew G Knepley   Input Parameters:
10999298eaa6SMatthew G Knepley + comm - The communicator
11009298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
11019298eaa6SMatthew G Knepley . numCells - The number of cells
11029298eaa6SMatthew G Knepley . numVertices - The number of vertices
11039298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
11049298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
11059298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
11069298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
11079298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
11089298eaa6SMatthew G Knepley 
11099298eaa6SMatthew G Knepley   Output Parameter:
11109298eaa6SMatthew G Knepley . dm - The DM
11119298eaa6SMatthew G Knepley 
11129298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
11139298eaa6SMatthew G Knepley $
11149298eaa6SMatthew G Knepley $        2
11159298eaa6SMatthew G Knepley $      / | \
11169298eaa6SMatthew G Knepley $     /  |  \
11179298eaa6SMatthew G Knepley $    /   |   \
11189298eaa6SMatthew G Knepley $   0  0 | 1  3
11199298eaa6SMatthew G Knepley $    \   |   /
11209298eaa6SMatthew G Knepley $     \  |  /
11219298eaa6SMatthew G Knepley $      \ | /
11229298eaa6SMatthew G Knepley $        1
11239298eaa6SMatthew G Knepley would have input
11249298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
11259298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
11269298eaa6SMatthew G Knepley $
11279298eaa6SMatthew G Knepley which would result in the DMPlex
11289298eaa6SMatthew G Knepley $
11299298eaa6SMatthew G Knepley $        4
11309298eaa6SMatthew G Knepley $      / | \
11319298eaa6SMatthew G Knepley $     /  |  \
11329298eaa6SMatthew G Knepley $    /   |   \
11339298eaa6SMatthew G Knepley $   2  0 | 1  5
11349298eaa6SMatthew G Knepley $    \   |   /
11359298eaa6SMatthew G Knepley $     \  |  /
11369298eaa6SMatthew G Knepley $      \ | /
11379298eaa6SMatthew G Knepley $        3
11389298eaa6SMatthew G Knepley 
11399298eaa6SMatthew G Knepley   Level: beginner
11409298eaa6SMatthew G Knepley 
1141939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
11429298eaa6SMatthew G Knepley @*/
11439298eaa6SMatthew 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)
11449298eaa6SMatthew G Knepley {
11459298eaa6SMatthew G Knepley   PetscErrorCode ierr;
11469298eaa6SMatthew G Knepley 
11479298eaa6SMatthew G Knepley   PetscFunctionBegin;
11489298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
11499298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
11509298eaa6SMatthew G Knepley   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
11519298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
11529298eaa6SMatthew G Knepley   if (interpolate) {
1153*e56d480eSMatthew G. Knepley     DM idm = NULL;
11549298eaa6SMatthew G Knepley 
11559298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
11569298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
11579298eaa6SMatthew G Knepley     *dm  = idm;
11589298eaa6SMatthew G Knepley   }
11599298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
11609298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
11619298eaa6SMatthew G Knepley }
11629298eaa6SMatthew G Knepley 
11639298eaa6SMatthew G Knepley #undef __FUNCT__
11649298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG"
1165939f6067SMatthew G. Knepley /*@
1166939f6067SMatthew 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
1167939f6067SMatthew G. Knepley 
1168939f6067SMatthew G. Knepley   Input Parameters:
1169939f6067SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension()
1170939f6067SMatthew G. Knepley . depth - The depth of the DAG
1171939f6067SMatthew G. Knepley . numPoints - The number of points at each depth
1172939f6067SMatthew G. Knepley . coneSize - The cone size of each point
1173939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1174939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point
1175939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1176939f6067SMatthew G. Knepley 
1177939f6067SMatthew G. Knepley   Output Parameter:
1178939f6067SMatthew G. Knepley . dm - The DM
1179939f6067SMatthew G. Knepley 
1180939f6067SMatthew G. Knepley   Note: Two triangles sharing a face would have input
1181939f6067SMatthew G. Knepley $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1182939f6067SMatthew G. Knepley $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1183939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1184939f6067SMatthew G. Knepley $
1185939f6067SMatthew G. Knepley which would result in the DMPlex
1186939f6067SMatthew G. Knepley $
1187939f6067SMatthew G. Knepley $        4
1188939f6067SMatthew G. Knepley $      / | \
1189939f6067SMatthew G. Knepley $     /  |  \
1190939f6067SMatthew G. Knepley $    /   |   \
1191939f6067SMatthew G. Knepley $   2  0 | 1  5
1192939f6067SMatthew G. Knepley $    \   |   /
1193939f6067SMatthew G. Knepley $     \  |  /
1194939f6067SMatthew G. Knepley $      \ | /
1195939f6067SMatthew G. Knepley $        3
1196939f6067SMatthew G. Knepley $
1197939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1198939f6067SMatthew G. Knepley 
1199939f6067SMatthew G. Knepley   Level: advanced
1200939f6067SMatthew G. Knepley 
1201939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1202939f6067SMatthew G. Knepley @*/
12039298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
12049298eaa6SMatthew G Knepley {
12059298eaa6SMatthew G Knepley   Vec            coordinates;
12069298eaa6SMatthew G Knepley   PetscSection   coordSection;
12079298eaa6SMatthew G Knepley   PetscScalar    *coords;
12089298eaa6SMatthew G Knepley   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
12099298eaa6SMatthew G Knepley   PetscErrorCode ierr;
12109298eaa6SMatthew G Knepley 
12119298eaa6SMatthew G Knepley   PetscFunctionBegin;
12129298eaa6SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
12139298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
12149298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
12159298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
12169298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
12179298eaa6SMatthew G Knepley   }
12189298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
12199298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
12209298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
12219298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
12229298eaa6SMatthew G Knepley   }
12239298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
12249298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
12259298eaa6SMatthew G Knepley   /* Build coordinates */
1226c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
12279298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
12289298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
12299298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
12309298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
12319298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
12326f8cbbeeSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
12339298eaa6SMatthew G Knepley   }
12349298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
12359298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
12369298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
12376f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
12389298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
12392eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
12409298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
12419298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
12429298eaa6SMatthew G Knepley     PetscInt off;
12439298eaa6SMatthew G Knepley 
12449298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
12459298eaa6SMatthew G Knepley     for (d = 0; d < dim; ++d) {
12469298eaa6SMatthew G Knepley       coords[off+d] = vertexCoords[v*dim+d];
12479298eaa6SMatthew G Knepley     }
12489298eaa6SMatthew G Knepley   }
12499298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
12509298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
12519298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
12529298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
12539298eaa6SMatthew G Knepley }
1254