Actual source code: plexcreate.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  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: . refinementUniform - Flag for uniform parallel refinement
 19: - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell

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

 24:   Level: beginner

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

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

 66:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 67:         for (p = 0; p < 4; ++p) {DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);}
 68:       } else {
 69:         PetscInt    numPoints[2]        = {6, 2};
 70:         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
 71:         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
 72:         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
 73:         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};

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

 87:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 88:         for (p = 0; p < 5; ++p) {DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);}
 89:       } else {
 90:         PetscInt    numPoints[2]         = {12, 2};
 91:         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
 92:         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
 93:         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
 94:         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,
 95:                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
 96:                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};

 98:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
 99:       }
100:       break;
101:     default:
102:       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %D", dim);
103:     }
104:   }
105:   *newdm = dm;
106:   if (refinementLimit > 0.0) {
107:     DM rdm;
108:     const char *name;

110:     DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);
111:     DMPlexSetRefinementLimit(*newdm, refinementLimit);
112:     DMRefine(*newdm, comm, &rdm);
113:     PetscObjectGetName((PetscObject) *newdm, &name);
114:     PetscObjectSetName((PetscObject)    rdm,  name);
115:     DMDestroy(newdm);
116:     *newdm = rdm;
117:   }
118:   if (interpolate) {
119:     DM idm;

121:     DMPlexInterpolate(*newdm, &idm);
122:     DMDestroy(newdm);
123:     *newdm = idm;
124:   }
125:   {
126:     DM refinedMesh     = NULL;
127:     DM distributedMesh = NULL;

129:     /* Distribute mesh over processes */
130:     DMPlexDistribute(*newdm, 0, NULL, &distributedMesh);
131:     if (distributedMesh) {
132:       DMDestroy(newdm);
133:       *newdm = distributedMesh;
134:     }
135:     if (refinementUniform) {
136:       DMPlexSetRefinementUniform(*newdm, refinementUniform);
137:       DMRefine(*newdm, comm, &refinedMesh);
138:       if (refinedMesh) {
139:         DMDestroy(newdm);
140:         *newdm = refinedMesh;
141:       }
142:     }
143:   }
144:   return(0);
145: }

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

150:   Collective

152:   Input Parameters:
153: + comm  - The communicator for the DM object
154: . lower - The lower left corner coordinates
155: . upper - The upper right corner coordinates
156: - edges - The number of cells in each direction

158:   Output Parameter:
159: . dm  - The DM object

161:   Note: Here is the numbering returned for 2 cells in each direction:
162: $ 18--5-17--4--16
163: $  |     |     |
164: $  6    10     3
165: $  |     |     |
166: $ 19-11-20--9--15
167: $  |     |     |
168: $  7     8     2
169: $  |     |     |
170: $ 12--0-13--1--14

172:   Level: beginner

174: .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
175: @*/
176: PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
177: {
178:   const PetscInt numVertices    = (edges[0]+1)*(edges[1]+1);
179:   const PetscInt numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
180:   PetscInt       markerTop      = 1;
181:   PetscInt       markerBottom   = 1;
182:   PetscInt       markerRight    = 1;
183:   PetscInt       markerLeft     = 1;
184:   PetscBool      markerSeparate = PETSC_FALSE;
185:   Vec            coordinates;
186:   PetscSection   coordSection;
187:   PetscScalar    *coords;
188:   PetscInt       coordSize;
189:   PetscMPIInt    rank;
190:   PetscInt       v, vx, vy;

194:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
195:   if (markerSeparate) {
196:     markerTop    = 3;
197:     markerBottom = 1;
198:     markerRight  = 2;
199:     markerLeft   = 4;
200:   }
201:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
202:   if (!rank) {
203:     PetscInt e, ex, ey;

205:     DMPlexSetChart(dm, 0, numEdges+numVertices);
206:     for (e = 0; e < numEdges; ++e) {
207:       DMPlexSetConeSize(dm, e, 2);
208:     }
209:     DMSetUp(dm); /* Allocate space for cones */
210:     for (vx = 0; vx <= edges[0]; vx++) {
211:       for (ey = 0; ey < edges[1]; ey++) {
212:         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
213:         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
214:         PetscInt cone[2];

216:         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
217:         DMPlexSetCone(dm, edge, cone);
218:         if (vx == edges[0]) {
219:           DMSetLabelValue(dm, "marker", edge,    markerRight);
220:           DMSetLabelValue(dm, "marker", cone[0], markerRight);
221:           if (ey == edges[1]-1) {
222:             DMSetLabelValue(dm, "marker", cone[1], markerRight);
223:             DMSetLabelValue(dm, "Face Sets", cone[1], markerRight);
224:           }
225:         } else if (vx == 0) {
226:           DMSetLabelValue(dm, "marker", edge,    markerLeft);
227:           DMSetLabelValue(dm, "marker", cone[0], markerLeft);
228:           if (ey == edges[1]-1) {
229:             DMSetLabelValue(dm, "marker", cone[1], markerLeft);
230:             DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft);
231:           }
232:         }
233:       }
234:     }
235:     for (vy = 0; vy <= edges[1]; vy++) {
236:       for (ex = 0; ex < edges[0]; ex++) {
237:         PetscInt edge   = vy*edges[0]     + ex;
238:         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
239:         PetscInt cone[2];

241:         cone[0] = vertex; cone[1] = vertex+1;
242:         DMPlexSetCone(dm, edge, cone);
243:         if (vy == edges[1]) {
244:           DMSetLabelValue(dm, "marker", edge,    markerTop);
245:           DMSetLabelValue(dm, "marker", cone[0], markerTop);
246:           if (ex == edges[0]-1) {
247:             DMSetLabelValue(dm, "marker", cone[1], markerTop);
248:             DMSetLabelValue(dm, "Face Sets", cone[1], markerTop);
249:           }
250:         } else if (vy == 0) {
251:           DMSetLabelValue(dm, "marker", edge,    markerBottom);
252:           DMSetLabelValue(dm, "marker", cone[0], markerBottom);
253:           if (ex == edges[0]-1) {
254:             DMSetLabelValue(dm, "marker", cone[1], markerBottom);
255:             DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom);
256:           }
257:         }
258:       }
259:     }
260:   }
261:   DMPlexSymmetrize(dm);
262:   DMPlexStratify(dm);
263:   /* Build coordinates */
264:   DMSetCoordinateDim(dm, 2);
265:   DMGetCoordinateSection(dm, &coordSection);
266:   PetscSectionSetNumFields(coordSection, 1);
267:   PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);
268:   PetscSectionSetFieldComponents(coordSection, 0, 2);
269:   for (v = numEdges; v < numEdges+numVertices; ++v) {
270:     PetscSectionSetDof(coordSection, v, 2);
271:     PetscSectionSetFieldDof(coordSection, v, 0, 2);
272:   }
273:   PetscSectionSetUp(coordSection);
274:   PetscSectionGetStorageSize(coordSection, &coordSize);
275:   VecCreate(PETSC_COMM_SELF, &coordinates);
276:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
277:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
278:   VecSetBlockSize(coordinates, 2);
279:   VecSetType(coordinates,VECSTANDARD);
280:   VecGetArray(coordinates, &coords);
281:   for (vy = 0; vy <= edges[1]; ++vy) {
282:     for (vx = 0; vx <= edges[0]; ++vx) {
283:       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
284:       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
285:     }
286:   }
287:   VecRestoreArray(coordinates, &coords);
288:   DMSetCoordinatesLocal(dm, coordinates);
289:   VecDestroy(&coordinates);
290:   return(0);
291: }

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

296:   Collective

298:   Input Parameters:
299: + comm  - The communicator for the DM object
300: . lower - The lower left front corner coordinates
301: . upper - The upper right back corner coordinates
302: - faces - The number of faces in each direction (the same as the number of cells)

304:   Output Parameter:
305: . dm  - The DM object

307:   Level: beginner

309: .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate()
310: @*/
311: PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
312: {
313:   PetscInt       vertices[3], numVertices;
314:   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
315:   Vec            coordinates;
316:   PetscSection   coordSection;
317:   PetscScalar    *coords;
318:   PetscInt       coordSize;
319:   PetscMPIInt    rank;
320:   PetscInt       v, vx, vy, vz;
321:   PetscInt       voffset, iface=0, cone[4];

325:   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");
326:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
327:   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
328:   numVertices = vertices[0]*vertices[1]*vertices[2];
329:   if (!rank) {
330:     PetscInt f;

332:     DMPlexSetChart(dm, 0, numFaces+numVertices);
333:     for (f = 0; f < numFaces; ++f) {
334:       DMPlexSetConeSize(dm, f, 4);
335:     }
336:     DMSetUp(dm); /* Allocate space for cones */

338:     /* Side 0 (Top) */
339:     for (vy = 0; vy < faces[1]; vy++) {
340:       for (vx = 0; vx < faces[0]; vx++) {
341:         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
342:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
343:         DMPlexSetCone(dm, iface, cone);
344:         DMSetLabelValue(dm, "marker", iface, 1);
345:         DMSetLabelValue(dm, "marker", voffset+0, 1);
346:         DMSetLabelValue(dm, "marker", voffset+1, 1);
347:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
348:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+1, 1);
349:         iface++;
350:       }
351:     }

353:     /* Side 1 (Bottom) */
354:     for (vy = 0; vy < faces[1]; vy++) {
355:       for (vx = 0; vx < faces[0]; vx++) {
356:         voffset = numFaces + vy*(faces[0]+1) + vx;
357:         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
358:         DMPlexSetCone(dm, iface, cone);
359:         DMSetLabelValue(dm, "marker", iface, 1);
360:         DMSetLabelValue(dm, "marker", voffset+0, 1);
361:         DMSetLabelValue(dm, "marker", voffset+1, 1);
362:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
363:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+1, 1);
364:         iface++;
365:       }
366:     }

368:     /* Side 2 (Front) */
369:     for (vz = 0; vz < faces[2]; vz++) {
370:       for (vx = 0; vx < faces[0]; vx++) {
371:         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
372:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
373:         DMPlexSetCone(dm, iface, cone);
374:         DMSetLabelValue(dm, "marker", iface, 1);
375:         DMSetLabelValue(dm, "marker", voffset+0, 1);
376:         DMSetLabelValue(dm, "marker", voffset+1, 1);
377:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
378:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+1, 1);
379:         iface++;
380:       }
381:     }

383:     /* Side 3 (Back) */
384:     for (vz = 0; vz < faces[2]; vz++) {
385:       for (vx = 0; vx < faces[0]; vx++) {
386:         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
387:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
388:         cone[2] = voffset+1; cone[3] = voffset;
389:         DMPlexSetCone(dm, iface, cone);
390:         DMSetLabelValue(dm, "marker", iface, 1);
391:         DMSetLabelValue(dm, "marker", voffset+0, 1);
392:         DMSetLabelValue(dm, "marker", voffset+1, 1);
393:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
394:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+1, 1);
395:         iface++;
396:       }
397:     }

399:     /* Side 4 (Left) */
400:     for (vz = 0; vz < faces[2]; vz++) {
401:       for (vy = 0; vy < faces[1]; vy++) {
402:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
403:         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
404:         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
405:         DMPlexSetCone(dm, iface, cone);
406:         DMSetLabelValue(dm, "marker", iface, 1);
407:         DMSetLabelValue(dm, "marker", voffset+0, 1);
408:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
409:         DMSetLabelValue(dm, "marker", voffset+vertices[1]+0, 1);
410:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+vertices[0], 1);
411:         iface++;
412:       }
413:     }

415:     /* Side 5 (Right) */
416:     for (vz = 0; vz < faces[2]; vz++) {
417:       for (vy = 0; vy < faces[1]; vy++) {
418:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + faces[0];
419:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
420:         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
421:         DMPlexSetCone(dm, iface, cone);
422:         DMSetLabelValue(dm, "marker", iface, 1);
423:         DMSetLabelValue(dm, "marker", voffset+0, 1);
424:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
425:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
426:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+vertices[0], 1);
427:         iface++;
428:       }
429:     }
430:   }
431:   DMPlexSymmetrize(dm);
432:   DMPlexStratify(dm);
433:   /* Build coordinates */
434:   DMSetCoordinateDim(dm, 3);
435:   DMGetCoordinateSection(dm, &coordSection);
436:   PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);
437:   for (v = numFaces; v < numFaces+numVertices; ++v) {
438:     PetscSectionSetDof(coordSection, v, 3);
439:   }
440:   PetscSectionSetUp(coordSection);
441:   PetscSectionGetStorageSize(coordSection, &coordSize);
442:   VecCreate(PETSC_COMM_SELF, &coordinates);
443:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
444:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
445:   VecSetBlockSize(coordinates, 3);
446:   VecSetType(coordinates,VECSTANDARD);
447:   VecGetArray(coordinates, &coords);
448:   for (vz = 0; vz <= faces[2]; ++vz) {
449:     for (vy = 0; vy <= faces[1]; ++vy) {
450:       for (vx = 0; vx <= faces[0]; ++vx) {
451:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
452:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
453:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
454:       }
455:     }
456:   }
457:   VecRestoreArray(coordinates, &coords);
458:   DMSetCoordinatesLocal(dm, coordinates);
459:   VecDestroy(&coordinates);
460:   return(0);
461: }

463: static PetscErrorCode DMPlexCreateLineMesh_Internal(MPI_Comm comm,PetscInt segments,PetscReal lower,PetscReal upper,DMBoundaryType bd,DM *dm)
464: {
465:   PetscInt       i,fStart,fEnd,numCells = 0,numVerts = 0;
466:   PetscInt       numPoints[2],*coneSize,*cones,*coneOrientations;
467:   PetscScalar    *vertexCoords;
468:   PetscReal      L,maxCell;
469:   PetscBool      markerSeparate = PETSC_FALSE;
470:   PetscInt       markerLeft  = 1, faceMarkerLeft  = 1;
471:   PetscInt       markerRight = 1, faceMarkerRight = 2;
472:   PetscBool      wrap = (bd == DM_BOUNDARY_PERIODIC || bd == DM_BOUNDARY_TWIST) ? PETSC_TRUE : PETSC_FALSE;
473:   PetscMPIInt    rank;


479:   DMCreate(comm,dm);
480:   DMSetType(*dm,DMPLEX);
481:   DMSetDimension(*dm,1);
482:   DMCreateLabel(*dm,"marker");
483:   DMCreateLabel(*dm,"Face Sets");

485:   MPI_Comm_rank(comm,&rank);
486:   if (!rank) numCells = segments;
487:   if (!rank) numVerts = segments + (wrap ? 0 : 1);

489:   numPoints[0] = numVerts ; numPoints[1] = numCells;
490:   PetscMalloc4(numCells+numVerts,&coneSize,numCells*2,&cones,numCells+numVerts,&coneOrientations,numVerts,&vertexCoords);
491:   PetscArrayzero(coneOrientations,numCells+numVerts);
492:   for (i = 0; i < numCells; ++i) { coneSize[i] = 2; }
493:   for (i = 0; i < numVerts; ++i) { coneSize[numCells+i] = 0; }
494:   for (i = 0; i < numCells; ++i) { cones[2*i] = numCells + i%numVerts; cones[2*i+1] = numCells + (i+1)%numVerts; }
495:   for (i = 0; i < numVerts; ++i) { vertexCoords[i] = lower + (upper-lower)*((PetscReal)i/(PetscReal)numCells); }
496:   DMPlexCreateFromDAG(*dm,1,numPoints,coneSize,cones,coneOrientations,vertexCoords);
497:   PetscFree4(coneSize,cones,coneOrientations,vertexCoords);

499:   PetscOptionsGetBool(((PetscObject)*dm)->options,((PetscObject)*dm)->prefix,"-dm_plex_separate_marker",&markerSeparate,NULL);
500:   if (markerSeparate) { markerLeft = faceMarkerLeft; markerRight = faceMarkerRight;}
501:   if (!wrap && !rank) {
502:     DMPlexGetHeightStratum(*dm,1,&fStart,&fEnd);
503:     DMSetLabelValue(*dm,"marker",fStart,markerLeft);
504:     DMSetLabelValue(*dm,"marker",fEnd-1,markerRight);
505:     DMSetLabelValue(*dm,"Face Sets",fStart,faceMarkerLeft);
506:     DMSetLabelValue(*dm,"Face Sets",fEnd-1,faceMarkerRight);
507:   }
508:   if (wrap) {
509:     L       = upper - lower;
510:     maxCell = (PetscReal)1.1*(L/(PetscReal)PetscMax(1,segments));
511:     DMSetPeriodicity(*dm,PETSC_TRUE,&maxCell,&L,&bd);
512:   }
513:   return(0);
514: }

