Actual source code: plexcreate.c

  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>
  3: #include <petsc/private/hashseti.h>
  4: #include <petscsf.h>

  6: PetscLogEvent DMPLEX_CreateFromFile, DMPLEX_BuildFromCellList, DMPLEX_BuildCoordinatesFromCellList;

  8: /*@
  9:   DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.

 11:   Collective

 13:   Input Parameters:
 14: + comm - The communicator for the DM object
 15: . dim - The spatial dimension
 16: . simplex - Flag for simplicial cells, otherwise they are tensor product cells
 17: . interpolate - Flag to create intermediate mesh pieces (edges, faces)
 18: - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell

 20:   Output Parameter:
 21: . dm - The DM object

 23:   Level: beginner

 25: .seealso: DMSetType(), DMCreate()
 26: @*/
 27: PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscReal refinementLimit, DM *newdm)
 28: {
 29:   DM             dm;
 30:   PetscMPIInt    rank;

 34:   DMCreate(comm, &dm);
 35:   DMSetType(dm, DMPLEX);
 36:   DMSetDimension(dm, dim);
 37:   MPI_Comm_rank(comm, &rank);
 38:   switch (dim) {
 39:   case 2:
 40:     if (simplex) {PetscObjectSetName((PetscObject) dm, "triangular");}
 41:     else         {PetscObjectSetName((PetscObject) dm, "quadrilateral");}
 42:     break;
 43:   case 3:
 44:     if (simplex) {PetscObjectSetName((PetscObject) dm, "tetrahedral");}
 45:     else         {PetscObjectSetName((PetscObject) dm, "hexahedral");}
 46:     break;
 47:   default:
 48:     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %D", dim);
 49:   }
 50:   if (rank) {
 51:     PetscInt numPoints[2] = {0, 0};
 52:     DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);
 53:   } else {
 54:     switch (dim) {
 55:     case 2:
 56:       if (simplex) {
 57:         PetscInt    numPoints[2]        = {4, 2};
 58:         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
 59:         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
 60:         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
 61:         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};

 63:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 64:       } else {
 65:         PetscInt    numPoints[2]        = {6, 2};
 66:         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
 67:         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
 68:         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
 69:         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};

 71:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 72:       }
 73:       break;
 74:     case 3:
 75:       if (simplex) {
 76:         PetscInt    numPoints[2]        = {5, 2};
 77:         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
 78:         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
 79:         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
 80:         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};

 82:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 83:       } else {
 84:         PetscInt    numPoints[2]         = {12, 2};
 85:         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
 86:         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
 87:         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
 88:         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,
 89:                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
 90:                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};

 92:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 93:       }
 94:       break;
 95:     default:
 96:       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %D", dim);
 97:     }
 98:   }
 99:   *newdm = dm;
100:   if (refinementLimit > 0.0) {
101:     DM rdm;
102:     const char *name;

104:     DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);
105:     DMPlexSetRefinementLimit(*newdm, refinementLimit);
106:     DMRefine(*newdm, comm, &rdm);
107:     PetscObjectGetName((PetscObject) *newdm, &name);
108:     PetscObjectSetName((PetscObject)    rdm,  name);
109:     DMDestroy(newdm);
110:     *newdm = rdm;
111:   }
112:   if (interpolate) {
113:     DM idm;

115:     DMPlexInterpolate(*newdm, &idm);
116:     DMDestroy(newdm);
117:     *newdm = idm;
118:   }
119:   return(0);
120: }

122: /*@
123:   DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice.

125:   Collective

127:   Input Parameters:
128: + comm  - The communicator for the DM object
129: . lower - The lower left corner coordinates
130: . upper - The upper right corner coordinates
131: - edges - The number of cells in each direction

133:   Output Parameter:
134: . dm  - The DM object

136:   Note: Here is the numbering returned for 2 cells in each direction:
137: $ 18--5-17--4--16
138: $  |     |     |
139: $  6    10     3
140: $  |     |     |
141: $ 19-11-20--9--15
142: $  |     |     |
143: $  7     8     2
144: $  |     |     |
145: $ 12--0-13--1--14

147:   Level: beginner

149: .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
150: @*/
151: PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
152: {
153:   const PetscInt numVertices    = (edges[0]+1)*(edges[1]+1);
154:   const PetscInt numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
155:   PetscInt       markerTop      = 1;
156:   PetscInt       markerBottom   = 1;
157:   PetscInt       markerRight    = 1;
158:   PetscInt       markerLeft     = 1;
159:   PetscBool      markerSeparate = PETSC_FALSE;
160:   Vec            coordinates;
161:   PetscSection   coordSection;
162:   PetscScalar    *coords;
163:   PetscInt       coordSize;
164:   PetscMPIInt    rank;
165:   PetscInt       v, vx, vy;

169:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
170:   if (markerSeparate) {
171:     markerTop    = 3;
172:     markerBottom = 1;
173:     markerRight  = 2;
174:     markerLeft   = 4;
175:   }
176:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
177:   if (!rank) {
178:     PetscInt e, ex, ey;

180:     DMPlexSetChart(dm, 0, numEdges+numVertices);
181:     for (e = 0; e < numEdges; ++e) {
182:       DMPlexSetConeSize(dm, e, 2);
183:     }
184:     DMSetUp(dm); /* Allocate space for cones */
185:     for (vx = 0; vx <= edges[0]; vx++) {
186:       for (ey = 0; ey < edges[1]; ey++) {
187:         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
188:         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
189:         PetscInt cone[2];

191:         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
192:         DMPlexSetCone(dm, edge, cone);
193:         if (vx == edges[0]) {
194:           DMSetLabelValue(dm, "marker", edge,    markerRight);
195:           DMSetLabelValue(dm, "marker", cone[0], markerRight);
196:           if (ey == edges[1]-1) {
197:             DMSetLabelValue(dm, "marker", cone[1], markerRight);
198:             DMSetLabelValue(dm, "Face Sets", cone[1], markerRight);
199:           }
200:         } else if (vx == 0) {
201:           DMSetLabelValue(dm, "marker", edge,    markerLeft);
202:           DMSetLabelValue(dm, "marker", cone[0], markerLeft);
203:           if (ey == edges[1]-1) {
204:             DMSetLabelValue(dm, "marker", cone[1], markerLeft);
205:             DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft);
206:           }
207:         }
208:       }
209:     }
210:     for (vy = 0; vy <= edges[1]; vy++) {
211:       for (ex = 0; ex < edges[0]; ex++) {
212:         PetscInt edge   = vy*edges[0]     + ex;
213:         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
214:         PetscInt cone[2];

216:         cone[0] = vertex; cone[1] = vertex+1;
217:         DMPlexSetCone(dm, edge, cone);
218:         if (vy == edges[1]) {
219:           DMSetLabelValue(dm, "marker", edge,    markerTop);
220:           DMSetLabelValue(dm, "marker", cone[0], markerTop);
221:           if (ex == edges[0]-1) {
222:             DMSetLabelValue(dm, "marker", cone[1], markerTop);
223:             DMSetLabelValue(dm, "Face Sets", cone[1], markerTop);
224:           }
225:         } else if (vy == 0) {
226:           DMSetLabelValue(dm, "marker", edge,    markerBottom);
227:           DMSetLabelValue(dm, "marker", cone[0], markerBottom);
228:           if (ex == edges[0]-1) {
229:             DMSetLabelValue(dm, "marker", cone[1], markerBottom);
230:             DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom);
231:           }
232:         }
233:       }
234:     }
235:   }
236:   DMPlexSymmetrize(dm);
237:   DMPlexStratify(dm);
238:   /* Build coordinates */
239:   DMSetCoordinateDim(dm, 2);
240:   DMGetCoordinateSection(dm, &coordSection);
241:   PetscSectionSetNumFields(coordSection, 1);
242:   PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);
243:   PetscSectionSetFieldComponents(coordSection, 0, 2);
244:   for (v = numEdges; v < numEdges+numVertices; ++v) {
245:     PetscSectionSetDof(coordSection, v, 2);
246:     PetscSectionSetFieldDof(coordSection, v, 0, 2);
247:   }
248:   PetscSectionSetUp(coordSection);
249:   PetscSectionGetStorageSize(coordSection, &coordSize);
250:   VecCreate(PETSC_COMM_SELF, &coordinates);
251:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
252:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
253:   VecSetBlockSize(coordinates, 2);
254:   VecSetType(coordinates,VECSTANDARD);
255:   VecGetArray(coordinates, &coords);
256:   for (vy = 0; vy <= edges[1]; ++vy) {
257:     for (vx = 0; vx <= edges[0]; ++vx) {
258:       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
259:       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
260:     }
261:   }
262:   VecRestoreArray(coordinates, &coords);
263:   DMSetCoordinatesLocal(dm, coordinates);
264:   VecDestroy(&coordinates);
265:   return(0);
266: }

268: /*@
269:   DMPlexCreateCubeBoundary - Creates a 2D mesh that is the boundary of a cubic lattice.

271:   Collective

273:   Input Parameters:
274: + comm  - The communicator for the DM object
275: . lower - The lower left front corner coordinates
276: . upper - The upper right back corner coordinates
277: - faces - The number of faces in each direction (the same as the number of cells)

279:   Output Parameter:
280: . dm  - The DM object

282:   Level: beginner

284: .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
285: @*/
286: PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
287: {
288:   PetscInt       vertices[3], numVertices;
289:   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
290:   Vec            coordinates;
291:   PetscSection   coordSection;
292:   PetscScalar    *coords;
293:   PetscInt       coordSize;
294:   PetscMPIInt    rank;
295:   PetscInt       v, vx, vy, vz;
296:   PetscInt       voffset, iface=0, cone[4];

300:   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");
301:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
302:   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
303:   numVertices = vertices[0]*vertices[1]*vertices[2];
304:   if (!rank) {
305:     PetscInt f;

307:     DMPlexSetChart(dm, 0, numFaces+numVertices);
308:     for (f = 0; f < numFaces; ++f) {
309:       DMPlexSetConeSize(dm, f, 4);
310:     }
311:     DMSetUp(dm); /* Allocate space for cones */

313:     /* Side 0 (Top) */
314:     for (vy = 0; vy < faces[1]; vy++) {
315:       for (vx = 0; vx < faces[0]; vx++) {
316:         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
317:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
318:         DMPlexSetCone(dm, iface, cone);
319:         DMSetLabelValue(dm, "marker", iface, 1);
320:         DMSetLabelValue(dm, "marker", voffset+0, 1);
321:         DMSetLabelValue(dm, "marker", voffset+1, 1);
322:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
323:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+1, 1);
324:         iface++;
325:       }
326:     }

328:     /* Side 1 (Bottom) */
329:     for (vy = 0; vy < faces[1]; vy++) {
330:       for (vx = 0; vx < faces[0]; vx++) {
331:         voffset = numFaces + vy*(faces[0]+1) + vx;
332:         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
333:         DMPlexSetCone(dm, iface, cone);
334:         DMSetLabelValue(dm, "marker", iface, 1);
335:         DMSetLabelValue(dm, "marker", voffset+0, 1);
336:         DMSetLabelValue(dm, "marker", voffset+1, 1);
337:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
338:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+1, 1);
339:         iface++;
340:       }
341:     }

343:     /* Side 2 (Front) */
344:     for (vz = 0; vz < faces[2]; vz++) {
345:       for (vx = 0; vx < faces[0]; vx++) {
346:         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
347:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
348:         DMPlexSetCone(dm, iface, cone);
349:         DMSetLabelValue(dm, "marker", iface, 1);
350:         DMSetLabelValue(dm, "marker", voffset+0, 1);
351:         DMSetLabelValue(dm, "marker", voffset+1, 1);
352:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
353:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+1, 1);
354:         iface++;
355:       }
356:     }

358:     /* Side 3 (Back) */
359:     for (vz = 0; vz < faces[2]; vz++) {
360:       for (vx = 0; vx < faces[0]; vx++) {
361:         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
362:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
363:         cone[2] = voffset+1; cone[3] = voffset;
364:         DMPlexSetCone(dm, iface, cone);
365:         DMSetLabelValue(dm, "marker", iface, 1);
366:         DMSetLabelValue(dm, "marker", voffset+0, 1);
367:         DMSetLabelValue(dm, "marker", voffset+1, 1);
368:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
369:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+1, 1);
370:         iface++;
371:       }
372:     }

374:     /* Side 4 (Left) */
375:     for (vz = 0; vz < faces[2]; vz++) {
376:       for (vy = 0; vy < faces[1]; vy++) {
377:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
378:         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
379:         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
380:         DMPlexSetCone(dm, iface, cone);
381:         DMSetLabelValue(dm, "marker", iface, 1);
382:         DMSetLabelValue(dm, "marker", voffset+0, 1);
383:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
384:         DMSetLabelValue(dm, "marker", voffset+vertices[1]+0, 1);
385:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+vertices[0], 1);
386:         iface++;
387:       }
388:     }

390:     /* Side 5 (Right) */
391:     for (vz = 0; vz < faces[2]; vz++) {
392:       for (vy = 0; vy < faces[1]; vy++) {
393:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + faces[0];
394:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
395:         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
396:         DMPlexSetCone(dm, iface, cone);
397:         DMSetLabelValue(dm, "marker", iface, 1);
398:         DMSetLabelValue(dm, "marker", voffset+0, 1);
399:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
400:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
401:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+vertices[0], 1);
402:         iface++;
403:       }
404:     }
405:   }
406:   DMPlexSymmetrize(dm);
407:   DMPlexStratify(dm);
408:   /* Build coordinates */
409:   DMSetCoordinateDim(dm, 3);
410:   DMGetCoordinateSection(dm, &coordSection);
411:   PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);
412:   for (v = numFaces; v < numFaces+numVertices; ++v) {
413:     PetscSectionSetDof(coordSection, v, 3);
414:   }
415:   PetscSectionSetUp(coordSection);
416:   PetscSectionGetStorageSize(coordSection, &coordSize);
417:   VecCreate(PETSC_COMM_SELF, &coordinates);
418:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
419:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
420:   VecSetBlockSize(coordinates, 3);
421:   VecSetType(coordinates,VECSTANDARD);
422:   VecGetArray(coordinates, &coords);
423:   for (vz = 0; vz <= faces[2]; ++vz) {
424:     for (vy = 0; vy <= faces[1]; ++vy) {
425:       for (vx = 0; vx <= faces[0]; ++vx) {
426:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
427:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
428:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
429:       }
430:     }
431:   }
432:   VecRestoreArray(coordinates, &coords);
433:   DMSetCoordinatesLocal(dm, coordinates);
434:   VecDestroy(&coordinates);
435:   return(0);
436: }

438: static PetscErrorCode DMPlexCreateLineMesh_Internal(MPI_Comm comm,PetscInt segments,PetscReal lower,PetscReal upper,DMBoundaryType bd,DM *dm)
439: {
440:   PetscInt       i,fStart,fEnd,numCells = 0,numVerts = 0;
441:   PetscInt       numPoints[2],*coneSize,*cones,*coneOrientations;
442:   PetscScalar    *vertexCoords;
443:   PetscReal      L,maxCell;
444:   PetscBool      markerSeparate = PETSC_FALSE;
445:   PetscInt       markerLeft  = 1, faceMarkerLeft  = 1;
446:   PetscInt       markerRight = 1, faceMarkerRight = 2;
447:   PetscBool      wrap = (bd == DM_BOUNDARY_PERIODIC || bd == DM_BOUNDARY_TWIST) ? PETSC_TRUE : PETSC_FALSE;
448:   PetscMPIInt    rank;


454:   DMCreate(comm,dm);
455:   DMSetType(*dm,DMPLEX);
456:   DMSetDimension(*dm,1);
457:   DMCreateLabel(*dm,"marker");
458:   DMCreateLabel(*dm,"Face Sets");

460:   MPI_Comm_rank(comm,&rank);
461:   if (!rank) numCells = segments;
462:   if (!rank) numVerts = segments + (wrap ? 0 : 1);

464:   numPoints[0] = numVerts ; numPoints[1] = numCells;
465:   PetscMalloc4(numCells+numVerts,&coneSize,numCells*2,&cones,numCells+numVerts,&coneOrientations,numVerts,&vertexCoords);
466:   PetscArrayzero(coneOrientations,numCells+numVerts);
467:   for (i = 0; i < numCells; ++i) { coneSize[i] = 2; }
468:   for (i = 0; i < numVerts; ++i) { coneSize[numCells+i] = 0; }
469:   for (i = 0; i < numCells; ++i) { cones[2*i] = numCells + i%numVerts; cones[2*i+1] = numCells + (i+1)%numVerts; }
470:   for (i = 0; i < numVerts; ++i) { vertexCoords[i] = lower + (upper-lower)*((PetscReal)i/(PetscReal)numCells); }
471:   DMPlexCreateFromDAG(*dm,1,numPoints,coneSize,cones,coneOrientations,vertexCoords);
472:   PetscFree4(coneSize,cones,coneOrientations,vertexCoords);

474:   PetscOptionsGetBool(((PetscObject)*dm)->options,((PetscObject)*dm)->prefix,"-dm_plex_separate_marker",&markerSeparate,NULL);
475:   if (markerSeparate) { markerLeft = faceMarkerLeft; markerRight = faceMarkerRight;}
476:   if (!wrap && !rank) {
477:     DMPlexGetHeightStratum(*dm,1,&fStart,&fEnd);
478:     DMSetLabelValue(*dm,"marker",fStart,markerLeft);
479:     DMSetLabelValue(*dm,"marker",fEnd-1,markerRight);
480:     DMSetLabelValue(*dm,"Face Sets",fStart,faceMarkerLeft);
481:     DMSetLabelValue(*dm,"Face Sets",fEnd-1,faceMarkerRight);
482:   }
483:   if (wrap) {
484:     L       = upper - lower;
485:     maxCell = (PetscReal)1.1*(L/(PetscReal)PetscMax(1,segments));
486:     DMSetPeriodicity(*dm,PETSC_TRUE,&maxCell,&L,&bd);
487:   }
488:   DMPlexSetRefinementUniform(*dm, PETSC_TRUE);
489:   return(0);
490: }

