xref: /petsc/src/dm/impls/plex/plexcreate.c (revision fbf5b4835cfa97ba15f94625ed9d322bf40d9a86)
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);
61ce78fa2fSMatthew G. Knepley   switch (dim) {
62ce78fa2fSMatthew G. Knepley   case 2:
63ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
64ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
65ce78fa2fSMatthew G. Knepley     break;
66ce78fa2fSMatthew G. Knepley   case 3:
67ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
68ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
69ce78fa2fSMatthew G. Knepley     break;
70ce78fa2fSMatthew G. Knepley   default:
71ce78fa2fSMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
72ce78fa2fSMatthew G. Knepley   }
731df5d5c5SMatthew G. Knepley   if (rank) {
741df5d5c5SMatthew G. Knepley     PetscInt numPoints[2] = {0, 0};
751df5d5c5SMatthew G. Knepley     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
761df5d5c5SMatthew G. Knepley   } else {
771df5d5c5SMatthew G. Knepley     switch (dim) {
781df5d5c5SMatthew G. Knepley     case 2:
791df5d5c5SMatthew G. Knepley       if (simplex) {
801df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {4, 2};
811df5d5c5SMatthew G. Knepley         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
821df5d5c5SMatthew G. Knepley         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
831df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
841df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
851df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
861df5d5c5SMatthew G. Knepley 
871df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
881df5d5c5SMatthew G. Knepley         for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
891df5d5c5SMatthew G. Knepley       } else {
901df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {6, 2};
911df5d5c5SMatthew G. Knepley         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
921df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
931df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
941df5d5c5SMatthew 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};
951df5d5c5SMatthew G. Knepley 
961df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
971df5d5c5SMatthew G. Knepley       }
981df5d5c5SMatthew G. Knepley       break;
991df5d5c5SMatthew G. Knepley     case 3:
1001df5d5c5SMatthew G. Knepley       if (simplex) {
1011df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {5, 2};
1021df5d5c5SMatthew G. Knepley         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
1031df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
1041df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
1051df5d5c5SMatthew 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};
1061df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
1071df5d5c5SMatthew G. Knepley 
1081df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
1091df5d5c5SMatthew G. Knepley         for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
1101df5d5c5SMatthew G. Knepley       } else {
1111df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]         = {12, 2};
1121df5d5c5SMatthew G. Knepley         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
1131df5d5c5SMatthew G. Knepley         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
1141df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
1151df5d5c5SMatthew 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,
1161df5d5c5SMatthew 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,
1171df5d5c5SMatthew 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};
1181df5d5c5SMatthew G. Knepley 
1191df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
1201df5d5c5SMatthew G. Knepley       }
1211df5d5c5SMatthew G. Knepley       break;
1221df5d5c5SMatthew G. Knepley     default:
1231df5d5c5SMatthew G. Knepley       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
1241df5d5c5SMatthew G. Knepley     }
1251df5d5c5SMatthew G. Knepley   }
1261df5d5c5SMatthew G. Knepley   *newdm = dm;
1271df5d5c5SMatthew G. Knepley   if (refinementLimit > 0.0) {
1281df5d5c5SMatthew G. Knepley     DM rdm;
1291df5d5c5SMatthew G. Knepley     const char *name;
1301df5d5c5SMatthew G. Knepley 
1311df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
1321df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
1331df5d5c5SMatthew G. Knepley     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
1341df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1351df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
1361df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1371df5d5c5SMatthew G. Knepley     *newdm = rdm;
1381df5d5c5SMatthew G. Knepley   }
1391df5d5c5SMatthew G. Knepley   if (interpolate) {
1401df5d5c5SMatthew G. Knepley     DM idm;
1411df5d5c5SMatthew G. Knepley     const char *name;
1421df5d5c5SMatthew G. Knepley 
1431df5d5c5SMatthew G. Knepley     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
1441df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1451df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
1461df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
1471df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr);
1481df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1491df5d5c5SMatthew G. Knepley     *newdm = idm;
1501df5d5c5SMatthew G. Knepley   }
1511df5d5c5SMatthew G. Knepley   {
1521df5d5c5SMatthew G. Knepley     DM refinedMesh     = NULL;
1531df5d5c5SMatthew G. Knepley     DM distributedMesh = NULL;
1541df5d5c5SMatthew G. Knepley 
1551df5d5c5SMatthew G. Knepley     /* Distribute mesh over processes */
15607104863SJed Brown     ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr);
1571df5d5c5SMatthew G. Knepley     if (distributedMesh) {
1581df5d5c5SMatthew G. Knepley       ierr = DMDestroy(newdm);CHKERRQ(ierr);
1591df5d5c5SMatthew G. Knepley       *newdm = distributedMesh;
1601df5d5c5SMatthew G. Knepley     }
1611df5d5c5SMatthew G. Knepley     if (refinementUniform) {
1621df5d5c5SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
1631df5d5c5SMatthew G. Knepley       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
1641df5d5c5SMatthew G. Knepley       if (refinedMesh) {
1651df5d5c5SMatthew G. Knepley         ierr = DMDestroy(newdm);CHKERRQ(ierr);
1661df5d5c5SMatthew G. Knepley         *newdm = refinedMesh;
1671df5d5c5SMatthew G. Knepley       }
1681df5d5c5SMatthew G. Knepley     }
1691df5d5c5SMatthew G. Knepley   }
1701df5d5c5SMatthew G. Knepley   PetscFunctionReturn(0);
1711df5d5c5SMatthew G. Knepley }
1721df5d5c5SMatthew G. Knepley 
1731df5d5c5SMatthew G. Knepley #undef __FUNCT__
174552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary"
17526492d91SMatthew G. Knepley /*@
17626492d91SMatthew G. Knepley   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
177552f7358SJed Brown 
17826492d91SMatthew G. Knepley   Collective on MPI_Comm
17926492d91SMatthew G. Knepley 
18026492d91SMatthew G. Knepley   Input Parameters:
18126492d91SMatthew G. Knepley + comm  - The communicator for the DM object
18226492d91SMatthew G. Knepley . lower - The lower left corner coordinates
18326492d91SMatthew G. Knepley . upper - The upper right corner coordinates
18426492d91SMatthew G. Knepley - edges - The number of cells in each direction
18526492d91SMatthew G. Knepley 
18626492d91SMatthew G. Knepley   Output Parameter:
18726492d91SMatthew G. Knepley . dm  - The DM object
18826492d91SMatthew G. Knepley 
18926492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
19026492d91SMatthew G. Knepley $ 18--5-17--4--16
19126492d91SMatthew G. Knepley $  |     |     |
19226492d91SMatthew G. Knepley $  6    10     3
19326492d91SMatthew G. Knepley $  |     |     |
19426492d91SMatthew G. Knepley $ 19-11-20--9--15
19526492d91SMatthew G. Knepley $  |     |     |
19626492d91SMatthew G. Knepley $  7     8     2
19726492d91SMatthew G. Knepley $  |     |     |
19826492d91SMatthew G. Knepley $ 12--0-13--1--14
19926492d91SMatthew G. Knepley 
20026492d91SMatthew G. Knepley   Level: beginner
20126492d91SMatthew G. Knepley 
20226492d91SMatthew G. Knepley .keywords: DM, create
20326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
20426492d91SMatthew G. Knepley @*/
205552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
206552f7358SJed Brown {
207552f7358SJed Brown   PetscInt       numVertices    = (edges[0]+1)*(edges[1]+1);
208552f7358SJed Brown   PetscInt       numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
209552f7358SJed Brown   PetscInt       markerTop      = 1;
210552f7358SJed Brown   PetscInt       markerBottom   = 1;
211552f7358SJed Brown   PetscInt       markerRight    = 1;
212552f7358SJed Brown   PetscInt       markerLeft     = 1;
213552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
214552f7358SJed Brown   Vec            coordinates;
215552f7358SJed Brown   PetscSection   coordSection;
216552f7358SJed Brown   PetscScalar    *coords;
217552f7358SJed Brown   PetscInt       coordSize;
218552f7358SJed Brown   PetscMPIInt    rank;
219552f7358SJed Brown   PetscInt       v, vx, vy;
220552f7358SJed Brown   PetscErrorCode ierr;
221552f7358SJed Brown 
222552f7358SJed Brown   PetscFunctionBegin;
2230298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
224552f7358SJed Brown   if (markerSeparate) {
225552f7358SJed Brown     markerTop    = 1;
226552f7358SJed Brown     markerBottom = 0;
227552f7358SJed Brown     markerRight  = 0;
228552f7358SJed Brown     markerLeft   = 0;
229552f7358SJed Brown   }
23082f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
231552f7358SJed Brown   if (!rank) {
232552f7358SJed Brown     PetscInt e, ex, ey;
233552f7358SJed Brown 
234552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
235552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
236552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
237552f7358SJed Brown     }
238552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
239552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
240552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
241552f7358SJed Brown         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
242552f7358SJed Brown         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
243da80777bSKarl Rupp         PetscInt cone[2];
244552f7358SJed Brown 
245da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
246552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
247552f7358SJed Brown         if (vx == edges[0]) {
248552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
249552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
250552f7358SJed Brown           if (ey == edges[1]-1) {
251552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
252552f7358SJed Brown           }
253552f7358SJed Brown         } else if (vx == 0) {
254552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
255552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
256552f7358SJed Brown           if (ey == edges[1]-1) {
257552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
258552f7358SJed Brown           }
259552f7358SJed Brown         }
260552f7358SJed Brown       }
261552f7358SJed Brown     }
262552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
263552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
264552f7358SJed Brown         PetscInt edge   = vy*edges[0]     + ex;
265552f7358SJed Brown         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
266da80777bSKarl Rupp         PetscInt cone[2];
267552f7358SJed Brown 
268da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
269552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
270552f7358SJed Brown         if (vy == edges[1]) {
271552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
272552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
273552f7358SJed Brown           if (ex == edges[0]-1) {
274552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
275552f7358SJed Brown           }
276552f7358SJed Brown         } else if (vy == 0) {
277552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
278552f7358SJed Brown           ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
279552f7358SJed Brown           if (ex == edges[0]-1) {
280552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
281552f7358SJed Brown           }
282552f7358SJed Brown         }
283552f7358SJed Brown       }
284552f7358SJed Brown     }
285552f7358SJed Brown   }
286552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
287552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
288552f7358SJed Brown   /* Build coordinates */
289c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
290552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
291552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
292552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
293552f7358SJed Brown   }
294552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
295552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
29682f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
297552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2982eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
299552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
300552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
301552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
302552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
303552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
304552f7358SJed Brown     }
305552f7358SJed Brown   }
306552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
307552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
308552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
309552f7358SJed Brown   PetscFunctionReturn(0);
310552f7358SJed Brown }
311552f7358SJed Brown 
312552f7358SJed Brown #undef __FUNCT__
313552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary"
31426492d91SMatthew G. Knepley /*@
31526492d91SMatthew G. Knepley   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
316552f7358SJed Brown 
31726492d91SMatthew G. Knepley   Collective on MPI_Comm
31826492d91SMatthew G. Knepley 
31926492d91SMatthew G. Knepley   Input Parameters:
32026492d91SMatthew G. Knepley + comm  - The communicator for the DM object
32126492d91SMatthew G. Knepley . lower - The lower left front corner coordinates
32226492d91SMatthew G. Knepley . upper - The upper right back corner coordinates
32326492d91SMatthew G. Knepley - edges - The number of cells in each direction
32426492d91SMatthew G. Knepley 
32526492d91SMatthew G. Knepley   Output Parameter:
32626492d91SMatthew G. Knepley . dm  - The DM object
32726492d91SMatthew G. Knepley 
32826492d91SMatthew G. Knepley   Level: beginner
32926492d91SMatthew G. Knepley 
33026492d91SMatthew G. Knepley .keywords: DM, create
33126492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
33226492d91SMatthew G. Knepley @*/
333552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
334552f7358SJed Brown {
335552f7358SJed Brown   PetscInt       numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
336552f7358SJed Brown   PetscInt       numFaces    = 6;
337552f7358SJed Brown   Vec            coordinates;
338552f7358SJed Brown   PetscSection   coordSection;
339552f7358SJed Brown   PetscScalar    *coords;
340552f7358SJed Brown   PetscInt       coordSize;
341552f7358SJed Brown   PetscMPIInt    rank;
342552f7358SJed Brown   PetscInt       v, vx, vy, vz;
343552f7358SJed Brown   PetscErrorCode ierr;
344552f7358SJed Brown 
345552f7358SJed Brown   PetscFunctionBegin;
34682f516ccSBarry 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");
34782f516ccSBarry 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");
34882f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
349552f7358SJed Brown   if (!rank) {
350552f7358SJed Brown     PetscInt f;
351552f7358SJed Brown 
352552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
353552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
354552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
355552f7358SJed Brown     }
356552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
357552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
358552f7358SJed Brown       ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
359552f7358SJed Brown     }
360552f7358SJed Brown     { /* Side 0 (Front) */
361da80777bSKarl Rupp       PetscInt cone[4];
362da80777bSKarl Rupp       cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6;
363552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr);
364552f7358SJed Brown     }
365552f7358SJed Brown     { /* Side 1 (Back) */
366da80777bSKarl Rupp       PetscInt cone[4];
367da80777bSKarl Rupp       cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3;
368552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr);
369552f7358SJed Brown     }
370552f7358SJed Brown     { /* Side 2 (Bottom) */
371da80777bSKarl Rupp       PetscInt cone[4];
372da80777bSKarl Rupp       cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4;
373552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr);
374552f7358SJed Brown     }
375552f7358SJed Brown     { /* Side 3 (Top) */
376da80777bSKarl Rupp       PetscInt cone[4];
377da80777bSKarl Rupp       cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2;
378552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr);
379552f7358SJed Brown     }
380552f7358SJed Brown     { /* Side 4 (Left) */
381da80777bSKarl Rupp       PetscInt cone[4];
382da80777bSKarl Rupp       cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2;
383552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr);
384552f7358SJed Brown     }
385552f7358SJed Brown     { /* Side 5 (Right) */
386da80777bSKarl Rupp       PetscInt cone[4];
387da80777bSKarl Rupp       cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7;
388552f7358SJed Brown       ierr    = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr);
389552f7358SJed Brown     }
390552f7358SJed Brown   }
391552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
392552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
393552f7358SJed Brown   /* Build coordinates */
394c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
395552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
396552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
397552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
398552f7358SJed Brown   }
399552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
400552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
40182f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
402552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
403552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4042eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
405552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
406552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
407552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
408552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
409552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
410552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
411552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
412552f7358SJed Brown       }
413552f7358SJed Brown     }
414552f7358SJed Brown   }
415552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
416552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
417552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
418552f7358SJed Brown   PetscFunctionReturn(0);
419552f7358SJed Brown }
420552f7358SJed Brown 
421552f7358SJed Brown #undef __FUNCT__
422552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh"
42326492d91SMatthew G. Knepley /*@
42426492d91SMatthew G. Knepley   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.
425552f7358SJed Brown 
42626492d91SMatthew G. Knepley   Collective on MPI_Comm
42726492d91SMatthew G. Knepley 
42826492d91SMatthew G. Knepley   Input Parameters:
42926492d91SMatthew G. Knepley + comm  - The communicator for the DM object
43026492d91SMatthew G. Knepley . lower - The lower left corner coordinates
43126492d91SMatthew G. Knepley . upper - The upper right corner coordinates
432*fbf5b483SMatthew G. Knepley . edges - The number of cells in each direction
433*fbf5b483SMatthew G. Knepley . bdX   - The boundary type for the X direction
434*fbf5b483SMatthew G. Knepley - bdY   - The boundary type for the Y direction
43526492d91SMatthew G. Knepley 
43626492d91SMatthew G. Knepley   Output Parameter:
43726492d91SMatthew G. Knepley . dm  - The DM object
43826492d91SMatthew G. Knepley 
43926492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
44026492d91SMatthew G. Knepley $ 22--8-23--9--24
44126492d91SMatthew G. Knepley $  |     |     |
44226492d91SMatthew G. Knepley $ 13  2 14  3  15
44326492d91SMatthew G. Knepley $  |     |     |
44426492d91SMatthew G. Knepley $ 19--6-20--7--21
44526492d91SMatthew G. Knepley $  |     |     |
44626492d91SMatthew G. Knepley $ 10  0 11  1 12
44726492d91SMatthew G. Knepley $  |     |     |
44826492d91SMatthew G. Knepley $ 16--4-17--5--18
44926492d91SMatthew G. Knepley 
45026492d91SMatthew G. Knepley   Level: beginner
45126492d91SMatthew G. Knepley 
45226492d91SMatthew G. Knepley .keywords: DM, create
45326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
45426492d91SMatthew G. Knepley @*/
455bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
456552f7358SJed Brown {
457552f7358SJed Brown   PetscInt       markerTop      = 1;
458552f7358SJed Brown   PetscInt       markerBottom   = 1;
459552f7358SJed Brown   PetscInt       markerRight    = 1;
460552f7358SJed Brown   PetscInt       markerLeft     = 1;
461552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
462552f7358SJed Brown   PetscMPIInt    rank;
463552f7358SJed Brown   PetscErrorCode ierr;
464552f7358SJed Brown 
465552f7358SJed Brown   PetscFunctionBegin;
46682f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
4670298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
468552f7358SJed Brown   if (markerSeparate) {
469552f7358SJed Brown     markerTop    = 3;
470552f7358SJed Brown     markerBottom = 1;
471552f7358SJed Brown     markerRight  = 2;
472552f7358SJed Brown     markerLeft   = 4;
473552f7358SJed Brown   }
474552f7358SJed Brown   {
475552f7358SJed Brown     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
476552f7358SJed Brown     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
47746675b8eSMatthew G. Knepley     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
47846675b8eSMatthew G. Knepley     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
479552f7358SJed Brown     const PetscInt numTotXEdges = numXEdges*numYVertices;
480552f7358SJed Brown     const PetscInt numTotYEdges = numYEdges*numXVertices;
481552f7358SJed Brown     const PetscInt numVertices  = numXVertices*numYVertices;
482552f7358SJed Brown     const PetscInt numEdges     = numTotXEdges + numTotYEdges;
483552f7358SJed Brown     const PetscInt numFaces     = numXEdges*numYEdges;
484552f7358SJed Brown     const PetscInt firstVertex  = numFaces;
485552f7358SJed Brown     const PetscInt firstXEdge   = numFaces + numVertices;
486552f7358SJed Brown     const PetscInt firstYEdge   = numFaces + numVertices + numTotXEdges;
487552f7358SJed Brown     Vec            coordinates;
488552f7358SJed Brown     PetscSection   coordSection;
489552f7358SJed Brown     PetscScalar   *coords;
490552f7358SJed Brown     PetscInt       coordSize;
491552f7358SJed Brown     PetscInt       v, vx, vy;
492552f7358SJed Brown     PetscInt       f, fx, fy, e, ex, ey;
493552f7358SJed Brown 
494552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr);
495552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
496552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
497552f7358SJed Brown     }
498552f7358SJed Brown     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
499552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
500552f7358SJed Brown     }
501552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
502552f7358SJed Brown     /* Build faces */
503552f7358SJed Brown     for (fy = 0; fy < numYEdges; fy++) {
504552f7358SJed Brown       for (fx = 0; fx < numXEdges; fx++) {
505552f7358SJed Brown         const PetscInt face    = fy*numXEdges + fx;
506552f7358SJed Brown         const PetscInt edgeL   = firstYEdge +   fx                 *numYEdges + fy;
507becd5721SMatthew G. Knepley         const PetscInt edgeR   = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy;
508552f7358SJed Brown         const PetscInt edgeB   = firstXEdge +   fy                 *numXEdges + fx;
509becd5721SMatthew G. Knepley         const PetscInt edgeT   = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx;
510552f7358SJed Brown         const PetscInt ornt[4] = {0, 0, -2, -2};
511da80777bSKarl Rupp         PetscInt       cone[4];
512552f7358SJed Brown 
513becd5721SMatthew G. Knepley         cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
514552f7358SJed Brown         ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
515552f7358SJed Brown         ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
516552f7358SJed Brown       }
517552f7358SJed Brown     }
518552f7358SJed Brown     /* Build Y edges*/
519552f7358SJed Brown     for (vx = 0; vx < numXVertices; vx++) {
520552f7358SJed Brown       for (ey = 0; ey < numYEdges; ey++) {
52146675b8eSMatthew G. Knepley         const PetscInt nextv   = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx;
522552f7358SJed Brown         const PetscInt edge    = firstYEdge  + vx*numYEdges + ey;
523becd5721SMatthew G. Knepley         const PetscInt vertexB = firstVertex + ey*numXVertices + vx;
52446675b8eSMatthew G. Knepley         const PetscInt vertexT = firstVertex + nextv;
525da80777bSKarl Rupp         PetscInt       cone[2];
526552f7358SJed Brown 
527becd5721SMatthew G. Knepley         cone[0] = vertexB; cone[1] = vertexT;
528552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
529c4b39bc1SMatthew G. Knepley         if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
530552f7358SJed Brown           if (vx == numXVertices-1) {
531552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
532552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
533552f7358SJed Brown             if (ey == numYEdges-1) {
534552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
535552f7358SJed Brown             }
536552f7358SJed Brown           } else if (vx == 0) {
537552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
538552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
539552f7358SJed Brown             if (ey == numYEdges-1) {
540552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
541552f7358SJed Brown             }
542552f7358SJed Brown           }
543552f7358SJed Brown         }
544552f7358SJed Brown       }
545c4b39bc1SMatthew G. Knepley     }
546552f7358SJed Brown     /* Build X edges*/
547552f7358SJed Brown     for (vy = 0; vy < numYVertices; vy++) {
548552f7358SJed Brown       for (ex = 0; ex < numXEdges; ex++) {
54946675b8eSMatthew G. Knepley         const PetscInt nextv   = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices;
550552f7358SJed Brown         const PetscInt edge    = firstXEdge  + vy*numXEdges + ex;
551becd5721SMatthew G. Knepley         const PetscInt vertexL = firstVertex + vy*numXVertices + ex;
55246675b8eSMatthew G. Knepley         const PetscInt vertexR = firstVertex + nextv;
553da80777bSKarl Rupp         PetscInt       cone[2];
554552f7358SJed Brown 
555becd5721SMatthew G. Knepley         cone[0] = vertexL; cone[1] = vertexR;
556552f7358SJed Brown         ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
557c4b39bc1SMatthew G. Knepley         if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
558552f7358SJed Brown           if (vy == numYVertices-1) {
559552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
560552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
561552f7358SJed Brown             if (ex == numXEdges-1) {
562552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
563552f7358SJed Brown             }
564552f7358SJed Brown           } else if (vy == 0) {
565552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
566552f7358SJed Brown             ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
567552f7358SJed Brown             if (ex == numXEdges-1) {
568552f7358SJed Brown               ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
569552f7358SJed Brown             }
570552f7358SJed Brown           }
571552f7358SJed Brown         }
572552f7358SJed Brown       }
573c4b39bc1SMatthew G. Knepley     }
574552f7358SJed Brown     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
575552f7358SJed Brown     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
576552f7358SJed Brown     /* Build coordinates */
577c2166f76SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
578552f7358SJed Brown     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
579552f7358SJed Brown     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
580552f7358SJed Brown       ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
581552f7358SJed Brown     }
582552f7358SJed Brown     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
583552f7358SJed Brown     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
58482f516ccSBarry Smith     ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
585552f7358SJed Brown     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
5862eb5907fSJed Brown     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
587552f7358SJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
588552f7358SJed Brown     for (vy = 0; vy < numYVertices; ++vy) {
589552f7358SJed Brown       for (vx = 0; vx < numXVertices; ++vx) {
590552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
591552f7358SJed Brown         coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
592552f7358SJed Brown       }
593552f7358SJed Brown     }
594552f7358SJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
595552f7358SJed Brown     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
596552f7358SJed Brown     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
597552f7358SJed Brown   }
598552f7358SJed Brown   PetscFunctionReturn(0);
599552f7358SJed Brown }
600552f7358SJed Brown 
601552f7358SJed Brown #undef __FUNCT__
602552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh"
6036b75479cSMatthew G Knepley /*@
60426492d91SMatthew G. Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.
6056b75479cSMatthew G Knepley 
6066b75479cSMatthew G Knepley   Collective on MPI_Comm
6076b75479cSMatthew G Knepley 
6086b75479cSMatthew G Knepley   Input Parameters:
6096b75479cSMatthew G Knepley + comm - The communicator for the DM object
6106b75479cSMatthew G Knepley . dim - The spatial dimension
6116b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
6126b75479cSMatthew G Knepley 
6136b75479cSMatthew G Knepley   Output Parameter:
6146b75479cSMatthew G Knepley . dm  - The DM object
6156b75479cSMatthew G Knepley 
6166b75479cSMatthew G Knepley   Level: beginner
6176b75479cSMatthew G Knepley 
6186b75479cSMatthew G Knepley .keywords: DM, create
61926492d91SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
6206b75479cSMatthew G Knepley @*/
6210adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm)
6220adebc6cSBarry Smith {
623552f7358SJed Brown   DM             boundary;
624552f7358SJed Brown   PetscErrorCode ierr;
625552f7358SJed Brown 
626552f7358SJed Brown   PetscFunctionBegin;
627552f7358SJed Brown   PetscValidPointer(dm, 4);
628552f7358SJed Brown   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
629552f7358SJed Brown   PetscValidLogicalCollectiveInt(boundary,dim,2);
630552f7358SJed Brown   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
631552f7358SJed Brown   ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr);
632552f7358SJed Brown   switch (dim) {
633552f7358SJed Brown   case 2:
634552f7358SJed Brown   {
635552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
636552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
637552f7358SJed Brown     PetscInt  edges[2] = {2, 2};
638552f7358SJed Brown 
639552f7358SJed Brown     ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr);
640552f7358SJed Brown     break;
641552f7358SJed Brown   }
642552f7358SJed Brown   case 3:
643552f7358SJed Brown   {
644552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
645552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
646552f7358SJed Brown     PetscInt  faces[3] = {1, 1, 1};
647552f7358SJed Brown 
648552f7358SJed Brown     ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);
649552f7358SJed Brown     break;
650552f7358SJed Brown   }
651552f7358SJed Brown   default:
652552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
653552f7358SJed Brown   }
6540298fd71SBarry Smith   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
655552f7358SJed Brown   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
656552f7358SJed Brown   PetscFunctionReturn(0);
657552f7358SJed Brown }
658552f7358SJed Brown 
659552f7358SJed Brown #undef __FUNCT__
660552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh"
66126492d91SMatthew G. Knepley /*@
66226492d91SMatthew G. Knepley   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.
66326492d91SMatthew G. Knepley 
66426492d91SMatthew G. Knepley   Collective on MPI_Comm
66526492d91SMatthew G. Knepley 
66626492d91SMatthew G. Knepley   Input Parameters:
66726492d91SMatthew G. Knepley + comm  - The communicator for the DM object
66826492d91SMatthew G. Knepley . dim   - The spatial dimension
669*fbf5b483SMatthew G. Knepley . periodicX - The boundary type for the X direction
670*fbf5b483SMatthew G. Knepley . periodicY - The boundary type for the Y direction
671*fbf5b483SMatthew G. Knepley . periodicZ - The boundary type for the Z direction
67226492d91SMatthew G. Knepley - cells - The number of cells in each direction
67326492d91SMatthew G. Knepley 
67426492d91SMatthew G. Knepley   Output Parameter:
67526492d91SMatthew G. Knepley . dm  - The DM object
67626492d91SMatthew G. Knepley 
67726492d91SMatthew G. Knepley   Level: beginner
67826492d91SMatthew G. Knepley 
67926492d91SMatthew G. Knepley .keywords: DM, create
68026492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
68126492d91SMatthew G. Knepley @*/
682bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
683a6dfd86eSKarl Rupp {
684552f7358SJed Brown   PetscErrorCode ierr;
685552f7358SJed Brown 
686552f7358SJed Brown   PetscFunctionBegin;
687552f7358SJed Brown   PetscValidPointer(dm, 4);
688552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
689552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
690552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
691552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
692552f7358SJed Brown   switch (dim) {
693552f7358SJed Brown   case 2:
694552f7358SJed Brown   {
695552f7358SJed Brown     PetscReal lower[2] = {0.0, 0.0};
696552f7358SJed Brown     PetscReal upper[2] = {1.0, 1.0};
697552f7358SJed Brown 
698becd5721SMatthew G. Knepley     ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr);
699552f7358SJed Brown     break;
700552f7358SJed Brown   }
701552f7358SJed Brown #if 0
702552f7358SJed Brown   case 3:
703552f7358SJed Brown   {
704552f7358SJed Brown     PetscReal lower[3] = {0.0, 0.0, 0.0};
705552f7358SJed Brown     PetscReal upper[3] = {1.0, 1.0, 1.0};
706552f7358SJed Brown 
707becd5721SMatthew G. Knepley     ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr);
708552f7358SJed Brown     break;
709552f7358SJed Brown   }
710552f7358SJed Brown #endif
711552f7358SJed Brown   default:
712552f7358SJed Brown     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
713552f7358SJed Brown   }
714552f7358SJed Brown   PetscFunctionReturn(0);
715552f7358SJed Brown }
716552f7358SJed Brown 
717552f7358SJed Brown /* External function declarations here */
718552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
719b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
720552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
721552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
7229a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
723552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
724552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
725552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
726552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
727552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS);
728552f7358SJed Brown 
729552f7358SJed Brown #undef __FUNCT__
730552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex"
731552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
732552f7358SJed Brown {
733552f7358SJed Brown   PetscErrorCode ierr;
734552f7358SJed Brown 
735552f7358SJed Brown   PetscFunctionBegin;
736552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
737552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
738552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr);
739552f7358SJed Brown   PetscFunctionReturn(0);
740552f7358SJed Brown }
741552f7358SJed Brown 
742552f7358SJed Brown #undef __FUNCT__
743552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex"
744552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
745552f7358SJed Brown {
746552f7358SJed Brown   PetscErrorCode ierr;
747552f7358SJed Brown 
748552f7358SJed Brown   PetscFunctionBegin;
749552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
750552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
751552f7358SJed Brown   PetscFunctionReturn(0);
752552f7358SJed Brown }
753552f7358SJed Brown 
75438221697SMatthew G. Knepley #undef __FUNCT__
755552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex"
756552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm)
757552f7358SJed Brown {
758552f7358SJed Brown   PetscFunctionBegin;
759552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
760552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
76138221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
762552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
763552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
764552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
765184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
766184d77edSJed Brown   dm->ops->getlocaltoglobalmappingblock    = NULL;
7670298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
768552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
769552f7358SJed Brown   dm->ops->getcoloring                     = 0;
770552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
771552f7358SJed Brown   dm->ops->createinterpolation             = 0;
772552f7358SJed Brown   dm->ops->getaggregates                   = 0;
773552f7358SJed Brown   dm->ops->getinjection                    = 0;
774552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
775552f7358SJed Brown   dm->ops->coarsen                         = 0;
776552f7358SJed Brown   dm->ops->refinehierarchy                 = 0;
777552f7358SJed Brown   dm->ops->coarsenhierarchy                = 0;
7780298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
7790298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
7800298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
7810298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
782552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
783552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
784552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
785552f7358SJed Brown   PetscFunctionReturn(0);
786552f7358SJed Brown }
787552f7358SJed Brown 
78863a16f15SMatthew G. Knepley #undef __FUNCT__
78963a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex"
79063a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
79163a16f15SMatthew G. Knepley {
79263a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
79363a16f15SMatthew G. Knepley   PetscErrorCode ierr;
79463a16f15SMatthew G. Knepley 
79563a16f15SMatthew G. Knepley   PetscFunctionBegin;
79663a16f15SMatthew G. Knepley   mesh->refct++;
79763a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
79863a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
79963a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
80063a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
80163a16f15SMatthew G. Knepley }
80263a16f15SMatthew G. Knepley 
8038818961aSMatthew G Knepley /*MC
8048818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
8058818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
8068818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
8078818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
8088818961aSMatthew G Knepley 
8098818961aSMatthew G Knepley   Level: intermediate
8108818961aSMatthew G Knepley 
8118818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
8128818961aSMatthew G Knepley M*/
8138818961aSMatthew G Knepley 
814552f7358SJed Brown #undef __FUNCT__
815552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex"
8168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
817552f7358SJed Brown {
818552f7358SJed Brown   DM_Plex        *mesh;
819770b213bSMatthew G Knepley   PetscInt       unit, d;
820552f7358SJed Brown   PetscErrorCode ierr;
821552f7358SJed Brown 
822552f7358SJed Brown   PetscFunctionBegin;
823552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
824b00a9115SJed Brown   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
825552f7358SJed Brown   dm->data = mesh;
826552f7358SJed Brown 
827552f7358SJed Brown   mesh->refct             = 1;
828552f7358SJed Brown   mesh->dim               = 0;
82982f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
830552f7358SJed Brown   mesh->maxConeSize       = 0;
8310298fd71SBarry Smith   mesh->cones             = NULL;
8320298fd71SBarry Smith   mesh->coneOrientations  = NULL;
83382f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
834552f7358SJed Brown   mesh->maxSupportSize    = 0;
8350298fd71SBarry Smith   mesh->supports          = NULL;
836552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
837552f7358SJed Brown   mesh->refinementLimit   = -1.0;
838552f7358SJed Brown 
8390298fd71SBarry Smith   mesh->facesTmp = NULL;
840552f7358SJed Brown 
8410298fd71SBarry Smith   mesh->subpointMap = NULL;
842552f7358SJed Brown 
8438865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
844552f7358SJed Brown 
8450298fd71SBarry Smith   mesh->labels              = NULL;
8468af19771SMatthew G. Knepley   mesh->depthLabel          = NULL;
8470298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
8480298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
8498865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
850552f7358SJed Brown   mesh->vtkCellHeight     = 0;
851552f7358SJed Brown   mesh->preallocCenterDim = -1;
852552f7358SJed Brown 
853552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
854552f7358SJed Brown   mesh->printFEM       = 0;
8556113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
856552f7358SJed Brown 
857552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
858552f7358SJed Brown   PetscFunctionReturn(0);
859552f7358SJed Brown }
860552f7358SJed Brown 
861552f7358SJed Brown #undef __FUNCT__
862552f7358SJed Brown #define __FUNCT__ "DMPlexCreate"
863552f7358SJed Brown /*@
864552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
865552f7358SJed Brown 
866552f7358SJed Brown   Collective on MPI_Comm
867552f7358SJed Brown 
868552f7358SJed Brown   Input Parameter:
869552f7358SJed Brown . comm - The communicator for the DMPlex object
870552f7358SJed Brown 
871552f7358SJed Brown   Output Parameter:
872552f7358SJed Brown . mesh  - The DMPlex object
873552f7358SJed Brown 
874552f7358SJed Brown   Level: beginner
875552f7358SJed Brown 
876552f7358SJed Brown .keywords: DMPlex, create
877552f7358SJed Brown @*/
878552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
879552f7358SJed Brown {
880552f7358SJed Brown   PetscErrorCode ierr;
881552f7358SJed Brown 
882552f7358SJed Brown   PetscFunctionBegin;
883552f7358SJed Brown   PetscValidPointer(mesh,2);
884552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
885552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
886552f7358SJed Brown   PetscFunctionReturn(0);
887552f7358SJed Brown }
888552f7358SJed Brown 
889552f7358SJed Brown #undef __FUNCT__
8909298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private"
8919298eaa6SMatthew G Knepley /*
8929298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
8939298eaa6SMatthew G Knepley */
8949298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
8959298eaa6SMatthew G Knepley {
8969298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
8979298eaa6SMatthew G Knepley   PetscErrorCode ierr;
8989298eaa6SMatthew G Knepley 
8999298eaa6SMatthew G Knepley   PetscFunctionBegin;
9009298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
9019298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
9029298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
9039298eaa6SMatthew G Knepley   }
9049298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
9059298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
9069298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
9079298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
9089298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
9099298eaa6SMatthew G Knepley     }
9109298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
9119298eaa6SMatthew G Knepley   }
9129298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
9139298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
9149298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
9159298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
9169298eaa6SMatthew G Knepley }
9179298eaa6SMatthew G Knepley 
9189298eaa6SMatthew G Knepley #undef __FUNCT__
9199298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private"
9209298eaa6SMatthew G Knepley /*
9219298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
9229298eaa6SMatthew G Knepley */
9239298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
9249298eaa6SMatthew G Knepley {
9259298eaa6SMatthew G Knepley   PetscSection   coordSection;
9269298eaa6SMatthew G Knepley   Vec            coordinates;
9279298eaa6SMatthew G Knepley   PetscScalar   *coords;
9289298eaa6SMatthew G Knepley   PetscInt       coordSize, v, d;
9299298eaa6SMatthew G Knepley   PetscErrorCode ierr;
9309298eaa6SMatthew G Knepley 
9319298eaa6SMatthew G Knepley   PetscFunctionBegin;
932c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
9339298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
9349298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
9359298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
9369298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
9379298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
9389298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
9399298eaa6SMatthew G Knepley   }
9409298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
9419298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
9429298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
9439298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
9449298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
9452eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
9469298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
9479298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
9489298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
9499298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
9509298eaa6SMatthew G Knepley     }
9519298eaa6SMatthew G Knepley   }
9529298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
9539298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
9549298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
9559298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
9569298eaa6SMatthew G Knepley }
9579298eaa6SMatthew G Knepley 
9589298eaa6SMatthew G Knepley #undef __FUNCT__
9599298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList"
9609298eaa6SMatthew G Knepley /*@C
9619298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
9629298eaa6SMatthew G Knepley 
9639298eaa6SMatthew G Knepley   Input Parameters:
9649298eaa6SMatthew G Knepley + comm - The communicator
9659298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
9669298eaa6SMatthew G Knepley . numCells - The number of cells
9679298eaa6SMatthew G Knepley . numVertices - The number of vertices
9689298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
9699298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
9709298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
9719298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
9729298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
9739298eaa6SMatthew G Knepley 
9749298eaa6SMatthew G Knepley   Output Parameter:
9759298eaa6SMatthew G Knepley . dm - The DM
9769298eaa6SMatthew G Knepley 
9779298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
9789298eaa6SMatthew G Knepley $
9799298eaa6SMatthew G Knepley $        2
9809298eaa6SMatthew G Knepley $      / | \
9819298eaa6SMatthew G Knepley $     /  |  \
9829298eaa6SMatthew G Knepley $    /   |   \
9839298eaa6SMatthew G Knepley $   0  0 | 1  3
9849298eaa6SMatthew G Knepley $    \   |   /
9859298eaa6SMatthew G Knepley $     \  |  /
9869298eaa6SMatthew G Knepley $      \ | /
9879298eaa6SMatthew G Knepley $        1
9889298eaa6SMatthew G Knepley would have input
9899298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
9909298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
9919298eaa6SMatthew G Knepley $
9929298eaa6SMatthew G Knepley which would result in the DMPlex
9939298eaa6SMatthew G Knepley $
9949298eaa6SMatthew G Knepley $        4
9959298eaa6SMatthew G Knepley $      / | \
9969298eaa6SMatthew G Knepley $     /  |  \
9979298eaa6SMatthew G Knepley $    /   |   \
9989298eaa6SMatthew G Knepley $   2  0 | 1  5
9999298eaa6SMatthew G Knepley $    \   |   /
10009298eaa6SMatthew G Knepley $     \  |  /
10019298eaa6SMatthew G Knepley $      \ | /
10029298eaa6SMatthew G Knepley $        3
10039298eaa6SMatthew G Knepley 
10049298eaa6SMatthew G Knepley   Level: beginner
10059298eaa6SMatthew G Knepley 
1006939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
10079298eaa6SMatthew G Knepley @*/
10089298eaa6SMatthew 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)
10099298eaa6SMatthew G Knepley {
10109298eaa6SMatthew G Knepley   PetscErrorCode ierr;
10119298eaa6SMatthew G Knepley 
10129298eaa6SMatthew G Knepley   PetscFunctionBegin;
10139298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
10149298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
10159298eaa6SMatthew G Knepley   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
10169298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
10179298eaa6SMatthew G Knepley   if (interpolate) {
10189298eaa6SMatthew G Knepley     DM idm;
10199298eaa6SMatthew G Knepley 
10209298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
10219298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
10229298eaa6SMatthew G Knepley     *dm  = idm;
10239298eaa6SMatthew G Knepley   }
10249298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
10259298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
10269298eaa6SMatthew G Knepley }
10279298eaa6SMatthew G Knepley 
10289298eaa6SMatthew G Knepley #undef __FUNCT__
10299298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG"
1030939f6067SMatthew G. Knepley /*@
1031939f6067SMatthew 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
1032939f6067SMatthew G. Knepley 
1033939f6067SMatthew G. Knepley   Input Parameters:
1034939f6067SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension()
1035939f6067SMatthew G. Knepley . depth - The depth of the DAG
1036939f6067SMatthew G. Knepley . numPoints - The number of points at each depth
1037939f6067SMatthew G. Knepley . coneSize - The cone size of each point
1038939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1039939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point
1040939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
1041939f6067SMatthew G. Knepley 
1042939f6067SMatthew G. Knepley   Output Parameter:
1043939f6067SMatthew G. Knepley . dm - The DM
1044939f6067SMatthew G. Knepley 
1045939f6067SMatthew G. Knepley   Note: Two triangles sharing a face would have input
1046939f6067SMatthew G. Knepley $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1047939f6067SMatthew G. Knepley $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1048939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1049939f6067SMatthew G. Knepley $
1050939f6067SMatthew G. Knepley which would result in the DMPlex
1051939f6067SMatthew G. Knepley $
1052939f6067SMatthew G. Knepley $        4
1053939f6067SMatthew G. Knepley $      / | \
1054939f6067SMatthew G. Knepley $     /  |  \
1055939f6067SMatthew G. Knepley $    /   |   \
1056939f6067SMatthew G. Knepley $   2  0 | 1  5
1057939f6067SMatthew G. Knepley $    \   |   /
1058939f6067SMatthew G. Knepley $     \  |  /
1059939f6067SMatthew G. Knepley $      \ | /
1060939f6067SMatthew G. Knepley $        3
1061939f6067SMatthew G. Knepley $
1062939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
1063939f6067SMatthew G. Knepley 
1064939f6067SMatthew G. Knepley   Level: advanced
1065939f6067SMatthew G. Knepley 
1066939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1067939f6067SMatthew G. Knepley @*/
10689298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
10699298eaa6SMatthew G Knepley {
10709298eaa6SMatthew G Knepley   Vec            coordinates;
10719298eaa6SMatthew G Knepley   PetscSection   coordSection;
10729298eaa6SMatthew G Knepley   PetscScalar    *coords;
10739298eaa6SMatthew G Knepley   PetscInt       coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off;
10749298eaa6SMatthew G Knepley   PetscErrorCode ierr;
10759298eaa6SMatthew G Knepley 
10769298eaa6SMatthew G Knepley   PetscFunctionBegin;
10779298eaa6SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
10789298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
10799298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
10809298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
10819298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
10829298eaa6SMatthew G Knepley   }
10839298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
10849298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
10859298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
10869298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
10879298eaa6SMatthew G Knepley   }
10889298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
10899298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
10909298eaa6SMatthew G Knepley   /* Build coordinates */
1091c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
10929298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
10939298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
10949298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
10959298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
10969298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
10976f8cbbeeSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
10989298eaa6SMatthew G Knepley   }
10999298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
11009298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
11019298eaa6SMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
11026f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
11039298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
11042eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
11059298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
11069298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
11079298eaa6SMatthew G Knepley     PetscInt off;
11089298eaa6SMatthew G Knepley 
11099298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
11109298eaa6SMatthew G Knepley     for (d = 0; d < dim; ++d) {
11119298eaa6SMatthew G Knepley       coords[off+d] = vertexCoords[v*dim+d];
11129298eaa6SMatthew G Knepley     }
11139298eaa6SMatthew G Knepley   }
11149298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
11159298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
11169298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
11179298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
11189298eaa6SMatthew G Knepley }
1119