516: 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)
517: {
518:   DM             boundary;
519:   PetscInt       i;

524:   for (i = 0; i < dim; ++i) if (periodicity[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity is not supported for simplex meshes");
525:   DMCreate(comm, &boundary);
527:   DMSetType(boundary, DMPLEX);
528:   DMSetDimension(boundary, dim-1);
529:   DMSetCoordinateDim(boundary, dim);
530:   switch (dim) {
531:   case 2: DMPlexCreateSquareBoundary(boundary, lower, upper, faces);break;
532:   case 3: DMPlexCreateCubeBoundary(boundary, lower, upper, faces);break;
533:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %D", dim);
534:   }
535:   DMPlexGenerate(boundary, NULL, interpolate, dm);
536:   DMDestroy(&boundary);
537:   return(0);
538: }

540: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
541: {
542:   DMLabel        cutLabel = NULL;
543:   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
544:   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
545:   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
546:   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
547:   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
548:   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
549:   PetscInt       dim;
550:   PetscBool      markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
551:   PetscMPIInt    rank;

555:   DMGetDimension(dm,&dim);
556:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
557:   DMCreateLabel(dm,"marker");
558:   DMCreateLabel(dm,"Face Sets");
559:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);
560:   if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ||
561:       bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ||
562:       bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {

564:     if (cutMarker) {DMCreateLabel(dm, "periodic_cut"); DMGetLabel(dm, "periodic_cut", &cutLabel);}
565:   }
566:   switch (dim) {
567:   case 2:
568:     faceMarkerTop    = 3;
569:     faceMarkerBottom = 1;
570:     faceMarkerRight  = 2;
571:     faceMarkerLeft   = 4;
572:     break;
573:   case 3:
574:     faceMarkerBottom = 1;
575:     faceMarkerTop    = 2;
576:     faceMarkerFront  = 3;
577:     faceMarkerBack   = 4;
578:     faceMarkerRight  = 5;
579:     faceMarkerLeft   = 6;
580:     break;
581:   default:
582:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D not supported",dim);
583:   }
584:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
585:   if (markerSeparate) {
586:     markerBottom = faceMarkerBottom;
587:     markerTop    = faceMarkerTop;
588:     markerFront  = faceMarkerFront;
589:     markerBack   = faceMarkerBack;
590:     markerRight  = faceMarkerRight;
591:     markerLeft   = faceMarkerLeft;
592:   }
593:   {
594:     const PetscInt numXEdges    = !rank ? edges[0] : 0;
595:     const PetscInt numYEdges    = !rank ? edges[1] : 0;
596:     const PetscInt numZEdges    = !rank ? edges[2] : 0;
597:     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
598:     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
599:     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
600:     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
601:     const PetscInt numXFaces    = numYEdges*numZEdges;
602:     const PetscInt numYFaces    = numXEdges*numZEdges;
603:     const PetscInt numZFaces    = numXEdges*numYEdges;
604:     const PetscInt numTotXFaces = numXVertices*numXFaces;
605:     const PetscInt numTotYFaces = numYVertices*numYFaces;
606:     const PetscInt numTotZFaces = numZVertices*numZFaces;
607:     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
608:     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
609:     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
610:     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
611:     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
612:     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
613:     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
614:     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
615:     const PetscInt firstYFace   = firstXFace + numTotXFaces;
616:     const PetscInt firstZFace   = firstYFace + numTotYFaces;
617:     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
618:     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
619:     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
620:     Vec            coordinates;
621:     PetscSection   coordSection;
622:     PetscScalar   *coords;
623:     PetscInt       coordSize;
624:     PetscInt       v, vx, vy, vz;
625:     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;

627:     DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);
628:     for (c = 0; c < numCells; c++) {
629:       DMPlexSetConeSize(dm, c, 6);
630:     }
631:     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
632:       DMPlexSetConeSize(dm, f, 4);
633:     }
634:     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
635:       DMPlexSetConeSize(dm, e, 2);
636:     }
637:     DMSetUp(dm); /* Allocate space for cones */
638:     /* Build cells */
639:     for (fz = 0; fz < numZEdges; ++fz) {
640:       for (fy = 0; fy < numYEdges; ++fy) {
641:         for (fx = 0; fx < numXEdges; ++fx) {
642:           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
643:           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
644:           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
645:           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
646:           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
647:           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
648:           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
649:                             /* B,  T,  F,  K,  R,  L */
650:           PetscInt ornt[6] = {-4,  0,  0, -1,  0, -4}; /* ??? */
651:           PetscInt cone[6];

653:           /* no boundary twisting in 3D */
654:           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
655:           DMPlexSetCone(dm, cell, cone);
656:           DMPlexSetConeOrientation(dm, cell, ornt);
657:           if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
658:           if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
659:           if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
660:         }
661:       }
662:     }
663:     /* Build x faces */
664:     for (fz = 0; fz < numZEdges; ++fz) {
665:       for (fy = 0; fy < numYEdges; ++fy) {
666:         for (fx = 0; fx < numXVertices; ++fx) {
667:           PetscInt face    = firstXFace + (fz*numYEdges+fy)     *numXVertices+fx;
668:           PetscInt edgeL   = firstZEdge + (fy                   *numXVertices+fx)*numZEdges + fz;
669:           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
670:           PetscInt edgeB   = firstYEdge + (fz                   *numXVertices+fx)*numYEdges + fy;
671:           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
672:           PetscInt ornt[4] = {0, 0, -2, -2};
673:           PetscInt cone[4];

675:           if (dim == 3) {
676:             /* markers */
677:             if (bdX != DM_BOUNDARY_PERIODIC) {
678:               if (fx == numXVertices-1) {
679:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);
680:                 DMSetLabelValue(dm, "marker", face, markerRight);
681:               }
682:               else if (fx == 0) {
683:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);
684:                 DMSetLabelValue(dm, "marker", face, markerLeft);
685:               }
686:             }
687:           }
688:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
689:           DMPlexSetCone(dm, face, cone);
690:           DMPlexSetConeOrientation(dm, face, ornt);
691:         }
692:       }
693:     }
694:     /* Build y faces */
695:     for (fz = 0; fz < numZEdges; ++fz) {
696:       for (fx = 0; fx < numXEdges; ++fx) {
697:         for (fy = 0; fy < numYVertices; ++fy) {
698:           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
699:           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx)*numZEdges + fz;
700:           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
701:           PetscInt edgeB   = firstXEdge + (fz                   *numYVertices+fy)*numXEdges + fx;
702:           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
703:           PetscInt ornt[4] = {0, 0, -2, -2};
704:           PetscInt cone[4];

706:           if (dim == 3) {
707:             /* markers */
708:             if (bdY != DM_BOUNDARY_PERIODIC) {
709:               if (fy == numYVertices-1) {
710:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);
711:                 DMSetLabelValue(dm, "marker", face, markerBack);
712:               }
713:               else if (fy == 0) {
714:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);
715:                 DMSetLabelValue(dm, "marker", face, markerFront);
716:               }
717:             }
718:           }
719:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
720:           DMPlexSetCone(dm, face, cone);
721:           DMPlexSetConeOrientation(dm, face, ornt);
722:         }
723:       }
724:     }
725:     /* Build z faces */
726:     for (fy = 0; fy < numYEdges; ++fy) {
727:       for (fx = 0; fx < numXEdges; ++fx) {
728:         for (fz = 0; fz < numZVertices; fz++) {
729:           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
730:           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx)*numYEdges + fy;
731:           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
732:           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy)*numXEdges + fx;
733:           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
734:           PetscInt ornt[4] = {0, 0, -2, -2};
735:           PetscInt cone[4];

737:           if (dim == 2) {
738:             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
739:             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
740:             if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, face, 2);}
741:             if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, face, 2);}
742:           } else {
743:             /* markers */
744:             if (bdZ != DM_BOUNDARY_PERIODIC) {
745:               if (fz == numZVertices-1) {
746:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);
747:                 DMSetLabelValue(dm, "marker", face, markerTop);
748:               }
749:               else if (fz == 0) {
750:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);
751:                 DMSetLabelValue(dm, "marker", face, markerBottom);
752:               }
753:             }
754:           }
755:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
756:           DMPlexSetCone(dm, face, cone);
757:           DMPlexSetConeOrientation(dm, face, ornt);
758:         }
759:       }
760:     }
761:     /* Build Z edges*/
762:     for (vy = 0; vy < numYVertices; vy++) {
763:       for (vx = 0; vx < numXVertices; vx++) {
764:         for (ez = 0; ez < numZEdges; ez++) {
765:           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
766:           const PetscInt vertexB = firstVertex + (ez                   *numYVertices+vy)*numXVertices + vx;
767:           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
768:           PetscInt       cone[2];

770:           if (dim == 3) {
771:             if (bdX != DM_BOUNDARY_PERIODIC) {
772:               if (vx == numXVertices-1) {
773:                 DMSetLabelValue(dm, "marker", edge, markerRight);
774:               }
775:               else if (vx == 0) {
776:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
777:               }
778:             }
779:             if (bdY != DM_BOUNDARY_PERIODIC) {
780:               if (vy == numYVertices-1) {
781:                 DMSetLabelValue(dm, "marker", edge, markerBack);
782:               }
783:               else if (vy == 0) {
784:                 DMSetLabelValue(dm, "marker", edge, markerFront);
785:               }
786:             }
787:           }
788:           cone[0] = vertexB; cone[1] = vertexT;
789:           DMPlexSetCone(dm, edge, cone);
790:         }
791:       }
792:     }
793:     /* Build Y edges*/
794:     for (vz = 0; vz < numZVertices; vz++) {
795:       for (vx = 0; vx < numXVertices; vx++) {
796:         for (ey = 0; ey < numYEdges; ey++) {
797:           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
798:           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
799:           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
800:           const PetscInt vertexK = firstVertex + nextv;
801:           PetscInt       cone[2];

803:           cone[0] = vertexF; cone[1] = vertexK;
804:           DMPlexSetCone(dm, edge, cone);
805:           if (dim == 2) {
806:             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
807:               if (vx == numXVertices-1) {
808:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);
809:                 DMSetLabelValue(dm, "marker", edge,    markerRight);
810:                 DMSetLabelValue(dm, "marker", cone[0], markerRight);
811:                 if (ey == numYEdges-1) {
812:                   DMSetLabelValue(dm, "marker", cone[1], markerRight);
813:                 }
814:               } else if (vx == 0) {
815:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);
816:                 DMSetLabelValue(dm, "marker", edge,    markerLeft);
817:                 DMSetLabelValue(dm, "marker", cone[0], markerLeft);
818:                 if (ey == numYEdges-1) {
819:                   DMSetLabelValue(dm, "marker", cone[1], markerLeft);
820:                 }
821:               }
822:             } else {
823:               if (vx == 0 && cutLabel) {
824:                 DMLabelSetValue(cutLabel, edge,    1);
825:                 DMLabelSetValue(cutLabel, cone[0], 1);
826:                 if (ey == numYEdges-1) {
827:                   DMLabelSetValue(cutLabel, cone[1], 1);
828:                 }
829:               }
830:             }
831:           } else {
832:             if (bdX != DM_BOUNDARY_PERIODIC) {
833:               if (vx == numXVertices-1) {
834:                 DMSetLabelValue(dm, "marker", edge, markerRight);
835:               } else if (vx == 0) {
836:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
837:               }
838:             }
839:             if (bdZ != DM_BOUNDARY_PERIODIC) {
840:               if (vz == numZVertices-1) {
841:                 DMSetLabelValue(dm, "marker", edge, markerTop);
842:               } else if (vz == 0) {
843:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
844:               }
845:             }
846:           }
847:         }
848:       }
849:     }
850:     /* Build X edges*/
851:     for (vz = 0; vz < numZVertices; vz++) {
852:       for (vy = 0; vy < numYVertices; vy++) {
853:         for (ex = 0; ex < numXEdges; ex++) {
854:           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
855:           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
856:           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
857:           const PetscInt vertexR = firstVertex + nextv;
858:           PetscInt       cone[2];

860:           cone[0] = vertexL; cone[1] = vertexR;
861:           DMPlexSetCone(dm, edge, cone);
862:           if (dim == 2) {
863:             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
864:               if (vy == numYVertices-1) {
865:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);
866:                 DMSetLabelValue(dm, "marker", edge,    markerTop);
867:                 DMSetLabelValue(dm, "marker", cone[0], markerTop);
868:                 if (ex == numXEdges-1) {
869:                   DMSetLabelValue(dm, "marker", cone[1], markerTop);
870:                 }
871:               } else if (vy == 0) {
872:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);
873:                 DMSetLabelValue(dm, "marker", edge,    markerBottom);
874:                 DMSetLabelValue(dm, "marker", cone[0], markerBottom);
875:                 if (ex == numXEdges-1) {
876:                   DMSetLabelValue(dm, "marker", cone[1], markerBottom);
877:                 }
878:               }
879:             } else {
880:               if (vy == 0 && cutLabel) {
881:                 DMLabelSetValue(cutLabel, edge,    1);
882:                 DMLabelSetValue(cutLabel, cone[0], 1);
883:                 if (ex == numXEdges-1) {
884:                   DMLabelSetValue(cutLabel, cone[1], 1);
885:                 }
886:               }
887:             }
888:           } else {
889:             if (bdY != DM_BOUNDARY_PERIODIC) {
890:               if (vy == numYVertices-1) {
891:                 DMSetLabelValue(dm, "marker", edge, markerBack);
892:               }
893:               else if (vy == 0) {
894:                 DMSetLabelValue(dm, "marker", edge, markerFront);
895:               }
896:             }
897:             if (bdZ != DM_BOUNDARY_PERIODIC) {
898:               if (vz == numZVertices-1) {
899:                 DMSetLabelValue(dm, "marker", edge, markerTop);
900:               }
901:               else if (vz == 0) {
902:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
903:               }
904:             }
905:           }
906:         }
907:       }
908:     }
909:     DMPlexSymmetrize(dm);
910:     DMPlexStratify(dm);
911:     /* Build coordinates */
912:     DMGetCoordinateSection(dm, &coordSection);
913:     PetscSectionSetNumFields(coordSection, 1);
914:     PetscSectionSetFieldComponents(coordSection, 0, dim);
915:     PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);
916:     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
917:       PetscSectionSetDof(coordSection, v, dim);
918:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
919:     }
920:     PetscSectionSetUp(coordSection);
921:     PetscSectionGetStorageSize(coordSection, &coordSize);
922:     VecCreate(PETSC_COMM_SELF, &coordinates);
923:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
924:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
925:     VecSetBlockSize(coordinates, dim);
926:     VecSetType(coordinates,VECSTANDARD);
927:     VecGetArray(coordinates, &coords);
928:     for (vz = 0; vz < numZVertices; ++vz) {
929:       for (vy = 0; vy < numYVertices; ++vy) {
930:         for (vx = 0; vx < numXVertices; ++vx) {
931:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
932:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
933:           if (dim == 3) {
934:             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
935:           }
936:         }
937:       }
938:     }
939:     VecRestoreArray(coordinates, &coords);
940:     DMSetCoordinatesLocal(dm, coordinates);
941:     VecDestroy(&coordinates);
942:   }
943:   return(0);
944: }

946: 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)
947: {
948:   PetscInt       i;

953:   DMCreate(comm, dm);
955:   DMSetType(*dm, DMPLEX);
956:   DMSetDimension(*dm, dim);
957:   switch (dim) {
958:   case 2: {DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], DM_BOUNDARY_NONE);break;}
959:   case 3: {DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], periodicity[2]);break;}
960:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %D", dim);
961:   }
962:   if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST ||
963:       periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST ||
964:       (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) {
965:     PetscReal L[3];
966:     PetscReal maxCell[3];

968:     for (i = 0; i < dim; i++) {
969:       L[i]       = upper[i] - lower[i];
970:       maxCell[i] = 1.1 * (L[i] / PetscMax(1,faces[i]));
971:     }
972:     DMSetPeriodicity(*dm,PETSC_TRUE,maxCell,L,periodicity);
973:   }
974:   if (!interpolate) {
975:     DM udm;

977:     DMPlexUninterpolate(*dm, &udm);
978:     DMPlexCopyCoordinates(*dm, udm);
979:     DMDestroy(dm);
980:     *dm  = udm;
981:   }
982:   return(0);
983: }

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

988:   Collective

990:   Input Parameters:
991: + comm        - The communicator for the DM object
992: . dim         - The spatial dimension
993: . simplex     - PETSC_TRUE for simplices, PETSC_FALSE for tensor cells
994: . faces       - Number of faces per dimension, or NULL for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
995: . lower       - The lower left corner, or NULL for (0, 0, 0)
996: . upper       - The upper right corner, or NULL for (1, 1, 1)
997: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
998: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1000:   Output Parameter:
1001: . dm  - The DM object

1003:   Options Database Keys:
1004:   These options override the hard-wired input values.
1005: + -dm_plex_box_dim <dim>          - Set the topological dimension
1006: . -dm_plex_box_simplex <bool>     - PETSC_TRUE for simplex elements, PETSC_FALSE for tensor elements
1007: . -dm_plex_box_lower <x,y,z>      - Specify lower-left-bottom coordinates for the box
1008: . -dm_plex_box_upper <x,y,z>      - Specify upper-right-top coordinates for the box
1009: . -dm_plex_box_faces <m,n,p>      - Number of faces in each linear direction
1010: . -dm_plex_box_bd <bx,by,bz>      - Specify the DMBoundaryType for each direction
1011: - -dm_plex_box_interpolate <bool> - PETSC_TRUE turns on topological interpolation (creating edges and faces)

1013:   Notes:
1014:   The options database keys above take lists of length d in d dimensions.

1016:   Here is the numbering returned for 2 faces in each direction for tensor cells:
1017: $ 10---17---11---18----12
1018: $  |         |         |
1019: $  |         |         |
1020: $ 20    2   22    3    24
1021: $  |         |         |
1022: $  |         |         |
1023: $  7---15----8---16----9
1024: $  |         |         |
1025: $  |         |         |
1026: $ 19    0   21    1   23
1027: $  |         |         |
1028: $  |         |         |
1029: $  4---13----5---14----6

1031: and for simplicial cells

1033: $ 14----8---15----9----16
1034: $  |\     5  |\      7 |
1035: $  | \       | \       |
1036: $ 13   2    14    3    15
1037: $  | 4   \   | 6   \   |
1038: $  |       \ |       \ |
1039: $ 11----6---12----7----13
1040: $  |\        |\        |
1041: $  | \    1  | \     3 |
1042: $ 10   0    11    1    12
1043: $  | 0   \   | 2   \   |
1044: $  |       \ |       \ |
1045: $  8----4----9----5----10

1047:   Level: beginner

