xref: /petsc/src/dm/impls/plex/plexcreate.c (revision 768d5fcef7e65bc5716c5e46069d6a09414963c7)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <petscdmda.h>
40c312b8eSJed Brown #include <petscsf.h>
5552f7358SJed Brown 
61df5d5c5SMatthew G. Knepley /*@
71df5d5c5SMatthew G. Knepley   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
81df5d5c5SMatthew G. Knepley 
91df5d5c5SMatthew G. Knepley   Collective on MPI_Comm
101df5d5c5SMatthew G. Knepley 
111df5d5c5SMatthew G. Knepley   Input Parameters:
121df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object
131df5d5c5SMatthew G. Knepley . dim - The spatial dimension
141df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells
151df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces)
161df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement
171df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
181df5d5c5SMatthew G. Knepley 
191df5d5c5SMatthew G. Knepley   Output Parameter:
201df5d5c5SMatthew G. Knepley . dm  - The DM object
211df5d5c5SMatthew G. Knepley 
221df5d5c5SMatthew G. Knepley   Level: beginner
231df5d5c5SMatthew G. Knepley 
241df5d5c5SMatthew G. Knepley .keywords: DM, create
251df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate()
261df5d5c5SMatthew G. Knepley @*/
271df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm)
281df5d5c5SMatthew G. Knepley {
291df5d5c5SMatthew G. Knepley   DM             dm;
301df5d5c5SMatthew G. Knepley   PetscInt       p;
311df5d5c5SMatthew G. Knepley   PetscMPIInt    rank;
321df5d5c5SMatthew G. Knepley   PetscErrorCode ierr;
331df5d5c5SMatthew G. Knepley 
341df5d5c5SMatthew G. Knepley   PetscFunctionBegin;
351df5d5c5SMatthew G. Knepley   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
361df5d5c5SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
37c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
381df5d5c5SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
39ce78fa2fSMatthew G. Knepley   switch (dim) {
40ce78fa2fSMatthew G. Knepley   case 2:
41ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);}
42ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);}
43ce78fa2fSMatthew G. Knepley     break;
44ce78fa2fSMatthew G. Knepley   case 3:
45ce78fa2fSMatthew G. Knepley     if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);}
46ce78fa2fSMatthew G. Knepley     else         {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);}
47ce78fa2fSMatthew G. Knepley     break;
48ce78fa2fSMatthew G. Knepley   default:
49ce78fa2fSMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
50ce78fa2fSMatthew G. Knepley   }
511df5d5c5SMatthew G. Knepley   if (rank) {
521df5d5c5SMatthew G. Knepley     PetscInt numPoints[2] = {0, 0};
531df5d5c5SMatthew G. Knepley     ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
541df5d5c5SMatthew G. Knepley   } else {
551df5d5c5SMatthew G. Knepley     switch (dim) {
561df5d5c5SMatthew G. Knepley     case 2:
571df5d5c5SMatthew G. Knepley       if (simplex) {
581df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {4, 2};
591df5d5c5SMatthew G. Knepley         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
601df5d5c5SMatthew G. Knepley         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
611df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
621df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
631df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
641df5d5c5SMatthew G. Knepley 
651df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
66c58f1c22SToby Isaac         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
671df5d5c5SMatthew G. Knepley       } else {
681df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {6, 2};
691df5d5c5SMatthew G. Knepley         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
701df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
711df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
721df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[12]    = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  -1.0, 0.5,  1.0, -0.5,  1.0, 0.5};
731df5d5c5SMatthew G. Knepley 
741df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
751df5d5c5SMatthew G. Knepley       }
761df5d5c5SMatthew G. Knepley       break;
771df5d5c5SMatthew G. Knepley     case 3:
781df5d5c5SMatthew G. Knepley       if (simplex) {
791df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]        = {5, 2};
801df5d5c5SMatthew G. Knepley         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
811df5d5c5SMatthew G. Knepley         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
821df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
831df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  1.0, 0.0, 0.0};
841df5d5c5SMatthew G. Knepley         PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
851df5d5c5SMatthew G. Knepley 
861df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
87c58f1c22SToby Isaac         for (p = 0; p < 5; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
881df5d5c5SMatthew G. Knepley       } else {
891df5d5c5SMatthew G. Knepley         PetscInt    numPoints[2]         = {12, 2};
901df5d5c5SMatthew G. Knepley         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
911df5d5c5SMatthew G. Knepley         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
921df5d5c5SMatthew G. Knepley         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
931df5d5c5SMatthew G. Knepley         PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5,  -1.0,  0.5, -0.5,  0.0,  0.5, -0.5,   0.0, -0.5, -0.5,
941df5d5c5SMatthew G. Knepley                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
951df5d5c5SMatthew G. Knepley                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
961df5d5c5SMatthew G. Knepley 
971df5d5c5SMatthew G. Knepley         ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
981df5d5c5SMatthew G. Knepley       }
991df5d5c5SMatthew G. Knepley       break;
1001df5d5c5SMatthew G. Knepley     default:
1011df5d5c5SMatthew G. Knepley       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
1021df5d5c5SMatthew G. Knepley     }
1031df5d5c5SMatthew G. Knepley   }
1041df5d5c5SMatthew G. Knepley   *newdm = dm;
1051df5d5c5SMatthew G. Knepley   if (refinementLimit > 0.0) {
1061df5d5c5SMatthew G. Knepley     DM rdm;
1071df5d5c5SMatthew G. Knepley     const char *name;
1081df5d5c5SMatthew G. Knepley 
1091df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
1101df5d5c5SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr);
1111df5d5c5SMatthew G. Knepley     ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr);
1121df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1131df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
1141df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1151df5d5c5SMatthew G. Knepley     *newdm = rdm;
1161df5d5c5SMatthew G. Knepley   }
1171df5d5c5SMatthew G. Knepley   if (interpolate) {
118e56d480eSMatthew G. Knepley     DM idm = NULL;
1191df5d5c5SMatthew G. Knepley     const char *name;
1201df5d5c5SMatthew G. Knepley 
1211df5d5c5SMatthew G. Knepley     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
1221df5d5c5SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
1231df5d5c5SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject)    idm,  name);CHKERRQ(ierr);
1241df5d5c5SMatthew G. Knepley     ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr);
125c58f1c22SToby Isaac     ierr = DMCopyLabels(*newdm, idm);CHKERRQ(ierr);
1261df5d5c5SMatthew G. Knepley     ierr = DMDestroy(newdm);CHKERRQ(ierr);
1271df5d5c5SMatthew G. Knepley     *newdm = idm;
1281df5d5c5SMatthew G. Knepley   }
1291df5d5c5SMatthew G. Knepley   {
1301df5d5c5SMatthew G. Knepley     DM refinedMesh     = NULL;
1311df5d5c5SMatthew G. Knepley     DM distributedMesh = NULL;
1321df5d5c5SMatthew G. Knepley 
1331df5d5c5SMatthew G. Knepley     /* Distribute mesh over processes */
13480cf41d5SMatthew G. Knepley     ierr = DMPlexDistribute(*newdm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
1351df5d5c5SMatthew G. Knepley     if (distributedMesh) {
1361df5d5c5SMatthew G. Knepley       ierr = DMDestroy(newdm);CHKERRQ(ierr);
1371df5d5c5SMatthew G. Knepley       *newdm = distributedMesh;
1381df5d5c5SMatthew G. Knepley     }
1391df5d5c5SMatthew G. Knepley     if (refinementUniform) {
1401df5d5c5SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr);
1411df5d5c5SMatthew G. Knepley       ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr);
1421df5d5c5SMatthew G. Knepley       if (refinedMesh) {
1431df5d5c5SMatthew G. Knepley         ierr = DMDestroy(newdm);CHKERRQ(ierr);
1441df5d5c5SMatthew G. Knepley         *newdm = refinedMesh;
1451df5d5c5SMatthew G. Knepley       }
1461df5d5c5SMatthew G. Knepley     }
1471df5d5c5SMatthew G. Knepley   }
1481df5d5c5SMatthew G. Knepley   PetscFunctionReturn(0);
1491df5d5c5SMatthew G. Knepley }
1501df5d5c5SMatthew G. Knepley 
15126492d91SMatthew G. Knepley /*@
15226492d91SMatthew G. Knepley   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.
153552f7358SJed Brown 
15426492d91SMatthew G. Knepley   Collective on MPI_Comm
15526492d91SMatthew G. Knepley 
15626492d91SMatthew G. Knepley   Input Parameters:
15726492d91SMatthew G. Knepley + comm  - The communicator for the DM object
15826492d91SMatthew G. Knepley . lower - The lower left corner coordinates
15926492d91SMatthew G. Knepley . upper - The upper right corner coordinates
16026492d91SMatthew G. Knepley - edges - The number of cells in each direction
16126492d91SMatthew G. Knepley 
16226492d91SMatthew G. Knepley   Output Parameter:
16326492d91SMatthew G. Knepley . dm  - The DM object
16426492d91SMatthew G. Knepley 
16526492d91SMatthew G. Knepley   Note: Here is the numbering returned for 2 cells in each direction:
16626492d91SMatthew G. Knepley $ 18--5-17--4--16
16726492d91SMatthew G. Knepley $  |     |     |
16826492d91SMatthew G. Knepley $  6    10     3
16926492d91SMatthew G. Knepley $  |     |     |
17026492d91SMatthew G. Knepley $ 19-11-20--9--15
17126492d91SMatthew G. Knepley $  |     |     |
17226492d91SMatthew G. Knepley $  7     8     2
17326492d91SMatthew G. Knepley $  |     |     |
17426492d91SMatthew G. Knepley $ 12--0-13--1--14
17526492d91SMatthew G. Knepley 
17626492d91SMatthew G. Knepley   Level: beginner
17726492d91SMatthew G. Knepley 
17826492d91SMatthew G. Knepley .keywords: DM, create
17926492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
18026492d91SMatthew G. Knepley @*/
181552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
182552f7358SJed Brown {
1831df21d24SMatthew G. Knepley   const PetscInt numVertices    = (edges[0]+1)*(edges[1]+1);
1841df21d24SMatthew G. Knepley   const PetscInt numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
185552f7358SJed Brown   PetscInt       markerTop      = 1;
186552f7358SJed Brown   PetscInt       markerBottom   = 1;
187552f7358SJed Brown   PetscInt       markerRight    = 1;
188552f7358SJed Brown   PetscInt       markerLeft     = 1;
189552f7358SJed Brown   PetscBool      markerSeparate = PETSC_FALSE;
190552f7358SJed Brown   Vec            coordinates;
191552f7358SJed Brown   PetscSection   coordSection;
192552f7358SJed Brown   PetscScalar    *coords;
193552f7358SJed Brown   PetscInt       coordSize;
194552f7358SJed Brown   PetscMPIInt    rank;
195552f7358SJed Brown   PetscInt       v, vx, vy;
196552f7358SJed Brown   PetscErrorCode ierr;
197552f7358SJed Brown 
198552f7358SJed Brown   PetscFunctionBegin;
199c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
200552f7358SJed Brown   if (markerSeparate) {
2011df21d24SMatthew G. Knepley     markerTop    = 3;
2021df21d24SMatthew G. Knepley     markerBottom = 1;
2031df21d24SMatthew G. Knepley     markerRight  = 2;
2041df21d24SMatthew G. Knepley     markerLeft   = 4;
205552f7358SJed Brown   }
20682f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
207552f7358SJed Brown   if (!rank) {
208552f7358SJed Brown     PetscInt e, ex, ey;
209552f7358SJed Brown 
210552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr);
211552f7358SJed Brown     for (e = 0; e < numEdges; ++e) {
212552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
213552f7358SJed Brown     }
214552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
215552f7358SJed Brown     for (vx = 0; vx <= edges[0]; vx++) {
216552f7358SJed Brown       for (ey = 0; ey < edges[1]; ey++) {
217552f7358SJed Brown         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
218552f7358SJed Brown         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
219da80777bSKarl Rupp         PetscInt cone[2];
220552f7358SJed Brown 
221da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
222552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
223552f7358SJed Brown         if (vx == edges[0]) {
224c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
225c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
226552f7358SJed Brown           if (ey == edges[1]-1) {
227c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
228c668006fSMatthew G. Knepley             ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerRight);CHKERRQ(ierr);
229552f7358SJed Brown           }
230552f7358SJed Brown         } else if (vx == 0) {
231c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
232c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
233552f7358SJed Brown           if (ey == edges[1]-1) {
234c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
235c668006fSMatthew G. Knepley             ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft);CHKERRQ(ierr);
236552f7358SJed Brown           }
237552f7358SJed Brown         }
238552f7358SJed Brown       }
239552f7358SJed Brown     }
240552f7358SJed Brown     for (vy = 0; vy <= edges[1]; vy++) {
241552f7358SJed Brown       for (ex = 0; ex < edges[0]; ex++) {
242552f7358SJed Brown         PetscInt edge   = vy*edges[0]     + ex;
243552f7358SJed Brown         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
244da80777bSKarl Rupp         PetscInt cone[2];
245552f7358SJed Brown 
246da80777bSKarl Rupp         cone[0] = vertex; cone[1] = vertex+1;
247552f7358SJed Brown         ierr    = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
248552f7358SJed Brown         if (vy == edges[1]) {
249c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
250c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
251552f7358SJed Brown           if (ex == edges[0]-1) {
252c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
253c668006fSMatthew G. Knepley             ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerTop);CHKERRQ(ierr);
254552f7358SJed Brown           }
255552f7358SJed Brown         } else if (vy == 0) {
256c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
257c58f1c22SToby Isaac           ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
258552f7358SJed Brown           if (ex == edges[0]-1) {
259c58f1c22SToby Isaac             ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
260c668006fSMatthew G. Knepley             ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom);CHKERRQ(ierr);
261552f7358SJed Brown           }
262552f7358SJed Brown         }
263552f7358SJed Brown       }
264552f7358SJed Brown     }
265552f7358SJed Brown   }
266552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
267552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
268552f7358SJed Brown   /* Build coordinates */
2699596c6baSMatthew G. Knepley   ierr = DMSetCoordinateDim(dm, 2);CHKERRQ(ierr);
270c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
271972bc18aSToby Isaac   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
272552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr);
273972bc18aSToby Isaac   ierr = PetscSectionSetFieldComponents(coordSection, 0, 2);CHKERRQ(ierr);
274552f7358SJed Brown   for (v = numEdges; v < numEdges+numVertices; ++v) {
275552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
276972bc18aSToby Isaac     ierr = PetscSectionSetFieldDof(coordSection, v, 0, 2);CHKERRQ(ierr);
277552f7358SJed Brown   }
278552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
279552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
2808b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
281da16285aSMichael Lange   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
282552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2838b9ced59SLisandro Dalcin   ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr);
2842eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
285552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
286552f7358SJed Brown   for (vy = 0; vy <= edges[1]; ++vy) {
287552f7358SJed Brown     for (vx = 0; vx <= edges[0]; ++vx) {
288552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
289552f7358SJed Brown       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
290552f7358SJed Brown     }
291552f7358SJed Brown   }
292552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
293552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
294552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
295552f7358SJed Brown   PetscFunctionReturn(0);
296552f7358SJed Brown }
297552f7358SJed Brown 
29826492d91SMatthew G. Knepley /*@
29926492d91SMatthew G. Knepley   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.
300552f7358SJed Brown 
30126492d91SMatthew G. Knepley   Collective on MPI_Comm
30226492d91SMatthew G. Knepley 
30326492d91SMatthew G. Knepley   Input Parameters:
30426492d91SMatthew G. Knepley + comm  - The communicator for the DM object
30526492d91SMatthew G. Knepley . lower - The lower left front corner coordinates
30626492d91SMatthew G. Knepley . upper - The upper right back corner coordinates
30726492d91SMatthew G. Knepley - edges - The number of cells in each direction
30826492d91SMatthew G. Knepley 
30926492d91SMatthew G. Knepley   Output Parameter:
31026492d91SMatthew G. Knepley . dm  - The DM object
31126492d91SMatthew G. Knepley 
31226492d91SMatthew G. Knepley   Level: beginner
31326492d91SMatthew G. Knepley 
31426492d91SMatthew G. Knepley .keywords: DM, create
31526492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
31626492d91SMatthew G. Knepley @*/
317552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
318552f7358SJed Brown {
3199e8abbc3SMichael Lange   PetscInt       vertices[3], numVertices;
3207b59f5a9SMichael Lange   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
321552f7358SJed Brown   Vec            coordinates;
322552f7358SJed Brown   PetscSection   coordSection;
323552f7358SJed Brown   PetscScalar    *coords;
324552f7358SJed Brown   PetscInt       coordSize;
325552f7358SJed Brown   PetscMPIInt    rank;
326552f7358SJed Brown   PetscInt       v, vx, vy, vz;
3277b59f5a9SMichael Lange   PetscInt       voffset, iface=0, cone[4];
328552f7358SJed Brown   PetscErrorCode ierr;
329552f7358SJed Brown 
330552f7358SJed Brown   PetscFunctionBegin;
33182f516ccSBarry 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");
33282f516ccSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
3339e8abbc3SMichael Lange   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
3349e8abbc3SMichael Lange   numVertices = vertices[0]*vertices[1]*vertices[2];
335552f7358SJed Brown   if (!rank) {
336552f7358SJed Brown     PetscInt f;
337552f7358SJed Brown 
338552f7358SJed Brown     ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr);
339552f7358SJed Brown     for (f = 0; f < numFaces; ++f) {
340552f7358SJed Brown       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
341552f7358SJed Brown     }
342552f7358SJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
343552f7358SJed Brown     for (v = 0; v < numFaces+numVertices; ++v) {
344c58f1c22SToby Isaac       ierr = DMSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr);
345552f7358SJed Brown     }
3467b59f5a9SMichael Lange 
3477b59f5a9SMichael Lange     /* Side 0 (Top) */
3487b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3497b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3507b59f5a9SMichael Lange         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
3517b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
3527b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3537b59f5a9SMichael Lange         iface++;
354552f7358SJed Brown       }
355552f7358SJed Brown     }
3567b59f5a9SMichael Lange 
3577b59f5a9SMichael Lange     /* Side 1 (Bottom) */
3587b59f5a9SMichael Lange     for (vy = 0; vy < faces[1]; vy++) {
3597b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3607b59f5a9SMichael Lange         voffset = numFaces + vy*(faces[0]+1) + vx;
3617b59f5a9SMichael Lange         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
3627b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3637b59f5a9SMichael Lange         iface++;
364552f7358SJed Brown       }
365552f7358SJed Brown     }
3667b59f5a9SMichael Lange 
3677b59f5a9SMichael Lange     /* Side 2 (Front) */
3687b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3697b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3707b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
3717b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
3727b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3737b59f5a9SMichael Lange         iface++;
374552f7358SJed Brown       }
3757b59f5a9SMichael Lange     }
3767b59f5a9SMichael Lange 
3777b59f5a9SMichael Lange     /* Side 3 (Back) */
3787b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3797b59f5a9SMichael Lange       for (vx = 0; vx < faces[0]; vx++) {
3807b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
3817b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
3827b59f5a9SMichael Lange         cone[2] = voffset+1; cone[3] = voffset;
3837b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3847b59f5a9SMichael Lange         iface++;
3857b59f5a9SMichael Lange       }
3867b59f5a9SMichael Lange     }
3877b59f5a9SMichael Lange 
3887b59f5a9SMichael Lange     /* Side 4 (Left) */
3897b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
3907b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
3917b59f5a9SMichael Lange         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
3927b59f5a9SMichael Lange         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
3937b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
3947b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
3957b59f5a9SMichael Lange         iface++;
3967b59f5a9SMichael Lange       }
3977b59f5a9SMichael Lange     }
3987b59f5a9SMichael Lange 
3997b59f5a9SMichael Lange     /* Side 5 (Right) */
4007b59f5a9SMichael Lange     for (vz = 0; vz < faces[2]; vz++) {
4017b59f5a9SMichael Lange       for (vy = 0; vy < faces[1]; vy++) {
402aab5bcd8SJed Brown         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + faces[0];
4037b59f5a9SMichael Lange         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
4047b59f5a9SMichael Lange         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
4057b59f5a9SMichael Lange         ierr    = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr);
4067b59f5a9SMichael Lange         iface++;
4077b59f5a9SMichael Lange       }
408552f7358SJed Brown     }
409552f7358SJed Brown   }
410552f7358SJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
411552f7358SJed Brown   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
412552f7358SJed Brown   /* Build coordinates */
4139596c6baSMatthew G. Knepley   ierr = DMSetCoordinateDim(dm, 3);CHKERRQ(ierr);
414c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
415552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr);
416552f7358SJed Brown   for (v = numFaces; v < numFaces+numVertices; ++v) {
417552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr);
418552f7358SJed Brown   }
419552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
420552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
4218b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
422552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
423552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4248b9ced59SLisandro Dalcin   ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr);
4252eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
426552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
427552f7358SJed Brown   for (vz = 0; vz <= faces[2]; ++vz) {
428552f7358SJed Brown     for (vy = 0; vy <= faces[1]; ++vy) {
429552f7358SJed Brown       for (vx = 0; vx <= faces[0]; ++vx) {
430552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
431552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
432552f7358SJed Brown         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
433552f7358SJed Brown       }
434552f7358SJed Brown     }
435552f7358SJed Brown   }
436552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
437552f7358SJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
438552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
439552f7358SJed Brown   PetscFunctionReturn(0);
440552f7358SJed Brown }
441552f7358SJed Brown 
442*768d5fceSMatthew G. Knepley static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
443d6218766SMatthew G. Knepley {
444d6218766SMatthew G. Knepley   DM             boundary;
445*768d5fceSMatthew G. Knepley   PetscInt       i;
446d6218766SMatthew G. Knepley   PetscErrorCode ierr;
447d6218766SMatthew G. Knepley 
448d6218766SMatthew G. Knepley   PetscFunctionBegin;
449d6218766SMatthew G. Knepley   PetscValidPointer(dm, 4);
450*768d5fceSMatthew G. Knepley   for (i = 0; i < dim; ++i) if (periodicity[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity is not supported for simplex meshes");
451d6218766SMatthew G. Knepley   ierr = DMCreate(comm, &boundary);CHKERRQ(ierr);
452d6218766SMatthew G. Knepley   PetscValidLogicalCollectiveInt(boundary,dim,2);
453d6218766SMatthew G. Knepley   ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr);
454d6218766SMatthew G. Knepley   ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr);
455d6218766SMatthew G. Knepley   ierr = DMSetCoordinateDim(boundary, dim);CHKERRQ(ierr);
456d6218766SMatthew G. Knepley   switch (dim) {
457*768d5fceSMatthew G. Knepley   case 2: ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);break;
458*768d5fceSMatthew G. Knepley   case 3: ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);break;
45942c1443fSMatthew G. Knepley   default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
460d6218766SMatthew G. Knepley   }
461d6218766SMatthew G. Knepley   ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
462d6218766SMatthew G. Knepley   ierr = DMDestroy(&boundary);CHKERRQ(ierr);
463d6218766SMatthew G. Knepley   PetscFunctionReturn(0);
464d6218766SMatthew G. Knepley }
465d6218766SMatthew G. Knepley 
4663dfda0b1SToby Isaac static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
4673dfda0b1SToby Isaac {
468ed0e4b50SMatthew G. Knepley   DMLabel        cutLabel = NULL;
469f4eb4c5dSMatthew G. Knepley   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
470f4eb4c5dSMatthew G. Knepley   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
471f4eb4c5dSMatthew G. Knepley   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
472f4eb4c5dSMatthew G. Knepley   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
473f4eb4c5dSMatthew G. Knepley   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
474f4eb4c5dSMatthew G. Knepley   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
4753dfda0b1SToby Isaac   PetscInt       dim;
476d8211ee3SMatthew G. Knepley   PetscBool      markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
4773dfda0b1SToby Isaac   PetscMPIInt    rank;
4783dfda0b1SToby Isaac   PetscErrorCode ierr;
4793dfda0b1SToby Isaac 
4803dfda0b1SToby Isaac   PetscFunctionBegin;
481f0226e14SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
4823dfda0b1SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
48350ae33c3SToby Isaac   ierr = DMCreateLabel(dm,"marker");CHKERRQ(ierr);
48450ae33c3SToby Isaac   ierr = DMCreateLabel(dm,"Face Sets");CHKERRQ(ierr);
485d8211ee3SMatthew G. Knepley   if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ||
486d8211ee3SMatthew G. Knepley       bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ||
487d8211ee3SMatthew G. Knepley       bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {
48846919275SMatthew G. Knepley     ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);CHKERRQ(ierr);
489d1c88043SMatthew G. Knepley     if (cutMarker) {ierr = DMCreateLabel(dm, "periodic_cut");CHKERRQ(ierr); ierr = DMGetLabel(dm, "periodic_cut", &cutLabel);CHKERRQ(ierr);}
490d8211ee3SMatthew G. Knepley   }
4913dfda0b1SToby Isaac   switch (dim) {
4923dfda0b1SToby Isaac   case 2:
493f4eb4c5dSMatthew G. Knepley     faceMarkerTop    = 3;
494f4eb4c5dSMatthew G. Knepley     faceMarkerBottom = 1;
495f4eb4c5dSMatthew G. Knepley     faceMarkerRight  = 2;
496f4eb4c5dSMatthew G. Knepley     faceMarkerLeft   = 4;
4973dfda0b1SToby Isaac     break;
4983dfda0b1SToby Isaac   case 3:
499f4eb4c5dSMatthew G. Knepley     faceMarkerBottom = 1;
500f4eb4c5dSMatthew G. Knepley     faceMarkerTop    = 2;
501f4eb4c5dSMatthew G. Knepley     faceMarkerFront  = 3;
502f4eb4c5dSMatthew G. Knepley     faceMarkerBack   = 4;
503f4eb4c5dSMatthew G. Knepley     faceMarkerRight  = 5;
504f4eb4c5dSMatthew G. Knepley     faceMarkerLeft   = 6;
5053dfda0b1SToby Isaac     break;
5063dfda0b1SToby Isaac   default:
5073dfda0b1SToby Isaac     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim);
5083dfda0b1SToby Isaac     break;
5093dfda0b1SToby Isaac   }
510c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr);
511f4eb4c5dSMatthew G. Knepley   if (markerSeparate) {
512f4eb4c5dSMatthew G. Knepley     markerBottom = faceMarkerBottom;
513f4eb4c5dSMatthew G. Knepley     markerTop    = faceMarkerTop;
514f4eb4c5dSMatthew G. Knepley     markerFront  = faceMarkerFront;
515f4eb4c5dSMatthew G. Knepley     markerBack   = faceMarkerBack;
516f4eb4c5dSMatthew G. Knepley     markerRight  = faceMarkerRight;
517f4eb4c5dSMatthew G. Knepley     markerLeft   = faceMarkerLeft;
5183dfda0b1SToby Isaac   }
5193dfda0b1SToby Isaac   {
5203dfda0b1SToby Isaac     const PetscInt numXEdges    = !rank ? edges[0] : 0;
5213dfda0b1SToby Isaac     const PetscInt numYEdges    = !rank ? edges[1] : 0;
5223dfda0b1SToby Isaac     const PetscInt numZEdges    = !rank ? edges[2] : 0;
5233dfda0b1SToby Isaac     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
5243dfda0b1SToby Isaac     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
52542206facSLisandro Dalcin     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
5263dfda0b1SToby Isaac     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
5273dfda0b1SToby Isaac     const PetscInt numXFaces    = numYEdges*numZEdges;
5283dfda0b1SToby Isaac     const PetscInt numYFaces    = numXEdges*numZEdges;
5293dfda0b1SToby Isaac     const PetscInt numZFaces    = numXEdges*numYEdges;
5303dfda0b1SToby Isaac     const PetscInt numTotXFaces = numXVertices*numXFaces;
5313dfda0b1SToby Isaac     const PetscInt numTotYFaces = numYVertices*numYFaces;
5323dfda0b1SToby Isaac     const PetscInt numTotZFaces = numZVertices*numZFaces;
5333dfda0b1SToby Isaac     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
5343dfda0b1SToby Isaac     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
5353dfda0b1SToby Isaac     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
5363dfda0b1SToby Isaac     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
5373dfda0b1SToby Isaac     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
5383dfda0b1SToby Isaac     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
5393dfda0b1SToby Isaac     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
5403dfda0b1SToby Isaac     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
5413dfda0b1SToby Isaac     const PetscInt firstYFace   = firstXFace + numTotXFaces;
5423dfda0b1SToby Isaac     const PetscInt firstZFace   = firstYFace + numTotYFaces;
5433dfda0b1SToby Isaac     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
5443dfda0b1SToby Isaac     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
5453dfda0b1SToby Isaac     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
5463dfda0b1SToby Isaac     Vec            coordinates;
5473dfda0b1SToby Isaac     PetscSection   coordSection;
5483dfda0b1SToby Isaac     PetscScalar   *coords;
5493dfda0b1SToby Isaac     PetscInt       coordSize;
5503dfda0b1SToby Isaac     PetscInt       v, vx, vy, vz;
5513dfda0b1SToby Isaac     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;
5523dfda0b1SToby Isaac 
5533dfda0b1SToby Isaac     ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr);
5543dfda0b1SToby Isaac     for (c = 0; c < numCells; c++) {
5553dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr);
5563dfda0b1SToby Isaac     }
5573dfda0b1SToby Isaac     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
5583dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr);
5593dfda0b1SToby Isaac     }
5603dfda0b1SToby Isaac     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
5613dfda0b1SToby Isaac       ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr);
5623dfda0b1SToby Isaac     }
5633dfda0b1SToby Isaac     ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
5643dfda0b1SToby Isaac     /* Build cells */
5653dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
5663dfda0b1SToby Isaac       for (fy = 0; fy < numYEdges; ++fy) {
5673dfda0b1SToby Isaac         for (fx = 0; fx < numXEdges; ++fx) {
5683dfda0b1SToby Isaac           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
5693dfda0b1SToby Isaac           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
5703dfda0b1SToby Isaac           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
5713dfda0b1SToby Isaac           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
5723dfda0b1SToby Isaac           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
5733dfda0b1SToby Isaac           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
5743dfda0b1SToby Isaac           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
5753dfda0b1SToby Isaac                             /* B,  T,  F,  K,  R,  L */
57642206facSLisandro Dalcin           PetscInt ornt[6] = {-4,  0,  0, -1,  0, -4}; /* ??? */
57742206facSLisandro Dalcin           PetscInt cone[6];
5783dfda0b1SToby Isaac 
5793dfda0b1SToby Isaac           /* no boundary twisting in 3D */
5803dfda0b1SToby Isaac           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
5813dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr);
5823dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr);
5838a5b437dSMatthew G. Knepley           if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);}
5848a5b437dSMatthew G. Knepley           if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);}
5858a5b437dSMatthew G. Knepley           if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);}
5863dfda0b1SToby Isaac         }
5873dfda0b1SToby Isaac       }
5883dfda0b1SToby Isaac     }
5893dfda0b1SToby Isaac     /* Build x faces */
5903dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
5913dfda0b1SToby Isaac       for (fy = 0; fy < numYEdges; ++fy) {
5923dfda0b1SToby Isaac         for (fx = 0; fx < numXVertices; ++fx) {
5933dfda0b1SToby Isaac           PetscInt face    = firstXFace + (fz*numYEdges+fy)*numXVertices + fx;
5943dfda0b1SToby Isaac           PetscInt edgeL   = firstZEdge + (  fy*                 numXVertices+fx)*numZEdges + fz;
5953dfda0b1SToby Isaac           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
5963dfda0b1SToby Isaac           PetscInt edgeB   = firstYEdge + (  fz*                 numXVertices+fx)*numYEdges + fy;
5973dfda0b1SToby Isaac           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
5983dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
5993dfda0b1SToby Isaac           PetscInt cone[4];
6003dfda0b1SToby Isaac 
6013dfda0b1SToby Isaac           if (dim == 3) {
6023dfda0b1SToby Isaac             /* markers */
6033dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
6043dfda0b1SToby Isaac               if (fx == numXVertices-1) {
605c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);CHKERRQ(ierr);
606c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr);
6073dfda0b1SToby Isaac               }
6083dfda0b1SToby Isaac               else if (fx == 0) {
609c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);CHKERRQ(ierr);
610c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr);
6113dfda0b1SToby Isaac               }
6123dfda0b1SToby Isaac             }
6133dfda0b1SToby Isaac           }
6143dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
6153dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
6163dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
6173dfda0b1SToby Isaac         }
6183dfda0b1SToby Isaac       }
6193dfda0b1SToby Isaac     }
6203dfda0b1SToby Isaac     /* Build y faces */
6213dfda0b1SToby Isaac     for (fz = 0; fz < numZEdges; ++fz) {
62242206facSLisandro Dalcin       for (fx = 0; fx < numXEdges; ++fx) {
6233dfda0b1SToby Isaac         for (fy = 0; fy < numYVertices; ++fy) {
6243dfda0b1SToby Isaac           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
6253dfda0b1SToby Isaac           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx                 )*numZEdges + fz;
6263dfda0b1SToby Isaac           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
6273dfda0b1SToby Isaac           PetscInt edgeB   = firstXEdge + (  fz                 *numYVertices+fy)*numXEdges + fx;
6283dfda0b1SToby Isaac           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
6293dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
6303dfda0b1SToby Isaac           PetscInt cone[4];
6313dfda0b1SToby Isaac 
6323dfda0b1SToby Isaac           if (dim == 3) {
6333dfda0b1SToby Isaac             /* markers */
6343dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
6353dfda0b1SToby Isaac               if (fy == numYVertices-1) {
636c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);CHKERRQ(ierr);
637c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr);
6383dfda0b1SToby Isaac               }
6393dfda0b1SToby Isaac               else if (fy == 0) {
640c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);CHKERRQ(ierr);
641c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr);
6423dfda0b1SToby Isaac               }
6433dfda0b1SToby Isaac             }
6443dfda0b1SToby Isaac           }
6453dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
6463dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
6473dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
6483dfda0b1SToby Isaac         }
6493dfda0b1SToby Isaac       }
6503dfda0b1SToby Isaac     }
6513dfda0b1SToby Isaac     /* Build z faces */
6523dfda0b1SToby Isaac     for (fy = 0; fy < numYEdges; ++fy) {
6533dfda0b1SToby Isaac       for (fx = 0; fx < numXEdges; ++fx) {
6543dfda0b1SToby Isaac         for (fz = 0; fz < numZVertices; fz++) {
6553dfda0b1SToby Isaac           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
6563dfda0b1SToby Isaac           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx                 )*numYEdges + fy;
6573dfda0b1SToby Isaac           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
6583dfda0b1SToby Isaac           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy                 )*numXEdges + fx;
6593dfda0b1SToby Isaac           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
6603dfda0b1SToby Isaac           PetscInt ornt[4] = {0, 0, -2, -2};
6613dfda0b1SToby Isaac           PetscInt cone[4];
6623dfda0b1SToby Isaac 
6633dfda0b1SToby Isaac           if (dim == 2) {
6643dfda0b1SToby Isaac             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
6653dfda0b1SToby Isaac             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
6668a5b437dSMatthew G. Knepley             if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, face, 2);CHKERRQ(ierr);}
6678a5b437dSMatthew G. Knepley             if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, face, 2);CHKERRQ(ierr);}
668d1c88043SMatthew G. Knepley           } else {
6693dfda0b1SToby Isaac             /* markers */
6703dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
6713dfda0b1SToby Isaac               if (fz == numZVertices-1) {
672c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);CHKERRQ(ierr);
673c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr);
6743dfda0b1SToby Isaac               }
6753dfda0b1SToby Isaac               else if (fz == 0) {
676c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);CHKERRQ(ierr);
677c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr);
6783dfda0b1SToby Isaac               }
6793dfda0b1SToby Isaac             }
6803dfda0b1SToby Isaac           }
6813dfda0b1SToby Isaac           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
6823dfda0b1SToby Isaac           ierr    = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr);
6833dfda0b1SToby Isaac           ierr    = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr);
6843dfda0b1SToby Isaac         }
6853dfda0b1SToby Isaac       }
6863dfda0b1SToby Isaac     }
6873dfda0b1SToby Isaac     /* Build Z edges*/
6883dfda0b1SToby Isaac     for (vy = 0; vy < numYVertices; vy++) {
6893dfda0b1SToby Isaac       for (vx = 0; vx < numXVertices; vx++) {
6903dfda0b1SToby Isaac         for (ez = 0; ez < numZEdges; ez++) {
6913dfda0b1SToby Isaac           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
6923dfda0b1SToby Isaac           const PetscInt vertexB = firstVertex + (  ez                 *numYVertices+vy)*numXVertices + vx;
6933dfda0b1SToby Isaac           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
6943dfda0b1SToby Isaac           PetscInt       cone[2];
6953dfda0b1SToby Isaac 
6963dfda0b1SToby Isaac           if (dim == 3) {
6973dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
6983dfda0b1SToby Isaac               if (vx == numXVertices-1) {
699c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
7003dfda0b1SToby Isaac               }
7013dfda0b1SToby Isaac               else if (vx == 0) {
702c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
7033dfda0b1SToby Isaac               }
7043dfda0b1SToby Isaac             }
7053dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
7063dfda0b1SToby Isaac               if (vy == numYVertices-1) {
707c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
7083dfda0b1SToby Isaac               }
7093dfda0b1SToby Isaac               else if (vy == 0) {
710c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
7113dfda0b1SToby Isaac               }
7123dfda0b1SToby Isaac             }
7133dfda0b1SToby Isaac           }
7143dfda0b1SToby Isaac           cone[0] = vertexB; cone[1] = vertexT;
7153dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
7163dfda0b1SToby Isaac         }
7173dfda0b1SToby Isaac       }
7183dfda0b1SToby Isaac     }
7193dfda0b1SToby Isaac     /* Build Y edges*/
7203dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; vz++) {
7213dfda0b1SToby Isaac       for (vx = 0; vx < numXVertices; vx++) {
7223dfda0b1SToby Isaac         for (ey = 0; ey < numYEdges; ey++) {
7233dfda0b1SToby Isaac           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
7243dfda0b1SToby Isaac           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
7253dfda0b1SToby Isaac           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
7263dfda0b1SToby Isaac           const PetscInt vertexK = firstVertex + nextv;
7273dfda0b1SToby Isaac           PetscInt       cone[2];
7283dfda0b1SToby Isaac 
7293dfda0b1SToby Isaac           cone[0] = vertexF; cone[1] = vertexK;
7303dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
7313dfda0b1SToby Isaac           if (dim == 2) {
7323dfda0b1SToby Isaac             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
7333dfda0b1SToby Isaac               if (vx == numXVertices-1) {
734c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);CHKERRQ(ierr);
735c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerRight);CHKERRQ(ierr);
736c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr);
7373dfda0b1SToby Isaac                 if (ey == numYEdges-1) {
738c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr);
7393dfda0b1SToby Isaac                 }
740d8211ee3SMatthew G. Knepley               } else if (vx == 0) {
741c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);CHKERRQ(ierr);
742c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerLeft);CHKERRQ(ierr);
743c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr);
7443dfda0b1SToby Isaac                 if (ey == numYEdges-1) {
745c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr);
7463dfda0b1SToby Isaac                 }
7473dfda0b1SToby Isaac               }
748d8211ee3SMatthew G. Knepley             } else {
749d8211ee3SMatthew G. Knepley               if (vx == 0 && cutMarker) {
750d1c88043SMatthew G. Knepley                 ierr = DMLabelSetValue(cutLabel, edge,    1);CHKERRQ(ierr);
751d1c88043SMatthew G. Knepley                 ierr = DMLabelSetValue(cutLabel, cone[0], 1);CHKERRQ(ierr);
752d8211ee3SMatthew G. Knepley                 if (ey == numYEdges-1) {
753d1c88043SMatthew G. Knepley                   ierr = DMLabelSetValue(cutLabel, cone[1], 1);CHKERRQ(ierr);
7543dfda0b1SToby Isaac                 }
7553dfda0b1SToby Isaac               }
756d8211ee3SMatthew G. Knepley             }
757d8211ee3SMatthew G. Knepley           } else {
7583dfda0b1SToby Isaac             if (bdX != DM_BOUNDARY_PERIODIC) {
7593dfda0b1SToby Isaac               if (vx == numXVertices-1) {
760c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr);
761d8211ee3SMatthew G. Knepley               } else if (vx == 0) {
762c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr);
7633dfda0b1SToby Isaac               }
7643dfda0b1SToby Isaac             }
7653dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
7663dfda0b1SToby Isaac               if (vz == numZVertices-1) {
767c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
768d8211ee3SMatthew G. Knepley               } else if (vz == 0) {
769c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
7703dfda0b1SToby Isaac               }
7713dfda0b1SToby Isaac             }
7723dfda0b1SToby Isaac           }
7733dfda0b1SToby Isaac         }
7743dfda0b1SToby Isaac       }
7753dfda0b1SToby Isaac     }
7763dfda0b1SToby Isaac     /* Build X edges*/
7773dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; vz++) {
7783dfda0b1SToby Isaac       for (vy = 0; vy < numYVertices; vy++) {
7793dfda0b1SToby Isaac         for (ex = 0; ex < numXEdges; ex++) {
7803dfda0b1SToby Isaac           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
7813dfda0b1SToby Isaac           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
7823dfda0b1SToby Isaac           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
7833dfda0b1SToby Isaac           const PetscInt vertexR = firstVertex + nextv;
7843dfda0b1SToby Isaac           PetscInt       cone[2];
7853dfda0b1SToby Isaac 
7863dfda0b1SToby Isaac           cone[0] = vertexL; cone[1] = vertexR;
7873dfda0b1SToby Isaac           ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr);
7883dfda0b1SToby Isaac           if (dim == 2) {
7893dfda0b1SToby Isaac             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
7903dfda0b1SToby Isaac               if (vy == numYVertices-1) {
791c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);CHKERRQ(ierr);
792c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerTop);CHKERRQ(ierr);
793c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr);
7943dfda0b1SToby Isaac                 if (ex == numXEdges-1) {
795c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr);
7963dfda0b1SToby Isaac                 }
797d8211ee3SMatthew G. Knepley               } else if (vy == 0) {
798c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);CHKERRQ(ierr);
799c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge,    markerBottom);CHKERRQ(ierr);
800c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr);
8013dfda0b1SToby Isaac                 if (ex == numXEdges-1) {
802c58f1c22SToby Isaac                   ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr);
8033dfda0b1SToby Isaac                 }
8043dfda0b1SToby Isaac               }
805d8211ee3SMatthew G. Knepley             } else {
806d8211ee3SMatthew G. Knepley               if (vy == 0 && cutMarker) {
807d1c88043SMatthew G. Knepley                 ierr = DMLabelSetValue(cutLabel, edge,    1);CHKERRQ(ierr);
808d1c88043SMatthew G. Knepley                 ierr = DMLabelSetValue(cutLabel, cone[0], 1);CHKERRQ(ierr);
809d8211ee3SMatthew G. Knepley                 if (ex == numXEdges-1) {
810d1c88043SMatthew G. Knepley                   ierr = DMLabelSetValue(cutLabel, cone[1], 1);CHKERRQ(ierr);
8113dfda0b1SToby Isaac                 }
8123dfda0b1SToby Isaac               }
813d8211ee3SMatthew G. Knepley             }
814d8211ee3SMatthew G. Knepley           } else {
8153dfda0b1SToby Isaac             if (bdY != DM_BOUNDARY_PERIODIC) {
8163dfda0b1SToby Isaac               if (vy == numYVertices-1) {
817c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr);
8183dfda0b1SToby Isaac               }
8193dfda0b1SToby Isaac               else if (vy == 0) {
820c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr);
8213dfda0b1SToby Isaac               }
8223dfda0b1SToby Isaac             }
8233dfda0b1SToby Isaac             if (bdZ != DM_BOUNDARY_PERIODIC) {
8243dfda0b1SToby Isaac               if (vz == numZVertices-1) {
825c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr);
8263dfda0b1SToby Isaac               }
8273dfda0b1SToby Isaac               else if (vz == 0) {
828c58f1c22SToby Isaac                 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr);
8293dfda0b1SToby Isaac               }
8303dfda0b1SToby Isaac             }
8313dfda0b1SToby Isaac           }
8323dfda0b1SToby Isaac         }
8333dfda0b1SToby Isaac       }
8343dfda0b1SToby Isaac     }
8353dfda0b1SToby Isaac     ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
8363dfda0b1SToby Isaac     ierr = DMPlexStratify(dm);CHKERRQ(ierr);
8373dfda0b1SToby Isaac     /* Build coordinates */
8383dfda0b1SToby Isaac     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
8393dfda0b1SToby Isaac     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
8403dfda0b1SToby Isaac     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
8413dfda0b1SToby Isaac     ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr);
8423dfda0b1SToby Isaac     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
8433dfda0b1SToby Isaac       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
8443dfda0b1SToby Isaac       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
8453dfda0b1SToby Isaac     }
8463dfda0b1SToby Isaac     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
8473dfda0b1SToby Isaac     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
8488b9ced59SLisandro Dalcin     ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
849da16285aSMichael Lange     ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
8503dfda0b1SToby Isaac     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
8518b9ced59SLisandro Dalcin     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
8523dfda0b1SToby Isaac     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
8533dfda0b1SToby Isaac     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
8543dfda0b1SToby Isaac     for (vz = 0; vz < numZVertices; ++vz) {
8553dfda0b1SToby Isaac       for (vy = 0; vy < numYVertices; ++vy) {
8563dfda0b1SToby Isaac         for (vx = 0; vx < numXVertices; ++vx) {
8573dfda0b1SToby Isaac           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
8583dfda0b1SToby Isaac           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
8593dfda0b1SToby Isaac           if (dim == 3) {
8603dfda0b1SToby Isaac             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
8613dfda0b1SToby Isaac           }
8623dfda0b1SToby Isaac         }
8633dfda0b1SToby Isaac       }
8643dfda0b1SToby Isaac     }
8653dfda0b1SToby Isaac     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
8663dfda0b1SToby Isaac     ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
8673dfda0b1SToby Isaac     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
8683dfda0b1SToby Isaac   }
8693dfda0b1SToby Isaac   PetscFunctionReturn(0);
8703dfda0b1SToby Isaac }
8713dfda0b1SToby Isaac 
872*768d5fceSMatthew G. Knepley static PetscErrorCode DMPlexCreateBoxMesh_Tensor_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
873a6dfd86eSKarl Rupp {
874c8c68bd8SToby Isaac   PetscInt       i;
875552f7358SJed Brown   PetscErrorCode ierr;
876552f7358SJed Brown 
877552f7358SJed Brown   PetscFunctionBegin;
878d67d2e29SLisandro Dalcin   PetscValidPointer(dm, 7);
879552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
880552f7358SJed Brown   PetscValidLogicalCollectiveInt(*dm,dim,2);
881552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
882c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
883552f7358SJed Brown   switch (dim) {
884*768d5fceSMatthew G. Knepley   case 2: {ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], DM_BOUNDARY_NONE);CHKERRQ(ierr);break;}
885*768d5fceSMatthew G. Knepley   case 3: {ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], periodicity[2]);CHKERRQ(ierr);break;}
886*768d5fceSMatthew G. Knepley   default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
887*768d5fceSMatthew G. Knepley   }
888*768d5fceSMatthew G. Knepley   if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST ||
889*768d5fceSMatthew G. Knepley       periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST ||
890*768d5fceSMatthew G. Knepley       (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) {
891*768d5fceSMatthew G. Knepley     PetscReal L[3];
892*768d5fceSMatthew G. Knepley     PetscReal maxCell[3];
893552f7358SJed Brown 
894c8c68bd8SToby Isaac     for (i = 0; i < dim; i++) {
895c8c68bd8SToby Isaac       L[i]       = upper[i] - lower[i];
896*768d5fceSMatthew G. Knepley       maxCell[i] = 1.1 * (L[i] / PetscMax(1,faces[i]));
897*768d5fceSMatthew G. Knepley     }
898*768d5fceSMatthew G. Knepley     ierr = DMSetPeriodicity(*dm,PETSC_TRUE,maxCell,L,periodicity);CHKERRQ(ierr);
899*768d5fceSMatthew G. Knepley   }
900*768d5fceSMatthew G. Knepley   if (!interpolate) {
901*768d5fceSMatthew G. Knepley     DM udm;
902*768d5fceSMatthew G. Knepley 
903*768d5fceSMatthew G. Knepley     ierr = DMPlexUninterpolate(*dm, &udm);CHKERRQ(ierr);
904*768d5fceSMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
905*768d5fceSMatthew G. Knepley     *dm  = udm;
906*768d5fceSMatthew G. Knepley   }
907*768d5fceSMatthew G. Knepley   PetscFunctionReturn(0);
908c8c68bd8SToby Isaac }
909c8c68bd8SToby Isaac 
910*768d5fceSMatthew G. Knepley /*@C
911*768d5fceSMatthew G. Knepley   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices or tensor cells (hexahedra).
912*768d5fceSMatthew G. Knepley 
913*768d5fceSMatthew G. Knepley   Collective on MPI_Comm
914*768d5fceSMatthew G. Knepley 
915*768d5fceSMatthew G. Knepley   Input Parameters:
916*768d5fceSMatthew G. Knepley + comm        - The communicator for the DM object
917*768d5fceSMatthew G. Knepley . dim         - The spatial dimension
918*768d5fceSMatthew G. Knepley . simplex     - PETSC_TRUE for simplices, PETSC_FALSE for tensor cells
919*768d5fceSMatthew G. Knepley . faces       - Number of faces per dimension, or NULL for (2, 2) in 2D and (1, 1, 1) in 3D
920*768d5fceSMatthew G. Knepley . lower       - The lower left corner, or NULL for (0, 0, 0)
921*768d5fceSMatthew G. Knepley . upper       - The upper right corner, or NULL for (1, 1, 1)
922*768d5fceSMatthew G. Knepley . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUDNARY_NONE
923*768d5fceSMatthew G. Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
924*768d5fceSMatthew G. Knepley 
925*768d5fceSMatthew G. Knepley   Output Parameter:
926*768d5fceSMatthew G. Knepley . dm  - The DM object
927*768d5fceSMatthew G. Knepley 
928*768d5fceSMatthew G. Knepley   Note: Here is the numbering returned for 2 faces in each direction for tensor cells:
929*768d5fceSMatthew G. Knepley $ 10---17---11---18----12
930*768d5fceSMatthew G. Knepley $  |         |         |
931*768d5fceSMatthew G. Knepley $  |         |         |
932*768d5fceSMatthew G. Knepley $ 20    2   22    3    24
933*768d5fceSMatthew G. Knepley $  |         |         |
934*768d5fceSMatthew G. Knepley $  |         |         |
935*768d5fceSMatthew G. Knepley $  7---15----8---16----9
936*768d5fceSMatthew G. Knepley $  |         |         |
937*768d5fceSMatthew G. Knepley $  |         |         |
938*768d5fceSMatthew G. Knepley $ 19    0   21    1   23
939*768d5fceSMatthew G. Knepley $  |         |         |
940*768d5fceSMatthew G. Knepley $  |         |         |
941*768d5fceSMatthew G. Knepley $  4---13----5---14----6
942*768d5fceSMatthew G. Knepley 
943*768d5fceSMatthew G. Knepley and for simplicial cells
944*768d5fceSMatthew G. Knepley 
945*768d5fceSMatthew G. Knepley $ 14----8---15----9----16
946*768d5fceSMatthew G. Knepley $  |\     5  |\      7 |
947*768d5fceSMatthew G. Knepley $  | \       | \       |
948*768d5fceSMatthew G. Knepley $ 13   2    14    3    15
949*768d5fceSMatthew G. Knepley $  | 4   \   | 6   \   |
950*768d5fceSMatthew G. Knepley $  |       \ |       \ |
951*768d5fceSMatthew G. Knepley $ 11----6---12----7----13
952*768d5fceSMatthew G. Knepley $  |\        |\        |
953*768d5fceSMatthew G. Knepley $  | \    1  | \     3 |
954*768d5fceSMatthew G. Knepley $ 10   0    11    1    12
955*768d5fceSMatthew G. Knepley $  | 0   \   | 2   \   |
956*768d5fceSMatthew G. Knepley $  |       \ |       \ |
957*768d5fceSMatthew G. Knepley $  8----4----9----5----10
958*768d5fceSMatthew G. Knepley 
959*768d5fceSMatthew G. Knepley   Level: beginner
960*768d5fceSMatthew G. Knepley 
961*768d5fceSMatthew G. Knepley .keywords: DM, create
962*768d5fceSMatthew G. Knepley .seealso: DMPlexCreateFromFile(), DMPlexCreateHexCylinderMesh(), DMSetType(), DMCreate()
963*768d5fceSMatthew G. Knepley @*/
964*768d5fceSMatthew G. Knepley PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
965552f7358SJed Brown {
966f56c5c00SMatthew G. Knepley   PetscReal      low[3] = {lower ? lower[0] : 0.0, lower ? lower[1] : 0.0, lower ? lower[2] : 0.0};
967f56c5c00SMatthew G. Knepley   PetscReal      upp[3] = {upper ? upper[0] : 1.0, upper ? upper[1] : 1.0, upper ? upper[2] : 1.0};
968*768d5fceSMatthew G. Knepley   PetscInt       fac[3] = {faces ? faces[0] : 4-dim, faces ? faces[1] : 4-dim, faces ? faces[2] : (dim > 2 ? 4-dim : 0)};
969*768d5fceSMatthew G. Knepley   DMBoundaryType bdt[3] = {periodicity ? periodicity[0] : DM_BOUNDARY_NONE, periodicity ? periodicity[1] : DM_BOUNDARY_NONE, periodicity ? periodicity[2] : DM_BOUNDARY_NONE};
970*768d5fceSMatthew G. Knepley   PetscErrorCode ierr;
971552f7358SJed Brown 
972*768d5fceSMatthew G. Knepley   PetscFunctionBegin;
973*768d5fceSMatthew G. Knepley   if (simplex) {ierr = DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);CHKERRQ(ierr);}
974*768d5fceSMatthew G. Knepley   else         {ierr = DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);CHKERRQ(ierr);}
975552f7358SJed Brown   PetscFunctionReturn(0);
976552f7358SJed Brown }
977552f7358SJed Brown 
978a9074c1eSMatthew G. Knepley /*@C
979a9074c1eSMatthew G. Knepley   DMPlexSetOptionsPrefix - Sets the prefix used for searching for all DM options in the database.
980a9074c1eSMatthew G. Knepley 
981a9074c1eSMatthew G. Knepley   Logically Collective on DM
982a9074c1eSMatthew G. Knepley 
983a9074c1eSMatthew G. Knepley   Input Parameters:
984a9074c1eSMatthew G. Knepley + dm - the DM context
985a9074c1eSMatthew G. Knepley - prefix - the prefix to prepend to all option names
986a9074c1eSMatthew G. Knepley 
987a9074c1eSMatthew G. Knepley   Notes:
988a9074c1eSMatthew G. Knepley   A hyphen (-) must NOT be given at the beginning of the prefix name.
989a9074c1eSMatthew G. Knepley   The first character of all runtime options is AUTOMATICALLY the hyphen.
990a9074c1eSMatthew G. Knepley 
991a9074c1eSMatthew G. Knepley   Level: advanced
992a9074c1eSMatthew G. Knepley 
993a9074c1eSMatthew G. Knepley .seealso: SNESSetFromOptions()
994a9074c1eSMatthew G. Knepley @*/
995a9074c1eSMatthew G. Knepley PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
996a9074c1eSMatthew G. Knepley {
997a9074c1eSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
998a9074c1eSMatthew G. Knepley   PetscErrorCode ierr;
999a9074c1eSMatthew G. Knepley 
1000a9074c1eSMatthew G. Knepley   PetscFunctionBegin;
1001a9074c1eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1002a9074c1eSMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);CHKERRQ(ierr);
1003a9074c1eSMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);CHKERRQ(ierr);
1004a9074c1eSMatthew G. Knepley   PetscFunctionReturn(0);
1005a9074c1eSMatthew G. Knepley }
1006a9074c1eSMatthew G. Knepley 
10070510c589SMatthew G. Knepley /*@
10080510c589SMatthew G. Knepley   DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra.
10090510c589SMatthew G. Knepley 
10100510c589SMatthew G. Knepley   Collective on MPI_Comm
10110510c589SMatthew G. Knepley 
10120510c589SMatthew G. Knepley   Input Parameters:
10130510c589SMatthew G. Knepley + comm      - The communicator for the DM object
1014dbc1dc17SMatthew G. Knepley . numRefine - The number of regular refinements to the basic 5 cell structure
10150510c589SMatthew G. Knepley - periodicZ - The boundary type for the Z direction
10160510c589SMatthew G. Knepley 
10170510c589SMatthew G. Knepley   Output Parameter:
10180510c589SMatthew G. Knepley . dm  - The DM object
10190510c589SMatthew G. Knepley 
10200510c589SMatthew G. Knepley   Note: Here is the output numbering looking from the bottom of the cylinder:
10210510c589SMatthew G. Knepley $       17-----14
10220510c589SMatthew G. Knepley $        |     |
10230510c589SMatthew G. Knepley $        |  2  |
10240510c589SMatthew G. Knepley $        |     |
10250510c589SMatthew G. Knepley $ 17-----8-----7-----14
10260510c589SMatthew G. Knepley $  |     |     |     |
10270510c589SMatthew G. Knepley $  |  3  |  0  |  1  |
10280510c589SMatthew G. Knepley $  |     |     |     |
10290510c589SMatthew G. Knepley $ 19-----5-----6-----13
10300510c589SMatthew G. Knepley $        |     |
10310510c589SMatthew G. Knepley $        |  4  |
10320510c589SMatthew G. Knepley $        |     |
10330510c589SMatthew G. Knepley $       19-----13
10340510c589SMatthew G. Knepley $
10350510c589SMatthew G. Knepley $ and up through the top
10360510c589SMatthew G. Knepley $
10370510c589SMatthew G. Knepley $       18-----16
10380510c589SMatthew G. Knepley $        |     |
10390510c589SMatthew G. Knepley $        |  2  |
10400510c589SMatthew G. Knepley $        |     |
10410510c589SMatthew G. Knepley $ 18----10----11-----16
10420510c589SMatthew G. Knepley $  |     |     |     |
10430510c589SMatthew G. Knepley $  |  3  |  0  |  1  |
10440510c589SMatthew G. Knepley $  |     |     |     |
10450510c589SMatthew G. Knepley $ 20-----9----12-----15
10460510c589SMatthew G. Knepley $        |     |
10470510c589SMatthew G. Knepley $        |  4  |
10480510c589SMatthew G. Knepley $        |     |
10490510c589SMatthew G. Knepley $       20-----15
10500510c589SMatthew G. Knepley 
10510510c589SMatthew G. Knepley   Level: beginner
10520510c589SMatthew G. Knepley 
10530510c589SMatthew G. Knepley .keywords: DM, create
1054*768d5fceSMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
10550510c589SMatthew G. Knepley @*/
1056dbc1dc17SMatthew G. Knepley PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
10570510c589SMatthew G. Knepley {
10580510c589SMatthew G. Knepley   const PetscInt dim = 3;
1059006a8963SMatthew G. Knepley   PetscInt       numCells, numVertices, r;
10600510c589SMatthew G. Knepley   PetscErrorCode ierr;
10610510c589SMatthew G. Knepley 
10620510c589SMatthew G. Knepley   PetscFunctionBegin;
10630510c589SMatthew G. Knepley   PetscValidPointer(dm, 4);
1064dbc1dc17SMatthew G. Knepley   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
10650510c589SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
10660510c589SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
10670510c589SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
10680510c589SMatthew G. Knepley   /* Create topology */
10690510c589SMatthew G. Knepley   {
10700510c589SMatthew G. Knepley     PetscInt cone[8], c;
10710510c589SMatthew G. Knepley 
1072006a8963SMatthew G. Knepley     numCells    = 5;
1073006a8963SMatthew G. Knepley     numVertices = 16;
1074006a8963SMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1075ae8bcbbbSMatthew G. Knepley       numCells   *= 3;
1076ae8bcbbbSMatthew G. Knepley       numVertices = 24;
1077006a8963SMatthew G. Knepley     }
10780510c589SMatthew G. Knepley     ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
10790510c589SMatthew G. Knepley     for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 8);CHKERRQ(ierr);}
10800510c589SMatthew G. Knepley     ierr = DMSetUp(*dm);CHKERRQ(ierr);
1081006a8963SMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1082ae8bcbbbSMatthew G. Knepley       cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1083ae8bcbbbSMatthew G. Knepley       cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1084006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr);
1085ae8bcbbbSMatthew G. Knepley       cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1086ae8bcbbbSMatthew G. Knepley       cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1087006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr);
1088ae8bcbbbSMatthew G. Knepley       cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1089ae8bcbbbSMatthew G. Knepley       cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1090006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr);
1091ae8bcbbbSMatthew G. Knepley       cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1092ae8bcbbbSMatthew G. Knepley       cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1093006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr);
1094ae8bcbbbSMatthew G. Knepley       cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1095ae8bcbbbSMatthew G. Knepley       cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1096006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr);
1097006a8963SMatthew G. Knepley 
1098ae8bcbbbSMatthew G. Knepley       cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1099ae8bcbbbSMatthew G. Knepley       cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1100006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr);
1101ae8bcbbbSMatthew G. Knepley       cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1102ae8bcbbbSMatthew G. Knepley       cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1103006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 6, cone);CHKERRQ(ierr);
1104ae8bcbbbSMatthew G. Knepley       cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1105ae8bcbbbSMatthew G. Knepley       cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1106006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 7, cone);CHKERRQ(ierr);
1107ae8bcbbbSMatthew G. Knepley       cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1108ae8bcbbbSMatthew G. Knepley       cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1109006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 8, cone);CHKERRQ(ierr);
1110ae8bcbbbSMatthew G. Knepley       cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1111ae8bcbbbSMatthew G. Knepley       cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1112006a8963SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 9, cone);CHKERRQ(ierr);
1113ae8bcbbbSMatthew G. Knepley 
1114ae8bcbbbSMatthew G. Knepley       cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1115ae8bcbbbSMatthew G. Knepley       cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1116ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 10, cone);CHKERRQ(ierr);
1117ae8bcbbbSMatthew G. Knepley       cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1118ae8bcbbbSMatthew G. Knepley       cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1119ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 11, cone);CHKERRQ(ierr);
1120ae8bcbbbSMatthew G. Knepley       cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1121ae8bcbbbSMatthew G. Knepley       cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1122ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 12, cone);CHKERRQ(ierr);
1123ae8bcbbbSMatthew G. Knepley       cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1124ae8bcbbbSMatthew G. Knepley       cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1125ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 13, cone);CHKERRQ(ierr);
1126ae8bcbbbSMatthew G. Knepley       cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1127ae8bcbbbSMatthew G. Knepley       cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1128ae8bcbbbSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr);
1129006a8963SMatthew G. Knepley     } else {
113010c6f908SMatthew G. Knepley       cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
113110c6f908SMatthew G. Knepley       cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
11320510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr);
113310c6f908SMatthew G. Knepley       cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
113410c6f908SMatthew G. Knepley       cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
11350510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr);
113610c6f908SMatthew G. Knepley       cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
113710c6f908SMatthew G. Knepley       cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
11380510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr);
113910c6f908SMatthew G. Knepley       cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
114010c6f908SMatthew G. Knepley       cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
11410510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr);
114210c6f908SMatthew G. Knepley       cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
114310c6f908SMatthew G. Knepley       cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
11440510c589SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr);
1145006a8963SMatthew G. Knepley     }
11460510c589SMatthew G. Knepley     ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
11470510c589SMatthew G. Knepley     ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
11480510c589SMatthew G. Knepley   }
11490510c589SMatthew G. Knepley   /* Interpolate */
11500510c589SMatthew G. Knepley   {
11510510c589SMatthew G. Knepley     DM idm = NULL;
11520510c589SMatthew G. Knepley 
11530510c589SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
11540510c589SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
11550510c589SMatthew G. Knepley     *dm  = idm;
11560510c589SMatthew G. Knepley   }
1157dbc1dc17SMatthew G. Knepley   /* Create cube geometry */
11580510c589SMatthew G. Knepley   {
11590510c589SMatthew G. Knepley     Vec             coordinates;
11600510c589SMatthew G. Knepley     PetscSection    coordSection;
11610510c589SMatthew G. Knepley     PetscScalar    *coords;
11620510c589SMatthew G. Knepley     PetscInt        coordSize, v;
11630510c589SMatthew G. Knepley     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
11640510c589SMatthew G. Knepley     const PetscReal ds2 = dis/2.0;
11650510c589SMatthew G. Knepley 
11660510c589SMatthew G. Knepley     /* Build coordinates */
11670510c589SMatthew G. Knepley     ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
11680510c589SMatthew G. Knepley     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
11690510c589SMatthew G. Knepley     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
11700510c589SMatthew G. Knepley     ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr);
11710510c589SMatthew G. Knepley     for (v = numCells; v < numCells+numVertices; ++v) {
11720510c589SMatthew G. Knepley       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
11730510c589SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
11740510c589SMatthew G. Knepley     }
11750510c589SMatthew G. Knepley     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
11760510c589SMatthew G. Knepley     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
11770510c589SMatthew G. Knepley     ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
11780510c589SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
11790510c589SMatthew G. Knepley     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
11800510c589SMatthew G. Knepley     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
11810510c589SMatthew G. Knepley     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
11820510c589SMatthew G. Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
11830510c589SMatthew G. Knepley     coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
11840510c589SMatthew G. Knepley     coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
11850510c589SMatthew G. Knepley     coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
11860510c589SMatthew G. Knepley     coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
11870510c589SMatthew G. Knepley     coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
11880510c589SMatthew G. Knepley     coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
11890510c589SMatthew G. Knepley     coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
11900510c589SMatthew G. Knepley     coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
11910510c589SMatthew G. Knepley     coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
11920510c589SMatthew G. Knepley     coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
11930510c589SMatthew G. Knepley     coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
11940510c589SMatthew G. Knepley     coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
11950510c589SMatthew G. Knepley     coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
11960510c589SMatthew G. Knepley     coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
11970510c589SMatthew G. Knepley     coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
11980510c589SMatthew G. Knepley     coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1199ae8bcbbbSMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1200ae8bcbbbSMatthew G. Knepley       /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1201ae8bcbbbSMatthew G. Knepley       /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1202ae8bcbbbSMatthew G. Knepley       /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1203ae8bcbbbSMatthew G. Knepley       /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1204ae8bcbbbSMatthew G. Knepley       /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1205ae8bcbbbSMatthew G. Knepley       /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1206ae8bcbbbSMatthew G. Knepley       /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1207ae8bcbbbSMatthew G. Knepley       /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1208ae8bcbbbSMatthew G. Knepley     }
12090510c589SMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
12100510c589SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
12110510c589SMatthew G. Knepley     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
12120510c589SMatthew G. Knepley   }
1213006a8963SMatthew G. Knepley   /* Create periodicity */
1214006a8963SMatthew G. Knepley   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1215006a8963SMatthew G. Knepley     PetscReal      L[3];
1216006a8963SMatthew G. Knepley     PetscReal      maxCell[3];
1217006a8963SMatthew G. Knepley     DMBoundaryType bdType[3];
1218006a8963SMatthew G. Knepley     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1219ae8bcbbbSMatthew G. Knepley     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1220ae8bcbbbSMatthew G. Knepley     PetscInt       i, numZCells = 3;
1221006a8963SMatthew G. Knepley 
1222006a8963SMatthew G. Knepley     bdType[0] = DM_BOUNDARY_NONE;
1223006a8963SMatthew G. Knepley     bdType[1] = DM_BOUNDARY_NONE;
1224006a8963SMatthew G. Knepley     bdType[2] = periodicZ;
1225006a8963SMatthew G. Knepley     for (i = 0; i < dim; i++) {
1226006a8963SMatthew G. Knepley       L[i]       = upper[i] - lower[i];
1227006a8963SMatthew G. Knepley       maxCell[i] = 1.1 * (L[i] / numZCells);
1228006a8963SMatthew G. Knepley     }
1229dd169d64SMatthew G. Knepley     ierr = DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);CHKERRQ(ierr);
1230006a8963SMatthew G. Knepley   }
1231dbc1dc17SMatthew G. Knepley   /* Refine topology */
1232dbc1dc17SMatthew G. Knepley   for (r = 0; r < numRefine; ++r) {
1233dbc1dc17SMatthew G. Knepley     DM rdm = NULL;
1234dbc1dc17SMatthew G. Knepley 
1235dbc1dc17SMatthew G. Knepley     ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr);
1236dbc1dc17SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
1237dbc1dc17SMatthew G. Knepley     *dm  = rdm;
1238dbc1dc17SMatthew G. Knepley   }
1239dbc1dc17SMatthew G. Knepley   /* Remap geometry to cylinder
1240dbc1dc17SMatthew G. Knepley        Interior square: Linear interpolation is correct
1241dbc1dc17SMatthew G. Knepley        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1242dbc1dc17SMatthew G. Knepley        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance
1243dbc1dc17SMatthew G. Knepley 
1244dbc1dc17SMatthew G. Knepley          phi     = arctan(y/x)
1245dbc1dc17SMatthew G. Knepley          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1246dbc1dc17SMatthew G. Knepley          d_far   = sqrt(1/2 + sin^2(phi))
1247dbc1dc17SMatthew G. Knepley 
1248dbc1dc17SMatthew G. Knepley        so we remap them using
1249dbc1dc17SMatthew G. Knepley 
1250dbc1dc17SMatthew G. Knepley          x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close)
1251dbc1dc17SMatthew G. Knepley          y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close)
1252dbc1dc17SMatthew G. Knepley 
1253dbc1dc17SMatthew G. Knepley        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1254dbc1dc17SMatthew G. Knepley   */
1255dbc1dc17SMatthew G. Knepley   {
1256dbc1dc17SMatthew G. Knepley     Vec           coordinates;
1257dbc1dc17SMatthew G. Knepley     PetscSection  coordSection;
1258dbc1dc17SMatthew G. Knepley     PetscScalar  *coords;
1259dbc1dc17SMatthew G. Knepley     PetscInt      vStart, vEnd, v;
1260dbc1dc17SMatthew G. Knepley     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1261dbc1dc17SMatthew G. Knepley     const PetscReal ds2 = 0.5*dis;
1262dbc1dc17SMatthew G. Knepley 
1263dbc1dc17SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1264dbc1dc17SMatthew G. Knepley     ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1265dbc1dc17SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr);
1266dbc1dc17SMatthew G. Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1267dbc1dc17SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
126871752167SMatthew G. Knepley       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1269dbc1dc17SMatthew G. Knepley       PetscInt  off;
1270dbc1dc17SMatthew G. Knepley 
1271dbc1dc17SMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
127271752167SMatthew G. Knepley       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
127371752167SMatthew G. Knepley       x    = PetscRealPart(coords[off]);
127471752167SMatthew G. Knepley       y    = PetscRealPart(coords[off+1]);
127571752167SMatthew G. Knepley       phi  = PetscAtan2Real(y, x);
1276dbc1dc17SMatthew G. Knepley       sinp = PetscSinReal(phi);
1277dbc1dc17SMatthew G. Knepley       cosp = PetscCosReal(phi);
1278dbc1dc17SMatthew G. Knepley       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
127971752167SMatthew G. Knepley         dc = PetscAbsReal(ds2/sinp);
128071752167SMatthew G. Knepley         df = PetscAbsReal(dis/sinp);
12819420fdb6SJose E. Roman         xc = ds2*x/PetscAbsReal(y);
128271752167SMatthew G. Knepley         yc = ds2*PetscSignReal(y);
1283dbc1dc17SMatthew G. Knepley       } else {
128471752167SMatthew G. Knepley         dc = PetscAbsReal(ds2/cosp);
128571752167SMatthew G. Knepley         df = PetscAbsReal(dis/cosp);
128671752167SMatthew G. Knepley         xc = ds2*PetscSignReal(x);
12879420fdb6SJose E. Roman         yc = ds2*y/PetscAbsReal(x);
1288dbc1dc17SMatthew G. Knepley       }
1289dbc1dc17SMatthew G. Knepley       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1290dbc1dc17SMatthew G. Knepley       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1291dbc1dc17SMatthew G. Knepley     }
1292dbc1dc17SMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
12930510c589SMatthew G. Knepley     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1294006a8963SMatthew G. Knepley       ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
12950510c589SMatthew G. Knepley     }
12960510c589SMatthew G. Knepley   }
12970510c589SMatthew G. Knepley   PetscFunctionReturn(0);
12980510c589SMatthew G. Knepley }
12990510c589SMatthew G. Knepley 
130024119c2aSMatthew G. Knepley /*@
130124119c2aSMatthew G. Knepley   DMPlexCreateWedgeCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using wedges.
130224119c2aSMatthew G. Knepley 
130324119c2aSMatthew G. Knepley   Collective on MPI_Comm
130424119c2aSMatthew G. Knepley 
130524119c2aSMatthew G. Knepley   Input Parameters:
130624119c2aSMatthew G. Knepley + comm - The communicator for the DM object
130724119c2aSMatthew G. Knepley . n    - The number of wedges around the origin
130824119c2aSMatthew G. Knepley - interpolate - Create edges and faces
130924119c2aSMatthew G. Knepley 
131024119c2aSMatthew G. Knepley   Output Parameter:
131124119c2aSMatthew G. Knepley . dm  - The DM object
131224119c2aSMatthew G. Knepley 
131324119c2aSMatthew G. Knepley   Level: beginner
131424119c2aSMatthew G. Knepley 
131524119c2aSMatthew G. Knepley .keywords: DM, create
1316*768d5fceSMatthew G. Knepley .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
131724119c2aSMatthew G. Knepley @*/
131824119c2aSMatthew G. Knepley PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
131924119c2aSMatthew G. Knepley {
132024119c2aSMatthew G. Knepley   const PetscInt dim = 3;
132124119c2aSMatthew G. Knepley   PetscInt       numCells, numVertices;
132224119c2aSMatthew G. Knepley   PetscErrorCode ierr;
132324119c2aSMatthew G. Knepley 
132424119c2aSMatthew G. Knepley   PetscFunctionBegin;
132524119c2aSMatthew G. Knepley   PetscValidPointer(dm, 3);
132624119c2aSMatthew G. Knepley   if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n);
132724119c2aSMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
132824119c2aSMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
132924119c2aSMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
133024119c2aSMatthew G. Knepley   /* Create topology */
133124119c2aSMatthew G. Knepley   {
133224119c2aSMatthew G. Knepley     PetscInt cone[6], c;
133324119c2aSMatthew G. Knepley 
133424119c2aSMatthew G. Knepley     numCells    = n;
133524119c2aSMatthew G. Knepley     numVertices = 2*(n+1);
133624119c2aSMatthew G. Knepley     ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
133724119c2aSMatthew G. Knepley     for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 6);CHKERRQ(ierr);}
133824119c2aSMatthew G. Knepley     ierr = DMSetUp(*dm);CHKERRQ(ierr);
133924119c2aSMatthew G. Knepley     for (c = 0; c < numCells; c++) {
134024119c2aSMatthew G. Knepley       cone[0] =  c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n;
134124119c2aSMatthew G. Knepley       cone[3] =  c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n;
134224119c2aSMatthew G. Knepley       ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
134324119c2aSMatthew G. Knepley     }
134424119c2aSMatthew G. Knepley     ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
134524119c2aSMatthew G. Knepley     ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
134624119c2aSMatthew G. Knepley   }
134724119c2aSMatthew G. Knepley   /* Interpolate */
134824119c2aSMatthew G. Knepley   if (interpolate) {
134924119c2aSMatthew G. Knepley     DM idm = NULL;
135024119c2aSMatthew G. Knepley 
135124119c2aSMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
135224119c2aSMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
135324119c2aSMatthew G. Knepley     *dm  = idm;
135424119c2aSMatthew G. Knepley   }
135524119c2aSMatthew G. Knepley   /* Create cylinder geometry */
135624119c2aSMatthew G. Knepley   {
135724119c2aSMatthew G. Knepley     Vec          coordinates;
135824119c2aSMatthew G. Knepley     PetscSection coordSection;
135924119c2aSMatthew G. Knepley     PetscScalar *coords;
136024119c2aSMatthew G. Knepley     PetscInt     coordSize, v, c;
136124119c2aSMatthew G. Knepley 
136224119c2aSMatthew G. Knepley     /* Build coordinates */
136324119c2aSMatthew G. Knepley     ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
136424119c2aSMatthew G. Knepley     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
136524119c2aSMatthew G. Knepley     ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
136624119c2aSMatthew G. Knepley     ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr);
136724119c2aSMatthew G. Knepley     for (v = numCells; v < numCells+numVertices; ++v) {
136824119c2aSMatthew G. Knepley       ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
136924119c2aSMatthew G. Knepley       ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
137024119c2aSMatthew G. Knepley     }
137124119c2aSMatthew G. Knepley     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
137224119c2aSMatthew G. Knepley     ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
137324119c2aSMatthew G. Knepley     ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
137424119c2aSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
137524119c2aSMatthew G. Knepley     ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
137624119c2aSMatthew G. Knepley     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
137724119c2aSMatthew G. Knepley     ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
137824119c2aSMatthew G. Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
137924119c2aSMatthew G. Knepley     for (c = 0; c < numCells; c++) {
138024119c2aSMatthew G. Knepley       coords[(c+0*n)*dim+0] = PetscCosReal(2.0*c*PETSC_PI/n); coords[(c+0*n)*dim+1] = PetscSinReal(2.0*c*PETSC_PI/n); coords[(c+0*n)*dim+2] = 1.0;
138124119c2aSMatthew G. Knepley       coords[(c+1*n)*dim+0] = PetscCosReal(2.0*c*PETSC_PI/n); coords[(c+1*n)*dim+1] = PetscSinReal(2.0*c*PETSC_PI/n); coords[(c+1*n)*dim+2] = 0.0;
138224119c2aSMatthew G. Knepley     }
138324119c2aSMatthew G. Knepley     coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0;
138424119c2aSMatthew G. Knepley     coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0;
138524119c2aSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
138624119c2aSMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
138724119c2aSMatthew G. Knepley     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
138824119c2aSMatthew G. Knepley   }
138924119c2aSMatthew G. Knepley   PetscFunctionReturn(0);
139024119c2aSMatthew G. Knepley }
139124119c2aSMatthew G. Knepley 
139265a81367SMatthew G. Knepley PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
139365a81367SMatthew G. Knepley {
139465a81367SMatthew G. Knepley   PetscReal prod = 0.0;
139565a81367SMatthew G. Knepley   PetscInt  i;
139665a81367SMatthew G. Knepley   for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
139765a81367SMatthew G. Knepley   return PetscSqrtReal(prod);
139865a81367SMatthew G. Knepley }
139965a81367SMatthew G. Knepley PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
140065a81367SMatthew G. Knepley {
140165a81367SMatthew G. Knepley   PetscReal prod = 0.0;
140265a81367SMatthew G. Knepley   PetscInt  i;
140365a81367SMatthew G. Knepley   for (i = 0; i < dim; ++i) prod += x[i]*y[i];
140465a81367SMatthew G. Knepley   return prod;
140565a81367SMatthew G. Knepley }
140665a81367SMatthew G. Knepley 
14072829fed8SMatthew G. Knepley /*@
140865a81367SMatthew G. Knepley   DMPlexCreateSphereMesh - Creates a mesh on the d-dimensional sphere, S^d.
14092829fed8SMatthew G. Knepley 
14102829fed8SMatthew G. Knepley   Collective on MPI_Comm
14112829fed8SMatthew G. Knepley 
14122829fed8SMatthew G. Knepley   Input Parameters:
14132829fed8SMatthew G. Knepley . comm  - The communicator for the DM object
141465a81367SMatthew G. Knepley . dim   - The dimension
141565a81367SMatthew G. Knepley - simplex - Use simplices, or tensor product cells
14162829fed8SMatthew G. Knepley 
14172829fed8SMatthew G. Knepley   Output Parameter:
14182829fed8SMatthew G. Knepley . dm  - The DM object
14192829fed8SMatthew G. Knepley 
14202829fed8SMatthew G. Knepley   Level: beginner
14212829fed8SMatthew G. Knepley 
14222829fed8SMatthew G. Knepley .keywords: DM, create
1423*768d5fceSMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
14242829fed8SMatthew G. Knepley @*/
142565a81367SMatthew G. Knepley PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *dm)
14262829fed8SMatthew G. Knepley {
142765a81367SMatthew G. Knepley   const PetscInt  embedDim = dim+1;
142865a81367SMatthew G. Knepley   PetscSection    coordSection;
142965a81367SMatthew G. Knepley   Vec             coordinates;
143065a81367SMatthew G. Knepley   PetscScalar    *coords;
143165a81367SMatthew G. Knepley   PetscReal      *coordsIn;
143265a81367SMatthew G. Knepley   PetscInt        numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e;
143365a81367SMatthew G. Knepley   PetscMPIInt     rank;
143465a81367SMatthew G. Knepley   PetscErrorCode  ierr;
143565a81367SMatthew G. Knepley 
143665a81367SMatthew G. Knepley   PetscFunctionBegin;
143765a81367SMatthew G. Knepley   PetscValidPointer(dm, 4);
143865a81367SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
143965a81367SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
144065a81367SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
144165a81367SMatthew G. Knepley   ierr = DMSetCoordinateDim(*dm, dim+1);CHKERRQ(ierr);
144265a81367SMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank);CHKERRQ(ierr);
144365a81367SMatthew G. Knepley   switch (dim) {
144465a81367SMatthew G. Knepley   case 2:
144565a81367SMatthew G. Knepley     if (simplex) {
1446116ded15SMatthew G. Knepley       DM              idm         = NULL;
1447116ded15SMatthew G. Knepley       const PetscReal edgeLen     = 2.0/(1.0 + PETSC_PHI);
1448116ded15SMatthew G. Knepley       const PetscReal vertex[3]   = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)};
144965a81367SMatthew G. Knepley       const PetscInt  degree      = 5;
145065a81367SMatthew G. Knepley       PetscInt        s[3]        = {1, 1, 1};
145165a81367SMatthew G. Knepley       PetscInt        cone[3];
145265a81367SMatthew G. Knepley       PetscInt       *graph, p, i, j, k;
145365a81367SMatthew G. Knepley 
145465a81367SMatthew G. Knepley       numCells    = !rank ? 20 : 0;
145565a81367SMatthew G. Knepley       numVerts    = !rank ? 12 : 0;
145665a81367SMatthew G. Knepley       firstVertex = numCells;
145765a81367SMatthew G. Knepley       /* Use icosahedron, which for a unit sphere has coordinates which are all cyclic permutations of
145865a81367SMatthew G. Knepley 
145965a81367SMatthew G. Knepley            (0, \pm 1/\phi+1, \pm \phi/\phi+1)
146065a81367SMatthew G. Knepley 
146165a81367SMatthew G. Knepley          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
146265a81367SMatthew G. Knepley          length is then given by 2/\phi = 2 * 2.73606 = 5.47214.
146365a81367SMatthew G. Knepley        */
146465a81367SMatthew G. Knepley       /* Construct vertices */
146565a81367SMatthew G. Knepley       ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr);
146665a81367SMatthew G. Knepley       for (p = 0, i = 0; p < embedDim; ++p) {
146765a81367SMatthew G. Knepley         for (s[1] = -1; s[1] < 2; s[1] += 2) {
146865a81367SMatthew G. Knepley           for (s[2] = -1; s[2] < 2; s[2] += 2) {
146965a81367SMatthew G. Knepley             for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim];
147065a81367SMatthew G. Knepley             ++i;
147165a81367SMatthew G. Knepley           }
147265a81367SMatthew G. Knepley         }
147365a81367SMatthew G. Knepley       }
147465a81367SMatthew G. Knepley       /* Construct graph */
147565a81367SMatthew G. Knepley       ierr = PetscCalloc1(numVerts * numVerts, &graph);CHKERRQ(ierr);
147665a81367SMatthew G. Knepley       for (i = 0; i < numVerts; ++i) {
147765a81367SMatthew G. Knepley         for (j = 0, k = 0; j < numVerts; ++j) {
147865a81367SMatthew G. Knepley           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
147965a81367SMatthew G. Knepley         }
148065a81367SMatthew G. Knepley         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid icosahedron, vertex %D degree %D != %D", i, k, degree);
148165a81367SMatthew G. Knepley       }
148265a81367SMatthew G. Knepley       /* Build Topology */
148365a81367SMatthew G. Knepley       ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr);
148465a81367SMatthew G. Knepley       for (c = 0; c < numCells; c++) {
1485116ded15SMatthew G. Knepley         ierr = DMPlexSetConeSize(*dm, c, embedDim);CHKERRQ(ierr);
148665a81367SMatthew G. Knepley       }
148765a81367SMatthew G. Knepley       ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */
148865a81367SMatthew G. Knepley       /* Cells */
148965a81367SMatthew G. Knepley       for (i = 0, c = 0; i < numVerts; ++i) {
149065a81367SMatthew G. Knepley         for (j = 0; j < i; ++j) {
149165a81367SMatthew G. Knepley           for (k = 0; k < j; ++k) {
149265a81367SMatthew G. Knepley             if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) {
149365a81367SMatthew G. Knepley               cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k;
149465a81367SMatthew G. Knepley               /* Check orientation */
149565a81367SMatthew G. Knepley               {
149665a81367SMatthew G. Knepley                 const PetscInt epsilon[3][3][3] = {{{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, {{0, 0, -1}, {0, 0, 0}, {1, 0, 0}}, {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0}}};
149765a81367SMatthew G. Knepley                 PetscReal normal[3];
149865a81367SMatthew G. Knepley                 PetscInt  e, f;
149965a81367SMatthew G. Knepley 
150065a81367SMatthew G. Knepley                 for (d = 0; d < embedDim; ++d) {
150165a81367SMatthew G. Knepley                   normal[d] = 0.0;
150265a81367SMatthew G. Knepley                   for (e = 0; e < embedDim; ++e) {
150365a81367SMatthew G. Knepley                     for (f = 0; f < embedDim; ++f) {
150465a81367SMatthew G. Knepley                       normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]);
150565a81367SMatthew G. Knepley                     }
150665a81367SMatthew G. Knepley                   }
150765a81367SMatthew G. Knepley                 }
150865a81367SMatthew G. Knepley                 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
150965a81367SMatthew G. Knepley               }
151065a81367SMatthew G. Knepley               ierr = DMPlexSetCone(*dm, c++, cone);CHKERRQ(ierr);
151165a81367SMatthew G. Knepley             }
151265a81367SMatthew G. Knepley           }
151365a81367SMatthew G. Knepley         }
151465a81367SMatthew G. Knepley       }
151565a81367SMatthew G. Knepley       ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
151665a81367SMatthew G. Knepley       ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
151765a81367SMatthew G. Knepley       ierr = PetscFree(graph);CHKERRQ(ierr);
151865a81367SMatthew G. Knepley       /* Interpolate mesh */
151965a81367SMatthew G. Knepley       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
152065a81367SMatthew G. Knepley       ierr = DMDestroy(dm);CHKERRQ(ierr);
152165a81367SMatthew G. Knepley       *dm  = idm;
152265a81367SMatthew G. Knepley     } else {
15232829fed8SMatthew G. Knepley       /*
15242829fed8SMatthew G. Knepley         12-21--13
15252829fed8SMatthew G. Knepley          |     |
15262829fed8SMatthew G. Knepley         25  4  24
15272829fed8SMatthew G. Knepley          |     |
15282829fed8SMatthew G. Knepley   12-25--9-16--8-24--13
15292829fed8SMatthew G. Knepley    |     |     |     |
15302829fed8SMatthew G. Knepley   23  5 17  0 15  3  22
15312829fed8SMatthew G. Knepley    |     |     |     |
15322829fed8SMatthew G. Knepley   10-20--6-14--7-19--11
15332829fed8SMatthew G. Knepley          |     |
15342829fed8SMatthew G. Knepley         20  1  19
15352829fed8SMatthew G. Knepley          |     |
15362829fed8SMatthew G. Knepley         10-18--11
15372829fed8SMatthew G. Knepley          |     |
15382829fed8SMatthew G. Knepley         23  2  22
15392829fed8SMatthew G. Knepley          |     |
15402829fed8SMatthew G. Knepley         12-21--13
15412829fed8SMatthew G. Knepley        */
154265a81367SMatthew G. Knepley       const PetscReal dist = 1.0/PetscSqrtReal(3.0);
15432829fed8SMatthew G. Knepley       PetscInt        cone[4], ornt[4];
15442829fed8SMatthew G. Knepley 
154565a81367SMatthew G. Knepley       numCells    = !rank ?  6 : 0;
154665a81367SMatthew G. Knepley       numEdges    = !rank ? 12 : 0;
154765a81367SMatthew G. Knepley       numVerts    = !rank ?  8 : 0;
154865a81367SMatthew G. Knepley       firstVertex = numCells;
154965a81367SMatthew G. Knepley       firstEdge   = numCells + numVerts;
15502829fed8SMatthew G. Knepley       /* Build Topology */
15512829fed8SMatthew G. Knepley       ierr = DMPlexSetChart(*dm, 0, numCells+numEdges+numVerts);CHKERRQ(ierr);
15522829fed8SMatthew G. Knepley       for (c = 0; c < numCells; c++) {
15532829fed8SMatthew G. Knepley         ierr = DMPlexSetConeSize(*dm, c, 4);CHKERRQ(ierr);
15542829fed8SMatthew G. Knepley       }
15552829fed8SMatthew G. Knepley       for (e = firstEdge; e < firstEdge+numEdges; ++e) {
15562829fed8SMatthew G. Knepley         ierr = DMPlexSetConeSize(*dm, e, 2);CHKERRQ(ierr);
15572829fed8SMatthew G. Knepley       }
15582829fed8SMatthew G. Knepley       ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */
15592829fed8SMatthew G. Knepley       /* Cell 0 */
15602829fed8SMatthew G. Knepley       cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17;
15612829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr);
15622829fed8SMatthew G. Knepley       ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0;
15632829fed8SMatthew G. Knepley       ierr = DMPlexSetConeOrientation(*dm, 0, ornt);CHKERRQ(ierr);
15642829fed8SMatthew G. Knepley       /* Cell 1 */
15652829fed8SMatthew G. Knepley       cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20;
15662829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr);
15672829fed8SMatthew G. Knepley       ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
15682829fed8SMatthew G. Knepley       ierr = DMPlexSetConeOrientation(*dm, 1, ornt);CHKERRQ(ierr);
15692829fed8SMatthew G. Knepley       /* Cell 2 */
15702829fed8SMatthew G. Knepley       cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23;
15712829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr);
15722829fed8SMatthew G. Knepley       ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
15732829fed8SMatthew G. Knepley       ierr = DMPlexSetConeOrientation(*dm, 2, ornt);CHKERRQ(ierr);
15742829fed8SMatthew G. Knepley       /* Cell 3 */
15752829fed8SMatthew G. Knepley       cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15;
15762829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr);
15772829fed8SMatthew G. Knepley       ornt[0] = -2; ornt[1] = -2; ornt[2] = 0; ornt[3] = -2;
15782829fed8SMatthew G. Knepley       ierr = DMPlexSetConeOrientation(*dm, 3, ornt);CHKERRQ(ierr);
15792829fed8SMatthew G. Knepley       /* Cell 4 */
15802829fed8SMatthew G. Knepley       cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25;
15812829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr);
15822829fed8SMatthew G. Knepley       ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = 0;
15832829fed8SMatthew G. Knepley       ierr = DMPlexSetConeOrientation(*dm, 4, ornt);CHKERRQ(ierr);
15842829fed8SMatthew G. Knepley       /* Cell 5 */
15852829fed8SMatthew G. Knepley       cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23;
15862829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr);
15872829fed8SMatthew G. Knepley       ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = -2;
15882829fed8SMatthew G. Knepley       ierr = DMPlexSetConeOrientation(*dm, 5, ornt);CHKERRQ(ierr);
15892829fed8SMatthew G. Knepley       /* Edges */
15902829fed8SMatthew G. Knepley       cone[0] =  6; cone[1] =  7;
15912829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr);
15922829fed8SMatthew G. Knepley       cone[0] =  7; cone[1] =  8;
15932829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 15, cone);CHKERRQ(ierr);
15942829fed8SMatthew G. Knepley       cone[0] =  8; cone[1] =  9;
15952829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 16, cone);CHKERRQ(ierr);
15962829fed8SMatthew G. Knepley       cone[0] =  9; cone[1] =  6;
15972829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 17, cone);CHKERRQ(ierr);
15982829fed8SMatthew G. Knepley       cone[0] = 10; cone[1] = 11;
15992829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 18, cone);CHKERRQ(ierr);
16002829fed8SMatthew G. Knepley       cone[0] = 11; cone[1] =  7;
16012829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 19, cone);CHKERRQ(ierr);
16022829fed8SMatthew G. Knepley       cone[0] =  6; cone[1] = 10;
16032829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 20, cone);CHKERRQ(ierr);
16042829fed8SMatthew G. Knepley       cone[0] = 12; cone[1] = 13;
16052829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 21, cone);CHKERRQ(ierr);
16062829fed8SMatthew G. Knepley       cone[0] = 13; cone[1] = 11;
16072829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 22, cone);CHKERRQ(ierr);
16082829fed8SMatthew G. Knepley       cone[0] = 10; cone[1] = 12;
16092829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 23, cone);CHKERRQ(ierr);
16102829fed8SMatthew G. Knepley       cone[0] = 13; cone[1] =  8;
16112829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 24, cone);CHKERRQ(ierr);
16122829fed8SMatthew G. Knepley       cone[0] = 12; cone[1] =  9;
16132829fed8SMatthew G. Knepley       ierr = DMPlexSetCone(*dm, 25, cone);CHKERRQ(ierr);
16142829fed8SMatthew G. Knepley       ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
16152829fed8SMatthew G. Knepley       ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
16162829fed8SMatthew G. Knepley       /* Build coordinates */
161765a81367SMatthew G. Knepley       ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr);
161865a81367SMatthew G. Knepley       coordsIn[0*embedDim+0] = -dist; coordsIn[0*embedDim+1] =  dist; coordsIn[0*embedDim+2] = -dist;
161965a81367SMatthew G. Knepley       coordsIn[1*embedDim+0] =  dist; coordsIn[1*embedDim+1] =  dist; coordsIn[1*embedDim+2] = -dist;
162065a81367SMatthew G. Knepley       coordsIn[2*embedDim+0] =  dist; coordsIn[2*embedDim+1] = -dist; coordsIn[2*embedDim+2] = -dist;
162165a81367SMatthew G. Knepley       coordsIn[3*embedDim+0] = -dist; coordsIn[3*embedDim+1] = -dist; coordsIn[3*embedDim+2] = -dist;
162265a81367SMatthew G. Knepley       coordsIn[4*embedDim+0] = -dist; coordsIn[4*embedDim+1] =  dist; coordsIn[4*embedDim+2] =  dist;
162365a81367SMatthew G. Knepley       coordsIn[5*embedDim+0] =  dist; coordsIn[5*embedDim+1] =  dist; coordsIn[5*embedDim+2] =  dist;
162465a81367SMatthew G. Knepley       coordsIn[6*embedDim+0] = -dist; coordsIn[6*embedDim+1] = -dist; coordsIn[6*embedDim+2] =  dist;
162565a81367SMatthew G. Knepley       coordsIn[7*embedDim+0] =  dist; coordsIn[7*embedDim+1] = -dist; coordsIn[7*embedDim+2] =  dist;
162665a81367SMatthew G. Knepley     }
162765a81367SMatthew G. Knepley     break;
162865a81367SMatthew G. Knepley   case 3:
1629116ded15SMatthew G. Knepley     if (simplex) {
1630116ded15SMatthew G. Knepley       DM              idm             = NULL;
1631116ded15SMatthew G. Knepley       const PetscReal edgeLen         = 1.0/PETSC_PHI;
1632116ded15SMatthew G. Knepley       const PetscReal vertexA[4]      = {0.5, 0.5, 0.5, 0.5};
1633116ded15SMatthew G. Knepley       const PetscReal vertexB[4]      = {1.0, 0.0, 0.0, 0.0};
1634116ded15SMatthew G. Knepley       const PetscReal vertexC[4]      = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0};
1635116ded15SMatthew G. Knepley       const PetscInt  degree          = 12;
1636116ded15SMatthew G. Knepley       PetscInt        s[4]            = {1, 1, 1};
1637116ded15SMatthew G. Knepley       PetscInt        evenPerm[12][4] = {{0, 1, 2, 3}, {0, 2, 3, 1}, {0, 3, 1, 2}, {1, 0, 3, 2}, {1, 2, 0, 3}, {1, 3, 2, 0},
1638116ded15SMatthew G. Knepley                                          {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}};
1639116ded15SMatthew G. Knepley       PetscInt        cone[4];
1640116ded15SMatthew G. Knepley       PetscInt       *graph, p, i, j, k, l;
1641116ded15SMatthew G. Knepley 
1642116ded15SMatthew G. Knepley       numCells    = !rank ? 600 : 0;
1643116ded15SMatthew G. Knepley       numVerts    = !rank ? 120 : 0;
1644116ded15SMatthew G. Knepley       firstVertex = numCells;
1645116ded15SMatthew G. Knepley       /* Use the 600-cell, which for a unit sphere has coordinates which are
1646116ded15SMatthew G. Knepley 
1647116ded15SMatthew G. Knepley            1/2 (\pm 1, \pm 1,    \pm 1, \pm 1)                          16
1648116ded15SMatthew G. Knepley                (\pm 1,    0,       0,      0)  all cyclic permutations   8
1649116ded15SMatthew G. Knepley            1/2 (\pm 1, \pm phi, \pm 1/phi, 0)  all even permutations    96
1650116ded15SMatthew G. Knepley 
1651116ded15SMatthew G. Knepley          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
1652116ded15SMatthew G. Knepley          length is then given by 1/\phi = 2.73606.
1653116ded15SMatthew G. Knepley 
1654116ded15SMatthew G. Knepley          http://buzzard.pugetsound.edu/sage-practice/ch03s03.html
1655116ded15SMatthew G. Knepley          http://mathworld.wolfram.com/600-Cell.html
1656116ded15SMatthew G. Knepley        */
1657116ded15SMatthew G. Knepley       /* Construct vertices */
1658116ded15SMatthew G. Knepley       ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr);
1659116ded15SMatthew G. Knepley       i    = 0;
1660116ded15SMatthew G. Knepley       for (s[0] = -1; s[0] < 2; s[0] += 2) {
1661116ded15SMatthew G. Knepley         for (s[1] = -1; s[1] < 2; s[1] += 2) {
1662116ded15SMatthew G. Knepley           for (s[2] = -1; s[2] < 2; s[2] += 2) {
1663116ded15SMatthew G. Knepley             for (s[3] = -1; s[3] < 2; s[3] += 2) {
1664116ded15SMatthew G. Knepley               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d];
1665116ded15SMatthew G. Knepley               ++i;
1666116ded15SMatthew G. Knepley             }
1667116ded15SMatthew G. Knepley           }
1668116ded15SMatthew G. Knepley         }
1669116ded15SMatthew G. Knepley       }
1670116ded15SMatthew G. Knepley       for (p = 0; p < embedDim; ++p) {
1671116ded15SMatthew G. Knepley         s[1] = s[2] = s[3] = 1;
1672116ded15SMatthew G. Knepley         for (s[0] = -1; s[0] < 2; s[0] += 2) {
1673116ded15SMatthew G. Knepley           for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim];
1674116ded15SMatthew G. Knepley           ++i;
1675116ded15SMatthew G. Knepley         }
1676116ded15SMatthew G. Knepley       }
1677116ded15SMatthew G. Knepley       for (p = 0; p < 12; ++p) {
1678116ded15SMatthew G. Knepley         s[3] = 1;
1679116ded15SMatthew G. Knepley         for (s[0] = -1; s[0] < 2; s[0] += 2) {
1680116ded15SMatthew G. Knepley           for (s[1] = -1; s[1] < 2; s[1] += 2) {
1681116ded15SMatthew G. Knepley             for (s[2] = -1; s[2] < 2; s[2] += 2) {
1682116ded15SMatthew G. Knepley               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]];
1683116ded15SMatthew G. Knepley               ++i;
1684116ded15SMatthew G. Knepley             }
1685116ded15SMatthew G. Knepley           }
1686116ded15SMatthew G. Knepley         }
1687116ded15SMatthew G. Knepley       }
1688116ded15SMatthew G. Knepley       if (i != numVerts) SETERRQ2(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertices %D != %D", i, numVerts);
1689116ded15SMatthew G. Knepley       /* Construct graph */
1690116ded15SMatthew G. Knepley       ierr = PetscCalloc1(numVerts * numVerts, &graph);CHKERRQ(ierr);
1691116ded15SMatthew G. Knepley       for (i = 0; i < numVerts; ++i) {
1692116ded15SMatthew G. Knepley         for (j = 0, k = 0; j < numVerts; ++j) {
1693116ded15SMatthew G. Knepley           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
1694116ded15SMatthew G. Knepley         }
1695116ded15SMatthew G. Knepley         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertex %D degree %D != %D", i, k, degree);
1696116ded15SMatthew G. Knepley       }
1697116ded15SMatthew G. Knepley       /* Build Topology */
1698116ded15SMatthew G. Knepley       ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr);
1699116ded15SMatthew G. Knepley       for (c = 0; c < numCells; c++) {
1700116ded15SMatthew G. Knepley         ierr = DMPlexSetConeSize(*dm, c, embedDim);CHKERRQ(ierr);
1701116ded15SMatthew G. Knepley       }
1702116ded15SMatthew G. Knepley       ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */
1703116ded15SMatthew G. Knepley       /* Cells */
1704116ded15SMatthew G. Knepley       for (i = 0, c = 0; i < numVerts; ++i) {
1705116ded15SMatthew G. Knepley         for (j = 0; j < i; ++j) {
1706116ded15SMatthew G. Knepley           for (k = 0; k < j; ++k) {
1707116ded15SMatthew G. Knepley             for (l = 0; l < k; ++l) {
1708116ded15SMatthew G. Knepley               if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] &&
1709116ded15SMatthew G. Knepley                   graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) {
1710116ded15SMatthew G. Knepley                 cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l;
1711116ded15SMatthew G. Knepley                 /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */
1712116ded15SMatthew G. Knepley                 {
1713116ded15SMatthew G. Knepley                   const PetscInt epsilon[4][4][4][4] = {{{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
1714116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  1}, { 0,  0, -1, 0}},
1715116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  0}, { 0, 0,  0, -1}, { 0,  0, 0,  0}, { 0,  1,  0, 0}},
1716116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  0}, { 0, 0,  1,  0}, { 0, -1, 0,  0}, { 0,  0,  0, 0}}},
1717116ded15SMatthew G. Knepley 
1718116ded15SMatthew G. Knepley                                                         {{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0, -1}, { 0,  0,  1, 0}},
1719116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
1720116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, {-1,  0,  0, 0}},
1721116ded15SMatthew G. Knepley                                                          {{0,  0, -1,  0}, { 0, 0,  0,  0}, { 1,  0, 0,  0}, { 0,  0,  0, 0}}},
1722116ded15SMatthew G. Knepley 
1723116ded15SMatthew G. Knepley                                                         {{{0,  0,  0,  0}, { 0, 0,  0,  1}, { 0,  0, 0,  0}, { 0, -1,  0, 0}},
1724116ded15SMatthew G. Knepley                                                          {{0,  0,  0, -1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 1,  0,  0, 0}},
1725116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
1726116ded15SMatthew G. Knepley                                                          {{0,  1,  0,  0}, {-1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}},
1727116ded15SMatthew G. Knepley 
1728116ded15SMatthew G. Knepley                                                         {{{0,  0,  0,  0}, { 0, 0, -1,  0}, { 0,  1, 0,  0}, { 0,  0,  0, 0}},
1729116ded15SMatthew G. Knepley                                                          {{0,  0,  1,  0}, { 0, 0,  0,  0}, {-1,  0, 0,  0}, { 0,  0,  0, 0}},
1730116ded15SMatthew G. Knepley                                                          {{0, -1,  0,  0}, { 1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
1731116ded15SMatthew G. Knepley                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}}};
1732116ded15SMatthew G. Knepley                   PetscReal normal[4];
1733116ded15SMatthew G. Knepley                   PetscInt  e, f, g;
1734116ded15SMatthew G. Knepley 
1735116ded15SMatthew G. Knepley                   for (d = 0; d < embedDim; ++d) {
1736116ded15SMatthew G. Knepley                     normal[d] = 0.0;
1737116ded15SMatthew G. Knepley                     for (e = 0; e < embedDim; ++e) {
1738116ded15SMatthew G. Knepley                       for (f = 0; f < embedDim; ++f) {
1739116ded15SMatthew G. Knepley                         for (g = 0; g < embedDim; ++g) {
1740116ded15SMatthew G. Knepley                           normal[d] += epsilon[d][e][f][g]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f])*(coordsIn[l*embedDim+f] - coordsIn[i*embedDim+f]);
1741116ded15SMatthew G. Knepley                         }
1742116ded15SMatthew G. Knepley                       }
1743116ded15SMatthew G. Knepley                     }
1744116ded15SMatthew G. Knepley                   }
1745116ded15SMatthew G. Knepley                   if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
1746116ded15SMatthew G. Knepley                 }
1747116ded15SMatthew G. Knepley                 ierr = DMPlexSetCone(*dm, c++, cone);CHKERRQ(ierr);
1748116ded15SMatthew G. Knepley               }
1749116ded15SMatthew G. Knepley             }
1750116ded15SMatthew G. Knepley           }
1751116ded15SMatthew G. Knepley         }
1752116ded15SMatthew G. Knepley       }
1753116ded15SMatthew G. Knepley       ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
1754116ded15SMatthew G. Knepley       ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
1755116ded15SMatthew G. Knepley       ierr = PetscFree(graph);CHKERRQ(ierr);
1756116ded15SMatthew G. Knepley       /* Interpolate mesh */
1757116ded15SMatthew G. Knepley       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1758116ded15SMatthew G. Knepley       ierr = DMDestroy(dm);CHKERRQ(ierr);
1759116ded15SMatthew G. Knepley       *dm  = idm;
1760116ded15SMatthew G. Knepley       break;
1761116ded15SMatthew G. Knepley     }
176265a81367SMatthew G. Knepley   default: SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim);
176365a81367SMatthew G. Knepley   }
176465a81367SMatthew G. Knepley   /* Create coordinates */
17652829fed8SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
17662829fed8SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
17672829fed8SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
17682829fed8SMatthew G. Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);CHKERRQ(ierr);
17692829fed8SMatthew G. Knepley   for (v = firstVertex; v < firstVertex+numVerts; ++v) {
17702829fed8SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
17712829fed8SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
17722829fed8SMatthew G. Knepley   }
17732829fed8SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
17742829fed8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
17752829fed8SMatthew G. Knepley   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
17762829fed8SMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
17772829fed8SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
17782829fed8SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
17792829fed8SMatthew G. Knepley   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
17802829fed8SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
178165a81367SMatthew G. Knepley   for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];}
17822829fed8SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
17832829fed8SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
17842829fed8SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
178565a81367SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
17862829fed8SMatthew G. Knepley   PetscFunctionReturn(0);
17872829fed8SMatthew G. Knepley }
17882829fed8SMatthew G. Knepley 
1789552f7358SJed Brown /* External function declarations here */
1790552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
17916dbf9973SLawrence Mitchell extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
1792fd59a867SMatthew G. Knepley extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
179366ad2231SToby Isaac extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
1794b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
1795552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
1796d9879d6bSMatthew G. Knepley PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
1797552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm);
1798552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm);
1799552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
18002c40f234SMatthew G. Knepley extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
1801552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
180229c8bc5fSMatthew G. Knepley static PetscErrorCode DMInitialize_Plex(DM dm);
1803552f7358SJed Brown 
18040a6ba040SMatthew G. Knepley /* Replace dm with the contents of dmNew
18050a6ba040SMatthew G. Knepley    - Share the DM_Plex structure
18060a6ba040SMatthew G. Knepley    - Share the coordinates
1807d7973521SMatthew G. Knepley    - Share the SF
18080a6ba040SMatthew G. Knepley */
18090a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
18100a6ba040SMatthew G. Knepley {
1811d7973521SMatthew G. Knepley   PetscSF               sf;
1812acdc6f61SToby Isaac   DM                    coordDM, coarseDM;
18130a6ba040SMatthew G. Knepley   Vec                   coords;
181490b157c4SStefano Zampini   PetscBool             isper;
181555fbe3e3SMatthew G. Knepley   const PetscReal      *maxCell, *L;
18165dc8c3f7SMatthew G. Knepley   const DMBoundaryType *bd;
18170a6ba040SMatthew G. Knepley   PetscErrorCode        ierr;
18180a6ba040SMatthew G. Knepley 
18190a6ba040SMatthew G. Knepley   PetscFunctionBegin;
1820d7973521SMatthew G. Knepley   ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr);
1821d7973521SMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
182255fbe3e3SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr);
18230a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr);
182455fbe3e3SMatthew G. Knepley   ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr);
18250a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
182690b157c4SStefano Zampini   ierr = DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);CHKERRQ(ierr);
182790b157c4SStefano Zampini   ierr = DMSetPeriodicity(dmNew, isper, maxCell, L, bd);CHKERRQ(ierr);
18280a6ba040SMatthew G. Knepley   ierr = DMDestroy_Plex(dm);CHKERRQ(ierr);
182929c8bc5fSMatthew G. Knepley   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
18300a6ba040SMatthew G. Knepley   dm->data = dmNew->data;
18310a6ba040SMatthew G. Knepley   ((DM_Plex *) dmNew->data)->refct++;
1832c58f1c22SToby Isaac   dmNew->labels->refct++;
1833c58f1c22SToby Isaac   if (!--(dm->labels->refct)) {
1834c58f1c22SToby Isaac     DMLabelLink next = dm->labels->next;
1835c58f1c22SToby Isaac 
1836c58f1c22SToby Isaac     /* destroy the labels */
1837c58f1c22SToby Isaac     while (next) {
1838c58f1c22SToby Isaac       DMLabelLink tmp = next->next;
1839c58f1c22SToby Isaac 
1840c58f1c22SToby Isaac       ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr);
1841c58f1c22SToby Isaac       ierr = PetscFree(next);CHKERRQ(ierr);
1842c58f1c22SToby Isaac       next = tmp;
1843c58f1c22SToby Isaac     }
1844c58f1c22SToby Isaac     ierr = PetscFree(dm->labels);CHKERRQ(ierr);
1845c58f1c22SToby Isaac   }
1846c58f1c22SToby Isaac   dm->labels = dmNew->labels;
1847c58f1c22SToby Isaac   dm->depthLabel = dmNew->depthLabel;
1848acdc6f61SToby Isaac   ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr);
1849acdc6f61SToby Isaac   ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr);
18500a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
18510a6ba040SMatthew G. Knepley }
18520a6ba040SMatthew G. Knepley 
18530a6ba040SMatthew G. Knepley /* Swap dm with the contents of dmNew
18540a6ba040SMatthew G. Knepley    - Swap the DM_Plex structure
18550a6ba040SMatthew G. Knepley    - Swap the coordinates
185679a015ccSMatthew G. Knepley    - Swap the point PetscSF
18570a6ba040SMatthew G. Knepley */
18580a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
18590a6ba040SMatthew G. Knepley {
18600a6ba040SMatthew G. Knepley   DM              coordDMA, coordDMB;
18610a6ba040SMatthew G. Knepley   Vec             coordsA,  coordsB;
186279a015ccSMatthew G. Knepley   PetscSF         sfA,      sfB;
18630a6ba040SMatthew G. Knepley   void            *tmp;
1864c58f1c22SToby Isaac   DMLabelLinkList listTmp;
1865c58f1c22SToby Isaac   DMLabel         depthTmp;
18660a6ba040SMatthew G. Knepley   PetscErrorCode  ierr;
18670a6ba040SMatthew G. Knepley 
18680a6ba040SMatthew G. Knepley   PetscFunctionBegin;
186979a015ccSMatthew G. Knepley   ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr);
187079a015ccSMatthew G. Knepley   ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr);
187179a015ccSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr);
187279a015ccSMatthew G. Knepley   ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr);
187379a015ccSMatthew G. Knepley   ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr);
187479a015ccSMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr);
187579a015ccSMatthew G. Knepley 
18760a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr);
18770a6ba040SMatthew G. Knepley   ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr);
18780a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr);
18790a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr);
18800a6ba040SMatthew G. Knepley   ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr);
18810a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr);
18820a6ba040SMatthew G. Knepley 
18830a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr);
18840a6ba040SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr);
18850a6ba040SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr);
18860a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr);
18870a6ba040SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr);
18880a6ba040SMatthew G. Knepley   ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr);
1889acdc6f61SToby Isaac 
18900a6ba040SMatthew G. Knepley   tmp       = dmA->data;
18910a6ba040SMatthew G. Knepley   dmA->data = dmB->data;
18920a6ba040SMatthew G. Knepley   dmB->data = tmp;
1893c58f1c22SToby Isaac   listTmp   = dmA->labels;
1894c58f1c22SToby Isaac   dmA->labels = dmB->labels;
1895c58f1c22SToby Isaac   dmB->labels = listTmp;
1896c58f1c22SToby Isaac   depthTmp  = dmA->depthLabel;
1897c58f1c22SToby Isaac   dmA->depthLabel = dmB->depthLabel;
1898c58f1c22SToby Isaac   dmB->depthLabel = depthTmp;
18990a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
19000a6ba040SMatthew G. Knepley }
19010a6ba040SMatthew G. Knepley 
190247920aaeSMatthew G. Knepley PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
19030a6ba040SMatthew G. Knepley {
19040a6ba040SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
19050a6ba040SMatthew G. Knepley   PetscErrorCode ierr;
19060a6ba040SMatthew G. Knepley 
19070a6ba040SMatthew G. Knepley   PetscFunctionBegin;
19080a6ba040SMatthew G. Knepley   /* Handle viewing */
19090a6ba040SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr);
19100a6ba040SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr);
191194ae4db5SBarry Smith   ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr);
1912953fc75cSMatthew G. Knepley   /* Point Location */
1913953fc75cSMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr);
19142e62ab5aSMatthew G. Knepley   /* Generation and remeshing */
19152e62ab5aSMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr);
1916b29cfa1cSToby Isaac   /* Projection behavior */
1917b29cfa1cSToby Isaac   ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr);
191864141f95SMatthew G. Knepley   ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr);
19194f3833eaSMatthew G. Knepley 
19204f3833eaSMatthew G. Knepley   ierr = PetscPartitionerSetFromOptions(mesh->partitioner);CHKERRQ(ierr);
192168d4fef7SMatthew G. Knepley   PetscFunctionReturn(0);
192268d4fef7SMatthew G. Knepley }
192368d4fef7SMatthew G. Knepley 
192446fa42a0SMatthew G. Knepley static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
192568d4fef7SMatthew G. Knepley {
1926b653a561SMatthew G. Knepley   PetscInt       refine = 0, coarsen = 0, r;
1927b653a561SMatthew G. Knepley   PetscBool      isHierarchy;
192868d4fef7SMatthew G. Knepley   PetscErrorCode ierr;
192968d4fef7SMatthew G. Knepley 
193068d4fef7SMatthew G. Knepley   PetscFunctionBegin;
193168d4fef7SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
19321a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr);
193368d4fef7SMatthew G. Knepley   /* Handle DMPlex refinement */
193468d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr);
193568d4fef7SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr);
1936b6a0289aSMatthew G. Knepley   if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);}
193768d4fef7SMatthew G. Knepley   if (refine && isHierarchy) {
1938acdc6f61SToby Isaac     DM *dms, coarseDM;
193968d4fef7SMatthew G. Knepley 
1940acdc6f61SToby Isaac     ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr);
1941acdc6f61SToby Isaac     ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
194268d4fef7SMatthew G. Knepley     ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr);
194368d4fef7SMatthew G. Knepley     ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr);
194468d4fef7SMatthew G. Knepley     /* Total hack since we do not pass in a pointer */
194568d4fef7SMatthew G. Knepley     ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr);
194668d4fef7SMatthew G. Knepley     if (refine == 1) {
1947a8fb8f29SToby Isaac       ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr);
19480aef6b92SMatthew G. Knepley       ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr);
194968d4fef7SMatthew G. Knepley     } else {
1950a8fb8f29SToby Isaac       ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr);
19510aef6b92SMatthew G. Knepley       ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr);
1952a8fb8f29SToby Isaac       ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr);
19530aef6b92SMatthew G. Knepley       ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr);
195468d4fef7SMatthew G. Knepley     }
1955acdc6f61SToby Isaac     ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr);
1956acdc6f61SToby Isaac     ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr);
195768d4fef7SMatthew G. Knepley     /* Free DMs */
195868d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
1959547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr);
196068d4fef7SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
196168d4fef7SMatthew G. Knepley     }
196268d4fef7SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
196368d4fef7SMatthew G. Knepley   } else {
196468d4fef7SMatthew G. Knepley     for (r = 0; r < refine; ++r) {
196568d4fef7SMatthew G. Knepley       DM refinedMesh;
196668d4fef7SMatthew G. Knepley 
19671a1499c8SBarry Smith       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
196868d4fef7SMatthew G. Knepley       ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr);
196968d4fef7SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
197068d4fef7SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr);
1971547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
197268d4fef7SMatthew G. Knepley       ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr);
197368d4fef7SMatthew G. Knepley     }
197468d4fef7SMatthew G. Knepley   }
19753cf6fe12SMatthew G. Knepley   /* Handle DMPlex coarsening */
1976b653a561SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr);
1977b653a561SMatthew G. Knepley   ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr);
1978b653a561SMatthew G. Knepley   if (coarsen && isHierarchy) {
1979b653a561SMatthew G. Knepley     DM *dms;
1980b653a561SMatthew G. Knepley 
1981b653a561SMatthew G. Knepley     ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr);
1982b653a561SMatthew G. Knepley     ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr);
1983b653a561SMatthew G. Knepley     /* Free DMs */
1984b653a561SMatthew G. Knepley     for (r = 0; r < coarsen; ++r) {
1985547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr);
1986b653a561SMatthew G. Knepley       ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);
1987b653a561SMatthew G. Knepley     }
1988b653a561SMatthew G. Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
1989b653a561SMatthew G. Knepley   } else {
1990b653a561SMatthew G. Knepley     for (r = 0; r < coarsen; ++r) {
19913cf6fe12SMatthew G. Knepley       DM coarseMesh;
19923cf6fe12SMatthew G. Knepley 
19933cf6fe12SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
19943cf6fe12SMatthew G. Knepley       ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr);
19953cf6fe12SMatthew G. Knepley       /* Total hack since we do not pass in a pointer */
19963cf6fe12SMatthew G. Knepley       ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr);
1997547f7119SMatthew G. Knepley       ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
19983cf6fe12SMatthew G. Knepley       ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr);
19993cf6fe12SMatthew G. Knepley     }
2000b653a561SMatthew G. Knepley   }
20013cf6fe12SMatthew G. Knepley   /* Handle */
20021a1499c8SBarry Smith   ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr);
20030a6ba040SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
20040a6ba040SMatthew G. Knepley   PetscFunctionReturn(0);
20050a6ba040SMatthew G. Knepley }
20060a6ba040SMatthew G. Knepley 
2007552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
2008552f7358SJed Brown {
2009552f7358SJed Brown   PetscErrorCode ierr;
2010552f7358SJed Brown 
2011552f7358SJed Brown   PetscFunctionBegin;
2012552f7358SJed Brown   ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr);
2013552f7358SJed Brown   /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
2014552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr);
2015d930f514SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr);
20162c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr);
2017d930f514SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr);
2018552f7358SJed Brown   PetscFunctionReturn(0);
2019552f7358SJed Brown }
2020552f7358SJed Brown 
2021552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
2022552f7358SJed Brown {
2023552f7358SJed Brown   PetscErrorCode ierr;
2024552f7358SJed Brown 
2025552f7358SJed Brown   PetscFunctionBegin;
2026552f7358SJed Brown   ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr);
2027552f7358SJed Brown   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr);
20282c40f234SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr);
2029552f7358SJed Brown   PetscFunctionReturn(0);
2030552f7358SJed Brown }
2031552f7358SJed Brown 
2032793f3fe5SMatthew G. Knepley static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
2033793f3fe5SMatthew G. Knepley {
2034793f3fe5SMatthew G. Knepley   PetscInt       depth, d;
2035793f3fe5SMatthew G. Knepley   PetscErrorCode ierr;
2036793f3fe5SMatthew G. Knepley 
2037793f3fe5SMatthew G. Knepley   PetscFunctionBegin;
2038793f3fe5SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2039793f3fe5SMatthew G. Knepley   if (depth == 1) {
2040793f3fe5SMatthew G. Knepley     ierr = DMGetDimension(dm, &d);CHKERRQ(ierr);
2041793f3fe5SMatthew G. Knepley     if (dim == 0)      {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);}
2042793f3fe5SMatthew G. Knepley     else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);}
2043793f3fe5SMatthew G. Knepley     else               {*pStart = 0; *pEnd = 0;}
2044793f3fe5SMatthew G. Knepley   } else {
2045793f3fe5SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);
2046793f3fe5SMatthew G. Knepley   }
2047793f3fe5SMatthew G. Knepley   PetscFunctionReturn(0);
2048793f3fe5SMatthew G. Knepley }
2049793f3fe5SMatthew G. Knepley 
20502f356facSMatthew G. Knepley static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
2051502a2867SDave May {
2052502a2867SDave May   PetscSF        sf;
20532f356facSMatthew G. Knepley   PetscErrorCode ierr;
2054502a2867SDave May 
20552f356facSMatthew G. Knepley   PetscFunctionBegin;
2056502a2867SDave May   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
2057502a2867SDave May   ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr);
2058502a2867SDave May   PetscFunctionReturn(0);
2059502a2867SDave May }
2060502a2867SDave May 
206146fa42a0SMatthew G. Knepley static PetscErrorCode DMInitialize_Plex(DM dm)
2062552f7358SJed Brown {
2063713918a9SToby Isaac   PetscErrorCode ierr;
2064713918a9SToby Isaac 
2065552f7358SJed Brown   PetscFunctionBegin;
2066552f7358SJed Brown   dm->ops->view                            = DMView_Plex;
20672c40f234SMatthew G. Knepley   dm->ops->load                            = DMLoad_Plex;
2068552f7358SJed Brown   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
206938221697SMatthew G. Knepley   dm->ops->clone                           = DMClone_Plex;
2070552f7358SJed Brown   dm->ops->setup                           = DMSetUp_Plex;
2071fd59a867SMatthew G. Knepley   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
207266ad2231SToby Isaac   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
2073552f7358SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
2074552f7358SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
2075184d77edSJed Brown   dm->ops->getlocaltoglobalmapping         = NULL;
20760298fd71SBarry Smith   dm->ops->createfieldis                   = NULL;
2077552f7358SJed Brown   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
20780a6ba040SMatthew G. Knepley   dm->ops->getcoloring                     = NULL;
2079552f7358SJed Brown   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
2080bceba477SMatthew G. Knepley   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
2081bceba477SMatthew G. Knepley   dm->ops->getaggregates                   = NULL;
2082bceba477SMatthew G. Knepley   dm->ops->getinjection                    = DMCreateInjection_Plex;
2083552f7358SJed Brown   dm->ops->refine                          = DMRefine_Plex;
20840a6ba040SMatthew G. Knepley   dm->ops->coarsen                         = DMCoarsen_Plex;
20850a6ba040SMatthew G. Knepley   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
2086b653a561SMatthew G. Knepley   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
20870d1cd5e0SMatthew G. Knepley   dm->ops->adaptlabel                      = DMAdaptLabel_Plex;
20880d1cd5e0SMatthew G. Knepley   dm->ops->adaptmetric                     = DMAdaptMetric_Plex;
20890298fd71SBarry Smith   dm->ops->globaltolocalbegin              = NULL;
20900298fd71SBarry Smith   dm->ops->globaltolocalend                = NULL;
20910298fd71SBarry Smith   dm->ops->localtoglobalbegin              = NULL;
20920298fd71SBarry Smith   dm->ops->localtoglobalend                = NULL;
2093552f7358SJed Brown   dm->ops->destroy                         = DMDestroy_Plex;
2094552f7358SJed Brown   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
2095793f3fe5SMatthew G. Knepley   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
2096552f7358SJed Brown   dm->ops->locatepoints                    = DMLocatePoints_Plex;
20970709b2feSToby Isaac   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
20980709b2feSToby Isaac   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
2099bfc4295aSToby Isaac   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
21008c6c5593SMatthew G. Knepley   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
21010709b2feSToby Isaac   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
2102b698f381SToby Isaac   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
21032a16baeaSToby Isaac   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
2104502a2867SDave May   dm->ops->getneighbors                    = DMGetNeighors_Plex;
2105f1d73a7aSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);CHKERRQ(ierr);
21068135c375SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);CHKERRQ(ierr);
2107552f7358SJed Brown   PetscFunctionReturn(0);
2108552f7358SJed Brown }
2109552f7358SJed Brown 
211046fa42a0SMatthew G. Knepley PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
211163a16f15SMatthew G. Knepley {
211263a16f15SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
211363a16f15SMatthew G. Knepley   PetscErrorCode ierr;
211463a16f15SMatthew G. Knepley 
211563a16f15SMatthew G. Knepley   PetscFunctionBegin;
211663a16f15SMatthew G. Knepley   mesh->refct++;
211763a16f15SMatthew G. Knepley   (*newdm)->data = mesh;
211863a16f15SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr);
211963a16f15SMatthew G. Knepley   ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr);
212063a16f15SMatthew G. Knepley   PetscFunctionReturn(0);
212163a16f15SMatthew G. Knepley }
212263a16f15SMatthew G. Knepley 
21238818961aSMatthew G Knepley /*MC
21248818961aSMatthew G Knepley   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
21258818961aSMatthew G Knepley                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
21268818961aSMatthew G Knepley                     specified by a PetscSection object. Ownership in the global representation is determined by
21278818961aSMatthew G Knepley                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.
21288818961aSMatthew G Knepley 
21298818961aSMatthew G Knepley   Level: intermediate
21308818961aSMatthew G Knepley 
21318818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
21328818961aSMatthew G Knepley M*/
21338818961aSMatthew G Knepley 
21348cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
2135552f7358SJed Brown {
2136552f7358SJed Brown   DM_Plex        *mesh;
2137770b213bSMatthew G Knepley   PetscInt       unit, d;
2138552f7358SJed Brown   PetscErrorCode ierr;
2139552f7358SJed Brown 
2140552f7358SJed Brown   PetscFunctionBegin;
2141552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2142b00a9115SJed Brown   ierr     = PetscNewLog(dm,&mesh);CHKERRQ(ierr);
2143c73cfb54SMatthew G. Knepley   dm->dim  = 0;
2144552f7358SJed Brown   dm->data = mesh;
2145552f7358SJed Brown 
2146552f7358SJed Brown   mesh->refct             = 1;
214782f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr);
2148552f7358SJed Brown   mesh->maxConeSize       = 0;
21490298fd71SBarry Smith   mesh->cones             = NULL;
21500298fd71SBarry Smith   mesh->coneOrientations  = NULL;
215182f516ccSBarry Smith   ierr                    = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr);
2152552f7358SJed Brown   mesh->maxSupportSize    = 0;
21530298fd71SBarry Smith   mesh->supports          = NULL;
2154552f7358SJed Brown   mesh->refinementUniform = PETSC_TRUE;
2155552f7358SJed Brown   mesh->refinementLimit   = -1.0;
2156552f7358SJed Brown 
21570298fd71SBarry Smith   mesh->facesTmp = NULL;
2158552f7358SJed Brown 
2159d9deefdfSMatthew G. Knepley   mesh->tetgenOpts   = NULL;
2160d9deefdfSMatthew G. Knepley   mesh->triangleOpts = NULL;
216177623264SMatthew G. Knepley   ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr);
21622e62ab5aSMatthew G. Knepley   mesh->remeshBd     = PETSC_FALSE;
2163d9deefdfSMatthew G. Knepley 
21640298fd71SBarry Smith   mesh->subpointMap = NULL;
2165552f7358SJed Brown 
21668865f1eaSKarl Rupp   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
2167552f7358SJed Brown 
21680aef6b92SMatthew G. Knepley   mesh->regularRefinement   = PETSC_FALSE;
2169df0420ecSMatthew G. Knepley   mesh->depthState          = -1;
21700298fd71SBarry Smith   mesh->globalVertexNumbers = NULL;
21710298fd71SBarry Smith   mesh->globalCellNumbers   = NULL;
2172a68b90caSToby Isaac   mesh->anchorSection       = NULL;
2173a68b90caSToby Isaac   mesh->anchorIS            = NULL;
217441e6d900SToby Isaac   mesh->createanchors       = NULL;
2175fa73a4e1SToby Isaac   mesh->computeanchormatrix = NULL;
2176d961a43aSToby Isaac   mesh->parentSection       = NULL;
2177d961a43aSToby Isaac   mesh->parents             = NULL;
2178d961a43aSToby Isaac   mesh->childIDs            = NULL;
2179d961a43aSToby Isaac   mesh->childSection        = NULL;
2180d961a43aSToby Isaac   mesh->children            = NULL;
2181d6a7ad0dSToby Isaac   mesh->referenceTree       = NULL;
2182dcbd3bf7SToby Isaac   mesh->getchildsymmetry    = NULL;
21838865f1eaSKarl Rupp   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
2184552f7358SJed Brown   mesh->vtkCellHeight       = 0;
218570034214SMatthew G. Knepley   mesh->useCone             = PETSC_FALSE;
218670034214SMatthew G. Knepley   mesh->useClosure          = PETSC_TRUE;
2187e228b242SToby Isaac   mesh->useAnchors          = PETSC_FALSE;
2188552f7358SJed Brown 
2189b29cfa1cSToby Isaac   mesh->maxProjectionHeight = 0;
2190b29cfa1cSToby Isaac 
2191552f7358SJed Brown   mesh->printSetValues = PETSC_FALSE;
2192552f7358SJed Brown   mesh->printFEM       = 0;
21936113b454SMatthew G. Knepley   mesh->printTol       = 1.0e-10;
2194552f7358SJed Brown 
2195552f7358SJed Brown   ierr = DMInitialize_Plex(dm);CHKERRQ(ierr);
2196552f7358SJed Brown   PetscFunctionReturn(0);
2197552f7358SJed Brown }
2198552f7358SJed Brown 
2199552f7358SJed Brown /*@
2200552f7358SJed Brown   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
2201552f7358SJed Brown 
2202552f7358SJed Brown   Collective on MPI_Comm
2203552f7358SJed Brown 
2204552f7358SJed Brown   Input Parameter:
2205552f7358SJed Brown . comm - The communicator for the DMPlex object
2206552f7358SJed Brown 
2207552f7358SJed Brown   Output Parameter:
2208552f7358SJed Brown . mesh  - The DMPlex object
2209552f7358SJed Brown 
2210552f7358SJed Brown   Level: beginner
2211552f7358SJed Brown 
2212552f7358SJed Brown .keywords: DMPlex, create
2213552f7358SJed Brown @*/
2214552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
2215552f7358SJed Brown {
2216552f7358SJed Brown   PetscErrorCode ierr;
2217552f7358SJed Brown 
2218552f7358SJed Brown   PetscFunctionBegin;
2219552f7358SJed Brown   PetscValidPointer(mesh,2);
2220552f7358SJed Brown   ierr = DMCreate(comm, mesh);CHKERRQ(ierr);
2221552f7358SJed Brown   ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr);
2222552f7358SJed Brown   PetscFunctionReturn(0);
2223552f7358SJed Brown }
2224552f7358SJed Brown 
2225a47d0d45SMatthew G. Knepley /*
2226a47d0d45SMatthew G. Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell, but vertex numbers are global and an SF is built for them
2227a47d0d45SMatthew G. Knepley */
2228a47d0d45SMatthew G. Knepley static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert)
2229a47d0d45SMatthew G. Knepley {
2230a47d0d45SMatthew G. Knepley   PetscSF         sfPoint;
2231a47d0d45SMatthew G. Knepley   PetscLayout     vLayout;
2232a47d0d45SMatthew G. Knepley   PetscHashI      vhash;
2233a47d0d45SMatthew G. Knepley   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
2234a47d0d45SMatthew G. Knepley   const PetscInt *vrange;
2235a47d0d45SMatthew G. Knepley   PetscInt        numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
2236a47d0d45SMatthew G. Knepley   PETSC_UNUSED PetscHashIIter ret, iter;
22379852e123SBarry Smith   PetscMPIInt     rank, size;
2238a47d0d45SMatthew G. Knepley   PetscErrorCode  ierr;
2239a47d0d45SMatthew G. Knepley 
2240a47d0d45SMatthew G. Knepley   PetscFunctionBegin;
2241a47d0d45SMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
22429852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
2243a47d0d45SMatthew G. Knepley   /* Partition vertices */
2244a47d0d45SMatthew G. Knepley   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr);
2245a47d0d45SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr);
2246a47d0d45SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
2247a47d0d45SMatthew G. Knepley   ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
2248a47d0d45SMatthew G. Knepley   ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr);
2249a47d0d45SMatthew G. Knepley   /* Count vertices and map them to procs */
2250a47d0d45SMatthew G. Knepley   PetscHashICreate(vhash);
2251a47d0d45SMatthew G. Knepley   for (c = 0; c < numCells; ++c) {
2252a47d0d45SMatthew G. Knepley     for (p = 0; p < numCorners; ++p) {
2253a47d0d45SMatthew G. Knepley       PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter);
2254a47d0d45SMatthew G. Knepley     }
2255a47d0d45SMatthew G. Knepley   }
2256a47d0d45SMatthew G. Knepley   PetscHashISize(vhash, numVerticesAdj);
2257a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr);
2258324a3438SMichael Lange   off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr);
2259a47d0d45SMatthew G. Knepley   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
22608cd1fd5cSMatthew G. Knepley   ierr = PetscSortInt(numVerticesAdj, verticesAdj);CHKERRQ(ierr);
2261a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr);
2262a47d0d45SMatthew G. Knepley   for (v = 0; v < numVerticesAdj; ++v) {
2263a47d0d45SMatthew G. Knepley     const PetscInt gv = verticesAdj[v];
2264a47d0d45SMatthew G. Knepley     PetscInt       vrank;
2265a47d0d45SMatthew G. Knepley 
22669852e123SBarry Smith     ierr = PetscFindInt(gv, size+1, vrange, &vrank);CHKERRQ(ierr);
2267a47d0d45SMatthew G. Knepley     vrank = vrank < 0 ? -(vrank+2) : vrank;
2268a47d0d45SMatthew G. Knepley     remoteVerticesAdj[v].index = gv - vrange[vrank];
2269a47d0d45SMatthew G. Knepley     remoteVerticesAdj[v].rank  = vrank;
2270a47d0d45SMatthew G. Knepley   }
2271a47d0d45SMatthew G. Knepley   /* Create cones */
2272a47d0d45SMatthew G. Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr);
2273a47d0d45SMatthew G. Knepley   for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);}
2274a47d0d45SMatthew G. Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
2275a47d0d45SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
2276a47d0d45SMatthew G. Knepley   for (c = 0; c < numCells; ++c) {
2277a47d0d45SMatthew G. Knepley     for (p = 0; p < numCorners; ++p) {
2278a47d0d45SMatthew G. Knepley       const PetscInt gv = cells[c*numCorners+p];
2279a47d0d45SMatthew G. Knepley       PetscInt       lv;
2280a47d0d45SMatthew G. Knepley 
2281a47d0d45SMatthew G. Knepley       ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr);
2282a47d0d45SMatthew G. Knepley       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
2283a47d0d45SMatthew G. Knepley       cone[p] = lv+numCells;
2284a47d0d45SMatthew G. Knepley     }
2285a47d0d45SMatthew G. Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
2286a47d0d45SMatthew G. Knepley   }
2287a47d0d45SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
2288a47d0d45SMatthew G. Knepley   /* Create SF for vertices */
2289a47d0d45SMatthew G. Knepley   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr);
2290a47d0d45SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr);
2291a47d0d45SMatthew G. Knepley   ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr);
2292a47d0d45SMatthew G. Knepley   ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr);
2293a47d0d45SMatthew G. Knepley   ierr = PetscFree(verticesAdj);CHKERRQ(ierr);
2294a47d0d45SMatthew G. Knepley   /* Build pointSF */
2295a47d0d45SMatthew G. Knepley   ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr);
2296a47d0d45SMatthew G. Knepley   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
2297a47d0d45SMatthew G. Knepley   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
22985119f1e2SToby Isaac   ierr = PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr);
22995119f1e2SToby Isaac   ierr = PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr);
2300a47d0d45SMatthew G. Knepley   for (v = 0; v < numVertices;    ++v) if (vertexOwner[v].rank < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global vertex %d on rank %d was unclaimed", v + vrange[rank], rank);
23015119f1e2SToby Isaac   ierr = PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr);
23025119f1e2SToby Isaac   ierr = PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr);
2303a47d0d45SMatthew G. Knepley   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
2304a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr);
2305a47d0d45SMatthew G. Knepley   ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr);
2306a47d0d45SMatthew G. Knepley   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
2307a47d0d45SMatthew G. Knepley     if (vertexLocal[v].rank != rank) {
2308a47d0d45SMatthew G. Knepley       localVertex[g]        = v+numCells;
2309a47d0d45SMatthew G. Knepley       remoteVertex[g].index = vertexLocal[v].index;
2310a47d0d45SMatthew G. Knepley       remoteVertex[g].rank  = vertexLocal[v].rank;
2311a47d0d45SMatthew G. Knepley       ++g;
2312a47d0d45SMatthew G. Knepley     }
2313a47d0d45SMatthew G. Knepley   }
2314a47d0d45SMatthew G. Knepley   ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr);
2315a47d0d45SMatthew G. Knepley   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
2316a47d0d45SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
2317a47d0d45SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr);
2318a47d0d45SMatthew G. Knepley   ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr);
2319a47d0d45SMatthew G. Knepley   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
2320a47d0d45SMatthew G. Knepley   PetscHashIDestroy(vhash);
2321a47d0d45SMatthew G. Knepley   /* Fill in the rest of the topology structure */
2322a47d0d45SMatthew G. Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
2323a47d0d45SMatthew G. Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
2324a47d0d45SMatthew G. Knepley   PetscFunctionReturn(0);
2325a47d0d45SMatthew G. Knepley }
2326a47d0d45SMatthew G. Knepley 
2327a47d0d45SMatthew G. Knepley /*
2328a47d0d45SMatthew G. Knepley   This takes as input the coordinates for each owned vertex
2329a47d0d45SMatthew G. Knepley */
233021016a8bSBarry Smith static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[])
2331a47d0d45SMatthew G. Knepley {
2332a47d0d45SMatthew G. Knepley   PetscSection   coordSection;
2333a47d0d45SMatthew G. Knepley   Vec            coordinates;
2334a47d0d45SMatthew G. Knepley   PetscScalar   *coords;
2335c5e08f72SMatthew G. Knepley   PetscInt       numVertices, numVerticesAdj, coordSize, v;
2336a47d0d45SMatthew G. Knepley   PetscErrorCode ierr;
2337a47d0d45SMatthew G. Knepley 
2338a47d0d45SMatthew G. Knepley   PetscFunctionBegin;
23399596c6baSMatthew G. Knepley   ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr);
2340a47d0d45SMatthew G. Knepley   ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr);
2341a47d0d45SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2342a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
2343a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
2344a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr);
2345a47d0d45SMatthew G. Knepley   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
2346a47d0d45SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
2347a47d0d45SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
2348a47d0d45SMatthew G. Knepley   }
2349a47d0d45SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
2350a47d0d45SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
2351a47d0d45SMatthew G. Knepley   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr);
2352a47d0d45SMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
2353a47d0d45SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
2354a47d0d45SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2355a47d0d45SMatthew G. Knepley   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
2356a47d0d45SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
2357a47d0d45SMatthew G. Knepley   {
2358a47d0d45SMatthew G. Knepley     MPI_Datatype coordtype;
2359a47d0d45SMatthew G. Knepley 
2360a47d0d45SMatthew G. Knepley     /* Need a temp buffer for coords if we have complex/single */
236121016a8bSBarry Smith     ierr = MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);CHKERRQ(ierr);
2362a47d0d45SMatthew G. Knepley     ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr);
236321016a8bSBarry Smith #if defined(PETSC_USE_COMPLEX)
236421016a8bSBarry Smith     {
236521016a8bSBarry Smith     PetscScalar *svertexCoords;
236621016a8bSBarry Smith     PetscInt    i;
236721016a8bSBarry Smith     ierr = PetscMalloc1(numV*spaceDim,&svertexCoords);CHKERRQ(ierr);
236821016a8bSBarry Smith     for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
236921016a8bSBarry Smith     ierr = PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr);
237021016a8bSBarry Smith     ierr = PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr);
237121016a8bSBarry Smith     ierr = PetscFree(svertexCoords);CHKERRQ(ierr);
237221016a8bSBarry Smith     }
237321016a8bSBarry Smith #else
2374a47d0d45SMatthew G. Knepley     ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr);
2375a47d0d45SMatthew G. Knepley     ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr);
237621016a8bSBarry Smith #endif
2377a47d0d45SMatthew G. Knepley     ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr);
2378a47d0d45SMatthew G. Knepley   }
2379a47d0d45SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
2380a47d0d45SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
2381a47d0d45SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
2382a47d0d45SMatthew G. Knepley   PetscFunctionReturn(0);
2383a47d0d45SMatthew G. Knepley }
2384a47d0d45SMatthew G. Knepley 
2385a47d0d45SMatthew G. Knepley /*@C
2386a47d0d45SMatthew G. Knepley   DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
2387a47d0d45SMatthew G. Knepley 
2388a47d0d45SMatthew G. Knepley   Input Parameters:
2389a47d0d45SMatthew G. Knepley + comm - The communicator
2390a47d0d45SMatthew G. Knepley . dim - The topological dimension of the mesh
2391a47d0d45SMatthew G. Knepley . numCells - The number of cells owned by this process
2392a47d0d45SMatthew G. Knepley . numVertices - The number of vertices owned by this process
2393a47d0d45SMatthew G. Knepley . numCorners - The number of vertices for each cell
2394a47d0d45SMatthew G. Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
2395a47d0d45SMatthew G. Knepley . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
2396a47d0d45SMatthew G. Knepley . spaceDim - The spatial dimension used for coordinates
2397a47d0d45SMatthew G. Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
2398a47d0d45SMatthew G. Knepley 
2399a47d0d45SMatthew G. Knepley   Output Parameter:
240018d54ad4SMichael Lange + dm - The DM
240118d54ad4SMichael Lange - vertexSF - Optional, SF describing complete vertex ownership
2402a47d0d45SMatthew G. Knepley 
2403a47d0d45SMatthew G. Knepley   Note: Two triangles sharing a face
2404a47d0d45SMatthew G. Knepley $
2405a47d0d45SMatthew G. Knepley $        2
2406a47d0d45SMatthew G. Knepley $      / | \
2407a47d0d45SMatthew G. Knepley $     /  |  \
2408a47d0d45SMatthew G. Knepley $    /   |   \
2409a47d0d45SMatthew G. Knepley $   0  0 | 1  3
2410a47d0d45SMatthew G. Knepley $    \   |   /
2411a47d0d45SMatthew G. Knepley $     \  |  /
2412a47d0d45SMatthew G. Knepley $      \ | /
2413a47d0d45SMatthew G. Knepley $        1
2414a47d0d45SMatthew G. Knepley would have input
2415a47d0d45SMatthew G. Knepley $  numCells = 2, numVertices = 4
2416a47d0d45SMatthew G. Knepley $  cells = [0 1 2  1 3 2]
2417a47d0d45SMatthew G. Knepley $
2418a47d0d45SMatthew G. Knepley which would result in the DMPlex
2419a47d0d45SMatthew G. Knepley $
2420a47d0d45SMatthew G. Knepley $        4
2421a47d0d45SMatthew G. Knepley $      / | \
2422a47d0d45SMatthew G. Knepley $     /  |  \
2423a47d0d45SMatthew G. Knepley $    /   |   \
2424a47d0d45SMatthew G. Knepley $   2  0 | 1  5
2425a47d0d45SMatthew G. Knepley $    \   |   /
2426a47d0d45SMatthew G. Knepley $     \  |  /
2427a47d0d45SMatthew G. Knepley $      \ | /
2428a47d0d45SMatthew G. Knepley $        3
2429a47d0d45SMatthew G. Knepley 
2430a47d0d45SMatthew G. Knepley   Level: beginner
2431a47d0d45SMatthew G. Knepley 
2432a47d0d45SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
2433a47d0d45SMatthew G. Knepley @*/
2434fba955ccSBarry Smith PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, DM *dm)
2435a47d0d45SMatthew G. Knepley {
2436a47d0d45SMatthew G. Knepley   PetscSF        sfVert;
2437a47d0d45SMatthew G. Knepley   PetscErrorCode ierr;
2438a47d0d45SMatthew G. Knepley 
2439a47d0d45SMatthew G. Knepley   PetscFunctionBegin;
2440a47d0d45SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
2441a47d0d45SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2442a47d0d45SMatthew G. Knepley   PetscValidLogicalCollectiveInt(*dm, dim, 2);
2443a47d0d45SMatthew G. Knepley   PetscValidLogicalCollectiveInt(*dm, spaceDim, 8);
2444a47d0d45SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
2445a47d0d45SMatthew G. Knepley   ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr);
2446a47d0d45SMatthew G. Knepley   if (interpolate) {
2447a47d0d45SMatthew G. Knepley     DM idm = NULL;
2448a47d0d45SMatthew G. Knepley 
2449a47d0d45SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
2450a47d0d45SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
2451a47d0d45SMatthew G. Knepley     *dm  = idm;
2452a47d0d45SMatthew G. Knepley   }
245321016a8bSBarry Smith   ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);CHKERRQ(ierr);
245418d54ad4SMichael Lange   if (vertexSF) *vertexSF = sfVert;
2455fba955ccSBarry Smith   else {ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);}
2456a47d0d45SMatthew G. Knepley   PetscFunctionReturn(0);
2457a47d0d45SMatthew G. Knepley }
2458a47d0d45SMatthew G. Knepley 
24599298eaa6SMatthew G Knepley /*
24609298eaa6SMatthew G Knepley   This takes as input the common mesh generator output, a list of the vertices for each cell
24619298eaa6SMatthew G Knepley */
24628f767406SMatthew G. Knepley static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
24639298eaa6SMatthew G Knepley {
24649298eaa6SMatthew G Knepley   PetscInt      *cone, c, p;
24659298eaa6SMatthew G Knepley   PetscErrorCode ierr;
24669298eaa6SMatthew G Knepley 
24679298eaa6SMatthew G Knepley   PetscFunctionBegin;
24689298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);
24699298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
24709298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);
24719298eaa6SMatthew G Knepley   }
24729298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr);
24739298eaa6SMatthew G Knepley   ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
24749298eaa6SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
24759298eaa6SMatthew G Knepley     for (p = 0; p < numCorners; ++p) {
24769298eaa6SMatthew G Knepley       cone[p] = cells[c*numCorners+p]+numCells;
24779298eaa6SMatthew G Knepley     }
24789298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);
24799298eaa6SMatthew G Knepley   }
24809298eaa6SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);
24819298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
24829298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
24839298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
24849298eaa6SMatthew G Knepley }
24859298eaa6SMatthew G Knepley 
24869298eaa6SMatthew G Knepley /*
24879298eaa6SMatthew G Knepley   This takes as input the coordinates for each vertex
24889298eaa6SMatthew G Knepley */
24898f767406SMatthew G. Knepley static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
24909298eaa6SMatthew G Knepley {
24919298eaa6SMatthew G Knepley   PetscSection   coordSection;
24929298eaa6SMatthew G Knepley   Vec            coordinates;
24938b9ced59SLisandro Dalcin   DM             cdm;
24949298eaa6SMatthew G Knepley   PetscScalar   *coords;
24958b9ced59SLisandro Dalcin   PetscInt       v, d;
24969298eaa6SMatthew G Knepley   PetscErrorCode ierr;
24979298eaa6SMatthew G Knepley 
24989298eaa6SMatthew G Knepley   PetscFunctionBegin;
24999596c6baSMatthew G. Knepley   ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr);
2500c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
25019298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
25029298eaa6SMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
25039298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
25049298eaa6SMatthew G Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
25059298eaa6SMatthew G Knepley     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
25069298eaa6SMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
25079298eaa6SMatthew G Knepley   }
25089298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
25098b9ced59SLisandro Dalcin 
25108b9ced59SLisandro Dalcin   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
25118b9ced59SLisandro Dalcin   ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr);
25124b66c01dSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr);
25139298eaa6SMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
25149298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
25159298eaa6SMatthew G Knepley   for (v = 0; v < numVertices; ++v) {
25169298eaa6SMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
25179298eaa6SMatthew G Knepley       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
25189298eaa6SMatthew G Knepley     }
25199298eaa6SMatthew G Knepley   }
25209298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
25219298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
25229298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
25239298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
25249298eaa6SMatthew G Knepley }
25259298eaa6SMatthew G Knepley 
25269298eaa6SMatthew G Knepley /*@C
25279298eaa6SMatthew G Knepley   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM
25289298eaa6SMatthew G Knepley 
25299298eaa6SMatthew G Knepley   Input Parameters:
25309298eaa6SMatthew G Knepley + comm - The communicator
25319298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh
25329298eaa6SMatthew G Knepley . numCells - The number of cells
25339298eaa6SMatthew G Knepley . numVertices - The number of vertices
25349298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell
25359298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
25369298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell
25379298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates
25389298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
25399298eaa6SMatthew G Knepley 
25409298eaa6SMatthew G Knepley   Output Parameter:
25419298eaa6SMatthew G Knepley . dm - The DM
25429298eaa6SMatthew G Knepley 
25439298eaa6SMatthew G Knepley   Note: Two triangles sharing a face
25449298eaa6SMatthew G Knepley $
25459298eaa6SMatthew G Knepley $        2
25469298eaa6SMatthew G Knepley $      / | \
25479298eaa6SMatthew G Knepley $     /  |  \
25489298eaa6SMatthew G Knepley $    /   |   \
25499298eaa6SMatthew G Knepley $   0  0 | 1  3
25509298eaa6SMatthew G Knepley $    \   |   /
25519298eaa6SMatthew G Knepley $     \  |  /
25529298eaa6SMatthew G Knepley $      \ | /
25539298eaa6SMatthew G Knepley $        1
25549298eaa6SMatthew G Knepley would have input
25559298eaa6SMatthew G Knepley $  numCells = 2, numVertices = 4
25569298eaa6SMatthew G Knepley $  cells = [0 1 2  1 3 2]
25579298eaa6SMatthew G Knepley $
25589298eaa6SMatthew G Knepley which would result in the DMPlex
25599298eaa6SMatthew G Knepley $
25609298eaa6SMatthew G Knepley $        4
25619298eaa6SMatthew G Knepley $      / | \
25629298eaa6SMatthew G Knepley $     /  |  \
25639298eaa6SMatthew G Knepley $    /   |   \
25649298eaa6SMatthew G Knepley $   2  0 | 1  5
25659298eaa6SMatthew G Knepley $    \   |   /
25669298eaa6SMatthew G Knepley $     \  |  /
25679298eaa6SMatthew G Knepley $      \ | /
25689298eaa6SMatthew G Knepley $        3
25699298eaa6SMatthew G Knepley 
25709298eaa6SMatthew G Knepley   Level: beginner
25719298eaa6SMatthew G Knepley 
2572939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
25739298eaa6SMatthew G Knepley @*/
25749298eaa6SMatthew 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)
25759298eaa6SMatthew G Knepley {
25769298eaa6SMatthew G Knepley   PetscErrorCode ierr;
25779298eaa6SMatthew G Knepley 
25789298eaa6SMatthew G Knepley   PetscFunctionBegin;
25799298eaa6SMatthew G Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
25809298eaa6SMatthew G Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2581c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
25829298eaa6SMatthew G Knepley   ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr);
25839298eaa6SMatthew G Knepley   if (interpolate) {
2584e56d480eSMatthew G. Knepley     DM idm = NULL;
25859298eaa6SMatthew G Knepley 
25869298eaa6SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
25879298eaa6SMatthew G Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
25889298eaa6SMatthew G Knepley     *dm  = idm;
25899298eaa6SMatthew G Knepley   }
25909298eaa6SMatthew G Knepley   ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr);
25919298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
25929298eaa6SMatthew G Knepley }
25939298eaa6SMatthew G Knepley 
2594939f6067SMatthew G. Knepley /*@
2595939f6067SMatthew 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
2596939f6067SMatthew G. Knepley 
2597939f6067SMatthew G. Knepley   Input Parameters:
2598c73cfb54SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
2599939f6067SMatthew G. Knepley . depth - The depth of the DAG
2600939f6067SMatthew G. Knepley . numPoints - The number of points at each depth
2601939f6067SMatthew G. Knepley . coneSize - The cone size of each point
2602939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
2603939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point
2604939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex
2605939f6067SMatthew G. Knepley 
2606939f6067SMatthew G. Knepley   Output Parameter:
2607939f6067SMatthew G. Knepley . dm - The DM
2608939f6067SMatthew G. Knepley 
2609939f6067SMatthew G. Knepley   Note: Two triangles sharing a face would have input
2610939f6067SMatthew G. Knepley $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
2611939f6067SMatthew G. Knepley $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
2612939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
2613939f6067SMatthew G. Knepley $
2614939f6067SMatthew G. Knepley which would result in the DMPlex
2615939f6067SMatthew G. Knepley $
2616939f6067SMatthew G. Knepley $        4
2617939f6067SMatthew G. Knepley $      / | \
2618939f6067SMatthew G. Knepley $     /  |  \
2619939f6067SMatthew G. Knepley $    /   |   \
2620939f6067SMatthew G. Knepley $   2  0 | 1  5
2621939f6067SMatthew G. Knepley $    \   |   /
2622939f6067SMatthew G. Knepley $     \  |  /
2623939f6067SMatthew G. Knepley $      \ | /
2624939f6067SMatthew G. Knepley $        3
2625939f6067SMatthew G. Knepley $
2626939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()
2627939f6067SMatthew G. Knepley 
2628939f6067SMatthew G. Knepley   Level: advanced
2629939f6067SMatthew G. Knepley 
2630939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
2631939f6067SMatthew G. Knepley @*/
26329298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
26339298eaa6SMatthew G Knepley {
26349298eaa6SMatthew G Knepley   Vec            coordinates;
26359298eaa6SMatthew G Knepley   PetscSection   coordSection;
26369298eaa6SMatthew G Knepley   PetscScalar    *coords;
2637811e8653SToby Isaac   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;
26389298eaa6SMatthew G Knepley   PetscErrorCode ierr;
26399298eaa6SMatthew G Knepley 
26409298eaa6SMatthew G Knepley   PetscFunctionBegin;
2641c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2642811e8653SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
2643811e8653SToby Isaac   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
26449298eaa6SMatthew G Knepley   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
26459298eaa6SMatthew G Knepley   ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr);
26469298eaa6SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
26479298eaa6SMatthew G Knepley     ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr);
264897e052ccSToby Isaac     if (firstVertex < 0 && !coneSize[p - pStart]) {
264997e052ccSToby Isaac       firstVertex = p - pStart;
26509298eaa6SMatthew G Knepley     }
265197e052ccSToby Isaac   }
265297e052ccSToby Isaac   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
26539298eaa6SMatthew G Knepley   ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */
26549298eaa6SMatthew G Knepley   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
26559298eaa6SMatthew G Knepley     ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr);
26569298eaa6SMatthew G Knepley     ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr);
26579298eaa6SMatthew G Knepley   }
26589298eaa6SMatthew G Knepley   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
26599298eaa6SMatthew G Knepley   ierr = DMPlexStratify(dm);CHKERRQ(ierr);
26609298eaa6SMatthew G Knepley   /* Build coordinates */
2661c2166f76SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
26629298eaa6SMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
2663811e8653SToby Isaac   ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr);
26649298eaa6SMatthew G Knepley   ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr);
26659298eaa6SMatthew G Knepley   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
2666811e8653SToby Isaac     ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr);
2667811e8653SToby Isaac     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr);
26689298eaa6SMatthew G Knepley   }
26699298eaa6SMatthew G Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
26709298eaa6SMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
26718b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
26726f8cbbeeSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
26739298eaa6SMatthew G Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
26748b9ced59SLisandro Dalcin   ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr);
26752eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
26769298eaa6SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
26779298eaa6SMatthew G Knepley   for (v = 0; v < numPoints[0]; ++v) {
26789298eaa6SMatthew G Knepley     PetscInt off;
26799298eaa6SMatthew G Knepley 
26809298eaa6SMatthew G Knepley     ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr);
2681811e8653SToby Isaac     for (d = 0; d < dimEmbed; ++d) {
2682811e8653SToby Isaac       coords[off+d] = vertexCoords[v*dimEmbed+d];
26839298eaa6SMatthew G Knepley     }
26849298eaa6SMatthew G Knepley   }
26859298eaa6SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
26869298eaa6SMatthew G Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
26879298eaa6SMatthew G Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
26889298eaa6SMatthew G Knepley   PetscFunctionReturn(0);
26899298eaa6SMatthew G Knepley }
26908415267dSToby Isaac 
2691ca522641SMatthew G. Knepley /*@C
2692ca522641SMatthew G. Knepley   DMPlexCreateFromFile - This takes a filename and produces a DM
2693ca522641SMatthew G. Knepley 
2694ca522641SMatthew G. Knepley   Input Parameters:
2695ca522641SMatthew G. Knepley + comm - The communicator
2696ca522641SMatthew G. Knepley . filename - A file name
2697ca522641SMatthew G. Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces)
2698ca522641SMatthew G. Knepley 
2699ca522641SMatthew G. Knepley   Output Parameter:
2700ca522641SMatthew G. Knepley . dm - The DM
2701ca522641SMatthew G. Knepley 
2702ca522641SMatthew G. Knepley   Level: beginner
2703ca522641SMatthew G. Knepley 
2704ca522641SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
2705ca522641SMatthew G. Knepley @*/
2706ca522641SMatthew G. Knepley PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
2707ca522641SMatthew G. Knepley {
2708ca522641SMatthew G. Knepley   const char    *extGmsh    = ".msh";
2709ca522641SMatthew G. Knepley   const char    *extCGNS    = ".cgns";
2710ca522641SMatthew G. Knepley   const char    *extExodus  = ".exo";
271190c68965SMatthew G. Knepley   const char    *extGenesis = ".gen";
27122f0bd6dcSMichael Lange   const char    *extFluent  = ".cas";
2713cc2f8f65SMatthew G. Knepley   const char    *extHDF5    = ".h5";
2714707dd687SMichael Lange   const char    *extMed     = ".med";
2715f2801cd6SMatthew G. Knepley   const char    *extPLY     = ".ply";
2716ca522641SMatthew G. Knepley   size_t         len;
271790c68965SMatthew G. Knepley   PetscBool      isGmsh, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY;
2718ca522641SMatthew G. Knepley   PetscMPIInt    rank;
2719ca522641SMatthew G. Knepley   PetscErrorCode ierr;
2720ca522641SMatthew G. Knepley 
2721ca522641SMatthew G. Knepley   PetscFunctionBegin;
2722ca522641SMatthew G. Knepley   PetscValidPointer(filename, 2);
2723ca522641SMatthew G. Knepley   PetscValidPointer(dm, 4);
2724ca522641SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2725ca522641SMatthew G. Knepley   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
2726ca522641SMatthew G. Knepley   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
2727ca522641SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,    4, &isGmsh);CHKERRQ(ierr);
2728ca522641SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,    5, &isCGNS);CHKERRQ(ierr);
2729ca522641SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus,  4, &isExodus);CHKERRQ(ierr);
273090c68965SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);CHKERRQ(ierr);
27312f0bd6dcSMichael Lange   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent,  4, &isFluent);CHKERRQ(ierr);
2732cc2f8f65SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,    3, &isHDF5);CHKERRQ(ierr);
2733707dd687SMichael Lange   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,     4, &isMed);CHKERRQ(ierr);
2734f2801cd6SMatthew G. Knepley   ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY,     4, &isPLY);CHKERRQ(ierr);
2735ca522641SMatthew G. Knepley   if (isGmsh) {
27367d282ae0SMichael Lange     ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2737ca522641SMatthew G. Knepley   } else if (isCGNS) {
2738ca522641SMatthew G. Knepley     ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
273990c68965SMatthew G. Knepley   } else if (isExodus || isGenesis) {
2740ca522641SMatthew G. Knepley     ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
27412f0bd6dcSMichael Lange   } else if (isFluent) {
27422f0bd6dcSMichael Lange     ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2743cc2f8f65SMatthew G. Knepley   } else if (isHDF5) {
2744cc2f8f65SMatthew G. Knepley     PetscViewer viewer;
2745cc2f8f65SMatthew G. Knepley 
2746cc2f8f65SMatthew G. Knepley     ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
2747cc2f8f65SMatthew G. Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr);
2748cc2f8f65SMatthew G. Knepley     ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
2749cc2f8f65SMatthew G. Knepley     ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
2750cc2f8f65SMatthew G. Knepley     ierr = DMCreate(comm, dm);CHKERRQ(ierr);
2751cc2f8f65SMatthew G. Knepley     ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2752cc2f8f65SMatthew G. Knepley     ierr = DMLoad(*dm, viewer);CHKERRQ(ierr);
2753cc2f8f65SMatthew G. Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2754707dd687SMichael Lange   } else if (isMed) {
2755707dd687SMichael Lange     ierr = DMPlexCreateMedFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2756f2801cd6SMatthew G. Knepley   } else if (isPLY) {
2757f2801cd6SMatthew G. Knepley     ierr = DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
2758ca522641SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
2759ca522641SMatthew G. Knepley   PetscFunctionReturn(0);
2760ca522641SMatthew G. Knepley }
2761ca522641SMatthew G. Knepley 
27628415267dSToby Isaac /*@
27638415267dSToby Isaac   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
27648415267dSToby Isaac 
27658415267dSToby Isaac   Collective on comm
27668415267dSToby Isaac 
27678415267dSToby Isaac   Input Parameters:
27688415267dSToby Isaac + comm    - The communicator
27698415267dSToby Isaac . dim     - The spatial dimension
27708415267dSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
27718415267dSToby Isaac 
27728415267dSToby Isaac   Output Parameter:
27738415267dSToby Isaac . refdm - The reference cell
27748415267dSToby Isaac 
27758415267dSToby Isaac   Level: intermediate
27768415267dSToby Isaac 
27778415267dSToby Isaac .keywords: reference cell
27788415267dSToby Isaac .seealso:
27798415267dSToby Isaac @*/
27808415267dSToby Isaac PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
27818415267dSToby Isaac {
27828415267dSToby Isaac   DM             rdm;
2783cdf740e0SMatthew G. Knepley   Vec            coords;
27848415267dSToby Isaac   PetscErrorCode ierr;
27858415267dSToby Isaac 
27868415267dSToby Isaac   PetscFunctionBeginUser;
27878415267dSToby Isaac   ierr = DMCreate(comm, &rdm);CHKERRQ(ierr);
27888415267dSToby Isaac   ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr);
27898415267dSToby Isaac   ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr);
27908415267dSToby Isaac   switch (dim) {
27918415267dSToby Isaac   case 0:
27928415267dSToby Isaac   {
27938415267dSToby Isaac     PetscInt    numPoints[1]        = {1};
27948415267dSToby Isaac     PetscInt    coneSize[1]         = {0};
27958415267dSToby Isaac     PetscInt    cones[1]            = {0};
27968415267dSToby Isaac     PetscInt    coneOrientations[1] = {0};
27978415267dSToby Isaac     PetscScalar vertexCoords[1]     = {0.0};
27988415267dSToby Isaac 
27998415267dSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
28008415267dSToby Isaac   }
28018415267dSToby Isaac   break;
28028415267dSToby Isaac   case 1:
28038415267dSToby Isaac   {
28048415267dSToby Isaac     PetscInt    numPoints[2]        = {2, 1};
28058415267dSToby Isaac     PetscInt    coneSize[3]         = {2, 0, 0};
28068415267dSToby Isaac     PetscInt    cones[2]            = {1, 2};
28078415267dSToby Isaac     PetscInt    coneOrientations[2] = {0, 0};
28088415267dSToby Isaac     PetscScalar vertexCoords[2]     = {-1.0,  1.0};
28098415267dSToby Isaac 
28108415267dSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
28118415267dSToby Isaac   }
28128415267dSToby Isaac   break;
28138415267dSToby Isaac   case 2:
28148415267dSToby Isaac     if (simplex) {
28158415267dSToby Isaac       PetscInt    numPoints[2]        = {3, 1};
28168415267dSToby Isaac       PetscInt    coneSize[4]         = {3, 0, 0, 0};
28178415267dSToby Isaac       PetscInt    cones[3]            = {1, 2, 3};
28188415267dSToby Isaac       PetscInt    coneOrientations[3] = {0, 0, 0};
28198415267dSToby Isaac       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};
28208415267dSToby Isaac 
28218415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
28228415267dSToby Isaac     } else {
28238415267dSToby Isaac       PetscInt    numPoints[2]        = {4, 1};
28248415267dSToby Isaac       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
28258415267dSToby Isaac       PetscInt    cones[4]            = {1, 2, 3, 4};
28268415267dSToby Isaac       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
28278415267dSToby Isaac       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};
28288415267dSToby Isaac 
28298415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
28308415267dSToby Isaac     }
28318415267dSToby Isaac   break;
28328415267dSToby Isaac   case 3:
28338415267dSToby Isaac     if (simplex) {
28348415267dSToby Isaac       PetscInt    numPoints[2]        = {4, 1};
28358415267dSToby Isaac       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
28368415267dSToby Isaac       PetscInt    cones[4]            = {1, 3, 2, 4};
28378415267dSToby Isaac       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
28388415267dSToby Isaac       PetscScalar vertexCoords[12]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,  -1.0, -1.0, 1.0};
28398415267dSToby Isaac 
28408415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
28418415267dSToby Isaac     } else {
28428415267dSToby Isaac       PetscInt    numPoints[2]        = {8, 1};
28438415267dSToby Isaac       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
28448415267dSToby Isaac       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
28458415267dSToby Isaac       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
28468415267dSToby Isaac       PetscScalar vertexCoords[24]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  1.0, 1.0, -1.0,  -1.0, 1.0, -1.0,
28478415267dSToby Isaac                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};
28488415267dSToby Isaac 
28498415267dSToby Isaac       ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
28508415267dSToby Isaac     }
28518415267dSToby Isaac   break;
28528415267dSToby Isaac   default:
28538415267dSToby Isaac     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
28548415267dSToby Isaac   }
28558415267dSToby Isaac   *refdm = NULL;
28568415267dSToby Isaac   ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr);
2857cdf740e0SMatthew G. Knepley   if (rdm->coordinateDM) {
2858cdf740e0SMatthew G. Knepley     DM           ncdm;
2859cdf740e0SMatthew G. Knepley     PetscSection cs;
2860cdf740e0SMatthew G. Knepley     PetscInt     pEnd = -1;
2861cdf740e0SMatthew G. Knepley 
2862cdf740e0SMatthew G. Knepley     ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr);
2863cdf740e0SMatthew G. Knepley     if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);}
2864cdf740e0SMatthew G. Knepley     if (pEnd >= 0) {
2865cdf740e0SMatthew G. Knepley       ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr);
2866cdf740e0SMatthew G. Knepley       ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr);
2867a61e840bSMatthew G. Knepley       ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr);
2868cdf740e0SMatthew G. Knepley       ierr = DMDestroy(&ncdm);CHKERRQ(ierr);
2869cdf740e0SMatthew G. Knepley     }
2870cdf740e0SMatthew G. Knepley   }
2871cdf740e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr);
2872cdf740e0SMatthew G. Knepley   if (coords) {
2873cdf740e0SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr);
2874cdf740e0SMatthew G. Knepley   } else {
2875cdf740e0SMatthew G. Knepley     ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr);
2876cdf740e0SMatthew G. Knepley     if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);}
2877cdf740e0SMatthew G. Knepley   }
28788415267dSToby Isaac   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
28798415267dSToby Isaac   PetscFunctionReturn(0);
28808415267dSToby Isaac }
2881