492: 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)
493: {
494:   DM             boundary;
495:   PetscInt       i;

500:   for (i = 0; i < dim; ++i) if (periodicity[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity is not supported for simplex meshes");
501:   DMCreate(comm, &boundary);
503:   DMSetType(boundary, DMPLEX);
504:   DMSetDimension(boundary, dim-1);
505:   DMSetCoordinateDim(boundary, dim);
506:   switch (dim) {
507:   case 2: DMPlexCreateSquareBoundary(boundary, lower, upper, faces);break;
508:   case 3: DMPlexCreateCubeBoundary(boundary, lower, upper, faces);break;
509:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %D", dim);
510:   }
511:   DMPlexGenerate(boundary, NULL, interpolate, dm);
512:   DMDestroy(&boundary);
513:   return(0);
514: }

516: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
517: {
518:   DMLabel        cutLabel = NULL;
519:   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
520:   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
521:   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
522:   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
523:   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
524:   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
525:   PetscInt       dim;
526:   PetscBool      markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
527:   PetscMPIInt    rank;

531:   DMGetDimension(dm,&dim);
532:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
533:   DMCreateLabel(dm,"marker");
534:   DMCreateLabel(dm,"Face Sets");
535:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);
536:   if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ||
537:       bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ||
538:       bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {

540:     if (cutMarker) {DMCreateLabel(dm, "periodic_cut"); DMGetLabel(dm, "periodic_cut", &cutLabel);}
541:   }
542:   switch (dim) {
543:   case 2:
544:     faceMarkerTop    = 3;
545:     faceMarkerBottom = 1;
546:     faceMarkerRight  = 2;
547:     faceMarkerLeft   = 4;
548:     break;
549:   case 3:
550:     faceMarkerBottom = 1;
551:     faceMarkerTop    = 2;
552:     faceMarkerFront  = 3;
553:     faceMarkerBack   = 4;
554:     faceMarkerRight  = 5;
555:     faceMarkerLeft   = 6;
556:     break;
557:   default:
558:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D not supported",dim);
559:   }
560:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
561:   if (markerSeparate) {
562:     markerBottom = faceMarkerBottom;
563:     markerTop    = faceMarkerTop;
564:     markerFront  = faceMarkerFront;
565:     markerBack   = faceMarkerBack;
566:     markerRight  = faceMarkerRight;
567:     markerLeft   = faceMarkerLeft;
568:   }
569:   {
570:     const PetscInt numXEdges    = !rank ? edges[0] : 0;
571:     const PetscInt numYEdges    = !rank ? edges[1] : 0;
572:     const PetscInt numZEdges    = !rank ? edges[2] : 0;
573:     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
574:     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
575:     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
576:     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
577:     const PetscInt numXFaces    = numYEdges*numZEdges;
578:     const PetscInt numYFaces    = numXEdges*numZEdges;
579:     const PetscInt numZFaces    = numXEdges*numYEdges;
580:     const PetscInt numTotXFaces = numXVertices*numXFaces;
581:     const PetscInt numTotYFaces = numYVertices*numYFaces;
582:     const PetscInt numTotZFaces = numZVertices*numZFaces;
583:     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
584:     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
585:     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
586:     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
587:     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
588:     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
589:     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
590:     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
591:     const PetscInt firstYFace   = firstXFace + numTotXFaces;
592:     const PetscInt firstZFace   = firstYFace + numTotYFaces;
593:     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
594:     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
595:     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
596:     Vec            coordinates;
597:     PetscSection   coordSection;
598:     PetscScalar   *coords;
599:     PetscInt       coordSize;
600:     PetscInt       v, vx, vy, vz;
601:     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;

603:     DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);
604:     for (c = 0; c < numCells; c++) {
605:       DMPlexSetConeSize(dm, c, 6);
606:     }
607:     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
608:       DMPlexSetConeSize(dm, f, 4);
609:     }
610:     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
611:       DMPlexSetConeSize(dm, e, 2);
612:     }
613:     DMSetUp(dm); /* Allocate space for cones */
614:     /* Build cells */
615:     for (fz = 0; fz < numZEdges; ++fz) {
616:       for (fy = 0; fy < numYEdges; ++fy) {
617:         for (fx = 0; fx < numXEdges; ++fx) {
618:           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
619:           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
620:           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
621:           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
622:           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
623:           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
624:           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
625:                             /* B,  T,  F,  K,  R,  L */
626:           PetscInt ornt[6] = {-4,  0,  0, -1,  0, -4}; /* ??? */
627:           PetscInt cone[6];

629:           /* no boundary twisting in 3D */
630:           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
631:           DMPlexSetCone(dm, cell, cone);
632:           DMPlexSetConeOrientation(dm, cell, ornt);
633:           if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
634:           if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
635:           if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
636:         }
637:       }
638:     }
639:     /* Build x faces */
640:     for (fz = 0; fz < numZEdges; ++fz) {
641:       for (fy = 0; fy < numYEdges; ++fy) {
642:         for (fx = 0; fx < numXVertices; ++fx) {
643:           PetscInt face    = firstXFace + (fz*numYEdges+fy)     *numXVertices+fx;
644:           PetscInt edgeL   = firstZEdge + (fy                   *numXVertices+fx)*numZEdges + fz;
645:           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
646:           PetscInt edgeB   = firstYEdge + (fz                   *numXVertices+fx)*numYEdges + fy;
647:           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
648:           PetscInt ornt[4] = {0, 0, -2, -2};
649:           PetscInt cone[4];

651:           if (dim == 3) {
652:             /* markers */
653:             if (bdX != DM_BOUNDARY_PERIODIC) {
654:               if (fx == numXVertices-1) {
655:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);
656:                 DMSetLabelValue(dm, "marker", face, markerRight);
657:               }
658:               else if (fx == 0) {
659:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);
660:                 DMSetLabelValue(dm, "marker", face, markerLeft);
661:               }
662:             }
663:           }
664:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
665:           DMPlexSetCone(dm, face, cone);
666:           DMPlexSetConeOrientation(dm, face, ornt);
667:         }
668:       }
669:     }
670:     /* Build y faces */
671:     for (fz = 0; fz < numZEdges; ++fz) {
672:       for (fx = 0; fx < numXEdges; ++fx) {
673:         for (fy = 0; fy < numYVertices; ++fy) {
674:           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
675:           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx)*numZEdges + fz;
676:           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
677:           PetscInt edgeB   = firstXEdge + (fz                   *numYVertices+fy)*numXEdges + fx;
678:           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
679:           PetscInt ornt[4] = {0, 0, -2, -2};
680:           PetscInt cone[4];

682:           if (dim == 3) {
683:             /* markers */
684:             if (bdY != DM_BOUNDARY_PERIODIC) {
685:               if (fy == numYVertices-1) {
686:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);
687:                 DMSetLabelValue(dm, "marker", face, markerBack);
688:               }
689:               else if (fy == 0) {
690:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);
691:                 DMSetLabelValue(dm, "marker", face, markerFront);
692:               }
693:             }
694:           }
695:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
696:           DMPlexSetCone(dm, face, cone);
697:           DMPlexSetConeOrientation(dm, face, ornt);
698:         }
699:       }
700:     }
701:     /* Build z faces */
702:     for (fy = 0; fy < numYEdges; ++fy) {
703:       for (fx = 0; fx < numXEdges; ++fx) {
704:         for (fz = 0; fz < numZVertices; fz++) {
705:           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
706:           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx)*numYEdges + fy;
707:           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
708:           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy)*numXEdges + fx;
709:           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
710:           PetscInt ornt[4] = {0, 0, -2, -2};
711:           PetscInt cone[4];

713:           if (dim == 2) {
714:             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
715:             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
716:             if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, face, 2);}
717:             if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, face, 2);}
718:           } else {
719:             /* markers */
720:             if (bdZ != DM_BOUNDARY_PERIODIC) {
721:               if (fz == numZVertices-1) {
722:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);
723:                 DMSetLabelValue(dm, "marker", face, markerTop);
724:               }
725:               else if (fz == 0) {
726:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);
727:                 DMSetLabelValue(dm, "marker", face, markerBottom);
728:               }
729:             }
730:           }
731:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
732:           DMPlexSetCone(dm, face, cone);
733:           DMPlexSetConeOrientation(dm, face, ornt);
734:         }
735:       }
736:     }
737:     /* Build Z edges*/
738:     for (vy = 0; vy < numYVertices; vy++) {
739:       for (vx = 0; vx < numXVertices; vx++) {
740:         for (ez = 0; ez < numZEdges; ez++) {
741:           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
742:           const PetscInt vertexB = firstVertex + (ez                   *numYVertices+vy)*numXVertices + vx;
743:           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
744:           PetscInt       cone[2];

746:           if (dim == 3) {
747:             if (bdX != DM_BOUNDARY_PERIODIC) {
748:               if (vx == numXVertices-1) {
749:                 DMSetLabelValue(dm, "marker", edge, markerRight);
750:               }
751:               else if (vx == 0) {
752:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
753:               }
754:             }
755:             if (bdY != DM_BOUNDARY_PERIODIC) {
756:               if (vy == numYVertices-1) {
757:                 DMSetLabelValue(dm, "marker", edge, markerBack);
758:               }
759:               else if (vy == 0) {
760:                 DMSetLabelValue(dm, "marker", edge, markerFront);
761:               }
762:             }
763:           }
764:           cone[0] = vertexB; cone[1] = vertexT;
765:           DMPlexSetCone(dm, edge, cone);
766:         }
767:       }
768:     }
769:     /* Build Y edges*/
770:     for (vz = 0; vz < numZVertices; vz++) {
771:       for (vx = 0; vx < numXVertices; vx++) {
772:         for (ey = 0; ey < numYEdges; ey++) {
773:           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
774:           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
775:           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
776:           const PetscInt vertexK = firstVertex + nextv;
777:           PetscInt       cone[2];

779:           cone[0] = vertexF; cone[1] = vertexK;
780:           DMPlexSetCone(dm, edge, cone);
781:           if (dim == 2) {
782:             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
783:               if (vx == numXVertices-1) {
784:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);
785:                 DMSetLabelValue(dm, "marker", edge,    markerRight);
786:                 DMSetLabelValue(dm, "marker", cone[0], markerRight);
787:                 if (ey == numYEdges-1) {
788:                   DMSetLabelValue(dm, "marker", cone[1], markerRight);
789:                 }
790:               } else if (vx == 0) {
791:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);
792:                 DMSetLabelValue(dm, "marker", edge,    markerLeft);
793:                 DMSetLabelValue(dm, "marker", cone[0], markerLeft);
794:                 if (ey == numYEdges-1) {
795:                   DMSetLabelValue(dm, "marker", cone[1], markerLeft);
796:                 }
797:               }
798:             } else {
799:               if (vx == 0 && cutLabel) {
800:                 DMLabelSetValue(cutLabel, edge,    1);
801:                 DMLabelSetValue(cutLabel, cone[0], 1);
802:                 if (ey == numYEdges-1) {
803:                   DMLabelSetValue(cutLabel, cone[1], 1);
804:                 }
805:               }
806:             }
807:           } else {
808:             if (bdX != DM_BOUNDARY_PERIODIC) {
809:               if (vx == numXVertices-1) {
810:                 DMSetLabelValue(dm, "marker", edge, markerRight);
811:               } else if (vx == 0) {
812:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
813:               }
814:             }
815:             if (bdZ != DM_BOUNDARY_PERIODIC) {
816:               if (vz == numZVertices-1) {
817:                 DMSetLabelValue(dm, "marker", edge, markerTop);
818:               } else if (vz == 0) {
819:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
820:               }
821:             }
822:           }
823:         }
824:       }
825:     }
826:     /* Build X edges*/
827:     for (vz = 0; vz < numZVertices; vz++) {
828:       for (vy = 0; vy < numYVertices; vy++) {
829:         for (ex = 0; ex < numXEdges; ex++) {
830:           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
831:           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
832:           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
833:           const PetscInt vertexR = firstVertex + nextv;
834:           PetscInt       cone[2];

836:           cone[0] = vertexL; cone[1] = vertexR;
837:           DMPlexSetCone(dm, edge, cone);
838:           if (dim == 2) {
839:             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
840:               if (vy == numYVertices-1) {
841:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);
842:                 DMSetLabelValue(dm, "marker", edge,    markerTop);
843:                 DMSetLabelValue(dm, "marker", cone[0], markerTop);
844:                 if (ex == numXEdges-1) {
845:                   DMSetLabelValue(dm, "marker", cone[1], markerTop);
846:                 }
847:               } else if (vy == 0) {
848:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);
849:                 DMSetLabelValue(dm, "marker", edge,    markerBottom);
850:                 DMSetLabelValue(dm, "marker", cone[0], markerBottom);
851:                 if (ex == numXEdges-1) {
852:                   DMSetLabelValue(dm, "marker", cone[1], markerBottom);
853:                 }
854:               }
855:             } else {
856:               if (vy == 0 && cutLabel) {
857:                 DMLabelSetValue(cutLabel, edge,    1);
858:                 DMLabelSetValue(cutLabel, cone[0], 1);
859:                 if (ex == numXEdges-1) {
860:                   DMLabelSetValue(cutLabel, cone[1], 1);
861:                 }
862:               }
863:             }
864:           } else {
865:             if (bdY != DM_BOUNDARY_PERIODIC) {
866:               if (vy == numYVertices-1) {
867:                 DMSetLabelValue(dm, "marker", edge, markerBack);
868:               }
869:               else if (vy == 0) {
870:                 DMSetLabelValue(dm, "marker", edge, markerFront);
871:               }
872:             }
873:             if (bdZ != DM_BOUNDARY_PERIODIC) {
874:               if (vz == numZVertices-1) {
875:                 DMSetLabelValue(dm, "marker", edge, markerTop);
876:               }
877:               else if (vz == 0) {
878:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
879:               }
880:             }
881:           }
882:         }
883:       }
884:     }
885:     DMPlexSymmetrize(dm);
886:     DMPlexStratify(dm);
887:     /* Build coordinates */
888:     DMGetCoordinateSection(dm, &coordSection);
889:     PetscSectionSetNumFields(coordSection, 1);
890:     PetscSectionSetFieldComponents(coordSection, 0, dim);
891:     PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);
892:     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
893:       PetscSectionSetDof(coordSection, v, dim);
894:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
895:     }
896:     PetscSectionSetUp(coordSection);
897:     PetscSectionGetStorageSize(coordSection, &coordSize);
898:     VecCreate(PETSC_COMM_SELF, &coordinates);
899:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
900:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
901:     VecSetBlockSize(coordinates, dim);
902:     VecSetType(coordinates,VECSTANDARD);
903:     VecGetArray(coordinates, &coords);
904:     for (vz = 0; vz < numZVertices; ++vz) {
905:       for (vy = 0; vy < numYVertices; ++vy) {
906:         for (vx = 0; vx < numXVertices; ++vx) {
907:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
908:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
909:           if (dim == 3) {
910:             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
911:           }
912:         }
913:       }
914:     }
915:     VecRestoreArray(coordinates, &coords);
916:     DMSetCoordinatesLocal(dm, coordinates);
917:     VecDestroy(&coordinates);
918:   }
919:   return(0);
920: }

922: 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)
923: {
924:   PetscInt       i;

929:   DMCreate(comm, dm);
931:   DMSetType(*dm, DMPLEX);
932:   DMSetDimension(*dm, dim);
933:   switch (dim) {
934:   case 2: {DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], DM_BOUNDARY_NONE);break;}
935:   case 3: {DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], periodicity[2]);break;}
936:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %D", dim);
937:   }
938:   if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST ||
939:       periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST ||
940:       (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) {
941:     PetscReal L[3];
942:     PetscReal maxCell[3];

944:     for (i = 0; i < dim; i++) {
945:       L[i]       = upper[i] - lower[i];
946:       maxCell[i] = 1.1 * (L[i] / PetscMax(1,faces[i]));
947:     }
948:     DMSetPeriodicity(*dm,PETSC_TRUE,maxCell,L,periodicity);
949:   }
950:   if (!interpolate) {
951:     DM udm;

953:     DMPlexUninterpolate(*dm, &udm);
954:     DMPlexCopyCoordinates(*dm, udm);
955:     DMDestroy(dm);
956:     *dm  = udm;
957:   }
958:   DMPlexSetRefinementUniform(*dm, PETSC_TRUE);
959:   return(0);
960: }