1049: .seealso: DMPlexCreateFromFile(), DMPlexCreateHexCylinderMesh(), DMSetType(), DMCreate()
1050: @*/
1051: 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)
1052: {
1053:   PetscInt       fac[3] = {0, 0, 0};
1054:   PetscReal      low[3] = {0, 0, 0};
1055:   PetscReal      upp[3] = {1, 1, 1};
1056:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1057:   PetscInt       i, n;
1058:   PetscBool      flg;

1062:   PetscOptionsGetInt(NULL, NULL, "-dm_plex_box_dim", &dim, &flg);
1063:   if ((dim < 0) || (dim > 3)) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D should be in [1, 3]", dim);
1064:   PetscOptionsGetBool(NULL, NULL, "-dm_plex_box_simplex", &simplex, &flg);
1065:   n    = dim;
1066:   PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", fac, &n, &flg);
1067:   for (i = 0; i < dim; ++i) fac[i] = faces ? faces[i] : (flg && i < n ? fac[i] : (dim == 1 ? 1 : 4-dim));
1068:   if (lower) for (i = 0; i < dim; ++i) low[i] = lower[i];
1069:   if (upper) for (i = 0; i < dim; ++i) upp[i] = upper[i];
1070:   if (periodicity) for (i = 0; i < dim; ++i) bdt[i] = periodicity[i];
1071:   /* Allow bounds to be specified from the command line */
1072:   n    = 3;
1073:   PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_lower", low, &n, &flg);
1074:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Lower box point had %D values, should have been %D", n, dim);
1075:   n    = 3;
1076:   PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_upper", upp, &n, &flg);
1077:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Upper box point had %D values, should have been %D", n, dim);
1078:   n    = 3;
1079:   PetscOptionsGetEnumArray(NULL, NULL, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, &flg);
1080:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Box boundary types had %D values, should have been %D", n, dim);
1081:   PetscOptionsGetBool(NULL, NULL, "-dm_plex_box_interpolate", &interpolate, &flg);

1083:   if (dim == 1)      {DMPlexCreateLineMesh_Internal(comm, fac[0], low[0], upp[0], bdt[0], dm);}
1084:   else if (simplex)  {DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1085:   else               {DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1086:   return(0);
1087: }

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

1092:   Collective

1094:   Input Parameters:
1095: + comm        - The communicator for the DM object
1096: . faces       - Number of faces per dimension, or NULL for (1, 1, 1)
1097: . lower       - The lower left corner, or NULL for (0, 0, 0)
1098: . upper       - The upper right corner, or NULL for (1, 1, 1)
1099: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
1100: . orderHeight - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1101: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1103:   Output Parameter:
1104: . dm  - The DM object

1106:   Level: beginner

1108: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateWedgeCylinderMesh(), DMPlexExtrude(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1109: @*/
1110: PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool orderHeight, PetscBool interpolate, DM *dm)
1111: {
1112:   DM             bdm, botdm;
1113:   PetscInt       i;
1114:   PetscInt       fac[3] = {0, 0, 0};
1115:   PetscReal      low[3] = {0, 0, 0};
1116:   PetscReal      upp[3] = {1, 1, 1};
1117:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1121:   for (i = 0; i < 3; ++i) fac[i] = faces ? (faces[i] > 0 ? faces[i] : 1) : 1;
1122:   if (lower) for (i = 0; i < 3; ++i) low[i] = lower[i];
1123:   if (upper) for (i = 0; i < 3; ++i) upp[i] = upper[i];
1124:   if (periodicity) for (i = 0; i < 3; ++i) bdt[i] = periodicity[i];
1125:   for (i = 0; i < 3; ++i) if (bdt[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity not yet supported");

1127:   DMCreate(comm, &bdm);
1128:   DMSetType(bdm, DMPLEX);
1129:   DMSetDimension(bdm, 1);
1130:   DMSetCoordinateDim(bdm, 2);
1131:   DMPlexCreateSquareBoundary(bdm, low, upp, fac);
1132:   DMPlexGenerate(bdm, NULL, PETSC_FALSE, &botdm);
1133:   DMDestroy(&bdm);
1134:   DMPlexExtrude(botdm, fac[2], upp[2] - low[2], orderHeight, NULL, interpolate, dm);
1135:   if (low[2] != 0.0) {
1136:     Vec         v;
1137:     PetscScalar *x;
1138:     PetscInt    cDim, n;

1140:     DMGetCoordinatesLocal(*dm, &v);
1141:     VecGetBlockSize(v, &cDim);
1142:     VecGetLocalSize(v, &n);
1143:     VecGetArray(v, &x);
1144:     x   += cDim;
1145:     for (i=0; i<n; i+=cDim) x[i] += low[2];
1146:     VecRestoreArray(v,&x);
1147:     DMSetCoordinatesLocal(*dm, v);
1148:   }
1149:   DMDestroy(&botdm);
1150:   return(0);
1151: }

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

1156:   Collective on idm

1158:   Input Parameters:
1159: + idm         - The mesh to be extruded
1160: . layers      - The number of layers, or PETSC_DETERMINE to use the default
1161: . height      - The total height of the extrusion, or PETSC_DETERMINE to use the default
1162: . orderHeight - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1163: . extNormal   - The normal direction in which the mesh should be extruded, or NULL to extrude using the surface normal
1164: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1166:   Output Parameter:
1167: . dm  - The DM object

1169:   Notes:
1170:   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.

1172:   Options Database Keys:
1173: +   -dm_plex_extrude_layers <k> - Sets the number of layers k
1174: .   -dm_plex_extrude_height <h> - Sets the total height of the extrusion
1175: .   -dm_plex_extrude_heights <h0,h1,...> - Sets the height of each layer
1176: .   -dm_plex_extrude_order_height - If true, order cells by height first
1177: -   -dm_plex_extrude_normal <n0,...,nd> - Sets the normal vector along which to extrude

1179:   Level: advanced

1181: .seealso: DMPlexCreateWedgeCylinderMesh(), DMPlexCreateWedgeBoxMesh(), DMSetType(), DMCreate()
1182: @*/
1183: PetscErrorCode DMPlexExtrude(DM idm, PetscInt layers, PetscReal height, PetscBool orderHeight, const PetscReal extNormal[], PetscBool interpolate, DM* dm)
1184: {
1185:   PetscScalar       *coordsB;
1186:   const PetscScalar *coordsA;
1187:   PetscReal         *normals = NULL, *heights = NULL;
1188:   PetscReal         clNormal[3];
1189:   Vec               coordinatesA, coordinatesB;
1190:   PetscSection      coordSectionA, coordSectionB;
1191:   PetscInt          dim, cDim, cDimB, c, l, v, coordSize, *newCone, nl;
1192:   PetscInt          cStart, cEnd, vStart, vEnd, cellV, numCells, numVertices;
1193:   const char       *prefix;
1194:   PetscBool         haveCLNormal, flg;
1195:   PetscErrorCode    ierr;

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

1207:   PetscObjectGetOptionsPrefix((PetscObject) idm, &prefix);
1208:   if (layers < 0) layers = 1;
1209:   PetscOptionsGetInt(NULL, prefix, "-dm_plex_extrude_layers", &layers, NULL);
1210:   if (layers <= 0) SETERRQ1(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Number of layers %D must be positive", layers);
1211:   if (height < 0.) height = 1.;
1212:   PetscOptionsGetReal(NULL, prefix, "-dm_plex_extrude_height", &height, NULL);
1213:   if (height <= 0.) SETERRQ1(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Height of layers %g must be positive", (double) height);
1214:   PetscMalloc1(layers, &heights);
1215:   nl   = layers;
1216:   PetscOptionsGetRealArray(NULL, prefix, "-dm_plex_extrude_heights", heights, &nl, &flg);
1217:   if (flg) {
1218:     if (!nl) SETERRQ(PetscObjectComm((PetscObject) idm), PETSC_ERR_ARG_OUTOFRANGE, "Must give at least one height for -dm_plex_extrude_heights");
1219:     for (l = nl; l < layers; ++l) heights[l] = heights[l-1];
1220:     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);
1221:   } else {
1222:     for (l = 0; l < layers; ++l) heights[l] = height/layers;
1223:   }
1224:   PetscOptionsGetBool(NULL, prefix, "-dm_plex_extrude_order_height", &orderHeight, NULL);
1225:   c = 3;
1226:   PetscOptionsGetRealArray(NULL, prefix, "-dm_plex_extrude_normal", clNormal, &c, &haveCLNormal);
1227:   if (haveCLNormal && c != cDimB) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_ARG_SIZ, "Input normal has size %D != %D extruded coordinate dimension", c, cDimB);

1229:   DMPlexGetHeightStratum(idm, 0, &cStart, &cEnd);
1230:   DMPlexGetDepthStratum(idm, 0, &vStart, &vEnd);
1231:   numCells = (cEnd - cStart)*layers;
1232:   numVertices = (vEnd - vStart)*(layers+1);
1233:   DMCreate(PetscObjectComm((PetscObject)idm), dm);
1234:   DMSetType(*dm, DMPLEX);
1235:   DMSetDimension(*dm, dim+1);
1236:   DMPlexSetChart(*dm, 0, numCells+numVertices);
1237:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1238:   DMCreateLabel(*dm, "celltype");
1239:   for (c = cStart, cellV = 0; c < cEnd; ++c) {
1240:     DMPolytopeType ct, nct;
1241:     PetscInt      *closure = NULL;
1242:     PetscInt       closureSize, numCorners = 0;

1244:     DMPlexGetCellType(idm, c, &ct);
1245:     switch (ct) {
1246:       case DM_POLYTOPE_SEGMENT:       nct = DM_POLYTOPE_SEG_PRISM_TENSOR;break;
1247:       case DM_POLYTOPE_TRIANGLE:      nct = DM_POLYTOPE_TRI_PRISM_TENSOR;break;
1248:       case DM_POLYTOPE_QUADRILATERAL: nct = DM_POLYTOPE_QUAD_PRISM_TENSOR;break;
1249:       default: nct = DM_POLYTOPE_UNKNOWN;
1250:     }
1251:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1252:     for (v = 0; v < closureSize*2; v += 2) if ((closure[v] >= vStart) && (closure[v] < vEnd)) numCorners++;
1253:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1254:     for (l = 0; l < layers; ++l) {
1255:       const PetscInt cell = orderHeight ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart;

1257:       DMPlexSetConeSize(*dm, cell, 2*numCorners);
1258:       DMPlexSetCellType(*dm, cell, nct);
1259:     }
1260:     cellV = PetscMax(numCorners,cellV);
1261:   }
1262:   DMSetUp(*dm);

1264:   if (dim != cDim && !(extNormal || haveCLNormal)) {PetscCalloc1(cDim*(vEnd - vStart), &normals);}
1265:   PetscMalloc1(3*cellV,&newCone);
1266:   for (c = cStart; c < cEnd; ++c) {
1267:     PetscInt *closure = NULL;
1268:     PetscInt closureSize, numCorners = 0, l;
1269:     PetscReal normal[3] = {0, 0, 0};

1271:     if (normals) {DMPlexComputeCellGeometryFVM(idm, c, NULL, NULL, normal);}
1272:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1273:     for (v = 0; v < closureSize*2; v += 2) {
1274:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1275:         PetscInt d;

1277:         newCone[numCorners++] = closure[v] - vStart;
1278:         if (normals) {for (d = 0; d < cDim; ++d) normals[cDim*(closure[v]-vStart)+d] += normal[d];}
1279:       }
1280:     }
1281:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1282:     for (l = 0; l < layers; ++l) {
1283:       PetscInt i;

1285:       for (i = 0; i < numCorners; ++i) {
1286:         newCone[  numCorners + i] = orderHeight ? (layers+1)*newCone[i] + l     + numCells :     l*(vEnd - vStart) + newCone[i] + numCells;
1287:         newCone[2*numCorners + i] = orderHeight ? (layers+1)*newCone[i] + l + 1 + numCells : (l+1)*(vEnd - vStart) + newCone[i] + numCells;
1288:       }
1289:       DMPlexSetCone(*dm, orderHeight ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart, newCone + numCorners);
1290:     }
1291:   }
1292:   DMPlexSymmetrize(*dm);
1293:   DMPlexStratify(*dm);
1294:   PetscFree(newCone);

1296:   DMGetCoordinateSection(*dm, &coordSectionB);
1297:   PetscSectionSetNumFields(coordSectionB, 1);
1298:   PetscSectionSetFieldComponents(coordSectionB, 0, cDimB);
1299:   PetscSectionSetChart(coordSectionB, numCells, numCells+numVertices);
1300:   for (v = numCells; v < numCells+numVertices; ++v) {
1301:     PetscSectionSetDof(coordSectionB, v, cDimB);
1302:     PetscSectionSetFieldDof(coordSectionB, v, 0, cDimB);
1303:     DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);
1304:   }
1305:   PetscSectionSetUp(coordSectionB);
1306:   PetscSectionGetStorageSize(coordSectionB, &coordSize);
1307:   VecCreate(PETSC_COMM_SELF, &coordinatesB);
1308:   PetscObjectSetName((PetscObject) coordinatesB, "coordinates");
1309:   VecSetSizes(coordinatesB, coordSize, PETSC_DETERMINE);
1310:   VecSetBlockSize(coordinatesB, cDimB);
1311:   VecSetType(coordinatesB,VECSTANDARD);

1313:   DMGetCoordinateSection(idm, &coordSectionA);
1314:   DMGetCoordinatesLocal(idm, &coordinatesA);
1315:   VecGetArray(coordinatesB, &coordsB);
1316:   VecGetArrayRead(coordinatesA, &coordsA);
1317:   for (v = vStart; v < vEnd; ++v) {
1318:     const PetscScalar *cptr;
1319:     PetscReal         ones2[2] = { 0., 1.}, ones3[3] = { 0., 0., 1.};
1320:     PetscReal         normal[3];
1321:     PetscReal         norm;
1322:     PetscInt          offA, d, cDimA = cDim;

1324:     if (normals)           {for (d = 0; d < cDimB; ++d) normal[d] = normals[cDimB*(v - vStart)+d];}
1325:     else if (haveCLNormal) {for (d = 0; d < cDimB; ++d) normal[d] = clNormal[d];}
1326:     else if (extNormal)    {for (d = 0; d < cDimB; ++d) normal[d] = extNormal[d];}
1327:     else if (cDimB == 2)   {for (d = 0; d < cDimB; ++d) normal[d] = ones2[d];}
1328:     else if (cDimB == 3)   {for (d = 0; d < cDimB; ++d) normal[d] = ones3[d];}
1329:     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to determine normal for extrusion");
1330:     for (d = 0, norm = 0.0; d < cDimB; ++d) norm += normal[d]*normal[d];
1331:     for (d = 0; d < cDimB; ++d) normal[d] *= 1./PetscSqrtReal(norm);

1333:     PetscSectionGetOffset(coordSectionA, v, &offA);
1334:     cptr = coordsA + offA;
1335:     for (l = 0; l <= layers; ++l) {
1336:       PetscInt offB, d, newV;

1338:       newV = orderHeight ? (layers+1)*(v -vStart) + l + numCells : (vEnd -vStart)*l + (v -vStart) + numCells;
1339:       PetscSectionGetOffset(coordSectionB, newV, &offB);
1340:       for (d = 0; d < cDimA; ++d) { coordsB[offB+d]  = cptr[d]; }
1341:       for (d = 0; d < cDimB; ++d) { coordsB[offB+d] += l ? normal[d]*heights[l-1] : 0.0; }
1342:       cptr    = coordsB + offB;
1343:       cDimA   = cDimB;
1344:     }
1345:   }
1346:   VecRestoreArrayRead(coordinatesA, &coordsA);
1347:   VecRestoreArray(coordinatesB, &coordsB);
1348:   DMSetCoordinatesLocal(*dm, coordinatesB);
1349:   VecDestroy(&coordinatesB);
1350:   PetscFree(normals);
1351:   PetscFree(heights);
1352:   if (interpolate) {
1353:     DM idm;

1355:     DMPlexInterpolate(*dm, &idm);
1356:     DMPlexCopyCoordinates(*dm, idm);
1357:     DMDestroy(dm);
1358:     *dm  = idm;
1359:   }
1360:   return(0);
1361: }

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

1366:   Logically Collective on dm

1368:   Input Parameters:
1369: + dm - the DM context
1370: - prefix - the prefix to prepend to all option names

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

1376:   Level: advanced

1378: .seealso: SNESSetFromOptions()
1379: @*/
1380: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1381: {
1382:   DM_Plex       *mesh = (DM_Plex *) dm->data;

1387:   PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);
1388:   PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);
1389:   return(0);
1390: }

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

1395:   Collective

1397:   Input Parameters:
1398: + comm      - The communicator for the DM object
1399: . numRefine - The number of regular refinements to the basic 5 cell structure
1400: - periodicZ - The boundary type for the Z direction

1402:   Output Parameter:
1403: . dm  - The DM object

1405:   Options Database Keys:
1406: These options override the hard-wired input values.
1407: + -dm_plex_hex_cyl_refine <r> - Refine the sylinder r times
1408: - -dm_plex_hex_cyl_bd <bz>    - Specify the DMBoundaryType in the z-direction

1410:   Note:
1411:   Here is the output numbering looking from the bottom of the cylinder:
1412: $       17-----14
1413: $        |     |
1414: $        |  2  |
1415: $        |     |
1416: $ 17-----8-----7-----14
1417: $  |     |     |     |
1418: $  |  3  |  0  |  1  |
1419: $  |     |     |     |
1420: $ 19-----5-----6-----13
1421: $        |     |
1422: $        |  4  |
1423: $        |     |
1424: $       19-----13
1425: $
1426: $ and up through the top
1427: $
1428: $       18-----16
1429: $        |     |
1430: $        |  2  |
1431: $        |     |
1432: $ 18----10----11-----16
1433: $  |     |     |     |
1434: $  |  3  |  0  |  1  |
1435: $  |     |     |     |
1436: $ 20-----9----12-----15
1437: $        |     |
1438: $        |  4  |
1439: $        |     |
1440: $       20-----15

1442:   Level: beginner

