xref: /petsc/src/dm/impls/plex/plexcreate.c (revision 07104863d5d9dd167e823ec2293ba6df7a2d7e67)
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__
7552f7358SJed Brown #define __FUNCT__ "DMSetFromOptions_Plex"
8552f7358SJed Brown PetscErrorCode  DMSetFromOptions_Plex(DM dm)
9552f7358SJed Brown {
10552f7358SJed Brown   DM_Plex        *mesh = (DM_Plex*) dm->data;
11552f7358SJed Brown   PetscErrorCode ierr;
12552f7358SJed Brown 
13552f7358SJed Brown   PetscFunctionBegin;
14552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15552f7358SJed Brown   ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr);
16552f7358SJed Brown   /* Handle DMPlex refinement */
17552f7358SJed Brown   /* Handle associated vectors */
18552f7358SJed Brown   /* Handle viewing */
190298fd71SBarry Smith   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
200298fd71SBarry Smith   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
216113b454SMatthew G. Knepley   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr);
22552f7358SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
23552f7358SJed Brown   PetscFunctionReturn(0);
24552f7358SJed Brown }
25552f7358SJed Brown 
26552f7358SJed Brown #undef __FUNCT__
271df5d5c5SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateDoublet"
281df5d5c5SMatthew G. Knepley /*@
291df5d5c5SMatthew G. Knepley   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
301df5d5c5SMatthew G. Knepley 
311df5d5c5SMatthew G. Knepley   Collective on MPI_Comm
321df5d5c5SMatthew G. Knepley 
331df5d5c5SMatthew G. Knepley   Input Parameters:
341df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object
351df5d5c5SMatthew G. Knepley . dim - The spatial dimension
361df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells
371df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces)
381df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement
391df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
401df5d5c5SMatthew G. Knepley 
411df5d5c5SMatthew G. Knepley   Output Parameter:
421df5d5c5SMatthew G. Knepley . dm  - The DM object
431df5d5c5SMatthew G. Knepley 
441df5d5c5SMatthew G. Knepley   Level: beginner
451df5d5c5SMatthew G. Knepley 
461df5d5c5SMatthew G. Knepley .keywords: DM, create
471df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate()
481df5d5c5SMatthew G. Knepley @*/
491df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
501df5d5c5SMatthew G. Knepley {
511df5d5c5SMatthew G. Knepley   DM             dm;
521df5d5c5SMatthew G. Knepley   PetscInt       p;
531df5d5c5SMatthew G. Knepley   PetscMPIInt    rank;
541df5d5c5SMatthew G. Knepley   PetscErrorCode ierr;
551df5d5c5SMatthew G. Knepley 
561df5d5c5SMatthew G. Knepley   PetscFunctionBegin;
571df5d5c5SMatthew G. Knepley   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
581df5d5c5SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
591df5d5c5SMatthew G. Knepley   ierr = DMPlexSetDimension(dm, dim);CHKERRQ(ierr);
601df5d5c5SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
611df5d5c5SMatthew G. Knepley   if (rank) {
621df5d5c5SMatthew G. Knepley     PetscInt numPoints[2] = {0, 0};
631df5d5c5SMatthew G. Knepley     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
641df5d5c5SMatthew G. Knepley   } else {
651df5d5c5SMatthew G. Knepley     switch (dim) {
661df5d5c5SMatthew G. Knepley     case 2:
671df5d5c5SMatthew G. Knepley       if (simplex) {
681df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {4, 2};
691df5d5c5SMatthew G. Knepley         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
701df5d5c5SMatthew G. Knepley         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
711df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
721df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
731df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
741df5d5c5SMatthew G. Knepley 
751df5d5c5SMatthew G. Knepley         ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);
761df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
771df5d5c5SMatthew G. Knepley         for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
781df5d5c5SMatthew G. Knepley       } else {
791df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {6, 2};
801df5d5c5SMatthew G. Knepley         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
811df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
821df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
831df5d5c5SMatthew 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};
841df5d5c5SMatthew G. Knepley 
851df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
861df5d5c5SMatthew G. Knepley       }
871df5d5c5SMatthew G. Knepley       break;
881df5d5c5SMatthew G. Knepley     case 3:
891df5d5c5SMatthew G. Knepley       if (simplex) {
901df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {5, 2};
911df5d5c5SMatthew G. Knepley         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
921df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
931df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
941df5d5c5SMatthew 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};
951df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
961df5d5c5SMatthew G. Knepley 
971df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
981df5d5c5SMatthew G. Knepley         for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
991df5d5c5SMatthew G. Knepley       } else {
1001df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]         = {12, 2};
1011df5d5c5SMatthew G. Knepley         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
1021df5d5c5SMatthew G. Knepley         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
1031df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
1041df5d5c5SMatthew 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,
1051df5d5c5SMatthew 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,
1061df5d5c5SMatthew 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};
1071df5d5c5SMatthew G. Knepley 
1081df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
1091df5d5c5SMatthew G. Knepley       }
1101df5d5c5SMatthew G. Knepley       break;
1111df5d5c5SMatthew G. Knepley     default:
1121df5d5c5SMatthew G. Knepley       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
1131df5d5c5SMatthew G. Knepley     }
1141df5d5c5SMatthew G. Knepley   }
1151df5d5c5SMatthew G. Knepley   *newdm = dm;
1161df5d5c5SMatthew G. Knepley   if (refinementLimit > 0.0) {
1171df5d5c5SMatthew G. Knepley     DM rdm;
1181df5d5c5SMatthew G. Knepley     const char *name;
1191df5d5c5SMatthew G. Knepley 
1201df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
1211df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
1221df5d5c5SMatthew G. Knepley     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
1231df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1241df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
1251df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1261df5d5c5SMatthew G. Knepley     *newdm = rdm;
1271df5d5c5SMatthew G. Knepley   }
1281df5d5c5SMatthew G. Knepley   if (interpolate) {
1291df5d5c5SMatthew G. Knepley     DM idm;
1301df5d5c5SMatthew G. Knepley     const char *name;
1311df5d5c5SMatthew G. Knepley 
1321df5d5c5SMatthew G. Knepley     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
1331df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1341df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
1351df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
1361df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr);
1371df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1381df5d5c5SMatthew G. Knepley     *newdm = idm;
1391df5d5c5SMatthew G. Knepley   }
1401df5d5c5SMatthew G. Knepley   {
1411df5d5c5SMatthew G. Knepley     DM refinedMesh     = NULL;
1421df5d5c5SMatthew G. Knepley     DM distributedMesh = NULL;
1431df5d5c5SMatthew G. Knepley 
1441df5d5c5SMatthew G. Knepley     /* Distribute mesh over processes */
145*07104863SJed Brown     ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr);
1461df5d5c5SMatthew G. Knepley     if (distributedMesh) {
1471df5d5c5SMatthew G. Knepley       ierr = DMDestroy(newdm);CHKERRQ(ierr);
1481df5d5c5SMatthew G. Knepley       *newdm = distributedMesh;
1491df5d5c5SMatthew G. Knepley     }
1501df5d5c5SMatthew G. Knepley     if (refinementUniform) {
1511df5d5c5SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
1521df5d5c5SMatthew G. Knepley       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
1531df5d5c5SMatthew G. Knepley       if (refinedMesh) {
1541df5d5c5SMatthew G. Knepley         ierr = DMDestroy(newdm);CHKERRQ(ierr);
1551df5d5c5SMatthew G. Knepley         *newdm = refinedMesh;
1561df5d5c5SMatthew G. Knepley       }
1571df5d5c5SMatthew G. Knepley     }
1581df5d5c5SMatthew G. Knepley   }
1591df5d5c5SMatthew G. Knepley   PetscFunctionReturn(0);
1601df5d5c5SMatthew G. Knepley }
1611df5d5c5SMatthew G. Knepley 
1621df5d5c5SMatthew G. Knepley #undef __FUNCT__
163552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary"
164552f7358SJed Brown /*
165552f7358SJed Brown  Simple square boundary:
166552f7358SJed Brown 
167552f7358SJed Brown  18--5-17--4--16
168552f7358SJed Brown   |     |     |
169552f7358SJed Brown   6    10     3
170552f7358SJed Brown   |     |     |
171552f7358SJed Brown  19-11-20--9--15
172552f7358SJed Brown   |     |     |
173552f7358SJed Brown   7     8     2
174552f7358SJed Brown   |     |     |
175552f7358SJed Brown  12--0-13--1--14
176552f7358SJed Brown */
177552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
178552f7358SJed Brown {
179552f7358SJed Brown   PetscInt       numVertices    = (edges[0]+1)*(edges[1]+1);
180552f7358SJed Brown   PetscInt       numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
181552f7358SJed Brown   PetscInt       markerTop      = 1;
182552f7358SJed Brown   PetscInt       markerBottom   = 1;
183552f7358SJed Brown   PetscInt       markerRight    = 1;
184552f7358SJed Brown   PetscInt       markerLeft     = 1;
185552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
186552f7358SJed Brown   Vec            coordinates;
187552f7358SJed Brown   PetscSection   coordSection;
188552f7358SJed Brown   PetscScalar    *coords;
189552f7358SJed Brown   PetscInt       coordSize;
190552f7358SJed Brown   PetscMPIInt    rank;
191552f7358SJed Brown   PetscInt       v, vx, vy;
192552f7358SJed Brown   PetscErrorCode ierr;
193552f7358SJed Brown 
194552f7358SJed Brown   PetscFunctionBegin;
1950298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
196552f7358SJed Brown   if (markerSeparate) {
197552f7358SJed Brown     markerTop    = 1;
198552f7358SJed Brown     markerBottom = 0;
199552f7358SJed Brown     markerRight  = 0;
200552f7358SJed Brown     markerLeft   = 0;
201552f7358SJed Brown   }
20282f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
203552f7358SJed Brown   if (!rank) {
204552f7358SJed Brown     PetscInt e, ex, ey;
205552f7358SJed Brown 
206552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
207552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
208552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
209552f7358SJed Brown     }
210552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
211552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
212552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
213552f7358SJed Brown         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
214552f7358SJed Brown         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
215da80777bSKarl Rupp         PetscInt cone[2];
216552f7358SJed Brown 
217da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
218552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
219552f7358SJed Brown         if (vx == edges[0]) {
220552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
221552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
222552f7358SJed Brown           if (ey == edges[1]-1) {
223552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
224552f7358SJed Brown           }
225552f7358SJed Brown         } else if (vx == 0) {
226552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
227552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
228552f7358SJed Brown           if (ey == edges[1]-1) {
229552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
230552f7358SJed Brown           }
231552f7358SJed Brown         }
232552f7358SJed Brown       }
233552f7358SJed Brown     }
234552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
235552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
236552f7358SJed Brown         PetscInt edge   = vy*edges[0]     + ex;
237552f7358SJed Brown         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
238da80777bSKarl Rupp         PetscInt cone[2];
239552f7358SJed Brown 
240da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
241552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
242552f7358SJed Brown         if (vy == edges[1]) {
243552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
244552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
245552f7358SJed Brown           if (ex == edges[0]-1) {
246552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
247552f7358SJed Brown           }
248552f7358SJed Brown         } else if (vy == 0) {
249552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
250552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
251552f7358SJed Brown           if (ex == edges[0]-1) {
252552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
253552f7358SJed Brown           }
254552f7358SJed Brown         }
255552f7358SJed Brown       }
256552f7358SJed Brown     }
257552f7358SJed Brown   }
258552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
259552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
260552f7358SJed Brown   /* Build coordinates */
261552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
262552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
263552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
264552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
265552f7358SJed Brown   }
266552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
267552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
26882f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
269552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
270552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
271552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
272552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
273552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
274552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
275552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
276552f7358SJed Brown     }
277552f7358SJed Brown   }
278552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
279552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
280552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
281552f7358SJed Brown   PetscFunctionReturn(0);
282552f7358SJed Brown }
283552f7358SJed Brown 
284552f7358SJed Brown #undef __FUNCT__
285552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
286552f7358SJed Brown /*
287552f7358SJed Brown  Simple cubic boundary:
288552f7358SJed Brown 
289552f7358SJed Brown      2-------3
290552f7358SJed Brown     /|      /|
291552f7358SJed Brown    6-------7 |
292552f7358SJed Brown    | |     | |
293552f7358SJed Brown    | 0-----|-1
294552f7358SJed Brown    |/      |/
295552f7358SJed Brown    4-------5
296552f7358SJed Brown */
297552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
298552f7358SJed Brown {
299552f7358SJed Brown   PetscInt       numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
300552f7358SJed Brown   PetscInt       numFaces    = 6;
301552f7358SJed Brown   Vec            coordinates;
302552f7358SJed Brown   PetscSection   coordSection;
303552f7358SJed Brown   PetscScalar    *coords;
304552f7358SJed Brown   PetscInt       coordSize;
305552f7358SJed Brown   PetscMPIInt    rank;
306552f7358SJed Brown   PetscInt       v, vx, vy, vz;
307552f7358SJed Brown   PetscErrorCode ierr;
308552f7358SJed Brown 
309552f7358SJed Brown   PetscFunctionBegin;
31082f516ccSBarry 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");
31182f516ccSBarry 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");
312552f7358SJed Brown   ierr = PetscMalloc(numVertices*2 * sizeof(PetscReal), &coords);CHKERRQ(ierr);
31382f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
314552f7358SJed Brown   if (!rank) {
315552f7358SJed Brown     PetscInt f;
316552f7358SJed Brown 
317552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
318552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
319552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
320552f7358SJed Brown     }
321552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
322552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
323552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
324552f7358SJed Brown     }
325552f7358SJed Brown     { /* Side 0 (Front) */
326da80777bSKarl Rupp       PetscInt cone[4];
327da80777bSKarl Rupp       cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6;
328552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr);
329552f7358SJed Brown     }
330552f7358SJed Brown     { /* Side 1 (Back) */
331da80777bSKarl Rupp       PetscInt cone[4];
332da80777bSKarl Rupp       cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3;
333552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr);
334552f7358SJed Brown     }
335552f7358SJed Brown     { /* Side 2 (Bottom) */
336da80777bSKarl Rupp       PetscInt cone[4];
337da80777bSKarl Rupp       cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4;
338552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr);
339552f7358SJed Brown     }
340552f7358SJed Brown     { /* Side 3 (Top) */
341da80777bSKarl Rupp       PetscInt cone[4];
342da80777bSKarl Rupp       cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2;
343552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr);
344552f7358SJed Brown     }
345552f7358SJed Brown     { /* Side 4 (Left) */
346da80777bSKarl Rupp       PetscInt cone[4];
347da80777bSKarl Rupp       cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2;
348552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr);
349552f7358SJed Brown     }
350552f7358SJed Brown     { /* Side 5 (Right) */
351da80777bSKarl Rupp       PetscInt cone[4];
352da80777bSKarl Rupp       cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7;
353552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr);
354552f7358SJed Brown     }
355552f7358SJed Brown   }
356552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
357552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
358552f7358SJed Brown   /* Build coordinates */
359552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
360552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
361552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
362552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
363552f7358SJed Brown   }
364552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
365552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
36682f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
367552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
368552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
369552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
370552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
371552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
372552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
373552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
374552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
375552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
376552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
377552f7358SJed Brown       }
378552f7358SJed Brown     }
379552f7358SJed Brown   }
380552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
381552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
382552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
383552f7358SJed Brown   PetscFunctionReturn(0);
384552f7358SJed Brown }
385552f7358SJed Brown 
386552f7358SJed Brown #undef __FUNCT__
387552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
388552f7358SJed Brown /*
389552f7358SJed Brown  Simple square mesh:
390552f7358SJed Brown 
391552f7358SJed Brown  22--8-23--9--24
392552f7358SJed Brown   |     |     |
393552f7358SJed Brown  13  2 14  3  15
394552f7358SJed Brown   |     |     |
395552f7358SJed Brown  19--6-20--7--21
396552f7358SJed Brown   |     |     |
397552f7358SJed Brown  10  0 11  1 12
398552f7358SJed Brown   |     |     |
399552f7358SJed Brown  16--4-17--5--18
400552f7358SJed Brown */
401552f7358SJed Brown PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
402552f7358SJed Brown {
403552f7358SJed Brown   PetscInt       markerTop      = 1;
404552f7358SJed Brown   PetscInt       markerBottom   = 1;
405552f7358SJed Brown   PetscInt       markerRight    = 1;
406552f7358SJed Brown   PetscInt       markerLeft     = 1;
407552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
408552f7358SJed Brown   PetscMPIInt    rank;
409552f7358SJed Brown   PetscErrorCode ierr;
410552f7358SJed Brown 
411552f7358SJed Brown   PetscFunctionBegin;
41282f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
4130298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
414552f7358SJed Brown   if (markerSeparate) {
415552f7358SJed Brown     markerTop    = 3;
416552f7358SJed Brown     markerBottom = 1;
417552f7358SJed Brown     markerRight  = 2;
418552f7358SJed Brown     markerLeft   = 4;
419552f7358SJed Brown   }
420552f7358SJed Brown   {
421552f7358SJed Brown     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
422552f7358SJed Brown     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
423552f7358SJed Brown     const PetscInt numXVertices = !rank ? edges[0]+1 : 0;
424552f7358SJed Brown     const PetscInt numYVertices = !rank ? edges[1]+1 : 0;
425552f7358SJed Brown     const PetscInt numTotXEdges = numXEdges*numYVertices;
426552f7358SJed Brown     const PetscInt numTotYEdges = numYEdges*numXVertices;
427552f7358SJed Brown     const PetscInt numVertices  = numXVertices*numYVertices;
428552f7358SJed Brown     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
429552f7358SJed Brown     const PetscInt numFaces     = numXEdges*numYEdges;
430552f7358SJed Brown     const PetscInt firstVertex  = numFaces;
431552f7358SJed Brown     const PetscInt firstXEdge   = numFaces + numVertices;
432552f7358SJed Brown     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
433552f7358SJed Brown     Vec            coordinates;
434552f7358SJed Brown     PetscSection   coordSection;
435552f7358SJed Brown     PetscScalar    *coords;
436552f7358SJed Brown     PetscInt       coordSize;
437552f7358SJed Brown     PetscInt       v, vx, vy;
438552f7358SJed Brown     PetscInt       f, fx, fy, e, ex, ey;
439552f7358SJed Brown 
440552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
441552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
442552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
443552f7358SJed Brown     }
444552f7358SJed Brown     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
445552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
446552f7358SJed Brown     }
447552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
448552f7358SJed Brown     /* Build faces */
449552f7358SJed Brown     for (fy = 0; fy < numYEdges; fy++) {
450552f7358SJed Brown       for (fx = 0; fx < numXEdges; fx++) {
451552f7358SJed Brown         const PetscInt face    = fy*numXEdges + fx;
452552f7358SJed Brown         const PetscInt edgeL   = firstYEdge + fx*numYEdges + fy;
453552f7358SJed Brown         const PetscInt edgeB   = firstXEdge + fy*numXEdges + fx;
454552f7358SJed Brown         const PetscInt ornt[4] = {0,     0,               -2,              -2};
455da80777bSKarl Rupp         PetscInt       cone[4];
456552f7358SJed Brown 
457da80777bSKarl Rupp         cone[0] = edgeB; cone[1] = edgeL+numYEdges; cone[2] = edgeB+numXEdges; cone[3] = edgeL;
458552f7358SJed Brown         ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
459552f7358SJed Brown         ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
460552f7358SJed Brown       }
461552f7358SJed Brown     }
462552f7358SJed Brown     /* Build Y edges*/
463552f7358SJed Brown     for (vx = 0; vx < numXVertices; vx++) {
464552f7358SJed Brown       for (ey = 0; ey < numYEdges; ey++) {
465552f7358SJed Brown         const PetscInt edge   = firstYEdge  + vx*numYEdges + ey;
466552f7358SJed Brown         const PetscInt vertex = firstVertex + ey*numXVertices + vx;
467da80777bSKarl Rupp         PetscInt       cone[2];
468552f7358SJed Brown 
469da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+numXVertices;
470552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
471552f7358SJed Brown         if (vx == numXVertices-1) {
472552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
473552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
474552f7358SJed Brown           if (ey == numYEdges-1) {
475552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
476552f7358SJed Brown           }
477552f7358SJed Brown         } else if (vx == 0) {
478552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
479552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
480552f7358SJed Brown           if (ey == numYEdges-1) {
481552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
482552f7358SJed Brown           }
483552f7358SJed Brown         }
484552f7358SJed Brown       }
485552f7358SJed Brown     }
486552f7358SJed Brown     /* Build X edges*/
487552f7358SJed Brown     for (vy = 0; vy < numYVertices; vy++) {
488552f7358SJed Brown       for (ex = 0; ex < numXEdges; ex++) {
489552f7358SJed Brown         const PetscInt edge   = firstXEdge  + vy*numXEdges + ex;
490552f7358SJed Brown         const PetscInt vertex = firstVertex + vy*numXVertices + ex;
491da80777bSKarl Rupp         PetscInt       cone[2];
492552f7358SJed Brown 
493da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
494552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
495552f7358SJed Brown         if (vy == numYVertices-1) {
496552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
497552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
498552f7358SJed Brown           if (ex == numXEdges-1) {
499552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
500552f7358SJed Brown           }
501552f7358SJed Brown         } else if (vy == 0) {
502552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
503552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
504552f7358SJed Brown           if (ex == numXEdges-1) {
505552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
506552f7358SJed Brown           }
507552f7358SJed Brown         }
508552f7358SJed Brown       }
509552f7358SJed Brown     }
510552f7358SJed Brown     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
511552f7358SJed Brown     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
512552f7358SJed Brown     /* Build coordinates */
513552f7358SJed Brown     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
514552f7358SJed Brown     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
515552f7358SJed Brown     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
516552f7358SJed Brown       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
517552f7358SJed Brown     }
518552f7358SJed Brown     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
519552f7358SJed Brown     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
52082f516ccSBarry Smith     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
521552f7358SJed Brown     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
522552f7358SJed Brown     ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
523552f7358SJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
524552f7358SJed Brown     for (vy = 0; vy < numYVertices; ++vy) {
525552f7358SJed Brown       for (vx = 0; vx < numXVertices; ++vx) {
526552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
527552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
528552f7358SJed Brown       }
529552f7358SJed Brown     }
530552f7358SJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
531552f7358SJed Brown     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
532552f7358SJed Brown     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
533552f7358SJed Brown   }
534552f7358SJed Brown   PetscFunctionReturn(0);
535552f7358SJed Brown }
536552f7358SJed Brown 
537552f7358SJed Brown #undef __FUNCT__
538552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
5396b75479cSMatthew G Knepley /*@
5406b75479cSMatthew G Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box).
5416b75479cSMatthew G Knepley 
5426b75479cSMatthew G Knepley   Collective on MPI_Comm
5436b75479cSMatthew G Knepley 
5446b75479cSMatthew G Knepley   Input Parameters:
5456b75479cSMatthew G Knepley + comm - The communicator for the DM object
5466b75479cSMatthew G Knepley . dim - The spatial dimension
5476b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
5486b75479cSMatthew G Knepley 
5496b75479cSMatthew G Knepley   Output Parameter:
5506b75479cSMatthew G Knepley . dm  - The DM object
5516b75479cSMatthew G Knepley 
5526b75479cSMatthew G Knepley   Level: beginner
5536b75479cSMatthew G Knepley 
5546b75479cSMatthew G Knepley .keywords: DM, create
5556b75479cSMatthew G Knepley .seealso: DMSetType(), DMCreate()
5566b75479cSMatthew G Knepley @*/
5570adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
5580adebc6cSBarry Smith {
559552f7358SJed Brown   DM             boundary;
560552f7358SJed Brown   PetscErrorCode ierr;
561552f7358SJed Brown 
562552f7358SJed Brown   PetscFunctionBegin;
563552f7358SJed Brown   PetscValidPointer(dm, 4);
564552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
565552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
566552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
567552f7358SJed Brown   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
568552f7358SJed Brown   switch (dim) {
569552f7358SJed Brown   case 2:
570552f7358SJed Brown   {
571552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
572552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
573552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
574552f7358SJed Brown 
575552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
576552f7358SJed Brown     break;
577552f7358SJed Brown   }
578552f7358SJed Brown   case 3:
579552f7358SJed Brown   {
580552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
581552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
582552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
583552f7358SJed Brown 
584552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
585552f7358SJed Brown     break;
586552f7358SJed Brown   }
587552f7358SJed Brown   default:
588552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
589552f7358SJed Brown   }
5900298fd71SBarry Smith   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
591552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
592552f7358SJed Brown   PetscFunctionReturn(0);
593552f7358SJed Brown }
594552f7358SJed Brown 
595552f7358SJed Brown #undef __FUNCT__
596552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
597a6dfd86eSKarl Rupp PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DM *dm)
598a6dfd86eSKarl Rupp {
599552f7358SJed Brown   PetscErrorCode ierr;
600552f7358SJed Brown 
601552f7358SJed Brown   PetscFunctionBegin;
602552f7358SJed Brown   PetscValidPointer(dm, 4);
603552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
604552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
605552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
606552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
607552f7358SJed Brown   switch (dim) {
608552f7358SJed Brown   case 2:
609552f7358SJed Brown   {
610552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
611552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
612552f7358SJed Brown 
613552f7358SJed Brown     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells);CHKERRQ(ierr);
614552f7358SJed Brown     break;
615552f7358SJed Brown   }
616552f7358SJed Brown #if 0
617552f7358SJed Brown   case 3:
618552f7358SJed Brown   {
619552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
620552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
621552f7358SJed Brown 
622552f7358SJed Brown     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells);CHKERRQ(ierr);
623552f7358SJed Brown     break;
624552f7358SJed Brown   }
625552f7358SJed Brown #endif
626552f7358SJed Brown   default:
627552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
628552f7358SJed Brown   }
629552f7358SJed Brown   PetscFunctionReturn(0);
630552f7358SJed Brown }
631552f7358SJed Brown 
632552f7358SJed Brown /* External function declarations here */
633552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
634552f7358SJed Brown extern PetscErrorCode DMCreateMatrix_Plex(DM dm, MatType mtype, Mat *J);
635552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
636552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
6379a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
638552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
639552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
640552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
641552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
642552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
643552f7358SJed Brown 
644552f7358SJed Brown #undef __FUNCT__
645552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
646552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
647552f7358SJed Brown {
648552f7358SJed Brown   PetscErrorCode ierr;
649552f7358SJed Brown 
650552f7358SJed Brown   PetscFunctionBegin;
651552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
652552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
653552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr);
654552f7358SJed Brown   PetscFunctionReturn(0);
655552f7358SJed Brown }
656552f7358SJed Brown 
657552f7358SJed Brown #undef __FUNCT__
658552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
659552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
660552f7358SJed Brown {
661552f7358SJed Brown   PetscErrorCode ierr;
662552f7358SJed Brown 
663552f7358SJed Brown   PetscFunctionBegin;
664552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
665552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
666552f7358SJed Brown   PetscFunctionReturn(0);
667552f7358SJed Brown }
668552f7358SJed Brown 
669552f7358SJed Brown #undef __FUNCT__
670552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
671552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
672552f7358SJed Brown {
673552f7358SJed Brown   PetscErrorCode ierr;
674552f7358SJed Brown 
675552f7358SJed Brown   PetscFunctionBegin;
676552f7358SJed Brown   ierr = PetscStrallocpy(VECSTANDARD, (char**)&dm->vectype);CHKERRQ(ierr);
6778865f1eaSKarl Rupp 
678552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
679552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
68038221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
681552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
682552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
683552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
684184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
685184d77edSJed Brown   dm->ops->getlocaltoglobalmappingblock    = NULL;
6860298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
687552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
688552f7358SJed Brown   dm->ops->getcoloring                     = 0;
689552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
690552f7358SJed Brown   dm->ops->createinterpolation             = 0;
691552f7358SJed Brown   dm->ops->getaggregates                   = 0;
692552f7358SJed Brown   dm->ops->getinjection                    = 0;
693552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
694552f7358SJed Brown   dm->ops->coarsen                         = 0;
695552f7358SJed Brown   dm->ops->refinehierarchy                 = 0;
696552f7358SJed Brown   dm->ops->coarsenhierarchy                = 0;
6970298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
6980298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
6990298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
7000298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
701552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
702552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
703552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
704552f7358SJed Brown   PetscFunctionReturn(0);
705552f7358SJed Brown }
706552f7358SJed Brown 
70763a16f15SMatthew G. Knepley #undef __FUNCT__
70863a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex"
70963a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
71063a16f15SMatthew G. Knepley {
71163a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
71263a16f15SMatthew G. Knepley   PetscErrorCode ierr;
71363a16f15SMatthew G. Knepley 
71463a16f15SMatthew G. Knepley   PetscFunctionBegin;
71563a16f15SMatthew G. Knepley   mesh->refct++;
71663a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
71763a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
71863a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
71963a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
72063a16f15SMatthew G. Knepley }
72163a16f15SMatthew G. Knepley 
7228818961aSMatthew G Knepley /*MC
7238818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
7248818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
7258818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
7268818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
7278818961aSMatthew G Knepley 
7288818961aSMatthew G Knepley   Level: intermediate
7298818961aSMatthew G Knepley 
7308818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
7318818961aSMatthew G Knepley M*/
7328818961aSMatthew G Knepley 
733552f7358SJed Brown #undef __FUNCT__
734552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
7358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
736552f7358SJed Brown {
737552f7358SJed Brown   DM_Plex        *mesh;
738770b213bSMatthew G Knepley   PetscInt       unit, d;
739552f7358SJed Brown   PetscErrorCode ierr;
740552f7358SJed Brown 
741552f7358SJed Brown   PetscFunctionBegin;
742552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
743552f7358SJed Brown   ierr     = PetscNewLog(dm, DM_Plex, &mesh);CHKERRQ(ierr);
744552f7358SJed Brown   dm->data = mesh;
745552f7358SJed Brown 
746552f7358SJed Brown   mesh->refct             = 1;
747552f7358SJed Brown   mesh->dim               = 0;
74882f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
749552f7358SJed Brown   mesh->maxConeSize       = 0;
7500298fd71SBarry Smith   mesh->cones             = NULL;
7510298fd71SBarry Smith   mesh->coneOrientations  = NULL;
75282f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
753552f7358SJed Brown   mesh->maxSupportSize    = 0;
7540298fd71SBarry Smith   mesh->supports          = NULL;
755552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
756552f7358SJed Brown   mesh->refinementLimit   = -1.0;
757552f7358SJed Brown 
7580298fd71SBarry Smith   mesh->facesTmp = NULL;
759552f7358SJed Brown 
7600298fd71SBarry Smith   mesh->subpointMap = NULL;
761552f7358SJed Brown 
7628865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
763552f7358SJed Brown 
7640298fd71SBarry Smith   mesh->labels              = NULL;
7658af19771SMatthew G. Knepley   mesh->depthLabel          = NULL;
7660298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
7670298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
7688865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
769552f7358SJed Brown   mesh->vtkCellHeight     = 0;
770552f7358SJed Brown   mesh->preallocCenterDim = -1;
771552f7358SJed Brown 
772552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
773552f7358SJed Brown   mesh->printFEM       = 0;
7746113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
775552f7358SJed Brown 
776552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
777552f7358SJed Brown   PetscFunctionReturn(0);
778552f7358SJed Brown }
779552f7358SJed Brown 
780552f7358SJed Brown #undef __FUNCT__
781552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
782552f7358SJed Brown /*@
783552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
784552f7358SJed Brown 
785552f7358SJed Brown   Collective on MPI_Comm
786552f7358SJed Brown 
787552f7358SJed Brown   Input Parameter:
788552f7358SJed Brown . comm - The communicator for the DMPlex object
789552f7358SJed Brown 
790552f7358SJed Brown   Output Parameter:
791552f7358SJed Brown . mesh  - The DMPlex object
792552f7358SJed Brown 
793552f7358SJed Brown   Level: beginner
794552f7358SJed Brown 
795552f7358SJed Brown .keywords: DMPlex, create
796552f7358SJed Brown @*/
797552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
798552f7358SJed Brown {
799552f7358SJed Brown   PetscErrorCode ierr;
800552f7358SJed Brown 
801552f7358SJed Brown   PetscFunctionBegin;
802552f7358SJed Brown   PetscValidPointer(mesh,2);
803552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
804552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
805552f7358SJed Brown   PetscFunctionReturn(0);
806552f7358SJed Brown }
807552f7358SJed Brown 
808552f7358SJed Brown #undef __FUNCT__
8099298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private"
8109298eaa6SMatthew G Knepley /*
8119298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
8129298eaa6SMatthew G Knepley */
8139298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
8149298eaa6SMatthew G Knepley {
8159298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
8169298eaa6SMatthew G Knepley   PetscErrorCode ierr;
8179298eaa6SMatthew G Knepley 
8189298eaa6SMatthew G Knepley   PetscFunctionBegin;
8199298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
8209298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
8219298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
8229298eaa6SMatthew G Knepley   }
8239298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
8249298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
8259298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
8269298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
8279298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
8289298eaa6SMatthew G Knepley     }
8299298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
8309298eaa6SMatthew G Knepley   }
8319298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
8329298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
8339298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
8349298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
8359298eaa6SMatthew G Knepley }
8369298eaa6SMatthew G Knepley 
8379298eaa6SMatthew G Knepley #undef __FUNCT__
8389298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private"
8399298eaa6SMatthew G Knepley /*
8409298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
8419298eaa6SMatthew G Knepley */
8429298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
8439298eaa6SMatthew G Knepley {
8449298eaa6SMatthew G Knepley   PetscSection   coordSection;
8459298eaa6SMatthew G Knepley   Vec            coordinates;
8469298eaa6SMatthew G Knepley   PetscScalar   *coords;
8479298eaa6SMatthew G Knepley   PetscInt       coordSize, v, d;
8489298eaa6SMatthew G Knepley   PetscErrorCode ierr;
8499298eaa6SMatthew G Knepley 
8509298eaa6SMatthew G Knepley   PetscFunctionBegin;
8519298eaa6SMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
8529298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
8539298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
8549298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
8559298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
8569298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
8579298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
8589298eaa6SMatthew G Knepley   }
8599298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
8609298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
8619298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
8629298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
8639298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
8649298eaa6SMatthew G Knepley   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
8659298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
8669298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
8679298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
8689298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
8699298eaa6SMatthew G Knepley     }
8709298eaa6SMatthew G Knepley   }
8719298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
8729298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
8739298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
8749298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
8759298eaa6SMatthew G Knepley }
8769298eaa6SMatthew G Knepley 
8779298eaa6SMatthew G Knepley #undef __FUNCT__
8789298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList"
8799298eaa6SMatthew G Knepley /*@C
8809298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
8819298eaa6SMatthew G Knepley 
8829298eaa6SMatthew G Knepley   Input Parameters:
8839298eaa6SMatthew G Knepley + comm - The communicator
8849298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
8859298eaa6SMatthew G Knepley . numCells - The number of cells
8869298eaa6SMatthew G Knepley . numVertices - The number of vertices
8879298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
8889298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
8899298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
8909298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
8919298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
8929298eaa6SMatthew G Knepley 
8939298eaa6SMatthew G Knepley   Output Parameter:
8949298eaa6SMatthew G Knepley . dm - The DM
8959298eaa6SMatthew G Knepley 
8969298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
8979298eaa6SMatthew G Knepley $
8989298eaa6SMatthew G Knepley $        2
8999298eaa6SMatthew G Knepley $      / | \
9009298eaa6SMatthew G Knepley $     /  |  \
9019298eaa6SMatthew G Knepley $    /   |   \
9029298eaa6SMatthew G Knepley $   0  0 | 1  3
9039298eaa6SMatthew G Knepley $    \   |   /
9049298eaa6SMatthew G Knepley $     \  |  /
9059298eaa6SMatthew G Knepley $      \ | /
9069298eaa6SMatthew G Knepley $        1
9079298eaa6SMatthew G Knepley would have input
9089298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
9099298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
9109298eaa6SMatthew G Knepley $
9119298eaa6SMatthew G Knepley which would result in the DMPlex
9129298eaa6SMatthew G Knepley $
9139298eaa6SMatthew G Knepley $        4
9149298eaa6SMatthew G Knepley $      / | \
9159298eaa6SMatthew G Knepley $     /  |  \
9169298eaa6SMatthew G Knepley $    /   |   \
9179298eaa6SMatthew G Knepley $   2  0 | 1  5
9189298eaa6SMatthew G Knepley $    \   |   /
9199298eaa6SMatthew G Knepley $     \  |  /
9209298eaa6SMatthew G Knepley $      \ | /
9219298eaa6SMatthew G Knepley $        3
9229298eaa6SMatthew G Knepley 
9239298eaa6SMatthew G Knepley   Level: beginner
9249298eaa6SMatthew G Knepley 
9259298eaa6SMatthew G Knepley .seealso: DMPlexCreate()
9269298eaa6SMatthew G Knepley @*/
9279298eaa6SMatthew 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)
9289298eaa6SMatthew G Knepley {
9299298eaa6SMatthew G Knepley   PetscErrorCode ierr;
9309298eaa6SMatthew G Knepley 
9319298eaa6SMatthew G Knepley   PetscFunctionBegin;
9329298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
9339298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
9349298eaa6SMatthew G Knepley   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
9359298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
9369298eaa6SMatthew G Knepley   if (interpolate) {
9379298eaa6SMatthew G Knepley     DM idm;
9389298eaa6SMatthew G Knepley 
9399298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
9409298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
9419298eaa6SMatthew G Knepley     *dm  = idm;
9429298eaa6SMatthew G Knepley   }
9439298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
9449298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
9459298eaa6SMatthew G Knepley }
9469298eaa6SMatthew G Knepley 
9479298eaa6SMatthew G Knepley #undef __FUNCT__
9489298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG"
9499298eaa6SMatthew G Knepley /*
9509298eaa6SMatthew G Knepley   This takes as input the raw Hasse Diagram data
9519298eaa6SMatthew G Knepley */
9529298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
9539298eaa6SMatthew G Knepley {
9549298eaa6SMatthew G Knepley   Vec            coordinates;
9559298eaa6SMatthew G Knepley   PetscSection   coordSection;
9569298eaa6SMatthew G Knepley   PetscScalar    *coords;
9579298eaa6SMatthew G Knepley   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
9589298eaa6SMatthew G Knepley   PetscErrorCode ierr;
9599298eaa6SMatthew G Knepley 
9609298eaa6SMatthew G Knepley   PetscFunctionBegin;
9619298eaa6SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
9629298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
9639298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
9649298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
9659298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
9669298eaa6SMatthew G Knepley   }
9679298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
9689298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
9699298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
9709298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
9719298eaa6SMatthew G Knepley   }
9729298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
9739298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
9749298eaa6SMatthew G Knepley   /* Build coordinates */
9759298eaa6SMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
9769298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
9779298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
9789298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
9799298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
9809298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
9816f8cbbeeSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
9829298eaa6SMatthew G Knepley   }
9839298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
9849298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
9859298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
9866f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
9879298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
9889298eaa6SMatthew G Knepley   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
9899298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
9909298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
9919298eaa6SMatthew G Knepley     PetscInt off;
9929298eaa6SMatthew G Knepley 
9939298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
9949298eaa6SMatthew G Knepley     for (d = 0; d < dim; ++d) {
9959298eaa6SMatthew G Knepley       coords[off+d] = vertexCoords[v*dim+d];
9969298eaa6SMatthew G Knepley     }
9979298eaa6SMatthew G Knepley   }
9989298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
9999298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
10009298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
10019298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
10029298eaa6SMatthew G Knepley }
1003