962: /*@C
963:   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices or tensor cells (hexahedra).

965:   Collective

967:   Input Parameters:
968: + comm        - The communicator for the DM object
969: . dim         - The spatial dimension
970: . simplex     - PETSC_TRUE for simplices, PETSC_FALSE for tensor cells
971: . faces       - Number of faces per dimension, or NULL for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
972: . lower       - The lower left corner, or NULL for (0, 0, 0)
973: . upper       - The upper right corner, or NULL for (1, 1, 1)
974: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
975: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

977:   Output Parameter:
978: . dm  - The DM object

980:   Options Database Keys:
981:   These options override the hard-wired input values.
982: + -dm_plex_box_dim <dim>          - Set the topological dimension
983: . -dm_plex_box_simplex <bool>     - PETSC_TRUE for simplex elements, PETSC_FALSE for tensor elements
984: . -dm_plex_box_lower <x,y,z>      - Specify lower-left-bottom coordinates for the box
985: . -dm_plex_box_upper <x,y,z>      - Specify upper-right-top coordinates for the box
986: . -dm_plex_box_faces <m,n,p>      - Number of faces in each linear direction
987: . -dm_plex_box_bd <bx,by,bz>      - Specify the DMBoundaryType for each direction
988: - -dm_plex_box_interpolate <bool> - PETSC_TRUE turns on topological interpolation (creating edges and faces)

990:   Notes:
991:   The options database keys above take lists of length d in d dimensions.

993:   Here is the numbering returned for 2 faces in each direction for tensor cells:
994: $ 10---17---11---18----12
995: $  |         |         |
996: $  |         |         |
997: $ 20    2   22    3    24
998: $  |         |         |
999: $  |         |         |
1000: $  7---15----8---16----9
1001: $  |         |         |
1002: $  |         |         |
1003: $ 19    0   21    1   23
1004: $  |         |         |
1005: $  |         |         |
1006: $  4---13----5---14----6

1008: and for simplicial cells

1010: $ 14----8---15----9----16
1011: $  |\     5  |\      7 |
1012: $  | \       | \       |
1013: $ 13   2    14    3    15
1014: $  | 4   \   | 6   \   |
1015: $  |       \ |       \ |
1016: $ 11----6---12----7----13
1017: $  |\        |\        |
1018: $  | \    1  | \     3 |
1019: $ 10   0    11    1    12
1020: $  | 0   \   | 2   \   |
1021: $  |       \ |       \ |
1022: $  8----4----9----5----10

1024:   Level: beginner

1026: .seealso: DMPlexCreateFromFile(), DMPlexCreateHexCylinderMesh(), DMSetType(), DMCreate()
1027: @*/
1028: 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)
1029: {
1030:   PetscInt       fac[3] = {0, 0, 0};
1031:   PetscReal      low[3] = {0, 0, 0};
1032:   PetscReal      upp[3] = {1, 1, 1};
1033:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1034:   PetscInt       i, n;
1035:   PetscBool      flg;

1039:   PetscOptionsGetInt(NULL, NULL, "-dm_plex_box_dim", &dim, &flg);
1040:   if ((dim < 0) || (dim > 3)) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D should be in [1, 3]", dim);
1041:   PetscOptionsGetBool(NULL, NULL, "-dm_plex_box_simplex", &simplex, &flg);
1042:   n    = dim;
1043:   PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", fac, &n, &flg);
1044:   for (i = 0; i < dim; ++i) fac[i] = faces ? faces[i] : (flg && i < n ? fac[i] : (dim == 1 ? 1 : 4-dim));
1045:   if (lower) for (i = 0; i < dim; ++i) low[i] = lower[i];
1046:   if (upper) for (i = 0; i < dim; ++i) upp[i] = upper[i];
1047:   if (periodicity) for (i = 0; i < dim; ++i) bdt[i] = periodicity[i];
1048:   /* Allow bounds to be specified from the command line */
1049:   n    = 3;
1050:   PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_lower", low, &n, &flg);
1051:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Lower box point had %D values, should have been %D", n, dim);
1052:   n    = 3;
1053:   PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_upper", upp, &n, &flg);
1054:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Upper box point had %D values, should have been %D", n, dim);
1055:   n    = 3;
1056:   PetscOptionsGetEnumArray(NULL, NULL, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, &flg);
1057:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Box boundary types had %D values, should have been %D", n, dim);
1058:   PetscOptionsGetBool(NULL, NULL, "-dm_plex_box_interpolate", &interpolate, &flg);

1060:   if (dim == 1)      {DMPlexCreateLineMesh_Internal(comm, fac[0], low[0], upp[0], bdt[0], dm);}
1061:   else if (simplex)  {DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1062:   else               {DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1063:   return(0);
1064: }

1066: /*@
1067:   DMPlexCreateWedgeBoxMesh - Creates a 3-D mesh tesselating the (x,y) plane and extruding in the third direction using wedge cells.

1069:   Collective

1071:   Input Parameters:
1072: + comm        - The communicator for the DM object
1073: . faces       - Number of faces per dimension, or NULL for (1, 1, 1)
1074: . lower       - The lower left corner, or NULL for (0, 0, 0)
1075: . upper       - The upper right corner, or NULL for (1, 1, 1)
1076: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
1077: . orderHeight - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1078: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1080:   Output Parameter:
1081: . dm  - The DM object

1083:   Level: beginner

1085: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateWedgeCylinderMesh(), DMPlexExtrude(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1086: @*/
1087: PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool orderHeight, PetscBool interpolate, DM *dm)
1088: {
1089:   DM             bdm, botdm;
1090:   PetscInt       i;
1091:   PetscInt       fac[3] = {0, 0, 0};
1092:   PetscReal      low[3] = {0, 0, 0};
1093:   PetscReal      upp[3] = {1, 1, 1};
1094:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1098:   for (i = 0; i < 3; ++i) fac[i] = faces ? (faces[i] > 0 ? faces[i] : 1) : 1;
1099:   if (lower) for (i = 0; i < 3; ++i) low[i] = lower[i];
1100:   if (upper) for (i = 0; i < 3; ++i) upp[i] = upper[i];
1101:   if (periodicity) for (i = 0; i < 3; ++i) bdt[i] = periodicity[i];
1102:   for (i = 0; i < 3; ++i) if (bdt[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity not yet supported");

1104:   DMCreate(comm, &bdm);
1105:   DMSetType(bdm, DMPLEX);
1106:   DMSetDimension(bdm, 1);
1107:   DMSetCoordinateDim(bdm, 2);
1108:   DMPlexCreateSquareBoundary(bdm, low, upp, fac);
1109:   DMPlexGenerate(bdm, NULL, PETSC_FALSE, &botdm);
1110:   DMDestroy(&bdm);
1111:   DMPlexExtrude(botdm, fac[2], upp[2] - low[2], orderHeight, NULL, interpolate, dm);
1112:   if (low[2] != 0.0) {
1113:     Vec         v;
1114:     PetscScalar *x;
1115:     PetscInt    cDim, n;

1117:     DMGetCoordinatesLocal(*dm, &v);
1118:     VecGetBlockSize(v, &cDim);
1119:     VecGetLocalSize(v, &n);
1120:     VecGetArray(v, &x);
1121:     x   += cDim;
1122:     for (i=0; i<n; i+=cDim) x[i] += low[2];
1123:     VecRestoreArray(v,&x);
1124:     DMSetCoordinatesLocal(*dm, v);
1125:   }
1126:   DMDestroy(&botdm);
1127:   return(0);
1128: }

1130: /*@C
1131:   DMPlexExtrude - Creates a (d+1)-D mesh by extruding a d-D mesh in the normal direction using prismatic cells.

1133:   Collective on idm

1135:   Input Parameters:
1136: + idm         - The mesh to be extruded
1137: . layers      - The number of layers, or PETSC_DETERMINE to use the default
1138: . height      - The total height of the extrusion, or PETSC_DETERMINE to use the default
1139: . orderHeight - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1140: . extNormal   - The normal direction in which the mesh should be extruded, or NULL to extrude using the surface normal
1141: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1143:   Output Parameter:
1144: . dm  - The DM object

1146:   Notes:
1147:   The mesh created has prismatic cells, and the vertex ordering in the cone of the cell is that of the tensor prismatic cells. Not currently supported in Fortran.

1149:   Options Database Keys:
1150: +   -dm_plex_extrude_layers <k> - Sets the number of layers k
1151: .   -dm_plex_extrude_height <h> - Sets the total height of the extrusion
1152: .   -dm_plex_extrude_heights <h0,h1,...> - Sets the height of each layer
1153: .   -dm_plex_extrude_order_height - If true, order cells by height first
1154: -   -dm_plex_extrude_normal <n0,...,nd> - Sets the normal vector along which to extrude

1156:   Level: advanced

1158: .seealso: DMPlexCreateWedgeCylinderMesh(), DMPlexCreateWedgeBoxMesh(), DMSetType(), DMCreate()
1159: @*/
1160: PetscErrorCode DMPlexExtrude(DM idm, PetscInt layers, PetscReal height, PetscBool orderHeight, const PetscReal extNormal[], PetscBool interpolate, DM* dm)
1161: {
1162:   PetscScalar       *coordsB;
1163:   const PetscScalar *coordsA;
1164:   PetscReal         *normals = NULL, *heights = NULL;
1165:   PetscReal         clNormal[3];
1166:   Vec               coordinatesA, coordinatesB;
1167:   PetscSection      coordSectionA, coordSectionB;
1168:   PetscInt          dim, cDim, cDimB, c, l, v, coordSize, *newCone, nl;
1169:   PetscInt          cStart, cEnd, vStart, vEnd, cellV, numCells, numVertices;
1170:   const char       *prefix;
1171:   PetscBool         haveCLNormal, flg;
1172:   PetscErrorCode    ierr;

1179:   DMGetDimension(idm, &dim);
1180:   DMGetCoordinateDim(idm, &cDim);
1181:   cDimB = cDim == dim ? cDim+1 : cDim;
1182:   if (dim < 1 || dim > 3) SETERRQ1(PetscObjectComm((PetscObject)idm), PETSC_ERR_SUP, "Support for dimension %D not coded", dim);

1184:   PetscObjectGetOptionsPrefix((PetscObject) idm, &prefix);
1185:   if (layers < 0) layers = 1;
1186:   PetscOptionsGetInt(NULL, prefix, "-dm_plex_extrude_layers", &layers, NULL);
1187:   if (layers <= 0) SETERRQ1(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Number of layers %D must be positive", layers);
1188:   if (height < 0.) height = 1.;
1189:   PetscOptionsGetReal(NULL, prefix, "-dm_plex_extrude_height", &height, NULL);
1190:   if (height <= 0.) SETERRQ1(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Height of layers %g must be positive", (double) height);
1191:   PetscMalloc1(layers, &heights);
1192:   nl   = layers;
1193:   PetscOptionsGetRealArray(NULL, prefix, "-dm_plex_extrude_heights", heights, &nl, &flg);
1194:   if (flg) {
1195:     if (!nl) SETERRQ(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Must give at least one height for -dm_plex_extrude_heights");
1196:     for (l = nl; l < layers; ++l) heights[l] = heights[l-1];
1197:     for (l = 0; l < layers; ++l) if (heights[l] <= 0.) SETERRQ2(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Height %g of layers %D must be positive", (double) heights[l], l);
1198:   } else {
1199:     for (l = 0; l < layers; ++l) heights[l] = height/layers;
1200:   }
1201:   PetscOptionsGetBool(NULL, prefix, "-dm_plex_extrude_order_height", &orderHeight, NULL);
1202:   c = 3;
1203:   PetscOptionsGetRealArray(NULL, prefix, "-dm_plex_extrude_normal", clNormal, &c, &haveCLNormal);
1204:   if (haveCLNormal && c != cDimB) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_ARG_SIZ, "Input normal has size %D != %D extruded coordinate dimension", c, cDimB);

1206:   DMPlexGetHeightStratum(idm, 0, &cStart, &cEnd);
1207:   DMPlexGetDepthStratum(idm, 0, &vStart, &vEnd);
1208:   numCells = (cEnd - cStart)*layers;
1209:   numVertices = (vEnd - vStart)*(layers+1);
1210:   DMCreate(PetscObjectComm((PetscObject)idm), dm);
1211:   DMSetType(*dm, DMPLEX);
1212:   DMSetDimension(*dm, dim+1);
1213:   DMPlexSetChart(*dm, 0, numCells+numVertices);
1214:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1215:   DMCreateLabel(*dm, "celltype");
1216:   for (c = cStart, cellV = 0; c < cEnd; ++c) {
1217:     DMPolytopeType ct, nct;
1218:     PetscInt      *closure = NULL;
1219:     PetscInt       closureSize, numCorners = 0;

1221:     DMPlexGetCellType(idm, c, &ct);
1222:     switch (ct) {
1223:       case DM_POLYTOPE_SEGMENT:       nct = DM_POLYTOPE_SEG_PRISM_TENSOR;break;
1224:       case DM_POLYTOPE_TRIANGLE:      nct = DM_POLYTOPE_TRI_PRISM_TENSOR;break;
1225:       case DM_POLYTOPE_QUADRILATERAL: nct = DM_POLYTOPE_QUAD_PRISM_TENSOR;break;
1226:       default: nct = DM_POLYTOPE_UNKNOWN;
1227:     }
1228:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1229:     for (v = 0; v < closureSize*2; v += 2) if ((closure[v] >= vStart) && (closure[v] < vEnd)) numCorners++;
1230:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1231:     for (l = 0; l < layers; ++l) {
1232:       const PetscInt cell = orderHeight ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart;

1234:       DMPlexSetConeSize(*dm, cell, 2*numCorners);
1235:       DMPlexSetCellType(*dm, cell, nct);
1236:     }
1237:     cellV = PetscMax(numCorners,cellV);
1238:   }
1239:   DMSetUp(*dm);

1241:   if (dim != cDim && !(extNormal || haveCLNormal)) {PetscCalloc1(cDim*(vEnd - vStart), &normals);}
1242:   PetscMalloc1(3*cellV,&newCone);
1243:   for (c = cStart; c < cEnd; ++c) {
1244:     PetscInt *closure = NULL;
1245:     PetscInt closureSize, numCorners = 0, l;
1246:     PetscReal normal[3] = {0, 0, 0};

1248:     if (normals) {DMPlexComputeCellGeometryFVM(idm, c, NULL, NULL, normal);}
1249:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1250:     for (v = 0; v < closureSize*2; v += 2) {
1251:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1252:         PetscInt d;

1254:         newCone[numCorners++] = closure[v] - vStart;
1255:         if (normals) {for (d = 0; d < cDim; ++d) normals[cDim*(closure[v]-vStart)+d] += normal[d];}
1256:       }
1257:     }
1258:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1259:     for (l = 0; l < layers; ++l) {
1260:       PetscInt i;

1262:       for (i = 0; i < numCorners; ++i) {
1263:         newCone[  numCorners + i] = orderHeight ? (layers+1)*newCone[i] + l     + numCells :     l*(vEnd - vStart) + newCone[i] + numCells;
1264:         newCone[2*numCorners + i] = orderHeight ? (layers+1)*newCone[i] + l + 1 + numCells : (l+1)*(vEnd - vStart) + newCone[i] + numCells;
1265:       }
1266:       DMPlexSetCone(*dm, orderHeight ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart, newCone + numCorners);
1267:     }
1268:   }
1269:   DMPlexSymmetrize(*dm);
1270:   DMPlexStratify(*dm);
1271:   PetscFree(newCone);

1273:   DMGetCoordinateSection(*dm, &coordSectionB);
1274:   PetscSectionSetNumFields(coordSectionB, 1);
1275:   PetscSectionSetFieldComponents(coordSectionB, 0, cDimB);
1276:   PetscSectionSetChart(coordSectionB, numCells, numCells+numVertices);
1277:   for (v = numCells; v < numCells+numVertices; ++v) {
1278:     PetscSectionSetDof(coordSectionB, v, cDimB);
1279:     PetscSectionSetFieldDof(coordSectionB, v, 0, cDimB);
1280:     DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);
1281:   }
1282:   PetscSectionSetUp(coordSectionB);
1283:   PetscSectionGetStorageSize(coordSectionB, &coordSize);
1284:   VecCreate(PETSC_COMM_SELF, &coordinatesB);
1285:   PetscObjectSetName((PetscObject) coordinatesB, "coordinates");
1286:   VecSetSizes(coordinatesB, coordSize, PETSC_DETERMINE);
1287:   VecSetBlockSize(coordinatesB, cDimB);
1288:   VecSetType(coordinatesB,VECSTANDARD);

1290:   DMGetCoordinateSection(idm, &coordSectionA);
1291:   DMGetCoordinatesLocal(idm, &coordinatesA);
1292:   VecGetArray(coordinatesB, &coordsB);
1293:   VecGetArrayRead(coordinatesA, &coordsA);
1294:   for (v = vStart; v < vEnd; ++v) {
1295:     const PetscScalar *cptr;
1296:     PetscReal         ones2[2] = { 0., 1.}, ones3[3] = { 0., 0., 1.};
1297:     PetscReal         normal[3];
1298:     PetscReal         norm;
1299:     PetscInt          offA, d, cDimA = cDim;

1301:     if (normals)           {for (d = 0; d < cDimB; ++d) normal[d] = normals[cDimB*(v - vStart)+d];}
1302:     else if (haveCLNormal) {for (d = 0; d < cDimB; ++d) normal[d] = clNormal[d];}
1303:     else if (extNormal)    {for (d = 0; d < cDimB; ++d) normal[d] = extNormal[d];}
1304:     else if (cDimB == 2)   {for (d = 0; d < cDimB; ++d) normal[d] = ones2[d];}
1305:     else if (cDimB == 3)   {for (d = 0; d < cDimB; ++d) normal[d] = ones3[d];}
1306:     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to determine normal for extrusion");
1307:     for (d = 0, norm = 0.0; d < cDimB; ++d) norm += normal[d]*normal[d];
1308:     for (d = 0; d < cDimB; ++d) normal[d] *= 1./PetscSqrtReal(norm);

1310:     PetscSectionGetOffset(coordSectionA, v, &offA);
1311:     cptr = coordsA + offA;
1312:     for (l = 0; l <= layers; ++l) {
1313:       PetscInt offB, d, newV;

1315:       newV = orderHeight ? (layers+1)*(v -vStart) + l + numCells : (vEnd -vStart)*l + (v -vStart) + numCells;
1316:       PetscSectionGetOffset(coordSectionB, newV, &offB);
1317:       for (d = 0; d < cDimA; ++d) { coordsB[offB+d]  = cptr[d]; }
1318:       for (d = 0; d < cDimB; ++d) { coordsB[offB+d] += l ? normal[d]*heights[l-1] : 0.0; }
1319:       cptr    = coordsB + offB;
1320:       cDimA   = cDimB;
1321:     }
1322:   }
1323:   VecRestoreArrayRead(coordinatesA, &coordsA);
1324:   VecRestoreArray(coordinatesB, &coordsB);
1325:   DMSetCoordinatesLocal(*dm, coordinatesB);
1326:   VecDestroy(&coordinatesB);
1327:   PetscFree(normals);
1328:   PetscFree(heights);
1329:   if (interpolate) {
1330:     DM idm;

1332:     DMPlexInterpolate(*dm, &idm);
1333:     DMPlexCopyCoordinates(*dm, idm);
1334:     DMDestroy(dm);
1335:     *dm  = idm;
1336:   }
1337:   return(0);
1338: }

1340: /*@C
1341:   DMPlexSetOptionsPrefix - Sets the prefix used for searching for all DM options in the database.

1343:   Logically Collective on dm

1345:   Input Parameters:
1346: + dm - the DM context
1347: - prefix - the prefix to prepend to all option names

1349:   Notes:
1350:   A hyphen (-) must NOT be given at the beginning of the prefix name.
1351:   The first character of all runtime options is AUTOMATICALLY the hyphen.

1353:   Level: advanced

1355: .seealso: SNESSetFromOptions()
1356: @*/
1357: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1358: {
1359:   DM_Plex       *mesh = (DM_Plex *) dm->data;

1364:   PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);
1365:   PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);
1366:   return(0);
1367: }

1369: /*@
1370:   DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra.

1372:   Collective

1374:   Input Parameters:
1375: + comm      - The communicator for the DM object
1376: . numRefine - The number of regular refinements to the basic 5 cell structure
1377: - periodicZ - The boundary type for the Z direction

1379:   Output Parameter:
1380: . dm  - The DM object

1382:   Options Database Keys:
1383: These options override the hard-wired input values.
1384: + -dm_plex_hex_cyl_refine <r> - Refine the sylinder r times
1385: - -dm_plex_hex_cyl_bd <bz>    - Specify the DMBoundaryType in the z-direction

1387:   Note:
1388:   Here is the output numbering looking from the bottom of the cylinder:
1389: $       17-----14
1390: $        |     |
1391: $        |  2  |
1392: $        |     |
1393: $ 17-----8-----7-----14
1394: $  |     |     |     |
1395: $  |  3  |  0  |  1  |
1396: $  |     |     |     |
1397: $ 19-----5-----6-----13
1398: $        |     |
1399: $        |  4  |
1400: $        |     |
1401: $       19-----13
1402: $
1403: $ and up through the top
1404: $
1405: $       18-----16
1406: $        |     |
1407: $        |  2  |
1408: $        |     |
1409: $ 18----10----11-----16
1410: $  |     |     |     |
1411: $  |  3  |  0  |  1  |
1412: $  |     |     |     |
1413: $ 20-----9----12-----15
1414: $        |     |
1415: $        |  4  |
1416: $        |     |
1417: $       20-----15

1419:   Level: beginner

1421: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1422: @*/
1423: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
1424: {
1425:   const PetscInt dim = 3;
1426:   PetscInt       numCells, numVertices, r;
1427:   PetscMPIInt    rank;

1432:   MPI_Comm_rank(comm, &rank);
1433:   PetscOptionsGetInt(NULL, NULL, "-dm_plex_hex_cyl_refine", &numRefine, NULL);
1434:   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
1435:   PetscOptionsGetEnum(NULL, NULL, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) &periodicZ, NULL);
1436:   DMCreate(comm, dm);
1437:   DMSetType(*dm, DMPLEX);
1438:   DMSetDimension(*dm, dim);
1439:   /* Create topology */
1440:   {
1441:     PetscInt cone[8], c;

1443:     numCells    = !rank ?  5 : 0;
1444:     numVertices = !rank ? 16 : 0;
1445:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1446:       numCells   *= 3;
1447:       numVertices = !rank ? 24 : 0;
1448:     }
1449:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1450:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 8);}
1451:     DMSetUp(*dm);
1452:     if (!rank) {
1453:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1454:         cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1455:         cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1456:         DMPlexSetCone(*dm, 0, cone);
1457:         cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1458:         cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1459:         DMPlexSetCone(*dm, 1, cone);
1460:         cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1461:         cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1462:         DMPlexSetCone(*dm, 2, cone);
1463:         cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1464:         cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1465:         DMPlexSetCone(*dm, 3, cone);
1466:         cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1467:         cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1468:         DMPlexSetCone(*dm, 4, cone);

1470:         cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1471:         cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1472:         DMPlexSetCone(*dm, 5, cone);
1473:         cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1474:         cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1475:         DMPlexSetCone(*dm, 6, cone);
1476:         cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1477:         cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1478:         DMPlexSetCone(*dm, 7, cone);
1479:         cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1480:         cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1481:         DMPlexSetCone(*dm, 8, cone);
1482:         cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1483:         cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1484:         DMPlexSetCone(*dm, 9, cone);

1486:         cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1487:         cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1488:         DMPlexSetCone(*dm, 10, cone);
1489:         cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1490:         cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1491:         DMPlexSetCone(*dm, 11, cone);
1492:         cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1493:         cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1494:         DMPlexSetCone(*dm, 12, cone);
1495:         cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1496:         cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1497:         DMPlexSetCone(*dm, 13, cone);
1498:         cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1499:         cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1500:         DMPlexSetCone(*dm, 14, cone);
1501:       } else {
1502:         cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
1503:         cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
1504:         DMPlexSetCone(*dm, 0, cone);
1505:         cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
1506:         cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
1507:         DMPlexSetCone(*dm, 1, cone);
1508:         cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
1509:         cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
1510:         DMPlexSetCone(*dm, 2, cone);
1511:         cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
1512:         cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
1513:         DMPlexSetCone(*dm, 3, cone);
1514:         cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
1515:         cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
1516:         DMPlexSetCone(*dm, 4, cone);
1517:       }
1518:     }
1519:     DMPlexSymmetrize(*dm);
1520:     DMPlexStratify(*dm);
1521:   }
1522:   /* Interpolate */
1523:   {
1524:     DM idm;

1526:     DMPlexInterpolate(*dm, &idm);
1527:     DMDestroy(dm);
1528:     *dm  = idm;
1529:   }
1530:   /* Create cube geometry */
1531:   {
1532:     Vec             coordinates;
1533:     PetscSection    coordSection;
1534:     PetscScalar    *coords;
1535:     PetscInt        coordSize, v;
1536:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1537:     const PetscReal ds2 = dis/2.0;

1539:     /* Build coordinates */
1540:     DMGetCoordinateSection(*dm, &coordSection);
1541:     PetscSectionSetNumFields(coordSection, 1);
1542:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1543:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1544:     for (v = numCells; v < numCells+numVertices; ++v) {
1545:       PetscSectionSetDof(coordSection, v, dim);
1546:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1547:     }
1548:     PetscSectionSetUp(coordSection);
1549:     PetscSectionGetStorageSize(coordSection, &coordSize);
1550:     VecCreate(PETSC_COMM_SELF, &coordinates);
1551:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1552:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1553:     VecSetBlockSize(coordinates, dim);
1554:     VecSetType(coordinates,VECSTANDARD);
1555:     VecGetArray(coordinates, &coords);
1556:     if (!rank) {
1557:       coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
1558:       coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
1559:       coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
1560:       coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
1561:       coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
1562:       coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
1563:       coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
1564:       coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
1565:       coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
1566:       coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
1567:       coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
1568:       coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
1569:       coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
1570:       coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
1571:       coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
1572:       coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1573:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1574:         /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1575:         /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1576:         /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1577:         /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1578:         /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1579:         /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1580:         /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1581:         /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1582:       }
1583:     }
1584:     VecRestoreArray(coordinates, &coords);
1585:     DMSetCoordinatesLocal(*dm, coordinates);
1586:     VecDestroy(&coordinates);
1587:   }
1588:   /* Create periodicity */
1589:   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1590:     PetscReal      L[3];
1591:     PetscReal      maxCell[3];
1592:     DMBoundaryType bdType[3];
1593:     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1594:     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1595:     PetscInt       i, numZCells = 3;