1444: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1445: @*/
1446: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
1447: {
1448:   const PetscInt dim = 3;
1449:   PetscInt       numCells, numVertices, r;
1450:   PetscMPIInt    rank;

1455:   MPI_Comm_rank(comm, &rank);
1456:   PetscOptionsGetInt(NULL, NULL, "-dm_plex_hex_cyl_refine", &numRefine, NULL);
1457:   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
1458:   PetscOptionsGetEnum(NULL, NULL, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) &periodicZ, NULL);
1459:   DMCreate(comm, dm);
1460:   DMSetType(*dm, DMPLEX);
1461:   DMSetDimension(*dm, dim);
1462:   /* Create topology */
1463:   {
1464:     PetscInt cone[8], c;

1466:     numCells    = !rank ?  5 : 0;
1467:     numVertices = !rank ? 16 : 0;
1468:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1469:       numCells   *= 3;
1470:       numVertices = !rank ? 24 : 0;
1471:     }
1472:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1473:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 8);}
1474:     DMSetUp(*dm);
1475:     if (!rank) {
1476:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1477:         cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1478:         cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1479:         DMPlexSetCone(*dm, 0, cone);
1480:         cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1481:         cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1482:         DMPlexSetCone(*dm, 1, cone);
1483:         cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1484:         cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1485:         DMPlexSetCone(*dm, 2, cone);
1486:         cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1487:         cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1488:         DMPlexSetCone(*dm, 3, cone);
1489:         cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1490:         cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1491:         DMPlexSetCone(*dm, 4, cone);

1493:         cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1494:         cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1495:         DMPlexSetCone(*dm, 5, cone);
1496:         cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1497:         cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1498:         DMPlexSetCone(*dm, 6, cone);
1499:         cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1500:         cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1501:         DMPlexSetCone(*dm, 7, cone);
1502:         cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1503:         cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1504:         DMPlexSetCone(*dm, 8, cone);
1505:         cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1506:         cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1507:         DMPlexSetCone(*dm, 9, cone);

1509:         cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1510:         cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1511:         DMPlexSetCone(*dm, 10, cone);
1512:         cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1513:         cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1514:         DMPlexSetCone(*dm, 11, cone);
1515:         cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1516:         cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1517:         DMPlexSetCone(*dm, 12, cone);
1518:         cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1519:         cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1520:         DMPlexSetCone(*dm, 13, cone);
1521:         cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1522:         cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1523:         DMPlexSetCone(*dm, 14, cone);
1524:       } else {
1525:         cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
1526:         cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
1527:         DMPlexSetCone(*dm, 0, cone);
1528:         cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
1529:         cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
1530:         DMPlexSetCone(*dm, 1, cone);
1531:         cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
1532:         cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
1533:         DMPlexSetCone(*dm, 2, cone);
1534:         cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
1535:         cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
1536:         DMPlexSetCone(*dm, 3, cone);
1537:         cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
1538:         cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
1539:         DMPlexSetCone(*dm, 4, cone);
1540:       }
1541:     }
1542:     DMPlexSymmetrize(*dm);
1543:     DMPlexStratify(*dm);
1544:   }
1545:   /* Interpolate */
1546:   {
1547:     DM idm;

1549:     DMPlexInterpolate(*dm, &idm);
1550:     DMDestroy(dm);
1551:     *dm  = idm;
1552:   }
1553:   /* Create cube geometry */
1554:   {
1555:     Vec             coordinates;
1556:     PetscSection    coordSection;
1557:     PetscScalar    *coords;
1558:     PetscInt        coordSize, v;
1559:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1560:     const PetscReal ds2 = dis/2.0;

1562:     /* Build coordinates */
1563:     DMGetCoordinateSection(*dm, &coordSection);
1564:     PetscSectionSetNumFields(coordSection, 1);
1565:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1566:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1567:     for (v = numCells; v < numCells+numVertices; ++v) {
1568:       PetscSectionSetDof(coordSection, v, dim);
1569:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1570:     }
1571:     PetscSectionSetUp(coordSection);
1572:     PetscSectionGetStorageSize(coordSection, &coordSize);
1573:     VecCreate(PETSC_COMM_SELF, &coordinates);
1574:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1575:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1576:     VecSetBlockSize(coordinates, dim);
1577:     VecSetType(coordinates,VECSTANDARD);
1578:     VecGetArray(coordinates, &coords);
1579:     if (!rank) {
1580:       coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
1581:       coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
1582:       coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
1583:       coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
1584:       coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
1585:       coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
1586:       coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
1587:       coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
1588:       coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
1589:       coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
1590:       coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
1591:       coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
1592:       coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
1593:       coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
1594:       coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
1595:       coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1596:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1597:         /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1598:         /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1599:         /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1600:         /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1601:         /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1602:         /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1603:         /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1604:         /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1605:       }
1606:     }
1607:     VecRestoreArray(coordinates, &coords);
1608:     DMSetCoordinatesLocal(*dm, coordinates);
1609:     VecDestroy(&coordinates);
1610:   }
1611:   /* Create periodicity */
1612:   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1613:     PetscReal      L[3];
1614:     PetscReal      maxCell[3];
1615:     DMBoundaryType bdType[3];
1616:     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1617:     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1618:     PetscInt       i, numZCells = 3;

1620:     bdType[0] = DM_BOUNDARY_NONE;
1621:     bdType[1] = DM_BOUNDARY_NONE;
1622:     bdType[2] = periodicZ;
1623:     for (i = 0; i < dim; i++) {
1624:       L[i]       = upper[i] - lower[i];
1625:       maxCell[i] = 1.1 * (L[i] / numZCells);
1626:     }
1627:     DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);
1628:   }
1629:   /* Refine topology */
1630:   for (r = 0; r < numRefine; ++r) {
1631:     DM rdm = NULL;

1633:     DMRefine(*dm, comm, &rdm);
1634:     DMDestroy(dm);
1635:     *dm  = rdm;
1636:   }
1637:   /* Remap geometry to cylinder
1638:        Interior square: Linear interpolation is correct
1639:        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1640:        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance

1642:          phi     = arctan(y/x)
1643:          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1644:          d_far   = sqrt(1/2 + sin^2(phi))

1646:        so we remap them using

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

1651:        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1652:   */
1653:   {
1654:     Vec           coordinates;
1655:     PetscSection  coordSection;
1656:     PetscScalar  *coords;
1657:     PetscInt      vStart, vEnd, v;
1658:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1659:     const PetscReal ds2 = 0.5*dis;

1661:     DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);
1662:     DMGetCoordinateSection(*dm, &coordSection);
1663:     DMGetCoordinatesLocal(*dm, &coordinates);
1664:     VecGetArray(coordinates, &coords);
1665:     for (v = vStart; v < vEnd; ++v) {
1666:       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1667:       PetscInt  off;

1669:       PetscSectionGetOffset(coordSection, v, &off);
1670:       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
1671:       x    = PetscRealPart(coords[off]);
1672:       y    = PetscRealPart(coords[off+1]);
1673:       phi  = PetscAtan2Real(y, x);
1674:       sinp = PetscSinReal(phi);
1675:       cosp = PetscCosReal(phi);
1676:       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
1677:         dc = PetscAbsReal(ds2/sinp);
1678:         df = PetscAbsReal(dis/sinp);
1679:         xc = ds2*x/PetscAbsReal(y);
1680:         yc = ds2*PetscSignReal(y);
1681:       } else {
1682:         dc = PetscAbsReal(ds2/cosp);
1683:         df = PetscAbsReal(dis/cosp);
1684:         xc = ds2*PetscSignReal(x);
1685:         yc = ds2*y/PetscAbsReal(x);
1686:       }
1687:       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1688:       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1689:     }
1690:     VecRestoreArray(coordinates, &coords);
1691:     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1692:       DMLocalizeCoordinates(*dm);
1693:     }
1694:   }
1695:   return(0);
1696: }

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

1701:   Collective

1703:   Input Parameters:
1704: + comm - The communicator for the DM object
1705: . n    - The number of wedges around the origin
1706: - interpolate - Create edges and faces

1708:   Output Parameter:
1709: . dm  - The DM object

1711:   Level: beginner

1713: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1714: @*/
1715: PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
1716: {
1717:   const PetscInt dim = 3;
1718:   PetscInt       numCells, numVertices, v;
1719:   PetscMPIInt    rank;

1724:   MPI_Comm_rank(comm, &rank);
1725:   if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n);
1726:   DMCreate(comm, dm);
1727:   DMSetType(*dm, DMPLEX);
1728:   DMSetDimension(*dm, dim);
1729:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1730:   DMCreateLabel(*dm, "celltype");
1731:   /* Create topology */
1732:   {
1733:     PetscInt cone[6], c;

1735:     numCells    = !rank ?        n : 0;
1736:     numVertices = !rank ?  2*(n+1) : 0;
1737:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1738:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 6);}
1739:     DMSetUp(*dm);
1740:     for (c = 0; c < numCells; c++) {
1741:       cone[0] =  c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n;
1742:       cone[3] =  c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n;
1743:       DMPlexSetCone(*dm, c, cone);
1744:       DMPlexSetCellType(*dm, c, DM_POLYTOPE_TRI_PRISM_TENSOR);
1745:     }
1746:     DMPlexSymmetrize(*dm);
1747:     DMPlexStratify(*dm);
1748:   }
1749:   for (v = numCells; v < numCells+numVertices; ++v) {
1750:     DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);
1751:   }
1752:   /* Interpolate */
1753:   if (interpolate) {
1754:     DM idm;

1756:     DMPlexInterpolate(*dm, &idm);
1757:     DMDestroy(dm);
1758:     *dm  = idm;
1759:   }
1760:   /* Create cylinder geometry */
1761:   {
1762:     Vec          coordinates;
1763:     PetscSection coordSection;
1764:     PetscScalar *coords;
1765:     PetscInt     coordSize, c;

1767:     /* Build coordinates */
1768:     DMGetCoordinateSection(*dm, &coordSection);
1769:     PetscSectionSetNumFields(coordSection, 1);
1770:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1771:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1772:     for (v = numCells; v < numCells+numVertices; ++v) {
1773:       PetscSectionSetDof(coordSection, v, dim);
1774:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1775:     }
1776:     PetscSectionSetUp(coordSection);
1777:     PetscSectionGetStorageSize(coordSection, &coordSize);
1778:     VecCreate(PETSC_COMM_SELF, &coordinates);
1779:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1780:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1781:     VecSetBlockSize(coordinates, dim);
1782:     VecSetType(coordinates,VECSTANDARD);
1783:     VecGetArray(coordinates, &coords);
1784:     for (c = 0; c < numCells; c++) {
1785:       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;
1786:       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;
1787:     }
1788:     if (!rank) {
1789:       coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0;
1790:       coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0;
1791:     }
1792:     VecRestoreArray(coordinates, &coords);
1793:     DMSetCoordinatesLocal(*dm, coordinates);
1794:     VecDestroy(&coordinates);
1795:   }
1796:   return(0);
1797: }

1799: PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1800: {
1801:   PetscReal prod = 0.0;
1802:   PetscInt  i;
1803:   for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
1804:   return PetscSqrtReal(prod);
1805: }
1806: PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1807: {
1808:   PetscReal prod = 0.0;
1809:   PetscInt  i;
1810:   for (i = 0; i < dim; ++i) prod += x[i]*y[i];
1811:   return prod;
1812: }

1814: /* The first constant is the sphere radius */
1815: static void snapToSphere(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1816:                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1817:                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1818:                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1819: {
1820:   PetscReal r = PetscRealPart(constants[0]);
1821:   PetscReal norm2 = 0.0, fac;
1822:   PetscInt  n = uOff[1] - uOff[0], d;

1824:   for (d = 0; d < n; ++d) norm2 += PetscSqr(PetscRealPart(u[d]));
1825:   fac = r/PetscSqrtReal(norm2);
1826:   for (d = 0; d < n; ++d) f0[d] = u[d]*fac;
1827: }

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

1832:   Collective

1834:   Input Parameters:
1835: + comm    - The communicator for the DM object
1836: . dim     - The dimension
1837: . simplex - Use simplices, or tensor product cells
1838: - R       - The radius

1840:   Output Parameter:
1841: . dm  - The DM object

1843:   Level: beginner

1845: .seealso: DMPlexCreateBallMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1846: @*/
1847: PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscReal R, DM *dm)
1848: {
1849:   const PetscInt  embedDim = dim+1;
1850:   PetscSection    coordSection;
1851:   Vec             coordinates;
1852:   PetscScalar    *coords;
1853:   PetscReal      *coordsIn;
1854:   PetscInt        numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e;
1855:   PetscMPIInt     rank;
1856:   PetscErrorCode  ierr;

1860:   DMCreate(comm, dm);
1861:   DMSetType(*dm, DMPLEX);
1862:   DMSetDimension(*dm, dim);
1863:   DMSetCoordinateDim(*dm, dim+1);
1864:   MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank);
1865:   switch (dim) {
1866:   case 2:
1867:     if (simplex) {
1868:       DM              idm;
1869:       const PetscReal radius    = PetscSqrtReal(1 + PETSC_PHI*PETSC_PHI)/(1.0 + PETSC_PHI);
1870:       const PetscReal edgeLen   = 2.0/(1.0 + PETSC_PHI) * (R/radius);
1871:       const PetscInt  degree    = 5;
1872:       PetscReal       vertex[3] = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)};
1873:       PetscInt        s[3]      = {1, 1, 1};
1874:       PetscInt        cone[3];
1875:       PetscInt       *graph, p, i, j, k;

1877:       vertex[0] *= R/radius; vertex[1] *= R/radius; vertex[2] *= R/radius;
1878:       numCells    = !rank ? 20 : 0;
1879:       numVerts    = !rank ? 12 : 0;
1880:       firstVertex = numCells;
1881:       /* Use icosahedron, which for a R-sphere has coordinates which are all cyclic permutations of

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

1885:          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
1886:          length is then given by 2/(1+\phi) = 2 * 0.38197 = 0.76393.
1887:       */
1888:       /* Construct vertices */
1889:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1890:       if (!rank) {
1891:         for (p = 0, i = 0; p < embedDim; ++p) {
1892:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
1893:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
1894:               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim];
1895:               ++i;
1896:             }
1897:           }
1898:         }
1899:       }
1900:       /* Construct graph */
1901:       PetscCalloc1(numVerts * numVerts, &graph);
1902:       for (i = 0; i < numVerts; ++i) {
1903:         for (j = 0, k = 0; j < numVerts; ++j) {
1904:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
1905:         }
1906:         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid icosahedron, vertex %D degree %D != %D", i, k, degree);
1907:       }
1908:       /* Build Topology */
1909:       DMPlexSetChart(*dm, 0, numCells+numVerts);
1910:       for (c = 0; c < numCells; c++) {
1911:         DMPlexSetConeSize(*dm, c, embedDim);
1912:       }
1913:       DMSetUp(*dm); /* Allocate space for cones */
1914:       /* Cells */
1915:       for (i = 0, c = 0; i < numVerts; ++i) {
1916:         for (j = 0; j < i; ++j) {
1917:           for (k = 0; k < j; ++k) {
1918:             if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) {
1919:               cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k;
1920:               /* Check orientation */
1921:               {
1922:                 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}}};
1923:                 PetscReal normal[3];
1924:                 PetscInt  e, f;

1926:                 for (d = 0; d < embedDim; ++d) {
1927:                   normal[d] = 0.0;
1928:                   for (e = 0; e < embedDim; ++e) {
1929:                     for (f = 0; f < embedDim; ++f) {
1930:                       normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]);
1931:                     }
1932:                   }
1933:                 }
1934:                 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
1935:               }
1936:               DMPlexSetCone(*dm, c++, cone);
1937:             }
1938:           }
1939:         }
1940:       }
1941:       DMPlexSymmetrize(*dm);
1942:       DMPlexStratify(*dm);
1943:       PetscFree(graph);
1944:       /* Interpolate mesh */
1945:       DMPlexInterpolate(*dm, &idm);
1946:       DMDestroy(dm);
1947:       *dm  = idm;
1948:     } else {
1949:       /*
1950:         12-21--13
1951:          |     |
1952:         25  4  24
1953:          |     |
1954:   12-25--9-16--8-24--13
1955:    |     |     |     |
1956:   23  5 17  0 15  3  22
1957:    |     |     |     |
1958:   10-20--6-14--7-19--11
1959:          |     |
1960:         20  1  19
1961:          |     |
1962:         10-18--11
1963:          |     |
1964:         23  2  22
1965:          |     |
1966:         12-21--13
1967:        */
1968:       PetscInt cone[4], ornt[4];

1970:       numCells    = !rank ?  6 : 0;
1971:       numEdges    = !rank ? 12 : 0;
1972:       numVerts    = !rank ?  8 : 0;
1973:       firstVertex = numCells;
1974:       firstEdge   = numCells + numVerts;
1975:       /* Build Topology */
1976:       DMPlexSetChart(*dm, 0, numCells+numEdges+numVerts);
1977:       for (c = 0; c < numCells; c++) {
1978:         DMPlexSetConeSize(*dm, c, 4);
1979:       }
1980:       for (e = firstEdge; e < firstEdge+numEdges; ++e) {
1981:         DMPlexSetConeSize(*dm, e, 2);
1982:       }
1983:       DMSetUp(*dm); /* Allocate space for cones */
1984:       if (!rank) {
1985:         /* Cell 0 */
1986:         cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17;
1987:         DMPlexSetCone(*dm, 0, cone);
1988:         ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0;
1989:         DMPlexSetConeOrientation(*dm, 0, ornt);
1990:         /* Cell 1 */
1991:         cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20;
1992:         DMPlexSetCone(*dm, 1, cone);
1993:         ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
1994:         DMPlexSetConeOrientation(*dm, 1, ornt);
1995:         /* Cell 2 */
1996:         cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23;
1997:         DMPlexSetCone(*dm, 2, cone);
1998:         ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
1999:         DMPlexSetConeOrientation(*dm, 2, ornt);
2000:         /* Cell 3 */
2001:         cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15;
2002:         DMPlexSetCone(*dm, 3, cone);
2003:         ornt[0] = -2; ornt[1] = -2; ornt[2] = 0; ornt[3] = -2;
2004:         DMPlexSetConeOrientation(*dm, 3, ornt);
2005:         /* Cell 4 */
2006:         cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25;
2007:         DMPlexSetCone(*dm, 4, cone);
2008:         ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = 0;
2009:         DMPlexSetConeOrientation(*dm, 4, ornt);
2010:         /* Cell 5 */
2011:         cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23;
2012:         DMPlexSetCone(*dm, 5, cone);
2013:         ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = -2;
2014:         DMPlexSetConeOrientation(*dm, 5, ornt);
2015:         /* Edges */
2016:         cone[0] =  6; cone[1] =  7;
2017:         DMPlexSetCone(*dm, 14, cone);
2018:         cone[0] =  7; cone[1] =  8;
2019:         DMPlexSetCone(*dm, 15, cone);
2020:         cone[0] =  8; cone[1] =  9;
2021:         DMPlexSetCone(*dm, 16, cone);
2022:         cone[0] =  9; cone[1] =  6;
2023:         DMPlexSetCone(*dm, 17, cone);
2024:         cone[0] = 10; cone[1] = 11;
2025:         DMPlexSetCone(*dm, 18, cone);
2026:         cone[0] = 11; cone[1] =  7;
2027:         DMPlexSetCone(*dm, 19, cone);
2028:         cone[0] =  6; cone[1] = 10;
2029:         DMPlexSetCone(*dm, 20, cone);
2030:         cone[0] = 12; cone[1] = 13;
2031:         DMPlexSetCone(*dm, 21, cone);
2032:         cone[0] = 13; cone[1] = 11;
2033:         DMPlexSetCone(*dm, 22, cone);
2034:         cone[0] = 10; cone[1] = 12;
2035:         DMPlexSetCone(*dm, 23, cone);
2036:         cone[0] = 13; cone[1] =  8;
2037:         DMPlexSetCone(*dm, 24, cone);
2038:         cone[0] = 12; cone[1] =  9;
2039:         DMPlexSetCone(*dm, 25, cone);
2040:       }
2041:       DMPlexSymmetrize(*dm);
2042:       DMPlexStratify(*dm);
2043:       /* Build coordinates */
2044:       PetscCalloc1(numVerts * embedDim, &coordsIn);
2045:       if (!rank) {
2046:         coordsIn[0*embedDim+0] = -R; coordsIn[0*embedDim+1] =  R; coordsIn[0*embedDim+2] = -R;
2047:         coordsIn[1*embedDim+0] =  R; coordsIn[1*embedDim+1] =  R; coordsIn[1*embedDim+2] = -R;
2048:         coordsIn[2*embedDim+0] =  R; coordsIn[2*embedDim+1] = -R; coordsIn[2*embedDim+2] = -R;
2049:         coordsIn[3*embedDim+0] = -R; coordsIn[3*embedDim+1] = -R; coordsIn[3*embedDim+2] = -R;
2050:         coordsIn[4*embedDim+0] = -R; coordsIn[4*embedDim+1] =  R; coordsIn[4*embedDim+2] =  R;
2051:         coordsIn[5*embedDim+0] =  R; coordsIn[5*embedDim+1] =  R; coordsIn[5*embedDim+2] =  R;
2052:         coordsIn[6*embedDim+0] = -R; coordsIn[6*embedDim+1] = -R; coordsIn[6*embedDim+2] =  R;
2053:         coordsIn[7*embedDim+0] =  R; coordsIn[7*embedDim+1] = -R; coordsIn[7*embedDim+2] =  R;
2054:       }
2055:     }
2056:     break;
2057:   case 3:
2058:     if (simplex) {
2059:       DM              idm;
2060:       const PetscReal edgeLen         = 1.0/PETSC_PHI;
2061:       PetscReal       vertexA[4]      = {0.5, 0.5, 0.5, 0.5};
2062:       PetscReal       vertexB[4]      = {1.0, 0.0, 0.0, 0.0};
2063:       PetscReal       vertexC[4]      = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0};
2064:       const PetscInt  degree          = 12;
2065:       PetscInt        s[4]            = {1, 1, 1};
2066:       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},
2067:                                          {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}};
2068:       PetscInt        cone[4];
2069:       PetscInt       *graph, p, i, j, k, l;

2071:       vertexA[0] *= R; vertexA[1] *= R; vertexA[2] *= R; vertexA[3] *= R;
2072:       vertexB[0] *= R; vertexB[1] *= R; vertexB[2] *= R; vertexB[3] *= R;
2073:       vertexC[0] *= R; vertexC[1] *= R; vertexC[2] *= R; vertexC[3] *= R;
2074:       numCells    = !rank ? 600 : 0;
2075:       numVerts    = !rank ? 120 : 0;
2076:       firstVertex = numCells;
2077:       /* Use the 600-cell, which for a unit sphere has coordinates which are

2079:            1/2 (\pm 1, \pm 1,    \pm 1, \pm 1)                          16
2080:                (\pm 1,    0,       0,      0)  all cyclic permutations   8
2081:            1/2 (\pm 1, \pm phi, \pm 1/phi, 0)  all even permutations    96

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

2086:          http://buzzard.pugetsound.edu/sage-practice/ch03s03.html
2087:          http://mathworld.wolfram.com/600-Cell.html
2088:       */
2089:       /* Construct vertices */
2090:       PetscCalloc1(numVerts * embedDim, &coordsIn);
2091:       i    = 0;
2092:       if (!rank) {
2093:         for (s[0] = -1; s[0] < 2; s[0] += 2) {
2094:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
2095:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
2096:               for (s[3] = -1; s[3] < 2; s[3] += 2) {
2097:                 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d];
2098:                 ++i;
2099:               }
2100:             }
2101:           }
2102:         }
2103:         for (p = 0; p < embedDim; ++p) {
2104:           s[1] = s[2] = s[3] = 1;
2105:           for (s[0] = -1; s[0] < 2; s[0] += 2) {
2106:             for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim];
2107:             ++i;
2108:           }
2109:         }
2110:         for (p = 0; p < 12; ++p) {
2111:           s[3] = 1;
2112:           for (s[0] = -1; s[0] < 2; s[0] += 2) {
2113:             for (s[1] = -1; s[1] < 2; s[1] += 2) {
2114:               for (s[2] = -1; s[2] < 2; s[2] += 2) {
2115:                 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]];
2116:                 ++i;
2117:               }
2118:             }
2119:           }
2120:         }
2121:       }
2122:       if (i != numVerts) SETERRQ2(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertices %D != %D", i, numVerts);
2123:       /* Construct graph */
2124:       PetscCalloc1(numVerts * numVerts, &graph);
2125:       for (i = 0; i < numVerts; ++i) {
2126:         for (j = 0, k = 0; j < numVerts; ++j) {
2127:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
2128:         }
2129:         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertex %D degree %D != %D", i, k, degree);
2130:       }
2131:       /* Build Topology */
2132:       DMPlexSetChart(*dm, 0, numCells+numVerts);
2133:       for (c = 0; c < numCells; c++) {
2134:         DMPlexSetConeSize(*dm, c, embedDim);
2135:       }
2136:       DMSetUp(*dm); /* Allocate space for cones */
2137:       /* Cells */
2138:       if (!rank) {
2139:         for (i = 0, c = 0; i < numVerts; ++i) {
2140:           for (j = 0; j < i; ++j) {
2141:             for (k = 0; k < j; ++k) {
2142:               for (l = 0; l < k; ++l) {
2143:                 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] &&
2144:                     graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) {
2145:                   cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l;
2146:                   /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */
2147:                   {
2148:                     const PetscInt epsilon[4][4][4][4] = {{{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2149:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  1}, { 0,  0, -1, 0}},
2150:                                                            {{0,  0,  0,  0}, { 0, 0,  0, -1}, { 0,  0, 0,  0}, { 0,  1,  0, 0}},
2151:                                                            {{0,  0,  0,  0}, { 0, 0,  1,  0}, { 0, -1, 0,  0}, { 0,  0,  0, 0}}},

2153:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0, -1}, { 0,  0,  1, 0}},
2154:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2155:                                                            {{0,  0,  0,  1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, {-1,  0,  0, 0}},
2156:                                                            {{0,  0, -1,  0}, { 0, 0,  0,  0}, { 1,  0, 0,  0}, { 0,  0,  0, 0}}},

2158:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  1}, { 0,  0, 0,  0}, { 0, -1,  0, 0}},
2159:                                                            {{0,  0,  0, -1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 1,  0,  0, 0}},
2160:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2161:                                                            {{0,  1,  0,  0}, {-1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}},

2163:                                                           {{{0,  0,  0,  0}, { 0, 0, -1,  0}, { 0,  1, 0,  0}, { 0,  0,  0, 0}},
2164:                                                            {{0,  0,  1,  0}, { 0, 0,  0,  0}, {-1,  0, 0,  0}, { 0,  0,  0, 0}},
2165:                                                            {{0, -1,  0,  0}, { 1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2166:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}}};
2167:                     PetscReal normal[4];
2168:                     PetscInt  e, f, g;

2170:                     for (d = 0; d < embedDim; ++d) {
2171:                       normal[d] = 0.0;
2172:                       for (e = 0; e < embedDim; ++e) {
2173:                         for (f = 0; f < embedDim; ++f) {
2174:                           for (g = 0; g < embedDim; ++g) {
2175:                             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]);
2176:                           }
2177:                         }
2178:                       }
2179:                     }
2180:                     if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
2181:                   }
2182:                   DMPlexSetCone(*dm, c++, cone);
2183:                 }
2184:               }
2185:             }
2186:           }
2187:         }
2188:       }
2189:       DMPlexSymmetrize(*dm);
2190:       DMPlexStratify(*dm);
2191:       PetscFree(graph);
2192:       /* Interpolate mesh */
2193:       DMPlexInterpolate(*dm, &idm);
2194:       DMDestroy(dm);
2195:       *dm  = idm;
2196:       break;
2197:     }
2198:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim);
2199:   }
2200:   /* Create coordinates */
2201:   DMGetCoordinateSection(*dm, &coordSection);
2202:   PetscSectionSetNumFields(coordSection, 1);
2203:   PetscSectionSetFieldComponents(coordSection, 0, embedDim);
2204:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);
2205:   for (v = firstVertex; v < firstVertex+numVerts; ++v) {
2206:     PetscSectionSetDof(coordSection, v, embedDim);
2207:     PetscSectionSetFieldDof(coordSection, v, 0, embedDim);
2208:   }
2209:   PetscSectionSetUp(coordSection);
2210:   PetscSectionGetStorageSize(coordSection, &coordSize);
2211:   VecCreate(PETSC_COMM_SELF, &coordinates);
2212:   VecSetBlockSize(coordinates, embedDim);
2213:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2214:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2215:   VecSetType(coordinates,VECSTANDARD);
2216:   VecGetArray(coordinates, &coords);
2217:   for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];}
2218:   VecRestoreArray(coordinates, &coords);
2219:   DMSetCoordinatesLocal(*dm, coordinates);
2220:   VecDestroy(&coordinates);
2221:   PetscFree(coordsIn);
2222:   /* Create coordinate function space */
2223:   {
2224:     DM          cdm;
2225:     PetscDS     cds;
2226:     PetscFE     fe;
2227:     PetscScalar radius = R;
2228:     PetscInt    dT, dE;

2230:     DMGetCoordinateDM(*dm, &cdm);
2231:     DMGetDimension(*dm, &dT);
2232:     DMGetCoordinateDim(*dm, &dE);
2233:     PetscFECreateLagrange(PETSC_COMM_SELF, dT, dE, simplex, 1, -1, &fe);
2234:     DMSetField(cdm, 0, NULL, (PetscObject) fe);
2235:     PetscFEDestroy(&fe);
2236:     DMCreateDS(cdm);

2238:     DMGetDS(cdm, &cds);
2239:     PetscDSSetConstants(cds, 1, &radius);
2240:   }
2241:   ((DM_Plex *) (*dm)->data)->coordFunc = snapToSphere;
2242:   return(0);
2243: }

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

2248:   Collective

2250:   Input Parameters:
2251: + comm  - The communicator for the DM object
2252: . dim   - The dimension
2253: - R     - The radius

2255:   Output Parameter:
2256: . dm  - The DM object

2258:   Options Database Keys:
2259: - bd_dm_refine - This will refine the surface mesh preserving the sphere geometry

2261:   Level: beginner

2263: .seealso: DMPlexCreateSphereMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
2264: @*/
2265: PetscErrorCode DMPlexCreateBallMesh(MPI_Comm comm, PetscInt dim, PetscReal R, DM *dm)
2266: {
2267:   DM             sdm;
2268:   DMLabel        bdlabel;

2272:   DMPlexCreateSphereMesh(comm, dim-1, PETSC_TRUE, R, &sdm);
2273:   PetscObjectSetOptionsPrefix((PetscObject) sdm, "bd_");
2274:   DMSetFromOptions(sdm);
2275:   DMPlexGenerate(sdm, NULL, PETSC_TRUE, dm);
2276:   DMDestroy(&sdm);
2277:   DMCreateLabel(*dm, "marker");
2278:   DMGetLabel(*dm, "marker", &bdlabel);
2279:   DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
2280:   DMPlexLabelComplete(*dm, bdlabel);
2281:   return(0);
2282: }

2284: /* External function declarations here */
2285: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
2286: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2287: extern PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2288: extern PetscErrorCode DMCreateLocalSection_Plex(DM dm);
2289: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
2290: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
2291: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
2292: extern PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
2293: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
2294: extern PetscErrorCode DMSetUp_Plex(DM dm);
2295: extern PetscErrorCode DMDestroy_Plex(DM dm);
2296: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
2297: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
2298: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
2299: extern PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);
2300: static PetscErrorCode DMInitialize_Plex(DM dm);

2302: /* Replace dm with the contents of dmNew
2303:    - Share the DM_Plex structure
2304:    - Share the coordinates
2305:    - Share the SF
2306: */
2307: static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
2308: {
2309:   PetscSF               sf;
2310:   DM                    coordDM, coarseDM;
2311:   DMField               coordField;
2312:   Vec                   coords;
2313:   PetscBool             isper;
2314:   const PetscReal      *maxCell, *L;
2315:   const DMBoundaryType *bd;
2316:   PetscErrorCode        ierr;

2319:   DMGetPointSF(dmNew, &sf);
2320:   DMSetPointSF(dm, sf);
2321:   DMGetCoordinateDM(dmNew, &coordDM);
2322:   DMGetCoordinatesLocal(dmNew, &coords);
2323:   DMGetCoordinateField(dmNew, &coordField);
2324:   DMSetCoordinateDM(dm, coordDM);
2325:   DMSetCoordinatesLocal(dm, coords);
2326:   DMSetCoordinateField(dm, coordField);
2327:   DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);
2328:   DMSetPeriodicity(dmNew, isper, maxCell, L, bd);
2329:   DMDestroy_Plex(dm);
2330:   DMInitialize_Plex(dm);
2331:   dm->data = dmNew->data;
2332:   ((DM_Plex *) dmNew->data)->refct++;
2333:   DMDestroyLabelLinkList_Internal(dm);
2334:   DMCopyLabels(dmNew, dm, PETSC_OWN_POINTER, PETSC_TRUE);
2335:   DMGetCoarseDM(dmNew,&coarseDM);
2336:   DMSetCoarseDM(dm,coarseDM);
2337:   return(0);
2338: }

2340: /* Swap dm with the contents of dmNew
2341:    - Swap the DM_Plex structure
2342:    - Swap the coordinates
2343:    - Swap the point PetscSF
2344: */
2345: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
2346: {
2347:   DM              coordDMA, coordDMB;
2348:   Vec             coordsA,  coordsB;
2349:   PetscSF         sfA,      sfB;
2350:   void            *tmp;
2351:   DMLabelLink     listTmp;
2352:   DMLabel         depthTmp;
2353:   PetscInt        tmpI;
2354:   PetscErrorCode  ierr;

2357:   DMGetPointSF(dmA, &sfA);
2358:   DMGetPointSF(dmB, &sfB);
2359:   PetscObjectReference((PetscObject) sfA);
2360:   DMSetPointSF(dmA, sfB);
2361:   DMSetPointSF(dmB, sfA);
2362:   PetscObjectDereference((PetscObject) sfA);

2364:   DMGetCoordinateDM(dmA, &coordDMA);
2365:   DMGetCoordinateDM(dmB, &coordDMB);
2366:   PetscObjectReference((PetscObject) coordDMA);
2367:   DMSetCoordinateDM(dmA, coordDMB);
2368:   DMSetCoordinateDM(dmB, coordDMA);
2369:   PetscObjectDereference((PetscObject) coordDMA);

2371:   DMGetCoordinatesLocal(dmA, &coordsA);
2372:   DMGetCoordinatesLocal(dmB, &coordsB);
2373:   PetscObjectReference((PetscObject) coordsA);
2374:   DMSetCoordinatesLocal(dmA, coordsB);
2375:   DMSetCoordinatesLocal(dmB, coordsA);
2376:   PetscObjectDereference((PetscObject) coordsA);

2378:   tmp       = dmA->data;
2379:   dmA->data = dmB->data;
2380:   dmB->data = tmp;
2381:   listTmp   = dmA->labels;
2382:   dmA->labels = dmB->labels;
2383:   dmB->labels = listTmp;
2384:   depthTmp  = dmA->depthLabel;
2385:   dmA->depthLabel = dmB->depthLabel;
2386:   dmB->depthLabel = depthTmp;
2387:   depthTmp  = dmA->celltypeLabel;
2388:   dmA->celltypeLabel = dmB->celltypeLabel;
2389:   dmB->celltypeLabel = depthTmp;
2390:   tmpI         = dmA->levelup;
2391:   dmA->levelup = dmB->levelup;
2392:   dmB->levelup = tmpI;
2393:   return(0);
2394: }