1597:     bdType[0] = DM_BOUNDARY_NONE;
1598:     bdType[1] = DM_BOUNDARY_NONE;
1599:     bdType[2] = periodicZ;
1600:     for (i = 0; i < dim; i++) {
1601:       L[i]       = upper[i] - lower[i];
1602:       maxCell[i] = 1.1 * (L[i] / numZCells);
1603:     }
1604:     DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);
1605:   }
1606:   /* Refine topology */
1607:   for (r = 0; r < numRefine; ++r) {
1608:     DM rdm = NULL;

1610:     DMRefine(*dm, comm, &rdm);
1611:     DMDestroy(dm);
1612:     *dm  = rdm;
1613:   }
1614:   /* Remap geometry to cylinder
1615:        Interior square: Linear interpolation is correct
1616:        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1617:        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance

1619:          phi     = arctan(y/x)
1620:          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1621:          d_far   = sqrt(1/2 + sin^2(phi))

1623:        so we remap them using

1625:          x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close)
1626:          y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close)

1628:        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1629:   */
1630:   {
1631:     Vec           coordinates;
1632:     PetscSection  coordSection;
1633:     PetscScalar  *coords;
1634:     PetscInt      vStart, vEnd, v;
1635:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1636:     const PetscReal ds2 = 0.5*dis;

1638:     DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);
1639:     DMGetCoordinateSection(*dm, &coordSection);
1640:     DMGetCoordinatesLocal(*dm, &coordinates);
1641:     VecGetArray(coordinates, &coords);
1642:     for (v = vStart; v < vEnd; ++v) {
1643:       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1644:       PetscInt  off;

1646:       PetscSectionGetOffset(coordSection, v, &off);
1647:       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
1648:       x    = PetscRealPart(coords[off]);
1649:       y    = PetscRealPart(coords[off+1]);
1650:       phi  = PetscAtan2Real(y, x);
1651:       sinp = PetscSinReal(phi);
1652:       cosp = PetscCosReal(phi);
1653:       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
1654:         dc = PetscAbsReal(ds2/sinp);
1655:         df = PetscAbsReal(dis/sinp);
1656:         xc = ds2*x/PetscAbsReal(y);
1657:         yc = ds2*PetscSignReal(y);
1658:       } else {
1659:         dc = PetscAbsReal(ds2/cosp);
1660:         df = PetscAbsReal(dis/cosp);
1661:         xc = ds2*PetscSignReal(x);
1662:         yc = ds2*y/PetscAbsReal(x);
1663:       }
1664:       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1665:       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1666:     }
1667:     VecRestoreArray(coordinates, &coords);
1668:     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1669:       DMLocalizeCoordinates(*dm);
1670:     }
1671:   }
1672:   return(0);
1673: }

1675: /*@
1676:   DMPlexCreateWedgeCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using wedges.

1678:   Collective

1680:   Input Parameters:
1681: + comm - The communicator for the DM object
1682: . n    - The number of wedges around the origin
1683: - interpolate - Create edges and faces

1685:   Output Parameter:
1686: . dm  - The DM object

1688:   Level: beginner

1690: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1691: @*/
1692: PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
1693: {
1694:   const PetscInt dim = 3;
1695:   PetscInt       numCells, numVertices, v;
1696:   PetscMPIInt    rank;

1701:   MPI_Comm_rank(comm, &rank);
1702:   if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n);
1703:   DMCreate(comm, dm);
1704:   DMSetType(*dm, DMPLEX);
1705:   DMSetDimension(*dm, dim);
1706:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1707:   DMCreateLabel(*dm, "celltype");
1708:   /* Create topology */
1709:   {
1710:     PetscInt cone[6], c;

1712:     numCells    = !rank ?        n : 0;
1713:     numVertices = !rank ?  2*(n+1) : 0;
1714:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1715:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 6);}
1716:     DMSetUp(*dm);
1717:     for (c = 0; c < numCells; c++) {
1718:       cone[0] =  c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n;
1719:       cone[3] =  c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n;
1720:       DMPlexSetCone(*dm, c, cone);
1721:       DMPlexSetCellType(*dm, c, DM_POLYTOPE_TRI_PRISM_TENSOR);
1722:     }
1723:     DMPlexSymmetrize(*dm);
1724:     DMPlexStratify(*dm);
1725:   }
1726:   for (v = numCells; v < numCells+numVertices; ++v) {
1727:     DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);
1728:   }
1729:   /* Interpolate */
1730:   if (interpolate) {
1731:     DM idm;

1733:     DMPlexInterpolate(*dm, &idm);
1734:     DMDestroy(dm);
1735:     *dm  = idm;
1736:   }
1737:   /* Create cylinder geometry */
1738:   {
1739:     Vec          coordinates;
1740:     PetscSection coordSection;
1741:     PetscScalar *coords;
1742:     PetscInt     coordSize, c;

1744:     /* Build coordinates */
1745:     DMGetCoordinateSection(*dm, &coordSection);
1746:     PetscSectionSetNumFields(coordSection, 1);
1747:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1748:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1749:     for (v = numCells; v < numCells+numVertices; ++v) {
1750:       PetscSectionSetDof(coordSection, v, dim);
1751:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1752:     }
1753:     PetscSectionSetUp(coordSection);
1754:     PetscSectionGetStorageSize(coordSection, &coordSize);
1755:     VecCreate(PETSC_COMM_SELF, &coordinates);
1756:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1757:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1758:     VecSetBlockSize(coordinates, dim);
1759:     VecSetType(coordinates,VECSTANDARD);
1760:     VecGetArray(coordinates, &coords);
1761:     for (c = 0; c < numCells; c++) {
1762:       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;
1763:       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;
1764:     }
1765:     if (!rank) {
1766:       coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0;
1767:       coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0;
1768:     }
1769:     VecRestoreArray(coordinates, &coords);
1770:     DMSetCoordinatesLocal(*dm, coordinates);
1771:     VecDestroy(&coordinates);
1772:   }
1773:   return(0);
1774: }

1776: PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1777: {
1778:   PetscReal prod = 0.0;
1779:   PetscInt  i;
1780:   for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
1781:   return PetscSqrtReal(prod);
1782: }
1783: PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1784: {
1785:   PetscReal prod = 0.0;
1786:   PetscInt  i;
1787:   for (i = 0; i < dim; ++i) prod += x[i]*y[i];
1788:   return prod;
1789: }

1791: /* The first constant is the sphere radius */
1792: static void snapToSphere(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1793:                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1794:                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1795:                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1796: {
1797:   PetscReal r = PetscRealPart(constants[0]);
1798:   PetscReal norm2 = 0.0, fac;
1799:   PetscInt  n = uOff[1] - uOff[0], d;

1801:   for (d = 0; d < n; ++d) norm2 += PetscSqr(PetscRealPart(u[d]));
1802:   fac = r/PetscSqrtReal(norm2);
1803:   for (d = 0; d < n; ++d) f0[d] = u[d]*fac;
1804: }

1806: /*@
1807:   DMPlexCreateSphereMesh - Creates a mesh on the d-dimensional sphere, S^d.

1809:   Collective

1811:   Input Parameters:
1812: + comm    - The communicator for the DM object
1813: . dim     - The dimension
1814: . simplex - Use simplices, or tensor product cells
1815: - R       - The radius

1817:   Output Parameter:
1818: . dm  - The DM object

1820:   Level: beginner

1822: .seealso: DMPlexCreateBallMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1823: @*/
1824: PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscReal R, DM *dm)
1825: {
1826:   const PetscInt  embedDim = dim+1;
1827:   PetscSection    coordSection;
1828:   Vec             coordinates;
1829:   PetscScalar    *coords;
1830:   PetscReal      *coordsIn;
1831:   PetscInt        numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e;
1832:   PetscMPIInt     rank;
1833:   PetscErrorCode  ierr;

1837:   DMCreate(comm, dm);
1838:   DMSetType(*dm, DMPLEX);
1839:   DMSetDimension(*dm, dim);
1840:   DMSetCoordinateDim(*dm, dim+1);
1841:   MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank);
1842:   switch (dim) {
1843:   case 2:
1844:     if (simplex) {
1845:       DM              idm;
1846:       const PetscReal radius    = PetscSqrtReal(1 + PETSC_PHI*PETSC_PHI)/(1.0 + PETSC_PHI);
1847:       const PetscReal edgeLen   = 2.0/(1.0 + PETSC_PHI) * (R/radius);
1848:       const PetscInt  degree    = 5;
1849:       PetscReal       vertex[3] = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)};
1850:       PetscInt        s[3]      = {1, 1, 1};
1851:       PetscInt        cone[3];
1852:       PetscInt       *graph, p, i, j, k;

1854:       vertex[0] *= R/radius; vertex[1] *= R/radius; vertex[2] *= R/radius;
1855:       numCells    = !rank ? 20 : 0;
1856:       numVerts    = !rank ? 12 : 0;
1857:       firstVertex = numCells;
1858:       /* Use icosahedron, which for a R-sphere has coordinates which are all cyclic permutations of

1860:            (0, \pm 1/\phi+1, \pm \phi/\phi+1)

1862:          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
1863:          length is then given by 2/(1+\phi) = 2 * 0.38197 = 0.76393.
1864:       */
1865:       /* Construct vertices */
1866:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1867:       if (!rank) {
1868:         for (p = 0, i = 0; p < embedDim; ++p) {
1869:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
1870:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
1871:               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim];
1872:               ++i;
1873:             }
1874:           }
1875:         }
1876:       }
1877:       /* Construct graph */
1878:       PetscCalloc1(numVerts * numVerts, &graph);
1879:       for (i = 0; i < numVerts; ++i) {
1880:         for (j = 0, k = 0; j < numVerts; ++j) {
1881:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
1882:         }
1883:         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid icosahedron, vertex %D degree %D != %D", i, k, degree);
1884:       }
1885:       /* Build Topology */
1886:       DMPlexSetChart(*dm, 0, numCells+numVerts);
1887:       for (c = 0; c < numCells; c++) {
1888:         DMPlexSetConeSize(*dm, c, embedDim);
1889:       }
1890:       DMSetUp(*dm); /* Allocate space for cones */
1891:       /* Cells */
1892:       for (i = 0, c = 0; i < numVerts; ++i) {
1893:         for (j = 0; j < i; ++j) {
1894:           for (k = 0; k < j; ++k) {
1895:             if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) {
1896:               cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k;
1897:               /* Check orientation */
1898:               {
1899:                 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}}};
1900:                 PetscReal normal[3];
1901:                 PetscInt  e, f;

1903:                 for (d = 0; d < embedDim; ++d) {
1904:                   normal[d] = 0.0;
1905:                   for (e = 0; e < embedDim; ++e) {
1906:                     for (f = 0; f < embedDim; ++f) {
1907:                       normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]);
1908:                     }
1909:                   }
1910:                 }
1911:                 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
1912:               }
1913:               DMPlexSetCone(*dm, c++, cone);
1914:             }
1915:           }
1916:         }
1917:       }
1918:       DMPlexSymmetrize(*dm);
1919:       DMPlexStratify(*dm);
1920:       PetscFree(graph);
1921:       /* Interpolate mesh */
1922:       DMPlexInterpolate(*dm, &idm);
1923:       DMDestroy(dm);
1924:       *dm  = idm;
1925:     } else {
1926:       /*
1927:         12-21--13
1928:          |     |
1929:         25  4  24
1930:          |     |
1931:   12-25--9-16--8-24--13
1932:    |     |     |     |
1933:   23  5 17  0 15  3  22
1934:    |     |     |     |
1935:   10-20--6-14--7-19--11
1936:          |     |
1937:         20  1  19
1938:          |     |
1939:         10-18--11
1940:          |     |
1941:         23  2  22
1942:          |     |
1943:         12-21--13
1944:        */
1945:       PetscInt cone[4], ornt[4];

1947:       numCells    = !rank ?  6 : 0;
1948:       numEdges    = !rank ? 12 : 0;
1949:       numVerts    = !rank ?  8 : 0;
1950:       firstVertex = numCells;
1951:       firstEdge   = numCells + numVerts;
1952:       /* Build Topology */
1953:       DMPlexSetChart(*dm, 0, numCells+numEdges+numVerts);
1954:       for (c = 0; c < numCells; c++) {
1955:         DMPlexSetConeSize(*dm, c, 4);
1956:       }
1957:       for (e = firstEdge; e < firstEdge+numEdges; ++e) {
1958:         DMPlexSetConeSize(*dm, e, 2);
1959:       }
1960:       DMSetUp(*dm); /* Allocate space for cones */
1961:       if (!rank) {
1962:         /* Cell 0 */
1963:         cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17;
1964:         DMPlexSetCone(*dm, 0, cone);
1965:         ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0;
1966:         DMPlexSetConeOrientation(*dm, 0, ornt);
1967:         /* Cell 1 */
1968:         cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20;
1969:         DMPlexSetCone(*dm, 1, cone);
1970:         ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
1971:         DMPlexSetConeOrientation(*dm, 1, ornt);
1972:         /* Cell 2 */
1973:         cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23;
1974:         DMPlexSetCone(*dm, 2, cone);
1975:         ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
1976:         DMPlexSetConeOrientation(*dm, 2, ornt);
1977:         /* Cell 3 */
1978:         cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15;
1979:         DMPlexSetCone(*dm, 3, cone);
1980:         ornt[0] = -2; ornt[1] = -2; ornt[2] = 0; ornt[3] = -2;
1981:         DMPlexSetConeOrientation(*dm, 3, ornt);
1982:         /* Cell 4 */
1983:         cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25;
1984:         DMPlexSetCone(*dm, 4, cone);
1985:         ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = 0;
1986:         DMPlexSetConeOrientation(*dm, 4, ornt);
1987:         /* Cell 5 */
1988:         cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23;
1989:         DMPlexSetCone(*dm, 5, cone);
1990:         ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = -2;
1991:         DMPlexSetConeOrientation(*dm, 5, ornt);
1992:         /* Edges */
1993:         cone[0] =  6; cone[1] =  7;
1994:         DMPlexSetCone(*dm, 14, cone);
1995:         cone[0] =  7; cone[1] =  8;
1996:         DMPlexSetCone(*dm, 15, cone);
1997:         cone[0] =  8; cone[1] =  9;
1998:         DMPlexSetCone(*dm, 16, cone);
1999:         cone[0] =  9; cone[1] =  6;
2000:         DMPlexSetCone(*dm, 17, cone);
2001:         cone[0] = 10; cone[1] = 11;
2002:         DMPlexSetCone(*dm, 18, cone);
2003:         cone[0] = 11; cone[1] =  7;
2004:         DMPlexSetCone(*dm, 19, cone);
2005:         cone[0] =  6; cone[1] = 10;
2006:         DMPlexSetCone(*dm, 20, cone);
2007:         cone[0] = 12; cone[1] = 13;
2008:         DMPlexSetCone(*dm, 21, cone);
2009:         cone[0] = 13; cone[1] = 11;
2010:         DMPlexSetCone(*dm, 22, cone);
2011:         cone[0] = 10; cone[1] = 12;
2012:         DMPlexSetCone(*dm, 23, cone);
2013:         cone[0] = 13; cone[1] =  8;
2014:         DMPlexSetCone(*dm, 24, cone);
2015:         cone[0] = 12; cone[1] =  9;
2016:         DMPlexSetCone(*dm, 25, cone);
2017:       }
2018:       DMPlexSymmetrize(*dm);
2019:       DMPlexStratify(*dm);
2020:       /* Build coordinates */
2021:       PetscCalloc1(numVerts * embedDim, &coordsIn);
2022:       if (!rank) {
2023:         coordsIn[0*embedDim+0] = -R; coordsIn[0*embedDim+1] =  R; coordsIn[0*embedDim+2] = -R;
2024:         coordsIn[1*embedDim+0] =  R; coordsIn[1*embedDim+1] =  R; coordsIn[1*embedDim+2] = -R;
2025:         coordsIn[2*embedDim+0] =  R; coordsIn[2*embedDim+1] = -R; coordsIn[2*embedDim+2] = -R;
2026:         coordsIn[3*embedDim+0] = -R; coordsIn[3*embedDim+1] = -R; coordsIn[3*embedDim+2] = -R;
2027:         coordsIn[4*embedDim+0] = -R; coordsIn[4*embedDim+1] =  R; coordsIn[4*embedDim+2] =  R;
2028:         coordsIn[5*embedDim+0] =  R; coordsIn[5*embedDim+1] =  R; coordsIn[5*embedDim+2] =  R;
2029:         coordsIn[6*embedDim+0] = -R; coordsIn[6*embedDim+1] = -R; coordsIn[6*embedDim+2] =  R;
2030:         coordsIn[7*embedDim+0] =  R; coordsIn[7*embedDim+1] = -R; coordsIn[7*embedDim+2] =  R;
2031:       }
2032:     }
2033:     break;
2034:   case 3:
2035:     if (simplex) {
2036:       DM              idm;
2037:       const PetscReal edgeLen         = 1.0/PETSC_PHI;
2038:       PetscReal       vertexA[4]      = {0.5, 0.5, 0.5, 0.5};
2039:       PetscReal       vertexB[4]      = {1.0, 0.0, 0.0, 0.0};
2040:       PetscReal       vertexC[4]      = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0};
2041:       const PetscInt  degree          = 12;
2042:       PetscInt        s[4]            = {1, 1, 1};
2043:       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},
2044:                                          {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}};
2045:       PetscInt        cone[4];
2046:       PetscInt       *graph, p, i, j, k, l;