2396: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2397: {
2398:   DM_Plex       *mesh = (DM_Plex*) dm->data;
2399:   PetscBool      flg;

2403:   /* Handle viewing */
2404:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL);
2405:   PetscOptionsBoundedInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL,0);
2406:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL);
2407:   PetscOptionsBoundedInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL,0);
2408:   DMMonitorSetFromOptions(dm, "-dm_plex_monitor_throughput", "Monitor the simulation throughput", "DMPlexMonitorThroughput", DMPlexMonitorThroughput, NULL, &flg);
2409:   if (flg) {PetscLogDefaultBegin();}
2410:   /* Point Location */
2411:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL);
2412:   /* Partitioning and distribution */
2413:   PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL);
2414:   /* Generation and remeshing */
2415:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL);
2416:   /* Projection behavior */
2417:   PetscOptionsBoundedInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL,0);
2418:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);
2419:   PetscOptionsEnum("-dm_plex_cell_refiner", "Strategy for cell refinment", "ex40.c", DMPlexCellRefinerTypes, (PetscEnum) mesh->cellRefiner, (PetscEnum *) &mesh->cellRefiner, NULL);
2420:   /* Checking structure */
2421:   {
2422:     PetscBool   flg = PETSC_FALSE, flg2 = PETSC_FALSE, all = PETSC_FALSE;

2424:     PetscOptionsBool("-dm_plex_check_all", "Perform all checks", NULL, PETSC_FALSE, &all, &flg2);
2425:     PetscOptionsBool("-dm_plex_check_symmetry", "Check that the adjacency information in the mesh is symmetric", "DMPlexCheckSymmetry", PETSC_FALSE, &flg, &flg2);
2426:     if (all || (flg && flg2)) {DMPlexCheckSymmetry(dm);}
2427:     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);
2428:     if (all || (flg && flg2)) {DMPlexCheckSkeleton(dm, 0);}
2429:     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);
2430:     if (all || (flg && flg2)) {DMPlexCheckFaces(dm, 0);}
2431:     PetscOptionsBool("-dm_plex_check_geometry", "Check that cells have positive volume", "DMPlexCheckGeometry", PETSC_FALSE, &flg, &flg2);
2432:     if (all || (flg && flg2)) {DMPlexCheckGeometry(dm);}
2433:     PetscOptionsBool("-dm_plex_check_pointsf", "Check some necessary conditions for PointSF", "DMPlexCheckPointSF", PETSC_FALSE, &flg, &flg2);
2434:     if (all || (flg && flg2)) {DMPlexCheckPointSF(dm);}
2435:     PetscOptionsBool("-dm_plex_check_interface_cones", "Check points on inter-partition interfaces have conforming order of cone points", "DMPlexCheckInterfaceCones", PETSC_FALSE, &flg, &flg2);
2436:     if (all || (flg && flg2)) {DMPlexCheckInterfaceCones(dm);}
2437:     PetscOptionsBool("-dm_plex_check_cell_shape", "Check cell shape", "DMPlexCheckCellShape", PETSC_FALSE, &flg, &flg2);
2438:     if (flg && flg2) {DMPlexCheckCellShape(dm, PETSC_TRUE, PETSC_DETERMINE);}
2439:   }

2441:   PetscPartitionerSetFromOptions(mesh->partitioner);
2442:   return(0);
2443: }

2445: static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2446: {
2447:   PetscReal      volume = -1.0;
2448:   PetscInt       prerefine = 0, refine = 0, r, coarsen = 0, overlap = 0;
2449:   PetscBool      uniformOrig, uniform = PETSC_TRUE, distribute = PETSC_FALSE, isHierarchy, flg;

2454:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
2455:   /* Handle DMPlex refinement before distribution */
2456:   DMPlexGetRefinementUniform(dm, &uniformOrig);
2457:   PetscOptionsBool("-dm_refine_uniform_pre", "Flag for uniform refinement before distribution", "DMCreate", uniform, &uniform, &flg);
2458:   if (flg) {DMPlexSetRefinementUniform(dm, uniform);}
2459:   PetscOptionsReal("-dm_refine_volume_limit_pre", "The maximum cell volume after refinement before distribution", "DMCreate", volume, &volume, &flg);
2460:   if (flg) {DMPlexSetRefinementLimit(dm, volume);}
2461:   PetscOptionsBoundedInt("-dm_refine_pre", "The number of refinements before distribution", "DMCreate", prerefine, &prerefine, NULL,0);
2462:   for (r = 0; r < prerefine; ++r) {
2463:     DM             rdm;
2464:     PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

2466:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2467:     DMRefine(dm, PetscObjectComm((PetscObject) dm), &rdm);
2468:     /* Total hack since we do not pass in a pointer */
2469:     DMPlexReplace_Static(dm, rdm);
2470:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2471:     if (coordFunc) {
2472:       DMPlexRemapGeometry(dm, 0.0, coordFunc);
2473:       ((DM_Plex*) dm->data)->coordFunc = coordFunc;
2474:     }
2475:     DMDestroy(&rdm);
2476:   }
2477:   DMPlexSetRefinementUniform(dm, uniformOrig);
2478:   /* Handle DMPlex distribution */
2479:   PetscOptionsBool("-dm_distribute", "Flag to redistribute a mesh among processes", "DMCreate", distribute, &distribute, NULL);
2480:   PetscOptionsBoundedInt("-dm_distribute_overlap", "The size of the overlap halo", "DMCreate", overlap, &overlap, NULL, 0);
2481:   if (distribute) {
2482:     DM               pdm = NULL;
2483:     PetscPartitioner part;

2485:     DMPlexGetPartitioner(dm, &part);
2486:     PetscPartitionerSetFromOptions(part);
2487:     DMPlexDistribute(dm, overlap, NULL, &pdm);
2488:     if (pdm) {
2489:       DMPlexReplace_Static(dm, pdm);
2490:       DMDestroy(&pdm);
2491:     }
2492:   }
2493:   /* Handle DMPlex refinement */
2494:   PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL,0);
2495:   PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy,0);
2496:   if (refine) {DMPlexSetRefinementUniform(dm, PETSC_TRUE);}
2497:   if (refine && isHierarchy) {
2498:     DM *dms, coarseDM;

2500:     DMGetCoarseDM(dm, &coarseDM);
2501:     PetscObjectReference((PetscObject)coarseDM);
2502:     PetscMalloc1(refine,&dms);
2503:     DMRefineHierarchy(dm, refine, dms);
2504:     /* Total hack since we do not pass in a pointer */
2505:     DMPlexSwap_Static(dm, dms[refine-1]);
2506:     if (refine == 1) {
2507:       DMSetCoarseDM(dm, dms[0]);
2508:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2509:     } else {
2510:       DMSetCoarseDM(dm, dms[refine-2]);
2511:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2512:       DMSetCoarseDM(dms[0], dms[refine-1]);
2513:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
2514:     }
2515:     DMSetCoarseDM(dms[refine-1], coarseDM);
2516:     PetscObjectDereference((PetscObject)coarseDM);
2517:     /* Free DMs */
2518:     for (r = 0; r < refine; ++r) {
2519:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2520:       DMDestroy(&dms[r]);
2521:     }
2522:     PetscFree(dms);
2523:   } else {
2524:     for (r = 0; r < refine; ++r) {
2525:       DM refinedMesh;
2526:       PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

2528:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2529:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);
2530:       /* Total hack since we do not pass in a pointer */
2531:       DMPlexReplace_Static(dm, refinedMesh);
2532:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2533:       if (coordFunc) {
2534:         DMPlexRemapGeometry(dm, 0.0, coordFunc);
2535:         ((DM_Plex*) dm->data)->coordFunc = coordFunc;
2536:       }
2537:       DMDestroy(&refinedMesh);
2538:     }
2539:   }
2540:   /* Handle DMPlex coarsening */
2541:   PetscOptionsBoundedInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL,0);
2542:   PetscOptionsBoundedInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy,0);
2543:   if (coarsen && isHierarchy) {
2544:     DM *dms;

2546:     PetscMalloc1(coarsen, &dms);
2547:     DMCoarsenHierarchy(dm, coarsen, dms);
2548:     /* Free DMs */
2549:     for (r = 0; r < coarsen; ++r) {
2550:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2551:       DMDestroy(&dms[r]);
2552:     }
2553:     PetscFree(dms);
2554:   } else {
2555:     for (r = 0; r < coarsen; ++r) {
2556:       DM coarseMesh;

2558:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2559:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);
2560:       /* Total hack since we do not pass in a pointer */
2561:       DMPlexReplace_Static(dm, coarseMesh);
2562:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2563:       DMDestroy(&coarseMesh);
2564:     }
2565:   }
2566:   /* Handle */
2567:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2568:   PetscOptionsTail();
2569:   return(0);
2570: }

2572: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
2573: {

2577:   DMCreateGlobalVector_Section_Private(dm,vec);
2578:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
2579:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
2580:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
2581:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
2582:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
2583:   return(0);
2584: }

2586: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
2587: {

2591:   DMCreateLocalVector_Section_Private(dm,vec);
2592:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
2593:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
2594:   return(0);
2595: }

2597: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
2598: {
2599:   PetscInt       depth, d;

2603:   DMPlexGetDepth(dm, &depth);
2604:   if (depth == 1) {
2605:     DMGetDimension(dm, &d);
2606:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
2607:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
2608:     else               {*pStart = 0; *pEnd = 0;}
2609:   } else {
2610:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
2611:   }
2612:   return(0);
2613: }

2615: static PetscErrorCode DMGetNeighbors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
2616: {
2617:   PetscSF           sf;
2618:   PetscInt          niranks, njranks, n;
2619:   const PetscMPIInt *iranks, *jranks;
2620:   DM_Plex           *data = (DM_Plex*) dm->data;
2621:   PetscErrorCode    ierr;

2624:   DMGetPointSF(dm, &sf);
2625:   if (!data->neighbors) {
2626:     PetscSFSetUp(sf);
2627:     PetscSFGetRootRanks(sf, &njranks, &jranks, NULL, NULL, NULL);
2628:     PetscSFGetLeafRanks(sf, &niranks, &iranks, NULL, NULL);
2629:     PetscMalloc1(njranks + niranks + 1, &data->neighbors);
2630:     PetscArraycpy(data->neighbors + 1, jranks, njranks);
2631:     PetscArraycpy(data->neighbors + njranks + 1, iranks, niranks);
2632:     n = njranks + niranks;
2633:     PetscSortRemoveDupsMPIInt(&n, data->neighbors + 1);
2634:     /* The following cast should never fail: can't have more neighbors than PETSC_MPI_INT_MAX */
2635:     PetscMPIIntCast(n, data->neighbors);
2636:   }
2637:   if (nranks) *nranks = data->neighbors[0];
2638:   if (ranks) {
2639:     if (data->neighbors[0]) *ranks = data->neighbors + 1;
2640:     else                    *ranks = NULL;
2641:   }
2642:   return(0);
2643: }

2645: static PetscErrorCode DMInitialize_Plex(DM dm)
2646: {

2650:   dm->ops->view                            = DMView_Plex;
2651:   dm->ops->load                            = DMLoad_Plex;
2652:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
2653:   dm->ops->clone                           = DMClone_Plex;
2654:   dm->ops->setup                           = DMSetUp_Plex;
2655:   dm->ops->createlocalsection              = DMCreateLocalSection_Plex;
2656:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
2657:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
2658:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
2659:   dm->ops->getlocaltoglobalmapping         = NULL;
2660:   dm->ops->createfieldis                   = NULL;
2661:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
2662:   dm->ops->createcoordinatefield           = DMCreateCoordinateField_Plex;
2663:   dm->ops->getcoloring                     = NULL;
2664:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
2665:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
2666:   dm->ops->createmassmatrix                = DMCreateMassMatrix_Plex;
2667:   dm->ops->createinjection                 = DMCreateInjection_Plex;
2668:   dm->ops->refine                          = DMRefine_Plex;
2669:   dm->ops->coarsen                         = DMCoarsen_Plex;
2670:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
2671:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
2672:   dm->ops->adaptlabel                      = DMAdaptLabel_Plex;
2673:   dm->ops->adaptmetric                     = DMAdaptMetric_Plex;
2674:   dm->ops->globaltolocalbegin              = NULL;
2675:   dm->ops->globaltolocalend                = NULL;
2676:   dm->ops->localtoglobalbegin              = NULL;
2677:   dm->ops->localtoglobalend                = NULL;
2678:   dm->ops->destroy                         = DMDestroy_Plex;
2679:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
2680:   dm->ops->createsuperdm                   = DMCreateSuperDM_Plex;
2681:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
2682:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
2683:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
2684:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
2685:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
2686:   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
2687:   dm->ops->projectbdfieldlabellocal        = DMProjectBdFieldLabelLocal_Plex;
2688:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
2689:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
2690:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
2691:   dm->ops->getneighbors                    = DMGetNeighbors_Plex;
2692:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);
2693:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertTimeDerviativeBoundaryValues_C",DMPlexInsertTimeDerivativeBoundaryValues_Plex);
2694:   PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);
2695:   PetscObjectComposeFunction((PetscObject)dm,"DMCreateNeumannOverlap_C",DMCreateNeumannOverlap_Plex);
2696:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexGetOverlap_C",DMPlexGetOverlap_Plex);
2697:   return(0);
2698: }

2700: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
2701: {
2702:   DM_Plex        *mesh = (DM_Plex *) dm->data;

2706:   mesh->refct++;
2707:   (*newdm)->data = mesh;
2708:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
2709:   DMInitialize_Plex(*newdm);
2710:   return(0);
2711: }

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

2719:   Options Database Keys:
2720: + -dm_refine_pre                     - Refine mesh before distribution
2721: + -dm_refine_uniform_pre             - Choose uniform or generator-based refinement
2722: + -dm_refine_volume_limit_pre        - Cell volume limit after pre-refinement using generator
2723: . -dm_distribute                     - Distribute mesh across processes
2724: . -dm_distribute_overlap             - Number of cells to overlap for distribution
2725: . -dm_refine                         - Refine mesh after distribution
2726: . -dm_plex_hash_location             - Use grid hashing for point location
2727: . -dm_plex_partition_balance         - Attempt to evenly divide points on partition boundary between processes
2728: . -dm_plex_remesh_bd                 - Allow changes to the boundary on remeshing
2729: . -dm_plex_max_projection_height     - Maxmimum mesh point height used to project locally
2730: . -dm_plex_regular_refinement        - Use special nested projection algorithm for regular refinement
2731: . -dm_plex_check_all                 - Perform all shecks below
2732: . -dm_plex_check_symmetry            - Check that the adjacency information in the mesh is symmetric
2733: . -dm_plex_check_skeleton <celltype> - Check that each cell has the correct number of vertices
2734: . -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
2735: . -dm_plex_check_geometry            - Check that cells have positive volume
2736: . -dm_view :mesh.tex:ascii_latex     - View the mesh in LaTeX/TikZ
2737: . -dm_plex_view_scale <num>          - Scale the TikZ
2738: - -dm_plex_print_fem <num>           - View FEM assembly information, such as element vectors and matrices


2741:   Level: intermediate

2743: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
2744: M*/

2746: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
2747: {
2748:   DM_Plex       *mesh;
2749:   PetscInt       unit;

2754:   PetscNewLog(dm,&mesh);
2755:   dm->dim  = 0;
2756:   dm->data = mesh;

2758:   mesh->refct             = 1;
2759:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
2760:   mesh->maxConeSize       = 0;
2761:   mesh->cones             = NULL;
2762:   mesh->coneOrientations  = NULL;
2763:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
2764:   mesh->maxSupportSize    = 0;
2765:   mesh->supports          = NULL;
2766:   mesh->refinementUniform = PETSC_TRUE;
2767:   mesh->refinementLimit   = -1.0;
2768:   mesh->interpolated      = DMPLEX_INTERPOLATED_INVALID;
2769:   mesh->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID;

2771:   mesh->facesTmp = NULL;

2773:   mesh->tetgenOpts   = NULL;
2774:   mesh->triangleOpts = NULL;
2775:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
2776:   mesh->remeshBd     = PETSC_FALSE;

2778:   mesh->subpointMap = NULL;

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

2782:   mesh->regularRefinement   = PETSC_FALSE;
2783:   mesh->depthState          = -1;
2784:   mesh->celltypeState       = -1;
2785:   mesh->globalVertexNumbers = NULL;
2786:   mesh->globalCellNumbers   = NULL;
2787:   mesh->anchorSection       = NULL;
2788:   mesh->anchorIS            = NULL;
2789:   mesh->createanchors       = NULL;
2790:   mesh->computeanchormatrix = NULL;
2791:   mesh->parentSection       = NULL;
2792:   mesh->parents             = NULL;
2793:   mesh->childIDs            = NULL;
2794:   mesh->childSection        = NULL;
2795:   mesh->children            = NULL;
2796:   mesh->referenceTree       = NULL;
2797:   mesh->getchildsymmetry    = NULL;
2798:   mesh->vtkCellHeight       = 0;
2799:   mesh->useAnchors          = PETSC_FALSE;

2801:   mesh->maxProjectionHeight = 0;

2803:   mesh->neighbors           = NULL;

2805:   mesh->printSetValues = PETSC_FALSE;
2806:   mesh->printFEM       = 0;
2807:   mesh->printTol       = 1.0e-10;

2809:   DMInitialize_Plex(dm);
2810:   return(0);
2811: }

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

2816:   Collective

2818:   Input Parameter:
2819: . comm - The communicator for the DMPlex object

2821:   Output Parameter:
2822: . mesh  - The DMPlex object

2824:   Level: beginner

2826: @*/
2827: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
2828: {

2833:   DMCreate(comm, mesh);
2834:   DMSetType(*mesh, DMPLEX);
2835:   return(0);
2836: }

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

2841:   Input Parameters:
2842: + dm - The DM
2843: . numCells - The number of cells owned by this process
2844: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
2845: . NVertices - The global number of vertices, or PETSC_DECIDE
2846: . numCorners - The number of vertices for each cell
2847: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell

2849:   Output Parameter:
2850: . vertexSF - (Optional) SF describing complete vertex ownership

2852:   Notes:
2853:   Two triangles sharing a face
2854: $
2855: $        2
2856: $      / | \
2857: $     /  |  \
2858: $    /   |   \
2859: $   0  0 | 1  3
2860: $    \   |   /
2861: $     \  |  /
2862: $      \ | /
2863: $        1
2864: would have input
2865: $  numCells = 2, numVertices = 4
2866: $  cells = [0 1 2  1 3 2]
2867: $
2868: which would result in the DMPlex
2869: $
2870: $        4
2871: $      / | \
2872: $     /  |  \
2873: $    /   |   \
2874: $   2  0 | 1  5
2875: $    \   |   /
2876: $     \  |  /
2877: $      \ | /
2878: $        3

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

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

2888:   Not currently supported in Fortran.

2890:   Level: advanced

2892: .seealso: DMPlexBuildFromCellList(), DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildCoordinatesFromCellListParallel()
2893: @*/
2894: PetscErrorCode DMPlexBuildFromCellListParallel(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, const PetscInt cells[], PetscSF *vertexSF)
2895: {
2896:   PetscSF         sfPoint, sfVert;
2897:   PetscLayout     vLayout;
2898:   PetscSFNode    *vertexLocal, *vertexOwner, *remoteVertex;
2899:   PetscInt        numVerticesAdj, *verticesAdj, numVerticesGhost = 0, *localVertex, *cones, c, p, v, g, dim;
2900:   PetscMPIInt     rank, size;
2901:   PetscErrorCode  ierr;

2905:   PetscLogEventBegin(DMPLEX_BuildFromCellList,dm,0,0,0);
2906:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
2907:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);
2908:   DMGetDimension(dm, &dim);
2909:   /* Get/check global number of vertices */
2910:   {
2911:     PetscInt NVerticesInCells, i;
2912:     const PetscInt len = numCells * numCorners;

2914:     /* NVerticesInCells = max(cells) + 1 */
2915:     NVerticesInCells = PETSC_MIN_INT;
2916:     for (i=0; i<len; i++) if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
2917:     ++NVerticesInCells;
2918:     MPI_Allreduce(MPI_IN_PLACE, &NVerticesInCells, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));

2920:     if (numVertices == PETSC_DECIDE && NVertices == PETSC_DECIDE) NVertices = NVerticesInCells;
2921:     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);
2922:   }
2923:   /* Partition vertices */
2924:   PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);
2925:   PetscLayoutSetLocalSize(vLayout, numVertices);
2926:   PetscLayoutSetSize(vLayout, NVertices);
2927:   PetscLayoutSetBlockSize(vLayout, 1);
2928:   PetscLayoutSetUp(vLayout);
2929:   /* Count locally unique vertices */
2930:   {
2931:     PetscHSetI vhash;
2932:     PetscInt off = 0;

2934:     PetscHSetICreate(&vhash);
2935:     for (c = 0; c < numCells; ++c) {
2936:       for (p = 0; p < numCorners; ++p) {
2937:         PetscHSetIAdd(vhash, cells[c*numCorners+p]);
2938:       }
2939:     }
2940:     PetscHSetIGetSize(vhash, &numVerticesAdj);
2941:     PetscMalloc1(numVerticesAdj, &verticesAdj);
2942:     PetscHSetIGetElems(vhash, &off, verticesAdj);
2943:     PetscHSetIDestroy(&vhash);
2944:     if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
2945:   }
2946:   PetscSortInt(numVerticesAdj, verticesAdj);
2947:   /* Create cones */
2948:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
2949:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
2950:   DMSetUp(dm);
2951:   DMPlexGetCones(dm,&cones);
2952:   for (c = 0; c < numCells; ++c) {
2953:     for (p = 0; p < numCorners; ++p) {
2954:       const PetscInt gv = cells[c*numCorners+p];
2955:       PetscInt       lv;

2957:       /* Positions within verticesAdj form 0-based local vertex numbering;
2958:          we need to shift it by numCells to get correct DAG points (cells go first) */
2959:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
2960:       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
2961:       cones[c*numCorners+p] = lv+numCells;
2962:     }
2963:   }
2964:   /* Create SF for vertices */
2965:   PetscSFCreate(PetscObjectComm((PetscObject)dm), &sfVert);
2966:   PetscObjectSetName((PetscObject) sfVert, "Vertex Ownership SF");
2967:   PetscSFSetFromOptions(sfVert);
2968:   PetscSFSetGraphLayout(sfVert, vLayout, numVerticesAdj, NULL, PETSC_OWN_POINTER, verticesAdj);
2969:   PetscFree(verticesAdj);
2970:   /* Build pointSF */
2971:   PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);
2972:   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
2973:   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
2974:   PetscSFReduceBegin(sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
2975:   PetscSFReduceEnd(sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
2976:   for (v = 0; v < numVertices;    ++v) if (vertexOwner[v].rank < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global vertex %D was unclaimed", v + vLayout->rstart);
2977:   PetscSFBcastBegin(sfVert, MPIU_2INT, vertexOwner, vertexLocal);
2978:   PetscSFBcastEnd(sfVert, MPIU_2INT, vertexOwner, vertexLocal);
2979:   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
2980:   PetscMalloc1(numVerticesGhost, &localVertex);
2981:   PetscMalloc1(numVerticesGhost, &remoteVertex);
2982:   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
2983:     if (vertexLocal[v].rank != rank) {
2984:       localVertex[g]        = v+numCells;
2985:       remoteVertex[g].index = vertexLocal[v].index;
2986:       remoteVertex[g].rank  = vertexLocal[v].rank;
2987:       ++g;
2988:     }
2989:   }
2990:   PetscFree2(vertexLocal, vertexOwner);
2991:   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
2992:   DMGetPointSF(dm, &sfPoint);
2993:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
2994:   PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);
2995:   PetscLayoutDestroy(&vLayout);
2996:   /* Fill in the rest of the topology structure */
2997:   DMPlexSymmetrize(dm);
2998:   DMPlexStratify(dm);
2999:   PetscLogEventEnd(DMPLEX_BuildFromCellList,dm,0,0,0);
3000:   if (vertexSF) *vertexSF = sfVert;
3001:   else {PetscSFDestroy(&sfVert);}
3002:   return(0);
3003: }

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

3008:   Input Parameters:
3009: + dm - The DM
3010: . spaceDim - The spatial dimension used for coordinates
3011: . sfVert - SF describing complete vertex ownership
3012: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3014:   Level: advanced

3016:   Notes:
3017:   Not currently supported in Fortran.

3019: .seealso: DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildFromCellListParallel()
3020: @*/
3021: PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM dm, PetscInt spaceDim, PetscSF sfVert, const PetscReal vertexCoords[])
3022: {
3023:   PetscSection   coordSection;
3024:   Vec            coordinates;
3025:   PetscScalar   *coords;
3026:   PetscInt       numVertices, numVerticesAdj, coordSize, v, vStart, vEnd;

3030:   PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3031:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
3032:   if (vStart < 0 || vEnd < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first.");
3033:   DMSetCoordinateDim(dm, spaceDim);
3034:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
3035:   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);
3036:   DMGetCoordinateSection(dm, &coordSection);
3037:   PetscSectionSetNumFields(coordSection, 1);
3038:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
3039:   PetscSectionSetChart(coordSection, vStart, vEnd);
3040:   for (v = vStart; v < vEnd; ++v) {
3041:     PetscSectionSetDof(coordSection, v, spaceDim);
3042:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
3043:   }
3044:   PetscSectionSetUp(coordSection);
3045:   PetscSectionGetStorageSize(coordSection, &coordSize);
3046:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
3047:   VecSetBlockSize(coordinates, spaceDim);
3048:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3049:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3050:   VecSetType(coordinates,VECSTANDARD);
3051:   VecGetArray(coordinates, &coords);
3052:   {
3053:     MPI_Datatype coordtype;

3055:     /* Need a temp buffer for coords if we have complex/single */
3056:     MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);
3057:     MPI_Type_commit(&coordtype);
3058: #if defined(PETSC_USE_COMPLEX)
3059:     {
3060:     PetscScalar *svertexCoords;
3061:     PetscInt    i;
3062:     PetscMalloc1(numVertices*spaceDim,&svertexCoords);
3063:     for (i=0; i<numVertices*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
3064:     PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);
3065:     PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);
3066:     PetscFree(svertexCoords);
3067:     }
3068: #else
3069:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);
3070:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);
3071: #endif
3072:     MPI_Type_free(&coordtype);
3073:   }
3074:   VecRestoreArray(coordinates, &coords);
3075:   DMSetCoordinatesLocal(dm, coordinates);
3076:   VecDestroy(&coordinates);
3077:   PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3078:   return(0);
3079: }

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

3084:   Input Parameters:
3085: + comm - The communicator
3086: . dim - The topological dimension of the mesh
3087: . numCells - The number of cells owned by this process
3088: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
3089: . NVertices - The global number of vertices, or PETSC_DECIDE
3090: . numCorners - The number of vertices for each cell
3091: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
3092: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
3093: . spaceDim - The spatial dimension used for coordinates
3094: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3096:   Output Parameter:
3097: + dm - The DM
3098: - vertexSF - (Optional) SF describing complete vertex ownership

3100:   Notes:
3101:   This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(),
3102:   DMPlexBuildFromCellListParallel(), DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellListParallel()

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

3107:   Level: intermediate

3109: .seealso: DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromDAG(), DMPlexCreate()
3110: @*/
3111: 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)
3112: {
3113:   PetscSF        sfVert;

3117:   DMCreate(comm, dm);
3118:   DMSetType(*dm, DMPLEX);
3121:   DMSetDimension(*dm, dim);
3122:   DMPlexBuildFromCellListParallel(*dm, numCells, numVertices, NVertices, numCorners, cells, &sfVert);
3123:   if (interpolate) {
3124:     DM idm;

3126:     DMPlexInterpolate(*dm, &idm);
3127:     DMDestroy(dm);
3128:     *dm  = idm;
3129:   }
3130:   DMPlexBuildCoordinatesFromCellListParallel(*dm, spaceDim, sfVert, vertexCoords);
3131:   if (vertexSF) *vertexSF = sfVert;
3132:   else {PetscSFDestroy(&sfVert);}
3133:   return(0);
3134: }


3137: /*@
3138:   DMPlexCreateFromCellListParallel - Deprecated, use DMPlexCreateFromCellListParallelPetsc()

3140:   Level: deprecated

3142: .seealso: DMPlexCreateFromCellListParallelPetsc()
3143: @*/
3144: 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)
3145: {
3147:   PetscInt       i;
3148:   PetscInt       *pintCells;

3151:   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));
3152:   if (sizeof(int) == sizeof(PetscInt)) {
3153:     pintCells = (PetscInt *) cells;
3154:   } else {
3155:     PetscMalloc1(numCells*numCorners, &pintCells);
3156:     for (i = 0; i < numCells*numCorners; i++) {
3157:       pintCells[i] = (PetscInt) cells[i];
3158:     }
3159:   }
3160:   DMPlexCreateFromCellListParallelPetsc(comm, dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, pintCells, spaceDim, vertexCoords, vertexSF, dm);
3161:   if (sizeof(int) != sizeof(PetscInt)) {
3162:     PetscFree(pintCells);
3163:   }
3164:   return(0);
3165: }

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

3170:   Input Parameters:
3171: + dm - The DM
3172: . numCells - The number of cells owned by this process
3173: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
3174: . numCorners - The number of vertices for each cell
3175: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell

3177:   Level: advanced

3179:   Notes:
3180:   Two triangles sharing a face
3181: $
3182: $        2
3183: $      / | \
3184: $     /  |  \
3185: $    /   |   \
3186: $   0  0 | 1  3
3187: $    \   |   /
3188: $     \  |  /
3189: $      \ | /
3190: $        1
3191: would have input
3192: $  numCells = 2, numVertices = 4
3193: $  cells = [0 1 2  1 3 2]
3194: $
3195: which would result in the DMPlex
3196: $
3197: $        4
3198: $      / | \
3199: $     /  |  \
3200: $    /   |   \
3201: $   2  0 | 1  5
3202: $    \   |   /
3203: $     \  |  /
3204: $      \ | /
3205: $        3

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

3209:   Not currently supported in Fortran.

3211: .seealso: DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromCellListPetsc()
3212: @*/
3213: PetscErrorCode DMPlexBuildFromCellList(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const PetscInt cells[])
3214: {
3215:   PetscInt      *cones, c, p, dim;

3219:   PetscLogEventBegin(DMPLEX_BuildFromCellList,dm,0,0,0);
3220:   DMGetDimension(dm, &dim);
3221:   /* Get/check global number of vertices */
3222:   {
3223:     PetscInt NVerticesInCells, i;
3224:     const PetscInt len = numCells * numCorners;

3226:     /* NVerticesInCells = max(cells) + 1 */
3227:     NVerticesInCells = PETSC_MIN_INT;
3228:     for (i=0; i<len; i++) if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
3229:     ++NVerticesInCells;

3231:     if (numVertices == PETSC_DECIDE) numVertices = NVerticesInCells;
3232:     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);
3233:   }
3234:   DMPlexSetChart(dm, 0, numCells+numVertices);
3235:   for (c = 0; c < numCells; ++c) {
3236:     DMPlexSetConeSize(dm, c, numCorners);
3237:   }
3238:   DMSetUp(dm);
3239:   DMPlexGetCones(dm,&cones);
3240:   for (c = 0; c < numCells; ++c) {
3241:     for (p = 0; p < numCorners; ++p) {
3242:       cones[c*numCorners+p] = cells[c*numCorners+p]+numCells;
3243:     }
3244:   }
3245:   DMPlexSymmetrize(dm);
3246:   DMPlexStratify(dm);
3247:   PetscLogEventEnd(DMPLEX_BuildFromCellList,dm,0,0,0);
3248:   return(0);
3249: }

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

3254:   Input Parameters:
3255: + dm - The DM
3256: . spaceDim - The spatial dimension used for coordinates
3257: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3259:   Level: advanced

3261:   Notes:
3262:   Not currently supported in Fortran.

3264: .seealso: DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellList()
3265: @*/
3266: PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM dm, PetscInt spaceDim, const PetscReal vertexCoords[])
3267: {
3268:   PetscSection   coordSection;
3269:   Vec            coordinates;
3270:   DM             cdm;
3271:   PetscScalar   *coords;
3272:   PetscInt       v, vStart, vEnd, d;

3276:   PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3277:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
3278:   if (vStart < 0 || vEnd < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first.");
3279:   DMSetCoordinateDim(dm, spaceDim);
3280:   DMGetCoordinateSection(dm, &coordSection);
3281:   PetscSectionSetNumFields(coordSection, 1);
3282:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
3283:   PetscSectionSetChart(coordSection, vStart, vEnd);
3284:   for (v = vStart; v < vEnd; ++v) {
3285:     PetscSectionSetDof(coordSection, v, spaceDim);
3286:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
3287:   }
3288:   PetscSectionSetUp(coordSection);

3290:   DMGetCoordinateDM(dm, &cdm);
3291:   DMCreateLocalVector(cdm, &coordinates);
3292:   VecSetBlockSize(coordinates, spaceDim);
3293:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3294:   VecGetArrayWrite(coordinates, &coords);
3295:   for (v = 0; v < vEnd-vStart; ++v) {
3296:     for (d = 0; d < spaceDim; ++d) {
3297:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
3298:     }
3299:   }
3300:   VecRestoreArrayWrite(coordinates, &coords);
3301:   DMSetCoordinatesLocal(dm, coordinates);
3302:   VecDestroy(&coordinates);
3303:   PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
3304:   return(0);
3305: }

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

3310:   Input Parameters:
3311: + comm - The communicator
3312: . dim - The topological dimension of the mesh
3313: . numCells - The number of cells
3314: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
3315: . numCorners - The number of vertices for each cell
3316: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
3317: . cells - An array of numCells*numCorners numbers, the vertices for each cell
3318: . spaceDim - The spatial dimension used for coordinates
3319: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

3321:   Output Parameter:
3322: . dm - The DM

3324:   Notes:
3325:   This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(), DMPlexBuildFromCellList(),
3326:   DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellList()

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

3331:   Level: intermediate

3333: .seealso: DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildFromCellList(), DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
3334: @*/
3335: 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)
3336: {

3340:   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.");
3341:   DMCreate(comm, dm);
3342:   DMSetType(*dm, DMPLEX);
3343:   DMSetDimension(*dm, dim);
3344:   DMPlexBuildFromCellList(*dm, numCells, numVertices, numCorners, cells);
3345:   if (interpolate) {
3346:     DM idm;

3348:     DMPlexInterpolate(*dm, &idm);
3349:     DMDestroy(dm);
3350:     *dm  = idm;
3351:   }
3352:   DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, vertexCoords);
3353:   return(0);
3354: }

3356: /*@
3357:   DMPlexCreateFromCellList - Deprecated, use DMPlexCreateFromCellListPetsc()

3359:   Level: deprecated

3361: .seealso: DMPlexCreateFromCellListPetsc()
3362: @*/
3363: 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)
3364: {
3366:   PetscInt       i;
3367:   PetscInt       *pintCells;
3368:   PetscReal      *prealVC;

3371:   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));
3372:   if (sizeof(int) == sizeof(PetscInt)) {
3373:     pintCells = (PetscInt *) cells;
3374:   } else {
3375:     PetscMalloc1(numCells*numCorners, &pintCells);
3376:     for (i = 0; i < numCells*numCorners; i++) {
3377:       pintCells[i] = (PetscInt) cells[i];
3378:     }
3379:   }
3380:   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));
3381:   if (sizeof(double) == sizeof(PetscReal)) {
3382:     prealVC = (PetscReal *) vertexCoords;
3383:   } else {
3384:     PetscMalloc1(numVertices*spaceDim, &prealVC);
3385:     for (i = 0; i < numVertices*spaceDim; i++) {
3386:       prealVC[i] = (PetscReal) vertexCoords[i];
3387:     }
3388:   }
3389:   DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, pintCells, spaceDim, prealVC, dm);
3390:   if (sizeof(int) != sizeof(PetscInt)) {
3391:     PetscFree(pintCells);
3392:   }
3393:   if (sizeof(double) != sizeof(PetscReal)) {
3394:     PetscFree(prealVC);
3395:   }
3396:   return(0);
3397: }

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

3402:   Input Parameters:
3403: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
3404: . depth - The depth of the DAG
3405: . numPoints - Array of size depth + 1 containing the number of points at each depth
3406: . coneSize - The cone size of each point
3407: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
3408: . coneOrientations - The orientation of each cone point
3409: - vertexCoords - An array of numPoints[0]*spacedim numbers representing the coordinates of each vertex, with spacedim the value set via DMSetCoordinateDim()

3411:   Output Parameter:
3412: . dm - The DM