2048:       vertexA[0] *= R; vertexA[1] *= R; vertexA[2] *= R; vertexA[3] *= R;
2049:       vertexB[0] *= R; vertexB[1] *= R; vertexB[2] *= R; vertexB[3] *= R;
2050:       vertexC[0] *= R; vertexC[1] *= R; vertexC[2] *= R; vertexC[3] *= R;
2051:       numCells    = !rank ? 600 : 0;
2052:       numVerts    = !rank ? 120 : 0;
2053:       firstVertex = numCells;
2054:       /* Use the 600-cell, which for a unit sphere has coordinates which are

2056:            1/2 (\pm 1, \pm 1,    \pm 1, \pm 1)                          16
2057:                (\pm 1,    0,       0,      0)  all cyclic permutations   8
2058:            1/2 (\pm 1, \pm phi, \pm 1/phi, 0)  all even permutations    96

2060:          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
2061:          length is then given by 1/\phi = 0.61803.

2063:          http://buzzard.pugetsound.edu/sage-practice/ch03s03.html
2064:          http://mathworld.wolfram.com/600-Cell.html
2065:       */
2066:       /* Construct vertices */
2067:       PetscCalloc1(numVerts * embedDim, &coordsIn);
2068:       i    = 0;
2069:       if (!rank) {
2070:         for (s[0] = -1; s[0] < 2; s[0] += 2) {
2071:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
2072:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
2073:               for (s[3] = -1; s[3] < 2; s[3] += 2) {
2074:                 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d];
2075:                 ++i;
2076:               }
2077:             }
2078:           }
2079:         }
2080:         for (p = 0; p < embedDim; ++p) {
2081:           s[1] = s[2] = s[3] = 1;
2082:           for (s[0] = -1; s[0] < 2; s[0] += 2) {
2083:             for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim];
2084:             ++i;
2085:           }
2086:         }
2087:         for (p = 0; p < 12; ++p) {
2088:           s[3] = 1;
2089:           for (s[0] = -1; s[0] < 2; s[0] += 2) {
2090:             for (s[1] = -1; s[1] < 2; s[1] += 2) {
2091:               for (s[2] = -1; s[2] < 2; s[2] += 2) {
2092:                 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]];
2093:                 ++i;
2094:               }
2095:             }
2096:           }
2097:         }
2098:       }
2099:       if (i != numVerts) SETERRQ2(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertices %D != %D", i, numVerts);
2100:       /* Construct graph */
2101:       PetscCalloc1(numVerts * numVerts, &graph);
2102:       for (i = 0; i < numVerts; ++i) {
2103:         for (j = 0, k = 0; j < numVerts; ++j) {
2104:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
2105:         }
2106:         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertex %D degree %D != %D", i, k, degree);
2107:       }
2108:       /* Build Topology */
2109:       DMPlexSetChart(*dm, 0, numCells+numVerts);
2110:       for (c = 0; c < numCells; c++) {
2111:         DMPlexSetConeSize(*dm, c, embedDim);
2112:       }
2113:       DMSetUp(*dm); /* Allocate space for cones */
2114:       /* Cells */
2115:       if (!rank) {
2116:         for (i = 0, c = 0; i < numVerts; ++i) {
2117:           for (j = 0; j < i; ++j) {
2118:             for (k = 0; k < j; ++k) {
2119:               for (l = 0; l < k; ++l) {
2120:                 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] &&
2121:                     graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) {
2122:                   cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l;
2123:                   /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */
2124:                   {
2125:                     const PetscInt epsilon[4][4][4][4] = {{{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2126:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  1}, { 0,  0, -1, 0}},
2127:                                                            {{0,  0,  0,  0}, { 0, 0,  0, -1}, { 0,  0, 0,  0}, { 0,  1,  0, 0}},
2128:                                                            {{0,  0,  0,  0}, { 0, 0,  1,  0}, { 0, -1, 0,  0}, { 0,  0,  0, 0}}},

2130:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0, -1}, { 0,  0,  1, 0}},
2131:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2132:                                                            {{0,  0,  0,  1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, {-1,  0,  0, 0}},
2133:                                                            {{0,  0, -1,  0}, { 0, 0,  0,  0}, { 1,  0, 0,  0}, { 0,  0,  0, 0}}},

2135:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  1}, { 0,  0, 0,  0}, { 0, -1,  0, 0}},
2136:                                                            {{0,  0,  0, -1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 1,  0,  0, 0}},
2137:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2138:                                                            {{0,  1,  0,  0}, {-1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}},

2140:                                                           {{{0,  0,  0,  0}, { 0, 0, -1,  0}, { 0,  1, 0,  0}, { 0,  0,  0, 0}},
2141:                                                            {{0,  0,  1,  0}, { 0, 0,  0,  0}, {-1,  0, 0,  0}, { 0,  0,  0, 0}},
2142:                                                            {{0, -1,  0,  0}, { 1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2143:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}}};
2144:                     PetscReal normal[4];
2145:                     PetscInt  e, f, g;

2147:                     for (d = 0; d < embedDim; ++d) {
2148:                       normal[d] = 0.0;
2149:                       for (e = 0; e < embedDim; ++e) {
2150:                         for (f = 0; f < embedDim; ++f) {
2151:                           for (g = 0; g < embedDim; ++g) {
2152:                             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]);
2153:                           }
2154:                         }
2155:                       }
2156:                     }
2157:                     if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
2158:                   }
2159:                   DMPlexSetCone(*dm, c++, cone);
2160:                 }
2161:               }
2162:             }
2163:           }
2164:         }
2165:       }
2166:       DMPlexSymmetrize(*dm);
2167:       DMPlexStratify(*dm);
2168:       PetscFree(graph);
2169:       /* Interpolate mesh */
2170:       DMPlexInterpolate(*dm, &idm);
2171:       DMDestroy(dm);
2172:       *dm  = idm;
2173:       break;
2174:     }
2175:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim);
2176:   }
2177:   /* Create coordinates */
2178:   DMGetCoordinateSection(*dm, &coordSection);
2179:   PetscSectionSetNumFields(coordSection, 1);
2180:   PetscSectionSetFieldComponents(coordSection, 0, embedDim);
2181:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);
2182:   for (v = firstVertex; v < firstVertex+numVerts; ++v) {
2183:     PetscSectionSetDof(coordSection, v, embedDim);
2184:     PetscSectionSetFieldDof(coordSection, v, 0, embedDim);
2185:   }
2186:   PetscSectionSetUp(coordSection);
2187:   PetscSectionGetStorageSize(coordSection, &coordSize);
2188:   VecCreate(PETSC_COMM_SELF, &coordinates);
2189:   VecSetBlockSize(coordinates, embedDim);
2190:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2191:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2192:   VecSetType(coordinates,VECSTANDARD);
2193:   VecGetArray(coordinates, &coords);
2194:   for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];}
2195:   VecRestoreArray(coordinates, &coords);
2196:   DMSetCoordinatesLocal(*dm, coordinates);
2197:   VecDestroy(&coordinates);
2198:   PetscFree(coordsIn);
2199:   /* Create coordinate function space */
2200:   {
2201:     DM          cdm;
2202:     PetscDS     cds;
2203:     PetscFE     fe;
2204:     PetscScalar radius = R;
2205:     PetscInt    dT, dE;

2207:     DMGetCoordinateDM(*dm, &cdm);
2208:     DMGetDimension(*dm, &dT);
2209:     DMGetCoordinateDim(*dm, &dE);
2210:     PetscFECreateLagrange(PETSC_COMM_SELF, dT, dE, simplex, 1, -1, &fe);
2211:     DMSetField(cdm, 0, NULL, (PetscObject) fe);
2212:     PetscFEDestroy(&fe);
2213:     DMCreateDS(cdm);

2215:     DMGetDS(cdm, &cds);
2216:     PetscDSSetConstants(cds, 1, &radius);
2217:   }
2218:   ((DM_Plex *) (*dm)->data)->coordFunc = snapToSphere;
2219:   return(0);
2220: }

2222: /*@
2223:   DMPlexCreateBallMesh - Creates a simplex mesh on the d-dimensional ball, B^d.

2225:   Collective

2227:   Input Parameters:
2228: + comm  - The communicator for the DM object
2229: . dim   - The dimension
2230: - R     - The radius

2232:   Output Parameter:
2233: . dm  - The DM object

2235:   Options Database Keys:
2236: - bd_dm_refine - This will refine the surface mesh preserving the sphere geometry

2238:   Level: beginner

2240: .seealso: DMPlexCreateSphereMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
2241: @*/
2242: PetscErrorCode DMPlexCreateBallMesh(MPI_Comm comm, PetscInt dim, PetscReal R, DM *dm)
2243: {
2244:   DM             sdm;
2245:   DMLabel        bdlabel;

2249:   DMPlexCreateSphereMesh(comm, dim-1, PETSC_TRUE, R, &sdm);
2250:   PetscObjectSetOptionsPrefix((PetscObject) sdm, "bd_");
2251:   DMSetFromOptions(sdm);
2252:   DMPlexGenerate(sdm, NULL, PETSC_TRUE, dm);
2253:   DMDestroy(&sdm);
2254:   DMCreateLabel(*dm, "marker");
2255:   DMGetLabel(*dm, "marker", &bdlabel);
2256:   DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
2257:   DMPlexLabelComplete(*dm, bdlabel);
2258:   return(0);
2259: }

2261: /* External function declarations here */
2262: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
2263: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2264: extern PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2265: extern PetscErrorCode DMCreateLocalSection_Plex(DM dm);
2266: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
2267: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
2268: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
2269: extern PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
2270: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
2271: extern PetscErrorCode DMSetUp_Plex(DM dm);
2272: extern PetscErrorCode DMDestroy_Plex(DM dm);
2273: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
2274: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
2275: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
2276: extern PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);
2277: static PetscErrorCode DMInitialize_Plex(DM dm);

2279: /* Replace dm with the contents of dmNew
2280:    - Share the DM_Plex structure
2281:    - Share the coordinates
2282:    - Share the SF
2283: */
2284: static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
2285: {
2286:   PetscSF               sf;
2287:   DM                    coordDM, coarseDM;
2288:   Vec                   coords;
2289:   PetscBool             isper;
2290:   const PetscReal      *maxCell, *L;
2291:   const DMBoundaryType *bd;
2292:   PetscErrorCode        ierr;

2295:   DMGetPointSF(dmNew, &sf);
2296:   DMSetPointSF(dm, sf);
2297:   DMGetCoordinateDM(dmNew, &coordDM);
2298:   DMGetCoordinatesLocal(dmNew, &coords);
2299:   DMSetCoordinateDM(dm, coordDM);
2300:   DMSetCoordinatesLocal(dm, coords);
2301:   /* Do not want to create the coordinate field if it does not already exist, so do not call DMGetCoordinateField() */
2302:   DMFieldDestroy(&dm->coordinateField);
2303:   dm->coordinateField = dmNew->coordinateField;
2304:   DMGetPeriodicity(dmNew, &isper, &maxCell, &L, &bd);
2305:   DMSetPeriodicity(dm, isper, maxCell, L, bd);
2306:   DMDestroy_Plex(dm);
2307:   DMInitialize_Plex(dm);
2308:   dm->data = dmNew->data;
2309:   ((DM_Plex *) dmNew->data)->refct++;
2310:   DMDestroyLabelLinkList_Internal(dm);
2311:   DMCopyLabels(dmNew, dm, PETSC_OWN_POINTER, PETSC_TRUE);
2312:   DMGetCoarseDM(dmNew,&coarseDM);
2313:   DMSetCoarseDM(dm,coarseDM);
2314:   return(0);
2315: }

2317: /* Swap dm with the contents of dmNew
2318:    - Swap the DM_Plex structure
2319:    - Swap the coordinates
2320:    - Swap the point PetscSF
2321: */
2322: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
2323: {
2324:   DM              coordDMA, coordDMB;
2325:   Vec             coordsA,  coordsB;
2326:   PetscSF         sfA,      sfB;
2327:   DMField         fieldTmp;
2328:   void            *tmp;
2329:   DMLabelLink     listTmp;
2330:   DMLabel         depthTmp;
2331:   PetscInt        tmpI;
2332:   PetscErrorCode  ierr;

2335:   DMGetPointSF(dmA, &sfA);
2336:   DMGetPointSF(dmB, &sfB);
2337:   PetscObjectReference((PetscObject) sfA);
2338:   DMSetPointSF(dmA, sfB);
2339:   DMSetPointSF(dmB, sfA);
2340:   PetscObjectDereference((PetscObject) sfA);

2342:   DMGetCoordinateDM(dmA, &coordDMA);
2343:   DMGetCoordinateDM(dmB, &coordDMB);
2344:   PetscObjectReference((PetscObject) coordDMA);
2345:   DMSetCoordinateDM(dmA, coordDMB);
2346:   DMSetCoordinateDM(dmB, coordDMA);
2347:   PetscObjectDereference((PetscObject) coordDMA);

2349:   DMGetCoordinatesLocal(dmA, &coordsA);
2350:   DMGetCoordinatesLocal(dmB, &coordsB);
2351:   PetscObjectReference((PetscObject) coordsA);
2352:   DMSetCoordinatesLocal(dmA, coordsB);
2353:   DMSetCoordinatesLocal(dmB, coordsA);
2354:   PetscObjectDereference((PetscObject) coordsA);

2356:   fieldTmp             = dmA->coordinateField;
2357:   dmA->coordinateField = dmB->coordinateField;
2358:   dmB->coordinateField = fieldTmp;
2359:   tmp       = dmA->data;
2360:   dmA->data = dmB->data;
2361:   dmB->data = tmp;
2362:   listTmp   = dmA->labels;
2363:   dmA->labels = dmB->labels;
2364:   dmB->labels = listTmp;
2365:   depthTmp  = dmA->depthLabel;
2366:   dmA->depthLabel = dmB->depthLabel;
2367:   dmB->depthLabel = depthTmp;
2368:   depthTmp  = dmA->celltypeLabel;
2369:   dmA->celltypeLabel = dmB->celltypeLabel;
2370:   dmB->celltypeLabel = depthTmp;
2371:   tmpI         = dmA->levelup;
2372:   dmA->levelup = dmB->levelup;
2373:   dmB->levelup = tmpI;
2374:   return(0);
2375: }

2377: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2378: {
2379:   DM_Plex       *mesh = (DM_Plex*) dm->data;
2380:   PetscBool      flg;

2384:   /* Handle viewing */
2385:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL);
2386:   PetscOptionsBoundedInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL,0);
2387:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL);
2388:   PetscOptionsBoundedInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL,0);
2389:   DMMonitorSetFromOptions(dm, "-dm_plex_monitor_throughput", "Monitor the simulation throughput", "DMPlexMonitorThroughput", DMPlexMonitorThroughput, NULL, &flg);
2390:   if (flg) {PetscLogDefaultBegin();}
2391:   /* Point Location */
2392:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL);
2393:   /* Partitioning and distribution */
2394:   PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL);
2395:   /* Generation and remeshing */
2396:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL);
2397:   /* Projection behavior */
2398:   PetscOptionsBoundedInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL,0);
2399:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);
2400:   PetscOptionsEnum("-dm_plex_cell_refiner", "Strategy for cell refinment", "ex40.c", DMPlexCellRefinerTypes, (PetscEnum) mesh->cellRefiner, (PetscEnum *) &mesh->cellRefiner, NULL);
2401:   /* Checking structure */
2402:   {
2403:     PetscBool   flg = PETSC_FALSE, flg2 = PETSC_FALSE, all = PETSC_FALSE;

2405:     PetscOptionsBool("-dm_plex_check_all", "Perform all checks", NULL, PETSC_FALSE, &all, &flg2);
2406:     PetscOptionsBool("-dm_plex_check_symmetry", "Check that the adjacency information in the mesh is symmetric", "DMPlexCheckSymmetry", PETSC_FALSE, &flg, &flg2);
2407:     if (all || (flg && flg2)) {DMPlexCheckSymmetry(dm);}
2408:     PetscOptionsBool("-dm_plex_check_skeleton", "Check that each cell has the correct number of vertices (only for homogeneous simplex or tensor meshes)", "DMPlexCheckSkeleton", PETSC_FALSE, &flg, &flg2);
2409:     if (all || (flg && flg2)) {DMPlexCheckSkeleton(dm, 0);}
2410:     PetscOptionsBool("-dm_plex_check_faces", "Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type", "DMPlexCheckFaces", PETSC_FALSE, &flg, &flg2);
2411:     if (all || (flg && flg2)) {DMPlexCheckFaces(dm, 0);}
2412:     PetscOptionsBool("-dm_plex_check_geometry", "Check that cells have positive volume", "DMPlexCheckGeometry", PETSC_FALSE, &flg, &flg2);
2413:     if (all || (flg && flg2)) {DMPlexCheckGeometry(dm);}
2414:     PetscOptionsBool("-dm_plex_check_pointsf", "Check some necessary conditions for PointSF", "DMPlexCheckPointSF", PETSC_FALSE, &flg, &flg2);
2415:     if (all || (flg && flg2)) {DMPlexCheckPointSF(dm);}
2416:     PetscOptionsBool("-dm_plex_check_interface_cones", "Check points on inter-partition interfaces have conforming order of cone points", "DMPlexCheckInterfaceCones", PETSC_FALSE, &flg, &flg2);
2417:     if (all || (flg && flg2)) {DMPlexCheckInterfaceCones(dm);}
2418:     PetscOptionsBool("-dm_plex_check_cell_shape", "Check cell shape", "DMPlexCheckCellShape", PETSC_FALSE, &flg, &flg2);
2419:     if (flg && flg2) {DMPlexCheckCellShape(dm, PETSC_TRUE, PETSC_DETERMINE);}
2420:   }

2422:   PetscPartitionerSetFromOptions(mesh->partitioner);
2423:   return(0);
2424: }

2426: static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2427: {
2428:   PetscReal      volume = -1.0;
2429:   PetscInt       prerefine = 0, refine = 0, r, coarsen = 0, overlap = 0;
2430:   PetscBool      uniformOrig, uniform = PETSC_TRUE, distribute = PETSC_FALSE, isHierarchy, flg;

2435:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
2436:   /* Handle DMPlex refinement before distribution */
2437:   DMPlexGetRefinementUniform(dm, &uniformOrig);
2438:   PetscOptionsBool("-dm_refine_uniform_pre", "Flag for uniform refinement before distribution", "DMCreate", uniform, &uniform, &flg);
2439:   if (flg) {DMPlexSetRefinementUniform(dm, uniform);}
2440:   PetscOptionsReal("-dm_refine_volume_limit_pre", "The maximum cell volume after refinement before distribution", "DMCreate", volume, &volume, &flg);
2441:   if (flg) {DMPlexSetRefinementLimit(dm, volume);}
2442:   PetscOptionsBoundedInt("-dm_refine_pre", "The number of refinements before distribution", "DMCreate", prerefine, &prerefine, NULL,0);
2443:   for (r = 0; r < prerefine; ++r) {
2444:     DM             rdm;
2445:     PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

2447:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2448:     DMRefine(dm, PetscObjectComm((PetscObject) dm), &rdm);
2449:     /* Total hack since we do not pass in a pointer */
2450:     DMPlexReplace_Static(dm, rdm);
2451:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2452:     if (coordFunc) {
2453:       DMPlexRemapGeometry(dm, 0.0, coordFunc);
2454:       ((DM_Plex*) dm->data)->coordFunc = coordFunc;
2455:     }
2456:     DMDestroy(&rdm);
2457:   }
2458:   DMPlexSetRefinementUniform(dm, uniformOrig);
2459:   /* Handle DMPlex distribution */
2460:   PetscOptionsBool("-dm_distribute", "Flag to redistribute a mesh among processes", "DMCreate", distribute, &distribute, NULL);
2461:   PetscOptionsBoundedInt("-dm_distribute_overlap", "The size of the overlap halo", "DMCreate", overlap, &overlap, NULL, 0);
2462:   if (distribute) {
2463:     DM               pdm = NULL;
2464:     PetscPartitioner part;

2466:     DMPlexGetPartitioner(dm, &part);
2467:     PetscPartitionerSetFromOptions(part);
2468:     DMPlexDistribute(dm, overlap, NULL, &pdm);
2469:     if (pdm) {
2470:       DMPlexReplace_Static(dm, pdm);
2471:       DMDestroy(&pdm);
2472:     }
2473:   }
2474:   /* Handle DMPlex refinement */
2475:   PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL,0);
2476:   PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy,0);
2477:   if (refine) {DMPlexSetRefinementUniform(dm, PETSC_TRUE);}
2478:   if (refine && isHierarchy) {
2479:     DM *dms, coarseDM;

2481:     DMGetCoarseDM(dm, &coarseDM);
2482:     PetscObjectReference((PetscObject)coarseDM);
2483:     PetscMalloc1(refine,&dms);
2484:     DMRefineHierarchy(dm, refine, dms);
2485:     /* Total hack since we do not pass in a pointer */
2486:     DMPlexSwap_Static(dm, dms[refine-1]);
2487:     if (refine == 1) {
2488:       DMSetCoarseDM(dm, dms[0]);
2489:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2490:     } else {
2491:       DMSetCoarseDM(dm, dms[refine-2]);
2492:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2493:       DMSetCoarseDM(dms[0], dms[refine-1]);
2494:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
2495:     }
2496:     DMSetCoarseDM(dms[refine-1], coarseDM);
2497:     PetscObjectDereference((PetscObject)coarseDM);
2498:     /* Free DMs */
2499:     for (r = 0; r < refine; ++r) {
2500:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2501:       DMDestroy(&dms[r]);
2502:     }
2503:     PetscFree(dms);
2504:   } else {
2505:     for (r = 0; r < refine; ++r) {
2506:       DM refinedMesh;
2507:       PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

2509:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2510:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);
2511:       /* Total hack since we do not pass in a pointer */
2512:       DMPlexReplace_Static(dm, refinedMesh);
2513:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2514:       if (coordFunc) {
2515:         DMPlexRemapGeometry(dm, 0.0, coordFunc);
2516:         ((DM_Plex*) dm->data)->coordFunc = coordFunc;
2517:       }
2518:       DMDestroy(&refinedMesh);
2519:     }
2520:   }
2521:   /* Handle DMPlex coarsening */
2522:   PetscOptionsBoundedInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL,0);
2523:   PetscOptionsBoundedInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy,0);
2524:   if (coarsen && isHierarchy) {
2525:     DM *dms;

2527:     PetscMalloc1(coarsen, &dms);
2528:     DMCoarsenHierarchy(dm, coarsen, dms);
2529:     /* Free DMs */
2530:     for (r = 0; r < coarsen; ++r) {
2531:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2532:       DMDestroy(&dms[r]);
2533:     }
2534:     PetscFree(dms);
2535:   } else {
2536:     for (r = 0; r < coarsen; ++r) {
2537:       DM coarseMesh;

2539:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2540:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);
2541:       /* Total hack since we do not pass in a pointer */
2542:       DMPlexReplace_Static(dm, coarseMesh);
2543:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2544:       DMDestroy(&coarseMesh);
2545:     }
2546:   }
2547:   /* Handle */
2548:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2549:   PetscOptionsTail();
2550:   return(0);
2551: }

2553: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
2554: {

2558:   DMCreateGlobalVector_Section_Private(dm,vec);
2559:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
2560:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
2561:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
2562:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
2563:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
2564:   return(0);
2565: }

2567: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
2568: {

2572:   DMCreateLocalVector_Section_Private(dm,vec);
2573:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
2574:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
2575:   return(0);
2576: }

2578: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
2579: {
2580:   PetscInt       depth, d;

2584:   DMPlexGetDepth(dm, &depth);
2585:   if (depth == 1) {
2586:     DMGetDimension(dm, &d);
2587:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
2588:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
2589:     else               {*pStart = 0; *pEnd = 0;}
2590:   } else {
2591:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
2592:   }
2593:   return(0);
2594: }

2596: static PetscErrorCode DMGetNeighbors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
2597: {
2598:   PetscSF           sf;
2599:   PetscInt          niranks, njranks, n;
2600:   const PetscMPIInt *iranks, *jranks;
2601:   DM_Plex           *data = (DM_Plex*) dm->data;
2602:   PetscErrorCode    ierr;

2605:   DMGetPointSF(dm, &sf);
2606:   if (!data->neighbors) {
2607:     PetscSFSetUp(sf);
2608:     PetscSFGetRootRanks(sf, &njranks, &jranks, NULL, NULL, NULL);
2609:     PetscSFGetLeafRanks(sf, &niranks, &iranks, NULL, NULL);
2610:     PetscMalloc1(njranks + niranks + 1, &data->neighbors);
2611:     PetscArraycpy(data->neighbors + 1, jranks, njranks);
2612:     PetscArraycpy(data->neighbors + njranks + 1, iranks, niranks);
2613:     n = njranks + niranks;
2614:     PetscSortRemoveDupsMPIInt(&n, data->neighbors + 1);
2615:     /* The following cast should never fail: can't have more neighbors than PETSC_MPI_INT_MAX */
2616:     PetscMPIIntCast(n, data->neighbors);
2617:   }
2618:   if (nranks) *nranks = data->neighbors[0];
2619:   if (ranks) {
2620:     if (data->neighbors[0]) *ranks = data->neighbors + 1;
2621:     else                    *ranks = NULL;
2622:   }
2623:   return(0);
2624: }

2626: PETSC_INTERN PetscErrorCode DMInterpolateSolution_Plex(DM, DM, Mat, Vec, Vec);

2628: static PetscErrorCode DMInitialize_Plex(DM dm)
2629: {

2633:   dm->ops->view                            = DMView_Plex;
2634:   dm->ops->load                            = DMLoad_Plex;
2635:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
2636:   dm->ops->clone                           = DMClone_Plex;
2637:   dm->ops->setup                           = DMSetUp_Plex;
2638:   dm->ops->createlocalsection              = DMCreateLocalSection_Plex;
2639:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
2640:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
2641:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
2642:   dm->ops->getlocaltoglobalmapping         = NULL;
2643:   dm->ops->createfieldis                   = NULL;
2644:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
2645:   dm->ops->createcoordinatefield           = DMCreateCoordinateField_Plex;
2646:   dm->ops->getcoloring                     = NULL;
2647:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
2648:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
2649:   dm->ops->createmassmatrix                = DMCreateMassMatrix_Plex;
2650:   dm->ops->createinjection                 = DMCreateInjection_Plex;
2651:   dm->ops->refine                          = DMRefine_Plex;
2652:   dm->ops->coarsen                         = DMCoarsen_Plex;
2653:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
2654:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
2655:   dm->ops->adaptlabel                      = DMAdaptLabel_Plex;
2656:   dm->ops->adaptmetric                     = DMAdaptMetric_Plex;
2657:   dm->ops->globaltolocalbegin              = NULL;
2658:   dm->ops->globaltolocalend                = NULL;
2659:   dm->ops->localtoglobalbegin              = NULL;
2660:   dm->ops->localtoglobalend                = NULL;
2661:   dm->ops->destroy                         = DMDestroy_Plex;
2662:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
2663:   dm->ops->createsuperdm                   = DMCreateSuperDM_Plex;
2664:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
2665:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
2666:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
2667:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
2668:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
2669:   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
2670:   dm->ops->projectbdfieldlabellocal        = DMProjectBdFieldLabelLocal_Plex;
2671:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
2672:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
2673:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
2674:   dm->ops->getneighbors                    = DMGetNeighbors_Plex;
2675:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);
2676:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertTimeDerviativeBoundaryValues_C",DMPlexInsertTimeDerivativeBoundaryValues_Plex);
2677:   PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);
2678:   PetscObjectComposeFunction((PetscObject)dm,"DMCreateNeumannOverlap_C",DMCreateNeumannOverlap_Plex);
2679:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexGetOverlap_C",DMPlexGetOverlap_Plex);
2680:   PetscObjectComposeFunction((PetscObject)dm,"DMInterpolateSolution_C",DMInterpolateSolution_Plex);
2681:   return(0);
2682: }

2684: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
2685: {
2686:   DM_Plex        *mesh = (DM_Plex *) dm->data;

2690:   mesh->refct++;
2691:   (*newdm)->data = mesh;
2692:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
2693:   DMInitialize_Plex(*newdm);
2694:   return(0);
2695: }

2697: /*MC
2698:   DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
2699:                     In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is
2700:                     specified by a PetscSection object. Ownership in the global representation is determined by
2701:                     ownership of the underlying DMPlex points. This is specified by another PetscSection object.

2703:   Options Database Keys:
2704: + -dm_refine_pre                     - Refine mesh before distribution
2705: + -dm_refine_uniform_pre             - Choose uniform or generator-based refinement
2706: + -dm_refine_volume_limit_pre        - Cell volume limit after pre-refinement using generator
2707: . -dm_distribute                     - Distribute mesh across processes
2708: . -dm_distribute_overlap             - Number of cells to overlap for distribution
2709: . -dm_refine                         - Refine mesh after distribution
2710: . -dm_plex_hash_location             - Use grid hashing for point location
2711: . -dm_plex_partition_balance         - Attempt to evenly divide points on partition boundary between processes
2712: . -dm_plex_remesh_bd                 - Allow changes to the boundary on remeshing
2713: . -dm_plex_max_projection_height     - Maxmimum mesh point height used to project locally
2714: . -dm_plex_regular_refinement        - Use special nested projection algorithm for regular refinement
2715: . -dm_plex_check_all                 - Perform all shecks below
2716: . -dm_plex_check_symmetry            - Check that the adjacency information in the mesh is symmetric
2717: . -dm_plex_check_skeleton <celltype> - Check that each cell has the correct number of vertices
2718: . -dm_plex_check_faces <celltype>    - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type
2719: . -dm_plex_check_geometry            - Check that cells have positive volume
2720: . -dm_view :mesh.tex:ascii_latex     - View the mesh in LaTeX/TikZ
2721: . -dm_plex_view_scale <num>          - Scale the TikZ
2722: - -dm_plex_print_fem <num>           - View FEM assembly information, such as element vectors and matrices


2725:   Level: intermediate

2727: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
2728: M*/

2730: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
2731: {
2732:   DM_Plex       *mesh;
2733:   PetscInt       unit;

2738:   PetscNewLog(dm,&mesh);
2739:   dm->dim  = 0;
2740:   dm->data = mesh;

2742:   mesh->refct             = 1;
2743:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
2744:   mesh->maxConeSize       = 0;
2745:   mesh->cones             = NULL;
2746:   mesh->coneOrientations  = NULL;
2747:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
2748:   mesh->maxSupportSize    = 0;
2749:   mesh->supports          = NULL;
2750:   mesh->refinementUniform = PETSC_TRUE;
2751:   mesh->refinementLimit   = -1.0;
2752:   mesh->interpolated      = DMPLEX_INTERPOLATED_INVALID;
2753:   mesh->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID;

2755:   mesh->facesTmp = NULL;

2757:   mesh->tetgenOpts   = NULL;
2758:   mesh->triangleOpts = NULL;
2759:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
2760:   mesh->remeshBd     = PETSC_FALSE;

2762:   mesh->subpointMap = NULL;

2764:   for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;

2766:   mesh->regularRefinement   = PETSC_FALSE;
2767:   mesh->depthState          = -1;
2768:   mesh->celltypeState       = -1;
2769:   mesh->globalVertexNumbers = NULL;
2770:   mesh->globalCellNumbers   = NULL;
2771:   mesh->anchorSection       = NULL;
2772:   mesh->anchorIS            = NULL;
2773:   mesh->createanchors       = NULL;
2774:   mesh->computeanchormatrix = NULL;
2775:   mesh->parentSection       = NULL;
2776:   mesh->parents             = NULL;
2777:   mesh->childIDs            = NULL;
2778:   mesh->childSection        = NULL;
2779:   mesh->children            = NULL;
2780:   mesh->referenceTree       = NULL;
2781:   mesh->getchildsymmetry    = NULL;
2782:   mesh->vtkCellHeight       = 0;
2783:   mesh->useAnchors          = PETSC_FALSE;

2785:   mesh->maxProjectionHeight = 0;

2787:   mesh->neighbors           = NULL;

2789:   mesh->printSetValues = PETSC_FALSE;
2790:   mesh->printFEM       = 0;
2791:   mesh->printTol       = 1.0e-10;

2793:   DMInitialize_Plex(dm);
2794:   return(0);
2795: }

2797: /*@
2798:   DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.

2800:   Collective

2802:   Input Parameter:
2803: . comm - The communicator for the DMPlex object

2805:   Output Parameter:
2806: . mesh  - The DMPlex object

2808:   Level: beginner

2810: @*/
2811: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
2812: {

2817:   DMCreate(comm, mesh);
2818:   DMSetType(*mesh, DMPLEX);
2819:   return(0);
2820: }

2822: /*@C
2823:   DMPlexBuildFromCellListParallel - Build distributed DMPLEX topology from a list of vertices for each cell (common mesh generator output)

2825:   Input Parameters:
2826: + dm - The DM
2827: . numCells - The number of cells owned by this process
2828: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
2829: . NVertices - The global number of vertices, or PETSC_DECIDE
2830: . numCorners - The number of vertices for each cell
2831: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell

2833:   Output Parameter:
2834: . vertexSF - (Optional) SF describing complete vertex ownership

2836:   Notes:
2837:   Two triangles sharing a face
2838: $
2839: $        2
2840: $      / | \
2841: $     /  |  \
2842: $    /   |   \
2843: $   0  0 | 1  3
2844: $    \   |   /
2845: $     \  |  /
2846: $      \ | /
2847: $        1
2848: would have input
2849: $  numCells = 2, numVertices = 4
2850: $  cells = [0 1 2  1 3 2]
2851: $
2852: which would result in the DMPlex
2853: $
2854: $        4
2855: $      / | \
2856: $     /  |  \
2857: $    /   |   \
2858: $   2  0 | 1  5
2859: $    \   |   /
2860: $     \  |  /
2861: $      \ | /
2862: $        3

2864:   Vertices are implicitly numbered consecutively 0,...,NVertices.
2865:   Each rank owns a chunk of numVertices consecutive vertices.
2866:   If numVertices is PETSC_DECIDE, PETSc will distribute them as evenly as possible using PetscLayout.
2867:   If both NVertices and numVertices are PETSC_DECIDE, NVertices is computed by PETSc as the maximum vertex index in cells + 1.
2868:   If only NVertices is PETSC_DECIDE, it is computed as the sum of numVertices over all ranks.

2870:   The cell distribution is arbitrary non-overlapping, independent of the vertex distribution.

2872:   Not currently supported in Fortran.

2874:   Level: advanced

2876: .seealso: DMPlexBuildFromCellList(), DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildCoordinatesFromCellListParallel()
2877: @*/
2878: PetscErrorCode DMPlexBuildFromCellListParallel(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, const PetscInt cells[], PetscSF *vertexSF)
2879: {
2880:   PetscSF         sfPoint;
2881:   PetscLayout     layout;
2882:   PetscInt        numVerticesAdj, *verticesAdj, *cones, c, p, dim;
2883:   PetscMPIInt     rank, size;
2884:   PetscErrorCode  ierr;

2888:   PetscLogEventBegin(DMPLEX_BuildFromCellList,dm,0,0,0);
2889:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
2890:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);
2891:   DMGetDimension(dm, &dim);
2892:   /* Get/check global number of vertices */
2893:   {
2894:     PetscInt NVerticesInCells, i;
2895:     const PetscInt len = numCells * numCorners;

2897:     /* NVerticesInCells = max(cells) + 1 */
2898:     NVerticesInCells = PETSC_MIN_INT;
2899:     for (i=0; i<len; i++) if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
2900:     ++NVerticesInCells;
2901:     MPI_Allreduce(MPI_IN_PLACE, &NVerticesInCells, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));