3414:   Note: Two triangles sharing a face would have input
3415: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
3416: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
3417: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
3418: $
3419: which would result in the DMPlex
3420: $
3421: $        4
3422: $      / | \
3423: $     /  |  \
3424: $    /   |   \
3425: $   2  0 | 1  5
3426: $    \   |   /
3427: $     \  |  /
3428: $      \ | /
3429: $        3
3430: $
3431: $ Notice that all points are numbered consecutively, unlike DMPlexCreateFromCellListPetsc()

3433:   Level: advanced

3435: .seealso: DMPlexCreateFromCellListPetsc(), DMPlexCreate()
3436: @*/
3437: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
3438: {
3439:   Vec            coordinates;
3440:   PetscSection   coordSection;
3441:   PetscScalar    *coords;
3442:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

3446:   DMGetDimension(dm, &dim);
3447:   DMGetCoordinateDim(dm, &dimEmbed);
3448:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %D cannot be less than intrinsic dimension %d",dimEmbed,dim);
3449:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
3450:   DMPlexSetChart(dm, pStart, pEnd);
3451:   for (p = pStart; p < pEnd; ++p) {
3452:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
3453:     if (firstVertex < 0 && !coneSize[p - pStart]) {
3454:       firstVertex = p - pStart;
3455:     }
3456:   }
3457:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %D vertices but could not find any", numPoints[0]);
3458:   DMSetUp(dm); /* Allocate space for cones */
3459:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
3460:     DMPlexSetCone(dm, p, &cones[off]);
3461:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
3462:   }
3463:   DMPlexSymmetrize(dm);
3464:   DMPlexStratify(dm);
3465:   /* Build coordinates */
3466:   DMGetCoordinateSection(dm, &coordSection);
3467:   PetscSectionSetNumFields(coordSection, 1);
3468:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
3469:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
3470:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
3471:     PetscSectionSetDof(coordSection, v, dimEmbed);
3472:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
3473:   }
3474:   PetscSectionSetUp(coordSection);
3475:   PetscSectionGetStorageSize(coordSection, &coordSize);
3476:   VecCreate(PETSC_COMM_SELF, &coordinates);
3477:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3478:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3479:   VecSetBlockSize(coordinates, dimEmbed);
3480:   VecSetType(coordinates,VECSTANDARD);
3481:   VecGetArray(coordinates, &coords);
3482:   for (v = 0; v < numPoints[0]; ++v) {
3483:     PetscInt off;

3485:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
3486:     for (d = 0; d < dimEmbed; ++d) {
3487:       coords[off+d] = vertexCoords[v*dimEmbed+d];
3488:     }
3489:   }
3490:   VecRestoreArray(coordinates, &coords);
3491:   DMSetCoordinatesLocal(dm, coordinates);
3492:   VecDestroy(&coordinates);
3493:   return(0);
3494: }

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

3499: + comm        - The MPI communicator
3500: . filename    - Name of the .dat file
3501: - interpolate - Create faces and edges in the mesh

3503:   Output Parameter:
3504: . dm  - The DM object representing the mesh

3506:   Note: The format is the simplest possible:
3507: $ Ne
3508: $ v0 v1 ... vk
3509: $ Nv
3510: $ x y z marker

3512:   Level: beginner

3514: .seealso: DMPlexCreateFromFile(), DMPlexCreateMedFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
3515: @*/
3516: PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3517: {
3518:   DMLabel         marker;
3519:   PetscViewer     viewer;
3520:   Vec             coordinates;
3521:   PetscSection    coordSection;
3522:   PetscScalar    *coords;
3523:   char            line[PETSC_MAX_PATH_LEN];
3524:   PetscInt        dim = 3, cdim = 3, coordSize, v, c, d;
3525:   PetscMPIInt     rank;
3526:   int             snum, Nv, Nc;
3527:   PetscErrorCode  ierr;

3530:   MPI_Comm_rank(comm, &rank);
3531:   PetscViewerCreate(comm, &viewer);
3532:   PetscViewerSetType(viewer, PETSCVIEWERASCII);
3533:   PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3534:   PetscViewerFileSetName(viewer, filename);
3535:   if (!rank) {
3536:     PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);
3537:     snum = sscanf(line, "%d %d", &Nc, &Nv);
3538:     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3539:   } else {
3540:     Nc = Nv = 0;
3541:   }
3542:   DMCreate(comm, dm);
3543:   DMSetType(*dm, DMPLEX);
3544:   DMPlexSetChart(*dm, 0, Nc+Nv);
3545:   DMSetDimension(*dm, dim);
3546:   DMSetCoordinateDim(*dm, cdim);
3547:   /* Read topology */
3548:   if (!rank) {
3549:     PetscInt cone[8], corners = 8;
3550:     int      vbuf[8], v;

3552:     for (c = 0; c < Nc; ++c) {DMPlexSetConeSize(*dm, c, corners);}
3553:     DMSetUp(*dm);
3554:     for (c = 0; c < Nc; ++c) {
3555:       PetscViewerRead(viewer, line, corners, NULL, PETSC_STRING);
3556:       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]);
3557:       if (snum != corners) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3558:       for (v = 0; v < corners; ++v) cone[v] = vbuf[v] + Nc;
3559:       /* Hexahedra are inverted */
3560:       {
3561:         PetscInt tmp = cone[1];
3562:         cone[1] = cone[3];
3563:         cone[3] = tmp;
3564:       }
3565:       DMPlexSetCone(*dm, c, cone);
3566:     }
3567:   }
3568:   DMPlexSymmetrize(*dm);
3569:   DMPlexStratify(*dm);
3570:   /* Read coordinates */
3571:   DMGetCoordinateSection(*dm, &coordSection);
3572:   PetscSectionSetNumFields(coordSection, 1);
3573:   PetscSectionSetFieldComponents(coordSection, 0, cdim);
3574:   PetscSectionSetChart(coordSection, Nc, Nc + Nv);
3575:   for (v = Nc; v < Nc+Nv; ++v) {
3576:     PetscSectionSetDof(coordSection, v, cdim);
3577:     PetscSectionSetFieldDof(coordSection, v, 0, cdim);
3578:   }
3579:   PetscSectionSetUp(coordSection);
3580:   PetscSectionGetStorageSize(coordSection, &coordSize);
3581:   VecCreate(PETSC_COMM_SELF, &coordinates);
3582:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3583:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3584:   VecSetBlockSize(coordinates, cdim);
3585:   VecSetType(coordinates, VECSTANDARD);
3586:   VecGetArray(coordinates, &coords);
3587:   if (!rank) {
3588:     double x[3];
3589:     int    val;

3591:     DMCreateLabel(*dm, "marker");
3592:     DMGetLabel(*dm, "marker", &marker);
3593:     for (v = 0; v < Nv; ++v) {
3594:       PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING);
3595:       snum = sscanf(line, "%lg %lg %lg %d", &x[0], &x[1], &x[2], &val);
3596:       if (snum != 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3597:       for (d = 0; d < cdim; ++d) coords[v*cdim+d] = x[d];
3598:       if (val) {DMLabelSetValue(marker, v+Nc, val);}
3599:     }
3600:   }
3601:   VecRestoreArray(coordinates, &coords);
3602:   DMSetCoordinatesLocal(*dm, coordinates);
3603:   VecDestroy(&coordinates);
3604:   PetscViewerDestroy(&viewer);
3605:   if (interpolate) {
3606:     DM      idm;
3607:     DMLabel bdlabel;

3609:     DMPlexInterpolate(*dm, &idm);
3610:     DMDestroy(dm);
3611:     *dm  = idm;

3613:     DMGetLabel(*dm, "marker", &bdlabel);
3614:     DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
3615:     DMPlexLabelComplete(*dm, bdlabel);
3616:   }
3617:   return(0);
3618: }

3620: /*@C
3621:   DMPlexCreateFromFile - This takes a filename and produces a DM

3623:   Input Parameters:
3624: + comm - The communicator
3625: . filename - A file name
3626: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

3628:   Output Parameter:
3629: . dm - The DM

3631:   Options Database Keys:
3632: . -dm_plex_create_from_hdf5_xdmf - use the PETSC_VIEWER_HDF5_XDMF format for reading HDF5

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

3637:   Level: beginner

3639: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellListPetsc(), DMPlexCreate()
3640: @*/
3641: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3642: {
3643:   const char    *extGmsh    = ".msh";
3644:   const char    *extGmsh2   = ".msh2";
3645:   const char    *extGmsh4   = ".msh4";
3646:   const char    *extCGNS    = ".cgns";
3647:   const char    *extExodus  = ".exo";
3648:   const char    *extGenesis = ".gen";
3649:   const char    *extFluent  = ".cas";
3650:   const char    *extHDF5    = ".h5";
3651:   const char    *extMed     = ".med";
3652:   const char    *extPLY     = ".ply";
3653:   const char    *extCV      = ".dat";
3654:   size_t         len;
3655:   PetscBool      isGmsh, isGmsh2, isGmsh4, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isCV;
3656:   PetscMPIInt    rank;

3662:   DMInitializePackage();
3663:   PetscLogEventBegin(DMPLEX_CreateFromFile,0,0,0,0);
3664:   MPI_Comm_rank(comm, &rank);
3665:   PetscStrlen(filename, &len);
3666:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
3667:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,    4, &isGmsh);
3668:   PetscStrncmp(&filename[PetscMax(0,len-5)], extGmsh2,   5, &isGmsh2);
3669:   PetscStrncmp(&filename[PetscMax(0,len-5)], extGmsh4,   5, &isGmsh4);
3670:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,    5, &isCGNS);
3671:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus,  4, &isExodus);
3672:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);
3673:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent,  4, &isFluent);
3674:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,    3, &isHDF5);
3675:   PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,     4, &isMed);
3676:   PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY,     4, &isPLY);
3677:   PetscStrncmp(&filename[PetscMax(0,len-4)], extCV,      4, &isCV);
3678:   if (isGmsh || isGmsh2 || isGmsh4) {
3679:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
3680:   } else if (isCGNS) {
3681:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
3682:   } else if (isExodus || isGenesis) {
3683:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
3684:   } else if (isFluent) {
3685:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
3686:   } else if (isHDF5) {
3687:     PetscBool      load_hdf5_xdmf = PETSC_FALSE;
3688:     PetscViewer viewer;

3690:     /* PETSC_VIEWER_HDF5_XDMF is used if the filename ends with .xdmf.h5, or if -dm_plex_create_from_hdf5_xdmf option is present */
3691:     PetscStrncmp(&filename[PetscMax(0,len-8)], ".xdmf",  5, &load_hdf5_xdmf);
3692:     PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &load_hdf5_xdmf, NULL);
3693:     PetscViewerCreate(comm, &viewer);
3694:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
3695:     PetscViewerSetOptionsPrefix(viewer, "dm_plex_create_");
3696:     PetscViewerSetFromOptions(viewer);
3697:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3698:     PetscViewerFileSetName(viewer, filename);
3699:     DMCreate(comm, dm);
3700:     DMSetType(*dm, DMPLEX);
3701:     if (load_hdf5_xdmf) {PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF);}
3702:     DMLoad(*dm, viewer);
3703:     if (load_hdf5_xdmf) {PetscViewerPopFormat(viewer);}
3704:     PetscViewerDestroy(&viewer);

3706:     if (interpolate) {
3707:       DM idm;

3709:       DMPlexInterpolate(*dm, &idm);
3710:       DMDestroy(dm);
3711:       *dm  = idm;
3712:     }
3713:   } else if (isMed) {
3714:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
3715:   } else if (isPLY) {
3716:     DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);
3717:   } else if (isCV) {
3718:     DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm);
3719:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
3720:   PetscLogEventEnd(DMPLEX_CreateFromFile,0,0,0,0);
3721:   return(0);
3722: }

3724: /*@
3725:   DMPlexCreateReferenceCellByType - Create a DMPLEX with the appropriate FEM reference cell

3727:   Collective

3729:   Input Parameters:
3730: + comm - The communicator
3731: - ct   - The cell type of the reference cell

3733:   Output Parameter:
3734: . refdm - The reference cell

3736:   Level: intermediate

3738: .seealso: DMPlexCreateReferenceCell(), DMPlexCreateBoxMesh()
3739: @*/
3740: PetscErrorCode DMPlexCreateReferenceCellByType(MPI_Comm comm, DMPolytopeType ct, DM *refdm)
3741: {
3742:   DM             rdm;
3743:   Vec            coords;

3747:   DMCreate(comm, &rdm);
3748:   DMSetType(rdm, DMPLEX);
3749:   switch (ct) {
3750:     case DM_POLYTOPE_POINT:
3751:     {
3752:       PetscInt    numPoints[1]        = {1};
3753:       PetscInt    coneSize[1]         = {0};
3754:       PetscInt    cones[1]            = {0};
3755:       PetscInt    coneOrientations[1] = {0};
3756:       PetscScalar vertexCoords[1]     = {0.0};

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

3770:       DMSetDimension(rdm, 1);
3771:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3772:     }
3773:     break;
3774:     case DM_POLYTOPE_TRIANGLE:
3775:     {
3776:       PetscInt    numPoints[2]        = {3, 1};
3777:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
3778:       PetscInt    cones[3]            = {1, 2, 3};
3779:       PetscInt    coneOrientations[3] = {0, 0, 0};
3780:       PetscScalar vertexCoords[6]     = {-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_QUADRILATERAL:
3787:     {
3788:       PetscInt    numPoints[2]        = {4, 1};
3789:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3790:       PetscInt    cones[4]            = {1, 2, 3, 4};
3791:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3792:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

3794:       DMSetDimension(rdm, 2);
3795:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3796:     }
3797:     break;
3798:     case DM_POLYTOPE_SEG_PRISM_TENSOR:
3799:     {
3800:       PetscInt    numPoints[2]        = {4, 1};
3801:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3802:       PetscInt    cones[4]            = {1, 2, 3, 4};
3803:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3804:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0,  1.0, 1.0};

3806:       DMSetDimension(rdm, 2);
3807:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3808:     }
3809:     break;
3810:     case DM_POLYTOPE_TETRAHEDRON:
3811:     {
3812:       PetscInt    numPoints[2]        = {4, 1};
3813:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3814:       PetscInt    cones[4]            = {1, 3, 2, 4};
3815:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3816:       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};

3818:       DMSetDimension(rdm, 3);
3819:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3820:     }
3821:     break;
3822:     case DM_POLYTOPE_HEXAHEDRON:
3823:     {
3824:       PetscInt    numPoints[2]        = {8, 1};
3825:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3826:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
3827:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3828:       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,
3829:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3831:       DMSetDimension(rdm, 3);
3832:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3833:     }
3834:     break;
3835:     case DM_POLYTOPE_TRI_PRISM:
3836:     {
3837:       PetscInt    numPoints[2]        = {6, 1};
3838:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3839:       PetscInt    cones[6]            = {1, 3, 2, 4, 5, 6};
3840:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3841:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3842:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3844:       DMSetDimension(rdm, 3);
3845:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3846:     }
3847:     break;
3848:     case DM_POLYTOPE_TRI_PRISM_TENSOR:
3849:     {
3850:       PetscInt    numPoints[2]        = {6, 1};
3851:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3852:       PetscInt    cones[6]            = {1, 2, 3, 4, 5, 6};
3853:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3854:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3855:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3857:       DMSetDimension(rdm, 3);
3858:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3859:     }
3860:     break;
3861:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
3862:     {
3863:       PetscInt    numPoints[2]        = {8, 1};
3864:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3865:       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
3866:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3867:       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,
3868:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3870:       DMSetDimension(rdm, 3);
3871:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3872:     }
3873:     break;
3874:     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for cell type %s", DMPolytopeTypes[ct]);
3875:   }
3876:   {
3877:     PetscInt Nv, v;

3879:     /* Must create the celltype label here so that we do not automatically try to compute the types */
3880:     DMCreateLabel(rdm, "celltype");
3881:     DMPlexSetCellType(rdm, 0, ct);
3882:     DMPlexGetChart(rdm, NULL, &Nv);
3883:     for (v = 1; v < Nv; ++v) {DMPlexSetCellType(rdm, v, DM_POLYTOPE_POINT);}
3884:   }
3885:   DMPlexInterpolate(rdm, refdm);
3886:   if (rdm->coordinateDM) {
3887:     DM           ncdm;
3888:     PetscSection cs;
3889:     PetscInt     pEnd = -1;

3891:     DMGetLocalSection(rdm->coordinateDM, &cs);
3892:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
3893:     if (pEnd >= 0) {
3894:       DMClone(*refdm, &ncdm);
3895:       DMCopyDisc(rdm->coordinateDM, ncdm);
3896:       DMSetLocalSection(ncdm, cs);
3897:       DMSetCoordinateDM(*refdm, ncdm);
3898:       DMDestroy(&ncdm);
3899:     }
3900:   }
3901:   DMGetCoordinatesLocal(rdm, &coords);
3902:   if (coords) {
3903:     DMSetCoordinatesLocal(*refdm, coords);
3904:   } else {
3905:     DMGetCoordinates(rdm, &coords);
3906:     if (coords) {DMSetCoordinates(*refdm, coords);}
3907:   }
3908:   DMDestroy(&rdm);
3909:   return(0);
3910: }

3912: /*@
3913:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

3915:   Collective

3917:   Input Parameters:
3918: + comm    - The communicator
3919: . dim     - The spatial dimension
3920: - simplex - Flag for simplex, otherwise use a tensor-product cell

3922:   Output Parameter:
3923: . refdm - The reference cell

3925:   Level: intermediate

3927: .seealso: DMPlexCreateReferenceCellByType(), DMPlexCreateBoxMesh()
3928: @*/
3929: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
3930: {

3934:   switch (dim) {
3935:   case 0: DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_POINT, refdm);break;
3936:   case 1: DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_SEGMENT, refdm);break;
3937:   case 2:
3938:     if (simplex) {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_TRIANGLE, refdm);}
3939:     else         {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_QUADRILATERAL, refdm);}
3940:     break;
3941:   case 3:
3942:     if (simplex) {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_TETRAHEDRON, refdm);}
3943:     else         {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_HEXAHEDRON, refdm);}
3944:     break;
3945:   default:
3946:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %D", dim);
3947:   }
3948:   return(0);
3949: }