2903:     if (numVertices == PETSC_DECIDE && NVertices == PETSC_DECIDE) NVertices = NVerticesInCells;
2904:     else if (NVertices != PETSC_DECIDE && NVertices < NVerticesInCells) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Specified global number of vertices %D must be greater than or equal to the number of vertices in cells %D",NVertices,NVerticesInCells);
2905:   }
2906:   /* Count locally unique vertices */
2907:   {
2908:     PetscHSetI vhash;
2909:     PetscInt off = 0;

2911:     PetscHSetICreate(&vhash);
2912:     for (c = 0; c < numCells; ++c) {
2913:       for (p = 0; p < numCorners; ++p) {
2914:         PetscHSetIAdd(vhash, cells[c*numCorners+p]);
2915:       }
2916:     }
2917:     PetscHSetIGetSize(vhash, &numVerticesAdj);
2918:     PetscMalloc1(numVerticesAdj, &verticesAdj);
2919:     PetscHSetIGetElems(vhash, &off, verticesAdj);
2920:     PetscHSetIDestroy(&vhash);
2921:     if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
2922:   }
2923:   PetscSortInt(numVerticesAdj, verticesAdj);
2924:   /* Create cones */
2925:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
2926:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
2927:   DMSetUp(dm);
2928:   DMPlexGetCones(dm,&cones);
2929:   for (c = 0; c < numCells; ++c) {
2930:     for (p = 0; p < numCorners; ++p) {
2931:       const PetscInt gv = cells[c*numCorners+p];
2932:       PetscInt       lv;

2934:       /* Positions within verticesAdj form 0-based local vertex numbering;
2935:          we need to shift it by numCells to get correct DAG points (cells go first) */
2936:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
2937:       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
2938:       cones[c*numCorners+p] = lv+numCells;
2939:     }
2940:   }
2941:   /* Build point sf */
2942:   PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &layout);
2943:   PetscLayoutSetSize(layout, NVertices);
2944:   PetscLayoutSetLocalSize(layout, numVertices);
2945:   PetscLayoutSetBlockSize(layout, 1);
2946:   PetscSFCreateByMatchingIndices(layout, numVerticesAdj, verticesAdj, NULL, numCells, numVerticesAdj, verticesAdj, NULL, numCells, vertexSF, &sfPoint);
2947:   PetscLayoutDestroy(&layout);
2948:   PetscFree(verticesAdj);
2949:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
2950:   if (dm->sf) {
2951:     const char *prefix;

2953:     PetscObjectGetOptionsPrefix((PetscObject)dm->sf, &prefix);
2954:     PetscObjectSetOptionsPrefix((PetscObject)sfPoint, prefix);
2955:   }
2956:   DMSetPointSF(dm, sfPoint);
2957:   PetscSFDestroy(&sfPoint);
2958:   if (vertexSF) {PetscObjectSetName((PetscObject)(*vertexSF), "Vertex Ownership SF");}
2959:   /* Fill in the rest of the topology structure */
2960:   DMPlexSymmetrize(dm);
2961:   DMPlexStratify(dm);
2962:   PetscLogEventEnd(DMPLEX_BuildFromCellList,dm,0,0,0);
2963:   return(0);
2964: }

2966: /*@C
2967:   DMPlexBuildCoordinatesFromCellListParallel - Build DM coordinates from a list of coordinates for each owned vertex (common mesh generator output)

2969:   Input Parameters:
2970: + dm - The DM
2971: . spaceDim - The spatial dimension used for coordinates
2972: . sfVert - SF describing complete vertex ownership
2973: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

2975:   Level: advanced

2977:   Notes:
2978:   Not currently supported in Fortran.

2980: .seealso: DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildFromCellListParallel()
2981: @*/
2982: PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM dm, PetscInt spaceDim, PetscSF sfVert, const PetscReal vertexCoords[])
2983: {
2984:   PetscSection   coordSection;
2985:   Vec            coordinates;
2986:   PetscScalar   *coords;
2987:   PetscInt       numVertices, numVerticesAdj, coordSize, v, vStart, vEnd;

2991:   PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
2992:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
2993:   if (vStart < 0 || vEnd < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first.");
2994:   DMSetCoordinateDim(dm, spaceDim);
2995:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
2996:   if (vEnd - vStart != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Supplied sfVert has wrong number of leaves = %D != %D = vEnd - vStart",numVerticesAdj,vEnd - vStart);
2997:   DMGetCoordinateSection(dm, &coordSection);
2998:   PetscSectionSetNumFields(coordSection, 1);
2999:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
3000:   PetscSectionSetChart(coordSection, vStart, vEnd);
3001:   for (v = vStart; v < vEnd; ++v) {
3002:     PetscSectionSetDof(coordSection, v, spaceDim);
3003:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
3004:   }
3005:   PetscSectionSetUp(coordSection);
3006:   PetscSectionGetStorageSize(coordSection, &coordSize);
3007:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
3008:   VecSetBlockSize(coordinates, spaceDim);
3009:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3010:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3011:   VecSetType(coordinates,VECSTANDARD);
3012:   VecGetArray(coordinates, &coords);
3013:   {
3014:     MPI_Datatype coordtype;

3016:     /* Need a temp buffer for coords if we have complex/single */
3017:     MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);
3018:     MPI_Type_commit(&coordtype);
3019: #if defined(PETSC_USE_COMPLEX)
3020:     {
3021:     PetscScalar *svertexCoords;
3022:     PetscInt    i;
3023:     PetscMalloc1(numVertices*spaceDim,&svertexCoords);
3024:     for (i=0; i<numVertices*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
3025:     PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords,MPI_REPLACE);
3026:     PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords,MPI_REPLACE);
3027:     PetscFree(svertexCoords);
3028:     }
3029: #else
3030:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords,MPI_REPLACE);
3031:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords,MPI_REPLACE);
3032: #endif
3033:     MPI_Type_free(&coordtype);
3034:   }
3035:   VecRestoreArray(coordinates, &coords);
3036:   DMSetCoordinatesLocal(dm, coordinates);
3037:   VecDestroy(&coordinates);
3038:   PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3039:   return(0);
3040: }

3042: /*@
3043:   DMPlexCreateFromCellListParallelPetsc - Create distributed DMPLEX from a list of vertices for each cell (common mesh generator output)

3045:   Input Parameters:
3046: + comm - The communicator
3047: . dim - The topological dimension of the mesh
3048: . numCells - The number of cells owned by this process
3049: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
3050: . NVertices - The global number of vertices, or PETSC_DECIDE
3051: . numCorners - The number of vertices for each cell
3052: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
3053: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
3054: . spaceDim - The spatial dimension used for coordinates
3055: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3057:   Output Parameter:
3058: + dm - The DM
3059: - vertexSF - (Optional) SF describing complete vertex ownership

3061:   Notes:
3062:   This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(),
3063:   DMPlexBuildFromCellListParallel(), DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellListParallel()

3065:   See DMPlexBuildFromCellListParallel() for an example and details about the topology-related parameters.
3066:   See DMPlexBuildCoordinatesFromCellListParallel() for details about the geometry-related parameters.

3068:   Level: intermediate

3070: .seealso: DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromDAG(), DMPlexCreate()
3071: @*/
3072: PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, DM *dm)
3073: {
3074:   PetscSF        sfVert;

3078:   DMCreate(comm, dm);
3079:   DMSetType(*dm, DMPLEX);
3082:   DMSetDimension(*dm, dim);
3083:   DMPlexBuildFromCellListParallel(*dm, numCells, numVertices, NVertices, numCorners, cells, &sfVert);
3084:   if (interpolate) {
3085:     DM idm;

3087:     DMPlexInterpolate(*dm, &idm);
3088:     DMDestroy(dm);
3089:     *dm  = idm;
3090:   }
3091:   DMPlexBuildCoordinatesFromCellListParallel(*dm, spaceDim, sfVert, vertexCoords);
3092:   if (vertexSF) *vertexSF = sfVert;
3093:   else {PetscSFDestroy(&sfVert);}
3094:   return(0);
3095: }


3098: /*@
3099:   DMPlexCreateFromCellListParallel - Deprecated, use DMPlexCreateFromCellListParallelPetsc()

3101:   Level: deprecated

3103: .seealso: DMPlexCreateFromCellListParallelPetsc()
3104: @*/
3105: 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)
3106: {
3108:   PetscInt       i;
3109:   PetscInt       *pintCells;

3112:   if (sizeof(int) > sizeof(PetscInt)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Size of int %zd greater than size of PetscInt %zd. Reconfigure PETSc --with-64-bit-indices=1", sizeof(int), sizeof(PetscInt));
3113:   if (sizeof(int) == sizeof(PetscInt)) {
3114:     pintCells = (PetscInt *) cells;
3115:   } else {
3116:     PetscMalloc1(numCells*numCorners, &pintCells);
3117:     for (i = 0; i < numCells*numCorners; i++) {
3118:       pintCells[i] = (PetscInt) cells[i];
3119:     }
3120:   }
3121:   DMPlexCreateFromCellListParallelPetsc(comm, dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, pintCells, spaceDim, vertexCoords, vertexSF, dm);
3122:   if (sizeof(int) != sizeof(PetscInt)) {
3123:     PetscFree(pintCells);
3124:   }
3125:   return(0);
3126: }

3128: /*@C
3129:   DMPlexBuildFromCellList - Build DMPLEX topology from a list of vertices for each cell (common mesh generator output)

3131:   Input Parameters:
3132: + dm - The DM
3133: . numCells - The number of cells owned by this process
3134: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
3135: . numCorners - The number of vertices for each cell
3136: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell

3138:   Level: advanced

3140:   Notes:
3141:   Two triangles sharing a face
3142: $
3143: $        2
3144: $      / | \
3145: $     /  |  \
3146: $    /   |   \
3147: $   0  0 | 1  3
3148: $    \   |   /
3149: $     \  |  /
3150: $      \ | /
3151: $        1
3152: would have input
3153: $  numCells = 2, numVertices = 4
3154: $  cells = [0 1 2  1 3 2]
3155: $
3156: which would result in the DMPlex
3157: $
3158: $        4
3159: $      / | \
3160: $     /  |  \
3161: $    /   |   \
3162: $   2  0 | 1  5
3163: $    \   |   /
3164: $     \  |  /
3165: $      \ | /
3166: $        3

3168:   If numVertices is PETSC_DECIDE, it is computed by PETSc as the maximum vertex index in cells + 1.

3170:   Not currently supported in Fortran.

3172: .seealso: DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromCellListPetsc()
3173: @*/
3174: PetscErrorCode DMPlexBuildFromCellList(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const PetscInt cells[])
3175: {
3176:   PetscInt      *cones, c, p, dim;

3180:   PetscLogEventBegin(DMPLEX_BuildFromCellList,dm,0,0,0);
3181:   DMGetDimension(dm, &dim);
3182:   /* Get/check global number of vertices */
3183:   {
3184:     PetscInt NVerticesInCells, i;
3185:     const PetscInt len = numCells * numCorners;

3187:     /* NVerticesInCells = max(cells) + 1 */
3188:     NVerticesInCells = PETSC_MIN_INT;
3189:     for (i=0; i<len; i++) if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
3190:     ++NVerticesInCells;

3192:     if (numVertices == PETSC_DECIDE) numVertices = NVerticesInCells;
3193:     else if (numVertices < NVerticesInCells) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Specified number of vertices %D must be greater than or equal to the number of vertices in cells %D",numVertices,NVerticesInCells);
3194:   }
3195:   DMPlexSetChart(dm, 0, numCells+numVertices);
3196:   for (c = 0; c < numCells; ++c) {
3197:     DMPlexSetConeSize(dm, c, numCorners);
3198:   }
3199:   DMSetUp(dm);
3200:   DMPlexGetCones(dm,&cones);
3201:   for (c = 0; c < numCells; ++c) {
3202:     for (p = 0; p < numCorners; ++p) {
3203:       cones[c*numCorners+p] = cells[c*numCorners+p]+numCells;
3204:     }
3205:   }
3206:   DMPlexSymmetrize(dm);
3207:   DMPlexStratify(dm);
3208:   PetscLogEventEnd(DMPLEX_BuildFromCellList,dm,0,0,0);
3209:   return(0);
3210: }

3212: /*@C
3213:   DMPlexBuildCoordinatesFromCellList - Build DM coordinates from a list of coordinates for each owned vertex (common mesh generator output)

3215:   Input Parameters:
3216: + dm - The DM
3217: . spaceDim - The spatial dimension used for coordinates
3218: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3220:   Level: advanced

3222:   Notes:
3223:   Not currently supported in Fortran.

3225: .seealso: DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellList()
3226: @*/
3227: PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM dm, PetscInt spaceDim, const PetscReal vertexCoords[])
3228: {
3229:   PetscSection   coordSection;
3230:   Vec            coordinates;
3231:   DM             cdm;
3232:   PetscScalar   *coords;
3233:   PetscInt       v, vStart, vEnd, d;

3237:   PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3238:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
3239:   if (vStart < 0 || vEnd < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first.");
3240:   DMSetCoordinateDim(dm, spaceDim);
3241:   DMGetCoordinateSection(dm, &coordSection);
3242:   PetscSectionSetNumFields(coordSection, 1);
3243:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
3244:   PetscSectionSetChart(coordSection, vStart, vEnd);
3245:   for (v = vStart; v < vEnd; ++v) {
3246:     PetscSectionSetDof(coordSection, v, spaceDim);
3247:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
3248:   }
3249:   PetscSectionSetUp(coordSection);

3251:   DMGetCoordinateDM(dm, &cdm);
3252:   DMCreateLocalVector(cdm, &coordinates);
3253:   VecSetBlockSize(coordinates, spaceDim);
3254:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3255:   VecGetArrayWrite(coordinates, &coords);
3256:   for (v = 0; v < vEnd-vStart; ++v) {
3257:     for (d = 0; d < spaceDim; ++d) {
3258:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
3259:     }
3260:   }
3261:   VecRestoreArrayWrite(coordinates, &coords);
3262:   DMSetCoordinatesLocal(dm, coordinates);
3263:   VecDestroy(&coordinates);
3264:   PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3265:   return(0);
3266: }

3268: /*@
3269:   DMPlexCreateFromCellListPetsc - Create DMPLEX from a list of vertices for each cell (common mesh generator output)

3271:   Input Parameters:
3272: + comm - The communicator
3273: . dim - The topological dimension of the mesh
3274: . numCells - The number of cells
3275: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
3276: . numCorners - The number of vertices for each cell
3277: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
3278: . cells - An array of numCells*numCorners numbers, the vertices for each cell
3279: . spaceDim - The spatial dimension used for coordinates
3280: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3282:   Output Parameter:
3283: . dm - The DM

3285:   Notes:
3286:   This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(), DMPlexBuildFromCellList(),
3287:   DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellList()

3289:   See DMPlexBuildFromCellList() for an example and details about the topology-related parameters.
3290:   See DMPlexBuildCoordinatesFromCellList() for details about the geometry-related parameters.

3292:   Level: intermediate

3294: .seealso: DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildFromCellList(), DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
3295: @*/
3296: PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], DM *dm)
3297: {

3301:   if (!dim) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "This is not appropriate for 0-dimensional meshes. Consider either creating the DM using DMPlexCreateFromDAG(), by hand, or using DMSwarm.");
3302:   DMCreate(comm, dm);
3303:   DMSetType(*dm, DMPLEX);
3304:   DMSetDimension(*dm, dim);
3305:   DMPlexBuildFromCellList(*dm, numCells, numVertices, numCorners, cells);
3306:   if (interpolate) {
3307:     DM idm;

3309:     DMPlexInterpolate(*dm, &idm);
3310:     DMDestroy(dm);
3311:     *dm  = idm;
3312:   }
3313:   DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, vertexCoords);
3314:   return(0);
3315: }

3317: /*@
3318:   DMPlexCreateFromCellList - Deprecated, use DMPlexCreateFromCellListPetsc()

3320:   Level: deprecated

3322: .seealso: DMPlexCreateFromCellListPetsc()
3323: @*/
3324: 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)
3325: {
3327:   PetscInt       i;
3328:   PetscInt       *pintCells;
3329:   PetscReal      *prealVC;

3332:   if (sizeof(int) > sizeof(PetscInt)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Size of int %zd greater than size of PetscInt %zd. Reconfigure PETSc --with-64-bit-indices=1", sizeof(int), sizeof(PetscInt));
3333:   if (sizeof(int) == sizeof(PetscInt)) {
3334:     pintCells = (PetscInt *) cells;
3335:   } else {
3336:     PetscMalloc1(numCells*numCorners, &pintCells);
3337:     for (i = 0; i < numCells*numCorners; i++) {
3338:       pintCells[i] = (PetscInt) cells[i];
3339:     }
3340:   }
3341:   if (sizeof(double) > sizeof(PetscReal)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Size of double %zd greater than size of PetscReal %zd. Reconfigure PETSc --with-precision=<higher precision>.", sizeof(double), sizeof(PetscReal));
3342:   if (sizeof(double) == sizeof(PetscReal)) {
3343:     prealVC = (PetscReal *) vertexCoords;
3344:   } else {
3345:     PetscMalloc1(numVertices*spaceDim, &prealVC);
3346:     for (i = 0; i < numVertices*spaceDim; i++) {
3347:       prealVC[i] = (PetscReal) vertexCoords[i];
3348:     }
3349:   }
3350:   DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, pintCells, spaceDim, prealVC, dm);
3351:   if (sizeof(int) != sizeof(PetscInt)) {
3352:     PetscFree(pintCells);
3353:   }
3354:   if (sizeof(double) != sizeof(PetscReal)) {
3355:     PetscFree(prealVC);
3356:   }
3357:   return(0);
3358: }

3360: /*@
3361:   DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM

3363:   Input Parameters:
3364: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
3365: . depth - The depth of the DAG
3366: . numPoints - Array of size depth + 1 containing the number of points at each depth
3367: . coneSize - The cone size of each point
3368: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
3369: . coneOrientations - The orientation of each cone point
3370: - vertexCoords - An array of numPoints[0]*spacedim numbers representing the coordinates of each vertex, with spacedim the value set via DMSetCoordinateDim()

3372:   Output Parameter:
3373: . dm - The DM

3375:   Note: Two triangles sharing a face would have input
3376: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
3377: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
3378: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
3379: $
3380: which would result in the DMPlex
3381: $
3382: $        4
3383: $      / | \
3384: $     /  |  \
3385: $    /   |   \
3386: $   2  0 | 1  5
3387: $    \   |   /
3388: $     \  |  /
3389: $      \ | /
3390: $        3
3391: $
3392: $ Notice that all points are numbered consecutively, unlike DMPlexCreateFromCellListPetsc()

3394:   Level: advanced

3396: .seealso: DMPlexCreateFromCellListPetsc(), DMPlexCreate()
3397: @*/
3398: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
3399: {
3400:   Vec            coordinates;
3401:   PetscSection   coordSection;
3402:   PetscScalar    *coords;
3403:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

3407:   DMGetDimension(dm, &dim);
3408:   DMGetCoordinateDim(dm, &dimEmbed);
3409:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %D cannot be less than intrinsic dimension %d",dimEmbed,dim);
3410:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
3411:   DMPlexSetChart(dm, pStart, pEnd);
3412:   for (p = pStart; p < pEnd; ++p) {
3413:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
3414:     if (firstVertex < 0 && !coneSize[p - pStart]) {
3415:       firstVertex = p - pStart;
3416:     }
3417:   }
3418:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %D vertices but could not find any", numPoints[0]);
3419:   DMSetUp(dm); /* Allocate space for cones */
3420:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
3421:     DMPlexSetCone(dm, p, &cones[off]);
3422:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
3423:   }
3424:   DMPlexSymmetrize(dm);
3425:   DMPlexStratify(dm);
3426:   /* Build coordinates */
3427:   DMGetCoordinateSection(dm, &coordSection);
3428:   PetscSectionSetNumFields(coordSection, 1);
3429:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
3430:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
3431:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
3432:     PetscSectionSetDof(coordSection, v, dimEmbed);
3433:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
3434:   }
3435:   PetscSectionSetUp(coordSection);
3436:   PetscSectionGetStorageSize(coordSection, &coordSize);
3437:   VecCreate(PETSC_COMM_SELF, &coordinates);
3438:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3439:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3440:   VecSetBlockSize(coordinates, dimEmbed);
3441:   VecSetType(coordinates,VECSTANDARD);
3442:   VecGetArray(coordinates, &coords);
3443:   for (v = 0; v < numPoints[0]; ++v) {
3444:     PetscInt off;

3446:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
3447:     for (d = 0; d < dimEmbed; ++d) {
3448:       coords[off+d] = vertexCoords[v*dimEmbed+d];
3449:     }
3450:   }
3451:   VecRestoreArray(coordinates, &coords);
3452:   DMSetCoordinatesLocal(dm, coordinates);
3453:   VecDestroy(&coordinates);
3454:   return(0);
3455: }

3457: /*@C
3458:   DMPlexCreateCellVertexFromFile - Create a DMPlex mesh from a simple cell-vertex file.

3460: + comm        - The MPI communicator
3461: . filename    - Name of the .dat file
3462: - interpolate - Create faces and edges in the mesh

3464:   Output Parameter:
3465: . dm  - The DM object representing the mesh

3467:   Note: The format is the simplest possible:
3468: $ Ne
3469: $ v0 v1 ... vk
3470: $ Nv
3471: $ x y z marker

3473:   Level: beginner

3475: .seealso: DMPlexCreateFromFile(), DMPlexCreateMedFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
3476: @*/
3477: PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3478: {
3479:   DMLabel         marker;
3480:   PetscViewer     viewer;
3481:   Vec             coordinates;
3482:   PetscSection    coordSection;
3483:   PetscScalar    *coords;
3484:   char            line[PETSC_MAX_PATH_LEN];
3485:   PetscInt        dim = 3, cdim = 3, coordSize, v, c, d;
3486:   PetscMPIInt     rank;
3487:   int             snum, Nv, Nc;
3488:   PetscErrorCode  ierr;

3491:   MPI_Comm_rank(comm, &rank);
3492:   PetscViewerCreate(comm, &viewer);
3493:   PetscViewerSetType(viewer, PETSCVIEWERASCII);
3494:   PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3495:   PetscViewerFileSetName(viewer, filename);
3496:   if (!rank) {
3497:     PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);
3498:     snum = sscanf(line, "%d %d", &Nc, &Nv);
3499:     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3500:   } else {
3501:     Nc = Nv = 0;
3502:   }
3503:   DMCreate(comm, dm);
3504:   DMSetType(*dm, DMPLEX);
3505:   DMPlexSetChart(*dm, 0, Nc+Nv);
3506:   DMSetDimension(*dm, dim);
3507:   DMSetCoordinateDim(*dm, cdim);
3508:   /* Read topology */
3509:   if (!rank) {
3510:     PetscInt cone[8], corners = 8;
3511:     int      vbuf[8], v;

3513:     for (c = 0; c < Nc; ++c) {DMPlexSetConeSize(*dm, c, corners);}
3514:     DMSetUp(*dm);
3515:     for (c = 0; c < Nc; ++c) {
3516:       PetscViewerRead(viewer, line, corners, NULL, PETSC_STRING);
3517:       snum = sscanf(line, "%d %d %d %d %d %d %d %d", &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5], &vbuf[6], &vbuf[7]);
3518:       if (snum != corners) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3519:       for (v = 0; v < corners; ++v) cone[v] = vbuf[v] + Nc;
3520:       /* Hexahedra are inverted */
3521:       {
3522:         PetscInt tmp = cone[1];
3523:         cone[1] = cone[3];
3524:         cone[3] = tmp;
3525:       }
3526:       DMPlexSetCone(*dm, c, cone);
3527:     }
3528:   }
3529:   DMPlexSymmetrize(*dm);
3530:   DMPlexStratify(*dm);
3531:   /* Read coordinates */
3532:   DMGetCoordinateSection(*dm, &coordSection);
3533:   PetscSectionSetNumFields(coordSection, 1);
3534:   PetscSectionSetFieldComponents(coordSection, 0, cdim);
3535:   PetscSectionSetChart(coordSection, Nc, Nc + Nv);
3536:   for (v = Nc; v < Nc+Nv; ++v) {
3537:     PetscSectionSetDof(coordSection, v, cdim);
3538:     PetscSectionSetFieldDof(coordSection, v, 0, cdim);
3539:   }
3540:   PetscSectionSetUp(coordSection);
3541:   PetscSectionGetStorageSize(coordSection, &coordSize);
3542:   VecCreate(PETSC_COMM_SELF, &coordinates);
3543:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3544:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3545:   VecSetBlockSize(coordinates, cdim);
3546:   VecSetType(coordinates, VECSTANDARD);
3547:   VecGetArray(coordinates, &coords);
3548:   if (!rank) {
3549:     double x[3];
3550:     int    val;

3552:     DMCreateLabel(*dm, "marker");
3553:     DMGetLabel(*dm, "marker", &marker);
3554:     for (v = 0; v < Nv; ++v) {
3555:       PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING);
3556:       snum = sscanf(line, "%lg %lg %lg %d", &x[0], &x[1], &x[2], &val);
3557:       if (snum != 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3558:       for (d = 0; d < cdim; ++d) coords[v*cdim+d] = x[d];
3559:       if (val) {DMLabelSetValue(marker, v+Nc, val);}
3560:     }
3561:   }
3562:   VecRestoreArray(coordinates, &coords);
3563:   DMSetCoordinatesLocal(*dm, coordinates);
3564:   VecDestroy(&coordinates);
3565:   PetscViewerDestroy(&viewer);
3566:   if (interpolate) {
3567:     DM      idm;
3568:     DMLabel bdlabel;

3570:     DMPlexInterpolate(*dm, &idm);
3571:     DMDestroy(dm);
3572:     *dm  = idm;

3574:     DMGetLabel(*dm, "marker", &bdlabel);
3575:     DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
3576:     DMPlexLabelComplete(*dm, bdlabel);
3577:   }
3578:   return(0);
3579: }

3581: /*@C
3582:   DMPlexCreateFromFile - This takes a filename and produces a DM

3584:   Input Parameters:
3585: + comm - The communicator
3586: . filename - A file name
3587: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

3589:   Output Parameter:
3590: . dm - The DM

3592:   Options Database Keys:
3593: . -dm_plex_create_from_hdf5_xdmf - use the PETSC_VIEWER_HDF5_XDMF format for reading HDF5

3595:   Use -dm_plex_create_ prefix to pass options to the internal PetscViewer, e.g.
3596: $ -dm_plex_create_viewer_hdf5_collective

3598:   Level: beginner

3600: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellListPetsc(), DMPlexCreate()
3601: @*/
3602: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3603: {
3604:   const char    *extGmsh    = ".msh";
3605:   const char    *extGmsh2   = ".msh2";
3606:   const char    *extGmsh4   = ".msh4";
3607:   const char    *extCGNS    = ".cgns";
3608:   const char    *extExodus  = ".exo";
3609:   const char    *extGenesis = ".gen";
3610:   const char    *extFluent  = ".cas";
3611:   const char    *extHDF5    = ".h5";
3612:   const char    *extMed     = ".med";
3613:   const char    *extPLY     = ".ply";
3614:   const char    *extEGADS   = ".egadslite";
3615:   const char    *extCV      = ".dat";
3616:   size_t         len;
3617:   PetscBool      isGmsh, isGmsh2, isGmsh4, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isEGADS, isCV;
3618:   PetscMPIInt    rank;

3624:   DMInitializePackage();
3625:   PetscLogEventBegin(DMPLEX_CreateFromFile,0,0,0,0);
3626:   MPI_Comm_rank(comm, &rank);
3627:   PetscStrlen(filename, &len);
3628:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
3629:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,    4, &isGmsh);
3630:   PetscStrncmp(&filename[PetscMax(0,len-5)], extGmsh2,   5, &isGmsh2);
3631:   PetscStrncmp(&filename[PetscMax(0,len-5)], extGmsh4,   5, &isGmsh4);
3632:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,    5, &isCGNS);
3633:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus,  4, &isExodus);
3634:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);
3635:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent,  4, &isFluent);
3636:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,    3, &isHDF5);
3637:   PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,     4, &isMed);
3638:   PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY,     4, &isPLY);
3639:   PetscStrncmp(&filename[PetscMax(0,len-10)], extEGADS,   9, &isEGADS);
3640:   PetscStrncmp(&filename[PetscMax(0,len-4)], extCV,      4, &isCV);
3641:   if (isGmsh || isGmsh2 || isGmsh4) {
3642:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
3643:   } else if (isCGNS) {
3644:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
3645:   } else if (isExodus || isGenesis) {
3646:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
3647:   } else if (isFluent) {
3648:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
3649:   } else if (isHDF5) {
3650:     PetscBool      load_hdf5_xdmf = PETSC_FALSE;
3651:     PetscViewer viewer;

3653:     /* PETSC_VIEWER_HDF5_XDMF is used if the filename ends with .xdmf.h5, or if -dm_plex_create_from_hdf5_xdmf option is present */
3654:     PetscStrncmp(&filename[PetscMax(0,len-8)], ".xdmf",  5, &load_hdf5_xdmf);
3655:     PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &load_hdf5_xdmf, NULL);
3656:     PetscViewerCreate(comm, &viewer);
3657:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
3658:     PetscViewerSetOptionsPrefix(viewer, "dm_plex_create_");
3659:     PetscViewerSetFromOptions(viewer);
3660:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3661:     PetscViewerFileSetName(viewer, filename);
3662:     DMCreate(comm, dm);
3663:     DMSetType(*dm, DMPLEX);
3664:     if (load_hdf5_xdmf) {PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF);}
3665:     DMLoad(*dm, viewer);
3666:     if (load_hdf5_xdmf) {PetscViewerPopFormat(viewer);}
3667:     PetscViewerDestroy(&viewer);

3669:     if (interpolate) {
3670:       DM idm;

3672:       DMPlexInterpolate(*dm, &idm);
3673:       DMDestroy(dm);
3674:       *dm  = idm;
3675:     }
3676:   } else if (isMed) {
3677:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
3678:   } else if (isPLY) {
3679:     DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);
3680:   } else if (isEGADS) {
3681:     DMPlexCreateEGADSFromFile(comm, filename, dm);
3682:     if (!interpolate) {
3683:       DM udm;

3685:       DMPlexUninterpolate(*dm, &udm);
3686:       DMDestroy(dm);
3687:       *dm  = udm;
3688:     }
3689:   } else if (isCV) {
3690:     DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm);
3691:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
3692:   PetscLogEventEnd(DMPLEX_CreateFromFile,0,0,0,0);
3693:   return(0);
3694: }

3696: /*@
3697:   DMPlexCreateReferenceCellByType - Create a DMPLEX with the appropriate FEM reference cell

3699:   Collective

3701:   Input Parameters:
3702: + comm - The communicator
3703: - ct   - The cell type of the reference cell

3705:   Output Parameter:
3706: . refdm - The reference cell

3708:   Options Database Keys:
3709: . -dm_plex_ref_type <ct> - Specify the celltyoe for the reference cell

3711:   Level: intermediate

3713: .seealso: DMPlexCreateReferenceCell(), DMPlexCreateBoxMesh()
3714: @*/
3715: PetscErrorCode DMPlexCreateReferenceCellByType(MPI_Comm comm, DMPolytopeType ct, DM *refdm)
3716: {
3717:   DM             rdm;
3718:   Vec            coords;

3722:   PetscOptionsGetEnum(NULL, NULL, "-dm_plex_ref_type", DMPolytopeTypes, (PetscEnum *) &ct, NULL);
3723:   DMCreate(comm, &rdm);
3724:   DMSetType(rdm, DMPLEX);
3725:   switch (ct) {
3726:     case DM_POLYTOPE_POINT:
3727:     {
3728:       PetscInt    numPoints[1]        = {1};
3729:       PetscInt    coneSize[1]         = {0};
3730:       PetscInt    cones[1]            = {0};
3731:       PetscInt    coneOrientations[1] = {0};
3732:       PetscScalar vertexCoords[1]     = {0.0};

3734:       DMSetDimension(rdm, 0);
3735:       DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3736:     }
3737:     break;
3738:     case DM_POLYTOPE_SEGMENT:
3739:     {
3740:       PetscInt    numPoints[2]        = {2, 1};
3741:       PetscInt    coneSize[3]         = {2, 0, 0};
3742:       PetscInt    cones[2]            = {1, 2};
3743:       PetscInt    coneOrientations[2] = {0, 0};
3744:       PetscScalar vertexCoords[2]     = {-1.0,  1.0};

3746:       DMSetDimension(rdm, 1);
3747:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3748:     }
3749:     break;
3750:     case DM_POLYTOPE_TRIANGLE:
3751:     {
3752:       PetscInt    numPoints[2]        = {3, 1};
3753:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
3754:       PetscInt    cones[3]            = {1, 2, 3};
3755:       PetscInt    coneOrientations[3] = {0, 0, 0};
3756:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

3758:       DMSetDimension(rdm, 2);
3759:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3760:     }
3761:     break;
3762:     case DM_POLYTOPE_QUADRILATERAL:
3763:     {
3764:       PetscInt    numPoints[2]        = {4, 1};
3765:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3766:       PetscInt    cones[4]            = {1, 2, 3, 4};
3767:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3768:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

3770:       DMSetDimension(rdm, 2);
3771:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3772:     }
3773:     break;
3774:     case DM_POLYTOPE_SEG_PRISM_TENSOR:
3775:     {
3776:       PetscInt    numPoints[2]        = {4, 1};
3777:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3778:       PetscInt    cones[4]            = {1, 2, 3, 4};
3779:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3780:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0,  1.0, 1.0};

3782:       DMSetDimension(rdm, 2);
3783:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3784:     }
3785:     break;
3786:     case DM_POLYTOPE_TETRAHEDRON:
3787:     {
3788:       PetscInt    numPoints[2]        = {4, 1};
3789:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3790:       PetscInt    cones[4]            = {1, 3, 2, 4};
3791:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3792:       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};

3794:       DMSetDimension(rdm, 3);
3795:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3796:     }
3797:     break;
3798:     case DM_POLYTOPE_HEXAHEDRON:
3799:     {
3800:       PetscInt    numPoints[2]        = {8, 1};
3801:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3802:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
3803:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3804:       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,
3805:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3807:       DMSetDimension(rdm, 3);
3808:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3809:     }
3810:     break;
3811:     case DM_POLYTOPE_TRI_PRISM:
3812:     {
3813:       PetscInt    numPoints[2]        = {6, 1};
3814:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3815:       PetscInt    cones[6]            = {1, 3, 2, 4, 5, 6};
3816:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3817:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3818:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3820:       DMSetDimension(rdm, 3);
3821:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3822:     }
3823:     break;
3824:     case DM_POLYTOPE_TRI_PRISM_TENSOR:
3825:     {
3826:       PetscInt    numPoints[2]        = {6, 1};
3827:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3828:       PetscInt    cones[6]            = {1, 2, 3, 4, 5, 6};
3829:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3830:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3831:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3833:       DMSetDimension(rdm, 3);
3834:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3835:     }
3836:     break;
3837:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
3838:     {
3839:       PetscInt    numPoints[2]        = {8, 1};
3840:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3841:       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
3842:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3843:       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,
3844:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3846:       DMSetDimension(rdm, 3);
3847:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3848:     }
3849:     break;
3850:     case DM_POLYTOPE_PYRAMID:
3851:     {
3852:       PetscInt    numPoints[2]        = {5, 1};
3853:       PetscInt    coneSize[6]         = {5, 0, 0, 0, 0, 0};
3854:       PetscInt    cones[5]            = {1, 4, 3, 2, 5};
3855:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3856:       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,
3857:                                           0.0,  0.0,  1.0};

3859:       DMSetDimension(rdm, 3);
3860:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3861:     }
3862:     break;
3863:     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for cell type %s", DMPolytopeTypes[ct]);
3864:   }
3865:   {
3866:     PetscInt Nv, v;

3868:     /* Must create the celltype label here so that we do not automatically try to compute the types */
3869:     DMCreateLabel(rdm, "celltype");
3870:     DMPlexSetCellType(rdm, 0, ct);
3871:     DMPlexGetChart(rdm, NULL, &Nv);
3872:     for (v = 1; v < Nv; ++v) {DMPlexSetCellType(rdm, v, DM_POLYTOPE_POINT);}
3873:   }
3874:   DMPlexInterpolate(rdm, refdm);
3875:   if (rdm->coordinateDM) {
3876:     DM           ncdm;
3877:     PetscSection cs;
3878:     PetscInt     pEnd = -1;

3880:     DMGetLocalSection(rdm->coordinateDM, &cs);
3881:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
3882:     if (pEnd >= 0) {
3883:       DMClone(*refdm, &ncdm);
3884:       DMCopyDisc(rdm->coordinateDM, ncdm);
3885:       DMSetLocalSection(ncdm, cs);
3886:       DMSetCoordinateDM(*refdm, ncdm);
3887:       DMDestroy(&ncdm);
3888:     }
3889:   }
3890:   DMGetCoordinatesLocal(rdm, &coords);
3891:   if (coords) {
3892:     DMSetCoordinatesLocal(*refdm, coords);
3893:   } else {
3894:     DMGetCoordinates(rdm, &coords);
3895:     if (coords) {DMSetCoordinates(*refdm, coords);}
3896:   }
3897:   DMDestroy(&rdm);
3898:   return(0);
3899: }

3901: /*@
3902:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

3904:   Collective

3906:   Input Parameters:
3907: + comm    - The communicator
3908: . dim     - The spatial dimension
3909: - simplex - Flag for simplex, otherwise use a tensor-product cell

3911:   Output Parameter:
3912: . refdm - The reference cell

3914:   Level: intermediate

3916: .seealso: DMPlexCreateReferenceCellByType(), DMPlexCreateBoxMesh()
3917: @*/
3918: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
3919: {

3923:   switch (dim) {
3924:   case 0: DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_POINT, refdm);break;
3925:   case 1: DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_SEGMENT, refdm);break;
3926:   case 2:
3927:     if (simplex) {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_TRIANGLE, refdm);}
3928:     else         {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_QUADRILATERAL, refdm);}
3929:     break;
3930:   case 3:
3931:     if (simplex) {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_TETRAHEDRON, refdm);}
3932:     else         {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_HEXAHEDRON, refdm);}
3933:     break;
3934:   default:
3935:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %D", dim);
3936:   }
3937:   return(0);
3938: }