Actual source code: plexcreate.c

  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>
  3: #include <petsc/private/hashseti.h>
  4: #include <petscsf.h>
  5: #include <petscdmplextransform.h>
  6: #include <petsc/private/kernels/blockmatmult.h>
  7: #include <petsc/private/kernels/blockinvert.h>

  9: PetscLogEvent DMPLEX_CreateFromFile, DMPLEX_BuildFromCellList, DMPLEX_BuildCoordinatesFromCellList;

 11: /* External function declarations here */
 12: static PetscErrorCode DMInitialize_Plex(DM dm);

 14: /* This copies internal things in the Plex structure that we generally want when making a new, related Plex */
 15: PetscErrorCode DMPlexCopy_Internal(DM dmin, PetscBool copyPeriodicity, DM dmout)
 16: {
 17:   const DMBoundaryType *bd;
 18:   const PetscReal      *maxCell, *L;
 19:   PetscBool             isper, dist;

 21:   if (copyPeriodicity) {
 22:     DMGetPeriodicity(dmin, &isper, &maxCell, &L, &bd);
 23:     DMSetPeriodicity(dmout, isper,  maxCell,  L,  bd);
 24:   }
 25:   DMPlexDistributeGetDefault(dmin, &dist);
 26:   DMPlexDistributeSetDefault(dmout, dist);
 27:   ((DM_Plex *) dmout->data)->useHashLocation = ((DM_Plex *) dmin->data)->useHashLocation;
 28:   return 0;
 29: }

 31: /* Replace dm with the contents of ndm, and then destroy ndm
 32:    - Share the DM_Plex structure
 33:    - Share the coordinates
 34:    - Share the SF
 35: */
 36: static PetscErrorCode DMPlexReplace_Static(DM dm, DM *ndm)
 37: {
 38:   PetscSF               sf;
 39:   DM                    dmNew = *ndm, coordDM, coarseDM;
 40:   Vec                   coords;
 41:   PetscBool             isper;
 42:   const PetscReal      *maxCell, *L;
 43:   const DMBoundaryType *bd;
 44:   PetscInt              dim, cdim;

 46:   if (dm == dmNew) {
 47:     DMDestroy(ndm);
 48:     return 0;
 49:   }
 50:   dm->setupcalled = dmNew->setupcalled;
 51:   DMGetDimension(dmNew, &dim);
 52:   DMSetDimension(dm, dim);
 53:   DMGetCoordinateDim(dmNew, &cdim);
 54:   DMSetCoordinateDim(dm, cdim);
 55:   DMGetPointSF(dmNew, &sf);
 56:   DMSetPointSF(dm, sf);
 57:   DMGetCoordinateDM(dmNew, &coordDM);
 58:   DMGetCoordinatesLocal(dmNew, &coords);
 59:   DMSetCoordinateDM(dm, coordDM);
 60:   DMSetCoordinatesLocal(dm, coords);
 61:   /* Do not want to create the coordinate field if it does not already exist, so do not call DMGetCoordinateField() */
 62:   DMFieldDestroy(&dm->coordinateField);
 63:   dm->coordinateField = dmNew->coordinateField;
 64:   ((DM_Plex *) dmNew->data)->coordFunc = ((DM_Plex *) dm->data)->coordFunc;
 65:   DMGetPeriodicity(dmNew, &isper, &maxCell, &L, &bd);
 66:   DMSetPeriodicity(dm, isper, maxCell, L, bd);
 67:   DMDestroy_Plex(dm);
 68:   DMInitialize_Plex(dm);
 69:   dm->data = dmNew->data;
 70:   ((DM_Plex *) dmNew->data)->refct++;
 71:   DMDestroyLabelLinkList_Internal(dm);
 72:   DMCopyLabels(dmNew, dm, PETSC_OWN_POINTER, PETSC_TRUE, DM_COPY_LABELS_FAIL);
 73:   DMGetCoarseDM(dmNew,&coarseDM);
 74:   DMSetCoarseDM(dm,coarseDM);
 75:   DMDestroy(ndm);
 76:   return 0;
 77: }

 79: /* Swap dm with the contents of dmNew
 80:    - Swap the DM_Plex structure
 81:    - Swap the coordinates
 82:    - Swap the point PetscSF
 83: */
 84: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
 85: {
 86:   DM              coordDMA, coordDMB;
 87:   Vec             coordsA,  coordsB;
 88:   PetscSF         sfA,      sfB;
 89:   DMField         fieldTmp;
 90:   void            *tmp;
 91:   DMLabelLink     listTmp;
 92:   DMLabel         depthTmp;
 93:   PetscInt        tmpI;

 95:   if (dmA == dmB) return 0;
 96:   DMGetPointSF(dmA, &sfA);
 97:   DMGetPointSF(dmB, &sfB);
 98:   PetscObjectReference((PetscObject) sfA);
 99:   DMSetPointSF(dmA, sfB);
100:   DMSetPointSF(dmB, sfA);
101:   PetscObjectDereference((PetscObject) sfA);

103:   DMGetCoordinateDM(dmA, &coordDMA);
104:   DMGetCoordinateDM(dmB, &coordDMB);
105:   PetscObjectReference((PetscObject) coordDMA);
106:   DMSetCoordinateDM(dmA, coordDMB);
107:   DMSetCoordinateDM(dmB, coordDMA);
108:   PetscObjectDereference((PetscObject) coordDMA);

110:   DMGetCoordinatesLocal(dmA, &coordsA);
111:   DMGetCoordinatesLocal(dmB, &coordsB);
112:   PetscObjectReference((PetscObject) coordsA);
113:   DMSetCoordinatesLocal(dmA, coordsB);
114:   DMSetCoordinatesLocal(dmB, coordsA);
115:   PetscObjectDereference((PetscObject) coordsA);

117:   fieldTmp             = dmA->coordinateField;
118:   dmA->coordinateField = dmB->coordinateField;
119:   dmB->coordinateField = fieldTmp;
120:   tmp       = dmA->data;
121:   dmA->data = dmB->data;
122:   dmB->data = tmp;
123:   listTmp   = dmA->labels;
124:   dmA->labels = dmB->labels;
125:   dmB->labels = listTmp;
126:   depthTmp  = dmA->depthLabel;
127:   dmA->depthLabel = dmB->depthLabel;
128:   dmB->depthLabel = depthTmp;
129:   depthTmp  = dmA->celltypeLabel;
130:   dmA->celltypeLabel = dmB->celltypeLabel;
131:   dmB->celltypeLabel = depthTmp;
132:   tmpI         = dmA->levelup;
133:   dmA->levelup = dmB->levelup;
134:   dmB->levelup = tmpI;
135:   return 0;
136: }

138: static PetscErrorCode DMPlexInterpolateInPlace_Internal(DM dm)
139: {
140:   DM             idm;

142:   DMPlexInterpolate(dm, &idm);
143:   DMPlexCopyCoordinates(dm, idm);
144:   DMPlexReplace_Static(dm, &idm);
145:   return 0;
146: }

148: /*@C
149:   DMPlexCreateCoordinateSpace - Creates a finite element space for the coordinates

151:   Collective

153:   Input Parameters:
154: + DM        - The DM
155: . degree    - The degree of the finite element or PETSC_DECIDE
156: - coordFunc - An optional function to map new points from refinement to the surface

158:   Level: advanced

160: .seealso: PetscFECreateLagrange(), DMGetCoordinateDM()
161: @*/
162: PetscErrorCode DMPlexCreateCoordinateSpace(DM dm, PetscInt degree, PetscPointFunc coordFunc)
163: {
164:   DM_Plex      *mesh = (DM_Plex *) dm->data;
165:   DM            cdm;
166:   PetscDS       cds;
167:   PetscFE       fe;
168:   PetscClassId  id;

170:   DMGetCoordinateDM(dm, &cdm);
171:   DMGetDS(cdm, &cds);
172:   PetscDSGetDiscretization(cds, 0, (PetscObject *) &fe);
173:   PetscObjectGetClassId((PetscObject) fe, &id);
174:   if (id != PETSCFE_CLASSID) {
175:     PetscBool      simplex;
176:     PetscInt       dim, dE, qorder;

179:     DMGetDimension(dm, &dim);
180:     DMGetCoordinateDim(dm, &dE);
181:     qorder = degree;
182:     PetscObjectOptionsBegin((PetscObject) cdm);
183:     PetscOptionsBoundedInt("-coord_dm_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "DMPlexCreateCoordinateSpace", qorder, &qorder, NULL, 0);
184:     PetscOptionsEnd();
185:     if (degree == PETSC_DECIDE) fe = NULL;
186:     else {
187:       DMPlexIsSimplex(dm, &simplex);
188:       PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, simplex, degree, qorder, &fe);
189:     }
190:     DMProjectCoordinates(dm, fe);
191:     PetscFEDestroy(&fe);
192:   }
193:   mesh->coordFunc = coordFunc;
194:   return 0;
195: }

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

200:   Collective

202:   Input Parameters:
203: + comm - The communicator for the DM object
204: . dim - The spatial dimension
205: . simplex - Flag for simplicial cells, otherwise they are tensor product cells
206: . interpolate - Flag to create intermediate mesh pieces (edges, faces)
207: - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell

209:   Output Parameter:
210: . dm - The DM object

212:   Level: beginner

214: .seealso: DMSetType(), DMCreate()
215: @*/
216: PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscReal refinementLimit, DM *newdm)
217: {
218:   DM             dm;
219:   PetscMPIInt    rank;

221:   DMCreate(comm, &dm);
222:   DMSetType(dm, DMPLEX);
223:   DMSetDimension(dm, dim);
224:   MPI_Comm_rank(comm, &rank);
225:   switch (dim) {
226:   case 2:
227:     if (simplex) PetscObjectSetName((PetscObject) dm, "triangular");
228:     else         PetscObjectSetName((PetscObject) dm, "quadrilateral");
229:     break;
230:   case 3:
231:     if (simplex) PetscObjectSetName((PetscObject) dm, "tetrahedral");
232:     else         PetscObjectSetName((PetscObject) dm, "hexahedral");
233:     break;
234:   default:
235:     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %D", dim);
236:   }
237:   if (rank) {
238:     PetscInt numPoints[2] = {0, 0};
239:     DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);
240:   } else {
241:     switch (dim) {
242:     case 2:
243:       if (simplex) {
244:         PetscInt    numPoints[2]        = {4, 2};
245:         PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
246:         PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
247:         PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
248:         PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};

250:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
251:       } else {
252:         PetscInt    numPoints[2]        = {6, 2};
253:         PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
254:         PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
255:         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
256:         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};

258:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
259:       }
260:       break;
261:     case 3:
262:       if (simplex) {
263:         PetscInt    numPoints[2]        = {5, 2};
264:         PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
265:         PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
266:         PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
267:         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};

269:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
270:       } else {
271:         PetscInt    numPoints[2]         = {12, 2};
272:         PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
273:         PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
274:         PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
275:         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,
276:                                             -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
277:                                              1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};

279:         DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
280:       }
281:       break;
282:     default:
283:       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %D", dim);
284:     }
285:   }
286:   *newdm = dm;
287:   if (refinementLimit > 0.0) {
288:     DM rdm;
289:     const char *name;

291:     DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);
292:     DMPlexSetRefinementLimit(*newdm, refinementLimit);
293:     DMRefine(*newdm, comm, &rdm);
294:     PetscObjectGetName((PetscObject) *newdm, &name);
295:     PetscObjectSetName((PetscObject)    rdm,  name);
296:     DMDestroy(newdm);
297:     *newdm = rdm;
298:   }
299:   if (interpolate) {
300:     DM idm;

302:     DMPlexInterpolate(*newdm, &idm);
303:     DMDestroy(newdm);
304:     *newdm = idm;
305:   }
306:   return 0;
307: }

309: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_1D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
310: {
311:   const PetscInt numVertices    = 2;
312:   PetscInt       markerRight    = 1;
313:   PetscInt       markerLeft     = 1;
314:   PetscBool      markerSeparate = PETSC_FALSE;
315:   Vec            coordinates;
316:   PetscSection   coordSection;
317:   PetscScalar   *coords;
318:   PetscInt       coordSize;
319:   PetscMPIInt    rank;
320:   PetscInt       cdim = 1, v;

322:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
323:   if (markerSeparate) {
324:     markerRight  = 2;
325:     markerLeft   = 1;
326:   }
327:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
328:   if (!rank) {
329:     DMPlexSetChart(dm, 0, numVertices);
330:     DMSetUp(dm); /* Allocate space for cones */
331:     DMSetLabelValue(dm, "marker", 0, markerLeft);
332:     DMSetLabelValue(dm, "marker", 1, markerRight);
333:   }
334:   DMPlexSymmetrize(dm);
335:   DMPlexStratify(dm);
336:   /* Build coordinates */
337:   DMSetCoordinateDim(dm, cdim);
338:   DMGetCoordinateSection(dm, &coordSection);
339:   PetscSectionSetNumFields(coordSection, 1);
340:   PetscSectionSetChart(coordSection, 0, numVertices);
341:   PetscSectionSetFieldComponents(coordSection, 0, cdim);
342:   for (v = 0; v < numVertices; ++v) {
343:     PetscSectionSetDof(coordSection, v, cdim);
344:     PetscSectionSetFieldDof(coordSection, v, 0, cdim);
345:   }
346:   PetscSectionSetUp(coordSection);
347:   PetscSectionGetStorageSize(coordSection, &coordSize);
348:   VecCreate(PETSC_COMM_SELF, &coordinates);
349:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
350:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
351:   VecSetBlockSize(coordinates, cdim);
352:   VecSetType(coordinates,VECSTANDARD);
353:   VecGetArray(coordinates, &coords);
354:   coords[0] = lower[0];
355:   coords[1] = upper[0];
356:   VecRestoreArray(coordinates, &coords);
357:   DMSetCoordinatesLocal(dm, coordinates);
358:   VecDestroy(&coordinates);
359:   return 0;
360: }

362: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_2D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
363: {
364:   const PetscInt numVertices    = (edges[0]+1)*(edges[1]+1);
365:   const PetscInt numEdges       = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
366:   PetscInt       markerTop      = 1;
367:   PetscInt       markerBottom   = 1;
368:   PetscInt       markerRight    = 1;
369:   PetscInt       markerLeft     = 1;
370:   PetscBool      markerSeparate = PETSC_FALSE;
371:   Vec            coordinates;
372:   PetscSection   coordSection;
373:   PetscScalar    *coords;
374:   PetscInt       coordSize;
375:   PetscMPIInt    rank;
376:   PetscInt       v, vx, vy;

378:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
379:   if (markerSeparate) {
380:     markerTop    = 3;
381:     markerBottom = 1;
382:     markerRight  = 2;
383:     markerLeft   = 4;
384:   }
385:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
386:   if (rank == 0) {
387:     PetscInt e, ex, ey;

389:     DMPlexSetChart(dm, 0, numEdges+numVertices);
390:     for (e = 0; e < numEdges; ++e) {
391:       DMPlexSetConeSize(dm, e, 2);
392:     }
393:     DMSetUp(dm); /* Allocate space for cones */
394:     for (vx = 0; vx <= edges[0]; vx++) {
395:       for (ey = 0; ey < edges[1]; ey++) {
396:         PetscInt edge   = vx*edges[1] + ey + edges[0]*(edges[1]+1);
397:         PetscInt vertex = ey*(edges[0]+1) + vx + numEdges;
398:         PetscInt cone[2];

400:         cone[0] = vertex; cone[1] = vertex+edges[0]+1;
401:         DMPlexSetCone(dm, edge, cone);
402:         if (vx == edges[0]) {
403:           DMSetLabelValue(dm, "marker", edge,    markerRight);
404:           DMSetLabelValue(dm, "marker", cone[0], markerRight);
405:           if (ey == edges[1]-1) {
406:             DMSetLabelValue(dm, "marker", cone[1], markerRight);
407:             DMSetLabelValue(dm, "Face Sets", cone[1], markerRight);
408:           }
409:         } else if (vx == 0) {
410:           DMSetLabelValue(dm, "marker", edge,    markerLeft);
411:           DMSetLabelValue(dm, "marker", cone[0], markerLeft);
412:           if (ey == edges[1]-1) {
413:             DMSetLabelValue(dm, "marker", cone[1], markerLeft);
414:             DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft);
415:           }
416:         }
417:       }
418:     }
419:     for (vy = 0; vy <= edges[1]; vy++) {
420:       for (ex = 0; ex < edges[0]; ex++) {
421:         PetscInt edge   = vy*edges[0]     + ex;
422:         PetscInt vertex = vy*(edges[0]+1) + ex + numEdges;
423:         PetscInt cone[2];

425:         cone[0] = vertex; cone[1] = vertex+1;
426:         DMPlexSetCone(dm, edge, cone);
427:         if (vy == edges[1]) {
428:           DMSetLabelValue(dm, "marker", edge,    markerTop);
429:           DMSetLabelValue(dm, "marker", cone[0], markerTop);
430:           if (ex == edges[0]-1) {
431:             DMSetLabelValue(dm, "marker", cone[1], markerTop);
432:             DMSetLabelValue(dm, "Face Sets", cone[1], markerTop);
433:           }
434:         } else if (vy == 0) {
435:           DMSetLabelValue(dm, "marker", edge,    markerBottom);
436:           DMSetLabelValue(dm, "marker", cone[0], markerBottom);
437:           if (ex == edges[0]-1) {
438:             DMSetLabelValue(dm, "marker", cone[1], markerBottom);
439:             DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom);
440:           }
441:         }
442:       }
443:     }
444:   }
445:   DMPlexSymmetrize(dm);
446:   DMPlexStratify(dm);
447:   /* Build coordinates */
448:   DMSetCoordinateDim(dm, 2);
449:   DMGetCoordinateSection(dm, &coordSection);
450:   PetscSectionSetNumFields(coordSection, 1);
451:   PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);
452:   PetscSectionSetFieldComponents(coordSection, 0, 2);
453:   for (v = numEdges; v < numEdges+numVertices; ++v) {
454:     PetscSectionSetDof(coordSection, v, 2);
455:     PetscSectionSetFieldDof(coordSection, v, 0, 2);
456:   }
457:   PetscSectionSetUp(coordSection);
458:   PetscSectionGetStorageSize(coordSection, &coordSize);
459:   VecCreate(PETSC_COMM_SELF, &coordinates);
460:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
461:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
462:   VecSetBlockSize(coordinates, 2);
463:   VecSetType(coordinates,VECSTANDARD);
464:   VecGetArray(coordinates, &coords);
465:   for (vy = 0; vy <= edges[1]; ++vy) {
466:     for (vx = 0; vx <= edges[0]; ++vx) {
467:       coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
468:       coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
469:     }
470:   }
471:   VecRestoreArray(coordinates, &coords);
472:   DMSetCoordinatesLocal(dm, coordinates);
473:   VecDestroy(&coordinates);
474:   return 0;
475: }

477: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_3D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
478: {
479:   PetscInt       vertices[3], numVertices;
480:   PetscInt       numFaces    = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2];
481:   Vec            coordinates;
482:   PetscSection   coordSection;
483:   PetscScalar    *coords;
484:   PetscInt       coordSize;
485:   PetscMPIInt    rank;
486:   PetscInt       v, vx, vy, vz;
487:   PetscInt       voffset, iface=0, cone[4];

490:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
491:   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
492:   numVertices = vertices[0]*vertices[1]*vertices[2];
493:   if (rank == 0) {
494:     PetscInt f;

496:     DMPlexSetChart(dm, 0, numFaces+numVertices);
497:     for (f = 0; f < numFaces; ++f) {
498:       DMPlexSetConeSize(dm, f, 4);
499:     }
500:     DMSetUp(dm); /* Allocate space for cones */

502:     /* Side 0 (Top) */
503:     for (vy = 0; vy < faces[1]; vy++) {
504:       for (vx = 0; vx < faces[0]; vx++) {
505:         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
506:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
507:         DMPlexSetCone(dm, iface, cone);
508:         DMSetLabelValue(dm, "marker", iface, 1);
509:         DMSetLabelValue(dm, "marker", voffset+0, 1);
510:         DMSetLabelValue(dm, "marker", voffset+1, 1);
511:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
512:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+1, 1);
513:         iface++;
514:       }
515:     }

517:     /* Side 1 (Bottom) */
518:     for (vy = 0; vy < faces[1]; vy++) {
519:       for (vx = 0; vx < faces[0]; vx++) {
520:         voffset = numFaces + vy*(faces[0]+1) + vx;
521:         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
522:         DMPlexSetCone(dm, iface, cone);
523:         DMSetLabelValue(dm, "marker", iface, 1);
524:         DMSetLabelValue(dm, "marker", voffset+0, 1);
525:         DMSetLabelValue(dm, "marker", voffset+1, 1);
526:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
527:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+1, 1);
528:         iface++;
529:       }
530:     }

532:     /* Side 2 (Front) */
533:     for (vz = 0; vz < faces[2]; vz++) {
534:       for (vx = 0; vx < faces[0]; vx++) {
535:         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
536:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
537:         DMPlexSetCone(dm, iface, cone);
538:         DMSetLabelValue(dm, "marker", iface, 1);
539:         DMSetLabelValue(dm, "marker", voffset+0, 1);
540:         DMSetLabelValue(dm, "marker", voffset+1, 1);
541:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
542:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+1, 1);
543:         iface++;
544:       }
545:     }

547:     /* Side 3 (Back) */
548:     for (vz = 0; vz < faces[2]; vz++) {
549:       for (vx = 0; vx < faces[0]; vx++) {
550:         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
551:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
552:         cone[2] = voffset+1; cone[3] = voffset;
553:         DMPlexSetCone(dm, iface, cone);
554:         DMSetLabelValue(dm, "marker", iface, 1);
555:         DMSetLabelValue(dm, "marker", voffset+0, 1);
556:         DMSetLabelValue(dm, "marker", voffset+1, 1);
557:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
558:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+1, 1);
559:         iface++;
560:       }
561:     }

563:     /* Side 4 (Left) */
564:     for (vz = 0; vz < faces[2]; vz++) {
565:       for (vy = 0; vy < faces[1]; vy++) {
566:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
567:         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
568:         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
569:         DMPlexSetCone(dm, iface, cone);
570:         DMSetLabelValue(dm, "marker", iface, 1);
571:         DMSetLabelValue(dm, "marker", voffset+0, 1);
572:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
573:         DMSetLabelValue(dm, "marker", voffset+vertices[1]+0, 1);
574:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+vertices[0], 1);
575:         iface++;
576:       }
577:     }

579:     /* Side 5 (Right) */
580:     for (vz = 0; vz < faces[2]; vz++) {
581:       for (vy = 0; vy < faces[1]; vy++) {
582:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + faces[0];
583:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset;
584:         cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0];
585:         DMPlexSetCone(dm, iface, cone);
586:         DMSetLabelValue(dm, "marker", iface, 1);
587:         DMSetLabelValue(dm, "marker", voffset+0, 1);
588:         DMSetLabelValue(dm, "marker", voffset+vertices[0]+0, 1);
589:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+0, 1);
590:         DMSetLabelValue(dm, "marker", voffset+vertices[0]*vertices[1]+vertices[0], 1);
591:         iface++;
592:       }
593:     }
594:   }
595:   DMPlexSymmetrize(dm);
596:   DMPlexStratify(dm);
597:   /* Build coordinates */
598:   DMSetCoordinateDim(dm, 3);
599:   DMGetCoordinateSection(dm, &coordSection);
600:   PetscSectionSetNumFields(coordSection, 1);
601:   PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);
602:   PetscSectionSetFieldComponents(coordSection, 0, 3);
603:   for (v = numFaces; v < numFaces+numVertices; ++v) {
604:     PetscSectionSetDof(coordSection, v, 3);
605:     PetscSectionSetFieldDof(coordSection, v, 0, 3);
606:   }
607:   PetscSectionSetUp(coordSection);
608:   PetscSectionGetStorageSize(coordSection, &coordSize);
609:   VecCreate(PETSC_COMM_SELF, &coordinates);
610:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
611:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
612:   VecSetBlockSize(coordinates, 3);
613:   VecSetType(coordinates,VECSTANDARD);
614:   VecGetArray(coordinates, &coords);
615:   for (vz = 0; vz <= faces[2]; ++vz) {
616:     for (vy = 0; vy <= faces[1]; ++vy) {
617:       for (vx = 0; vx <= faces[0]; ++vx) {
618:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
619:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
620:         coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
621:       }
622:     }
623:   }
624:   VecRestoreArray(coordinates, &coords);
625:   DMSetCoordinatesLocal(dm, coordinates);
626:   VecDestroy(&coordinates);
627:   return 0;
628: }

630: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], PetscBool interpolate)
631: {
633:   DMSetDimension(dm, dim-1);
634:   DMSetCoordinateDim(dm, dim);
635:   switch (dim) {
636:     case 1: DMPlexCreateBoxSurfaceMesh_Tensor_1D_Internal(dm, lower, upper, faces);break;
637:     case 2: DMPlexCreateBoxSurfaceMesh_Tensor_2D_Internal(dm, lower, upper, faces);break;
638:     case 3: DMPlexCreateBoxSurfaceMesh_Tensor_3D_Internal(dm, lower, upper, faces);break;
639:     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Dimension not supported: %D", dim);
640:   }
641:   if (interpolate) DMPlexInterpolateInPlace_Internal(dm);
642:   return 0;
643: }

645: /*@C
646:   DMPlexCreateBoxSurfaceMesh - Creates a mesh on the surface of the tensor product of unit intervals (box) using tensor cells (hexahedra).

648:   Collective

650:   Input Parameters:
651: + comm        - The communicator for the DM object
652: . dim         - The spatial dimension of the box, so the resulting mesh is has dimension dim-1
653: . faces       - Number of faces per dimension, or NULL for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
654: . lower       - The lower left corner, or NULL for (0, 0, 0)
655: . upper       - The upper right corner, or NULL for (1, 1, 1)
656: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

658:   Output Parameter:
659: . dm  - The DM object

661:   Level: beginner

663: .seealso: DMSetFromOptions(), DMPlexCreateBoxMesh(), DMPlexCreateFromFile(), DMSetType(), DMCreate()
664: @*/
665: PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], PetscBool interpolate, DM *dm)
666: {
667:   PetscInt       fac[3] = {1, 1, 1};
668:   PetscReal      low[3] = {0, 0, 0};
669:   PetscReal      upp[3] = {1, 1, 1};

671:   DMCreate(comm,dm);
672:   DMSetType(*dm,DMPLEX);
673:   DMPlexCreateBoxSurfaceMesh_Internal(*dm, dim, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, interpolate);
674:   return 0;
675: }

677: static PetscErrorCode DMPlexCreateLineMesh_Internal(DM dm,PetscInt segments,PetscReal lower,PetscReal upper,DMBoundaryType bd)
678: {
679:   PetscInt       i,fStart,fEnd,numCells = 0,numVerts = 0;
680:   PetscInt       numPoints[2],*coneSize,*cones,*coneOrientations;
681:   PetscScalar    *vertexCoords;
682:   PetscReal      L,maxCell;
683:   PetscBool      markerSeparate = PETSC_FALSE;
684:   PetscInt       markerLeft  = 1, faceMarkerLeft  = 1;
685:   PetscInt       markerRight = 1, faceMarkerRight = 2;
686:   PetscBool      wrap = (bd == DM_BOUNDARY_PERIODIC || bd == DM_BOUNDARY_TWIST) ? PETSC_TRUE : PETSC_FALSE;
687:   PetscMPIInt    rank;


691:   DMSetDimension(dm,1);
692:   DMCreateLabel(dm,"marker");
693:   DMCreateLabel(dm,"Face Sets");

695:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm),&rank);
696:   if (rank == 0) numCells = segments;
697:   if (rank == 0) numVerts = segments + (wrap ? 0 : 1);

699:   numPoints[0] = numVerts ; numPoints[1] = numCells;
700:   PetscMalloc4(numCells+numVerts,&coneSize,numCells*2,&cones,numCells+numVerts,&coneOrientations,numVerts,&vertexCoords);
701:   PetscArrayzero(coneOrientations,numCells+numVerts);
702:   for (i = 0; i < numCells; ++i) { coneSize[i] = 2; }
703:   for (i = 0; i < numVerts; ++i) { coneSize[numCells+i] = 0; }
704:   for (i = 0; i < numCells; ++i) { cones[2*i] = numCells + i%numVerts; cones[2*i+1] = numCells + (i+1)%numVerts; }
705:   for (i = 0; i < numVerts; ++i) { vertexCoords[i] = lower + (upper-lower)*((PetscReal)i/(PetscReal)numCells); }
706:   DMPlexCreateFromDAG(dm,1,numPoints,coneSize,cones,coneOrientations,vertexCoords);
707:   PetscFree4(coneSize,cones,coneOrientations,vertexCoords);

709:   PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dm_plex_separate_marker",&markerSeparate,NULL);
710:   if (markerSeparate) { markerLeft = faceMarkerLeft; markerRight = faceMarkerRight;}
711:   if (!wrap && rank == 0) {
712:     DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);
713:     DMSetLabelValue(dm,"marker",fStart,markerLeft);
714:     DMSetLabelValue(dm,"marker",fEnd-1,markerRight);
715:     DMSetLabelValue(dm,"Face Sets",fStart,faceMarkerLeft);
716:     DMSetLabelValue(dm,"Face Sets",fEnd-1,faceMarkerRight);
717:   }
718:   if (wrap) {
719:     L       = upper - lower;
720:     maxCell = (PetscReal)1.1*(L/(PetscReal)PetscMax(1,segments));
721:     DMSetPeriodicity(dm,PETSC_TRUE,&maxCell,&L,&bd);
722:   }
723:   DMPlexSetRefinementUniform(dm, PETSC_TRUE);
724:   return 0;
725: }

727: static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
728: {
729:   DM             boundary, vol;
730:   PetscInt       i;

734:   DMCreate(PetscObjectComm((PetscObject) dm), &boundary);
735:   DMSetType(boundary, DMPLEX);
736:   DMPlexCreateBoxSurfaceMesh_Internal(boundary, dim, faces, lower, upper, PETSC_FALSE);
737:   DMPlexGenerate(boundary, NULL, interpolate, &vol);
738:   DMPlexCopy_Internal(dm, PETSC_TRUE, vol);
739:   DMPlexReplace_Static(dm, &vol);
740:   DMDestroy(&boundary);
741:   return 0;
742: }

744: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
745: {
746:   DMLabel        cutLabel = NULL;
747:   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
748:   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
749:   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
750:   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
751:   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
752:   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
753:   PetscInt       dim;
754:   PetscBool      markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
755:   PetscMPIInt    rank;

757:   DMGetDimension(dm,&dim);
758:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
759:   DMCreateLabel(dm,"marker");
760:   DMCreateLabel(dm,"Face Sets");
761:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);
762:   if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ||
763:       bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ||
764:       bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {

766:     if (cutMarker) {DMCreateLabel(dm, "periodic_cut")); PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel);}
767:   }
768:   switch (dim) {
769:   case 2:
770:     faceMarkerTop    = 3;
771:     faceMarkerBottom = 1;
772:     faceMarkerRight  = 2;
773:     faceMarkerLeft   = 4;
774:     break;
775:   case 3:
776:     faceMarkerBottom = 1;
777:     faceMarkerTop    = 2;
778:     faceMarkerFront  = 3;
779:     faceMarkerBack   = 4;
780:     faceMarkerRight  = 5;
781:     faceMarkerLeft   = 6;
782:     break;
783:   default:
784:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D not supported",dim);
785:   }
786:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
787:   if (markerSeparate) {
788:     markerBottom = faceMarkerBottom;
789:     markerTop    = faceMarkerTop;
790:     markerFront  = faceMarkerFront;
791:     markerBack   = faceMarkerBack;
792:     markerRight  = faceMarkerRight;
793:     markerLeft   = faceMarkerLeft;
794:   }
795:   {
796:     const PetscInt numXEdges    = rank == 0 ? edges[0] : 0;
797:     const PetscInt numYEdges    = rank == 0 ? edges[1] : 0;
798:     const PetscInt numZEdges    = rank == 0 ? edges[2] : 0;
799:     const PetscInt numXVertices = rank == 0 ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
800:     const PetscInt numYVertices = rank == 0 ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
801:     const PetscInt numZVertices = rank == 0 ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
802:     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
803:     const PetscInt numXFaces    = numYEdges*numZEdges;
804:     const PetscInt numYFaces    = numXEdges*numZEdges;
805:     const PetscInt numZFaces    = numXEdges*numYEdges;
806:     const PetscInt numTotXFaces = numXVertices*numXFaces;
807:     const PetscInt numTotYFaces = numYVertices*numYFaces;
808:     const PetscInt numTotZFaces = numZVertices*numZFaces;
809:     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
810:     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
811:     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
812:     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
813:     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
814:     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
815:     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
816:     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
817:     const PetscInt firstYFace   = firstXFace + numTotXFaces;
818:     const PetscInt firstZFace   = firstYFace + numTotYFaces;
819:     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
820:     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
821:     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
822:     Vec            coordinates;
823:     PetscSection   coordSection;
824:     PetscScalar   *coords;
825:     PetscInt       coordSize;
826:     PetscInt       v, vx, vy, vz;
827:     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;

829:     DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);
830:     for (c = 0; c < numCells; c++) {
831:       DMPlexSetConeSize(dm, c, 6);
832:     }
833:     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
834:       DMPlexSetConeSize(dm, f, 4);
835:     }
836:     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
837:       DMPlexSetConeSize(dm, e, 2);
838:     }
839:     DMSetUp(dm); /* Allocate space for cones */
840:     /* Build cells */
841:     for (fz = 0; fz < numZEdges; ++fz) {
842:       for (fy = 0; fy < numYEdges; ++fy) {
843:         for (fx = 0; fx < numXEdges; ++fx) {
844:           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
845:           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
846:           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
847:           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
848:           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
849:           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
850:           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
851:                             /* B,  T,  F,  K,  R,  L */
852:           PetscInt ornt[6] = {-2,  0,  0, -3,  0, -2}; /* ??? */
853:           PetscInt cone[6];

855:           /* no boundary twisting in 3D */
856:           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
857:           DMPlexSetCone(dm, cell, cone);
858:           DMPlexSetConeOrientation(dm, cell, ornt);
859:           if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) DMLabelSetValue(cutLabel, cell, 2);
860:           if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) DMLabelSetValue(cutLabel, cell, 2);
861:           if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) DMLabelSetValue(cutLabel, cell, 2);
862:         }
863:       }
864:     }
865:     /* Build x faces */
866:     for (fz = 0; fz < numZEdges; ++fz) {
867:       for (fy = 0; fy < numYEdges; ++fy) {
868:         for (fx = 0; fx < numXVertices; ++fx) {
869:           PetscInt face    = firstXFace + (fz*numYEdges+fy)     *numXVertices+fx;
870:           PetscInt edgeL   = firstZEdge + (fy                   *numXVertices+fx)*numZEdges + fz;
871:           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
872:           PetscInt edgeB   = firstYEdge + (fz                   *numXVertices+fx)*numYEdges + fy;
873:           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
874:           PetscInt ornt[4] = {0, 0, -1, -1};
875:           PetscInt cone[4];

877:           if (dim == 3) {
878:             /* markers */
879:             if (bdX != DM_BOUNDARY_PERIODIC) {
880:               if (fx == numXVertices-1) {
881:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);
882:                 DMSetLabelValue(dm, "marker", face, markerRight);
883:               }
884:               else if (fx == 0) {
885:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);
886:                 DMSetLabelValue(dm, "marker", face, markerLeft);
887:               }
888:             }
889:           }
890:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
891:           DMPlexSetCone(dm, face, cone);
892:           DMPlexSetConeOrientation(dm, face, ornt);
893:         }
894:       }
895:     }
896:     /* Build y faces */
897:     for (fz = 0; fz < numZEdges; ++fz) {
898:       for (fx = 0; fx < numXEdges; ++fx) {
899:         for (fy = 0; fy < numYVertices; ++fy) {
900:           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
901:           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx)*numZEdges + fz;
902:           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
903:           PetscInt edgeB   = firstXEdge + (fz                   *numYVertices+fy)*numXEdges + fx;
904:           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
905:           PetscInt ornt[4] = {0, 0, -1, -1};
906:           PetscInt cone[4];

908:           if (dim == 3) {
909:             /* markers */
910:             if (bdY != DM_BOUNDARY_PERIODIC) {
911:               if (fy == numYVertices-1) {
912:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);
913:                 DMSetLabelValue(dm, "marker", face, markerBack);
914:               }
915:               else if (fy == 0) {
916:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);
917:                 DMSetLabelValue(dm, "marker", face, markerFront);
918:               }
919:             }
920:           }
921:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
922:           DMPlexSetCone(dm, face, cone);
923:           DMPlexSetConeOrientation(dm, face, ornt);
924:         }
925:       }
926:     }
927:     /* Build z faces */
928:     for (fy = 0; fy < numYEdges; ++fy) {
929:       for (fx = 0; fx < numXEdges; ++fx) {
930:         for (fz = 0; fz < numZVertices; fz++) {
931:           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
932:           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx)*numYEdges + fy;
933:           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
934:           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy)*numXEdges + fx;
935:           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
936:           PetscInt ornt[4] = {0, 0, -1, -1};
937:           PetscInt cone[4];

939:           if (dim == 2) {
940:             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -1;}
941:             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
942:             if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) DMLabelSetValue(cutLabel, face, 2);
943:             if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) DMLabelSetValue(cutLabel, face, 2);
944:           } else {
945:             /* markers */
946:             if (bdZ != DM_BOUNDARY_PERIODIC) {
947:               if (fz == numZVertices-1) {
948:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);
949:                 DMSetLabelValue(dm, "marker", face, markerTop);
950:               }
951:               else if (fz == 0) {
952:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);
953:                 DMSetLabelValue(dm, "marker", face, markerBottom);
954:               }
955:             }
956:           }
957:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
958:           DMPlexSetCone(dm, face, cone);
959:           DMPlexSetConeOrientation(dm, face, ornt);
960:         }
961:       }
962:     }
963:     /* Build Z edges*/
964:     for (vy = 0; vy < numYVertices; vy++) {
965:       for (vx = 0; vx < numXVertices; vx++) {
966:         for (ez = 0; ez < numZEdges; ez++) {
967:           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
968:           const PetscInt vertexB = firstVertex + (ez                   *numYVertices+vy)*numXVertices + vx;
969:           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
970:           PetscInt       cone[2];

972:           if (dim == 3) {
973:             if (bdX != DM_BOUNDARY_PERIODIC) {
974:               if (vx == numXVertices-1) {
975:                 DMSetLabelValue(dm, "marker", edge, markerRight);
976:               }
977:               else if (vx == 0) {
978:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
979:               }
980:             }
981:             if (bdY != DM_BOUNDARY_PERIODIC) {
982:               if (vy == numYVertices-1) {
983:                 DMSetLabelValue(dm, "marker", edge, markerBack);
984:               }
985:               else if (vy == 0) {
986:                 DMSetLabelValue(dm, "marker", edge, markerFront);
987:               }
988:             }
989:           }
990:           cone[0] = vertexB; cone[1] = vertexT;
991:           DMPlexSetCone(dm, edge, cone);
992:         }
993:       }
994:     }
995:     /* Build Y edges*/
996:     for (vz = 0; vz < numZVertices; vz++) {
997:       for (vx = 0; vx < numXVertices; vx++) {
998:         for (ey = 0; ey < numYEdges; ey++) {
999:           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
1000:           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
1001:           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
1002:           const PetscInt vertexK = firstVertex + nextv;
1003:           PetscInt       cone[2];

1005:           cone[0] = vertexF; cone[1] = vertexK;
1006:           DMPlexSetCone(dm, edge, cone);
1007:           if (dim == 2) {
1008:             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
1009:               if (vx == numXVertices-1) {
1010:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);
1011:                 DMSetLabelValue(dm, "marker", edge,    markerRight);
1012:                 DMSetLabelValue(dm, "marker", cone[0], markerRight);
1013:                 if (ey == numYEdges-1) {
1014:                   DMSetLabelValue(dm, "marker", cone[1], markerRight);
1015:                 }
1016:               } else if (vx == 0) {
1017:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);
1018:                 DMSetLabelValue(dm, "marker", edge,    markerLeft);
1019:                 DMSetLabelValue(dm, "marker", cone[0], markerLeft);
1020:                 if (ey == numYEdges-1) {
1021:                   DMSetLabelValue(dm, "marker", cone[1], markerLeft);
1022:                 }
1023:               }
1024:             } else {
1025:               if (vx == 0 && cutLabel) {
1026:                 DMLabelSetValue(cutLabel, edge,    1);
1027:                 DMLabelSetValue(cutLabel, cone[0], 1);
1028:                 if (ey == numYEdges-1) {
1029:                   DMLabelSetValue(cutLabel, cone[1], 1);
1030:                 }
1031:               }
1032:             }
1033:           } else {
1034:             if (bdX != DM_BOUNDARY_PERIODIC) {
1035:               if (vx == numXVertices-1) {
1036:                 DMSetLabelValue(dm, "marker", edge, markerRight);
1037:               } else if (vx == 0) {
1038:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
1039:               }
1040:             }
1041:             if (bdZ != DM_BOUNDARY_PERIODIC) {
1042:               if (vz == numZVertices-1) {
1043:                 DMSetLabelValue(dm, "marker", edge, markerTop);
1044:               } else if (vz == 0) {
1045:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
1046:               }
1047:             }
1048:           }
1049:         }
1050:       }
1051:     }
1052:     /* Build X edges*/
1053:     for (vz = 0; vz < numZVertices; vz++) {
1054:       for (vy = 0; vy < numYVertices; vy++) {
1055:         for (ex = 0; ex < numXEdges; ex++) {
1056:           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
1057:           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
1058:           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
1059:           const PetscInt vertexR = firstVertex + nextv;
1060:           PetscInt       cone[2];

1062:           cone[0] = vertexL; cone[1] = vertexR;
1063:           DMPlexSetCone(dm, edge, cone);
1064:           if (dim == 2) {
1065:             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
1066:               if (vy == numYVertices-1) {
1067:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);
1068:                 DMSetLabelValue(dm, "marker", edge,    markerTop);
1069:                 DMSetLabelValue(dm, "marker", cone[0], markerTop);
1070:                 if (ex == numXEdges-1) {
1071:                   DMSetLabelValue(dm, "marker", cone[1], markerTop);
1072:                 }
1073:               } else if (vy == 0) {
1074:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);
1075:                 DMSetLabelValue(dm, "marker", edge,    markerBottom);
1076:                 DMSetLabelValue(dm, "marker", cone[0], markerBottom);
1077:                 if (ex == numXEdges-1) {
1078:                   DMSetLabelValue(dm, "marker", cone[1], markerBottom);
1079:                 }
1080:               }
1081:             } else {
1082:               if (vy == 0 && cutLabel) {
1083:                 DMLabelSetValue(cutLabel, edge,    1);
1084:                 DMLabelSetValue(cutLabel, cone[0], 1);
1085:                 if (ex == numXEdges-1) {
1086:                   DMLabelSetValue(cutLabel, cone[1], 1);
1087:                 }
1088:               }
1089:             }
1090:           } else {
1091:             if (bdY != DM_BOUNDARY_PERIODIC) {
1092:               if (vy == numYVertices-1) {
1093:                 DMSetLabelValue(dm, "marker", edge, markerBack);
1094:               }
1095:               else if (vy == 0) {
1096:                 DMSetLabelValue(dm, "marker", edge, markerFront);
1097:               }
1098:             }
1099:             if (bdZ != DM_BOUNDARY_PERIODIC) {
1100:               if (vz == numZVertices-1) {
1101:                 DMSetLabelValue(dm, "marker", edge, markerTop);
1102:               }
1103:               else if (vz == 0) {
1104:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
1105:               }
1106:             }
1107:           }
1108:         }
1109:       }
1110:     }
1111:     DMPlexSymmetrize(dm);
1112:     DMPlexStratify(dm);
1113:     /* Build coordinates */
1114:     DMGetCoordinateSection(dm, &coordSection);
1115:     PetscSectionSetNumFields(coordSection, 1);
1116:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1117:     PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);
1118:     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
1119:       PetscSectionSetDof(coordSection, v, dim);
1120:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1121:     }
1122:     PetscSectionSetUp(coordSection);
1123:     PetscSectionGetStorageSize(coordSection, &coordSize);
1124:     VecCreate(PETSC_COMM_SELF, &coordinates);
1125:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1126:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1127:     VecSetBlockSize(coordinates, dim);
1128:     VecSetType(coordinates,VECSTANDARD);
1129:     VecGetArray(coordinates, &coords);
1130:     for (vz = 0; vz < numZVertices; ++vz) {
1131:       for (vy = 0; vy < numYVertices; ++vy) {
1132:         for (vx = 0; vx < numXVertices; ++vx) {
1133:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
1134:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
1135:           if (dim == 3) {
1136:             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
1137:           }
1138:         }
1139:       }
1140:     }
1141:     VecRestoreArray(coordinates, &coords);
1142:     DMSetCoordinatesLocal(dm, coordinates);
1143:     VecDestroy(&coordinates);
1144:   }
1145:   return 0;
1146: }

1148: static PetscErrorCode DMPlexCreateBoxMesh_Tensor_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[])
1149: {
1150:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1151:   PetscInt       fac[3] = {0, 0, 0}, d;

1155:   DMSetDimension(dm, dim);
1156:   for (d = 0; d < dim; ++d) {fac[d] = faces[d]; bdt[d] = periodicity[d];}
1157:   DMPlexCreateCubeMesh_Internal(dm, lower, upper, fac, bdt[0], bdt[1], bdt[2]);
1158:   if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST ||
1159:       periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST ||
1160:       (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) {
1161:     PetscReal L[3];
1162:     PetscReal maxCell[3];

1164:     for (d = 0; d < dim; ++d) {
1165:       L[d]       = upper[d] - lower[d];
1166:       maxCell[d] = 1.1 * (L[d] / PetscMax(1, faces[d]));
1167:     }
1168:     DMSetPeriodicity(dm, PETSC_TRUE, maxCell, L, periodicity);
1169:   }
1170:   DMPlexSetRefinementUniform(dm, PETSC_TRUE);
1171:   return 0;
1172: }

1174: static PetscErrorCode DMPlexCreateBoxMesh_Internal(DM dm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
1175: {
1176:   if (dim == 1)      DMPlexCreateLineMesh_Internal(dm, faces[0], lower[0], upper[0], periodicity[0]);
1177:   else if (simplex)  DMPlexCreateBoxMesh_Simplex_Internal(dm, dim, faces, lower, upper, periodicity, interpolate);
1178:   else               DMPlexCreateBoxMesh_Tensor_Internal(dm, dim, faces, lower, upper, periodicity);
1179:   if (!interpolate && dim > 1 && !simplex) {
1180:     DM udm;

1182:     DMPlexUninterpolate(dm, &udm);
1183:     DMPlexCopyCoordinates(dm, udm);
1184:     DMPlexReplace_Static(dm, &udm);
1185:   }
1186:   return 0;
1187: }

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

1192:   Collective

1194:   Input Parameters:
1195: + comm        - The communicator for the DM object
1196: . dim         - The spatial dimension
1197: . simplex     - PETSC_TRUE for simplices, PETSC_FALSE for tensor cells
1198: . faces       - Number of faces per dimension, or NULL for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
1199: . lower       - The lower left corner, or NULL for (0, 0, 0)
1200: . upper       - The upper right corner, or NULL for (1, 1, 1)
1201: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
1202: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1204:   Output Parameter:
1205: . dm  - The DM object

1207:   Note: If you want to customize this mesh using options, you just need to
1208: $  DMCreate(comm, &dm);
1209: $  DMSetType(dm, DMPLEX);
1210: $  DMSetFromOptions(dm);
1211: and use the options on the DMSetFromOptions() page.

1213:   Here is the numbering returned for 2 faces in each direction for tensor cells:
1214: $ 10---17---11---18----12
1215: $  |         |         |
1216: $  |         |         |
1217: $ 20    2   22    3    24
1218: $  |         |         |
1219: $  |         |         |
1220: $  7---15----8---16----9
1221: $  |         |         |
1222: $  |         |         |
1223: $ 19    0   21    1   23
1224: $  |         |         |
1225: $  |         |         |
1226: $  4---13----5---14----6

1228: and for simplicial cells

1230: $ 14----8---15----9----16
1231: $  |\     5  |\      7 |
1232: $  | \       | \       |
1233: $ 13   2    14    3    15
1234: $  | 4   \   | 6   \   |
1235: $  |       \ |       \ |
1236: $ 11----6---12----7----13
1237: $  |\        |\        |
1238: $  | \    1  | \     3 |
1239: $ 10   0    11    1    12
1240: $  | 0   \   | 2   \   |
1241: $  |       \ |       \ |
1242: $  8----4----9----5----10

1244:   Level: beginner

1246: .seealso: DMSetFromOptions(), DMPlexCreateFromFile(), DMPlexCreateHexCylinderMesh(), DMSetType(), DMCreate()
1247: @*/
1248: 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)
1249: {
1250:   PetscInt       fac[3] = {1, 1, 1};
1251:   PetscReal      low[3] = {0, 0, 0};
1252:   PetscReal      upp[3] = {1, 1, 1};
1253:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1255:   DMCreate(comm,dm);
1256:   DMSetType(*dm,DMPLEX);
1257:   DMPlexCreateBoxMesh_Internal(*dm, dim, simplex, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, periodicity ? periodicity : bdt, interpolate);
1258:   return 0;
1259: }

1261: static PetscErrorCode DMPlexCreateWedgeBoxMesh_Internal(DM dm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[])
1262: {
1263:   DM             bdm, vol;
1264:   PetscInt       i;

1267:   DMCreate(PetscObjectComm((PetscObject) dm), &bdm);
1268:   DMSetType(bdm, DMPLEX);
1269:   DMSetDimension(bdm, 2);
1270:   DMPlexCreateBoxMesh_Simplex_Internal(bdm, 2, faces, lower, upper, periodicity, PETSC_TRUE);
1271:   DMPlexExtrude(bdm, faces[2], upper[2] - lower[2], PETSC_TRUE, PETSC_FALSE, NULL, NULL, &vol);
1272:   DMDestroy(&bdm);
1273:   DMPlexReplace_Static(dm, &vol);
1274:   if (lower[2] != 0.0) {
1275:     Vec          v;
1276:     PetscScalar *x;
1277:     PetscInt     cDim, n;

1279:     DMGetCoordinatesLocal(dm, &v);
1280:     VecGetBlockSize(v, &cDim);
1281:     VecGetLocalSize(v, &n);
1282:     VecGetArray(v, &x);
1283:     x   += cDim;
1284:     for (i = 0; i < n; i += cDim) x[i] += lower[2];
1285:     VecRestoreArray(v,&x);
1286:     DMSetCoordinatesLocal(dm, v);
1287:   }
1288:   return 0;
1289: }

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

1294:   Collective

1296:   Input Parameters:
1297: + comm        - The communicator for the DM object
1298: . faces       - Number of faces per dimension, or NULL for (1, 1, 1)
1299: . lower       - The lower left corner, or NULL for (0, 0, 0)
1300: . upper       - The upper right corner, or NULL for (1, 1, 1)
1301: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
1302: . orderHeight - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1303: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1305:   Output Parameter:
1306: . dm  - The DM object

1308:   Level: beginner

1310: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateWedgeCylinderMesh(), DMExtrude(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1311: @*/
1312: PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool orderHeight, PetscBool interpolate, DM *dm)
1313: {
1314:   PetscInt       fac[3] = {1, 1, 1};
1315:   PetscReal      low[3] = {0, 0, 0};
1316:   PetscReal      upp[3] = {1, 1, 1};
1317:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1319:   DMCreate(comm,dm);
1320:   DMSetType(*dm,DMPLEX);
1321:   DMPlexCreateWedgeBoxMesh_Internal(*dm, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, periodicity ? periodicity : bdt);
1322:   if (!interpolate) {
1323:     DM udm;

1325:     DMPlexUninterpolate(*dm, &udm);
1326:     DMPlexReplace_Static(*dm, &udm);
1327:   }
1328:   return 0;
1329: }

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

1334:   Logically Collective on dm

1336:   Input Parameters:
1337: + dm - the DM context
1338: - prefix - the prefix to prepend to all option names

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

1344:   Level: advanced

1346: .seealso: SNESSetFromOptions()
1347: @*/
1348: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1349: {
1350:   DM_Plex       *mesh = (DM_Plex *) dm->data;

1353:   PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);
1354:   PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);
1355:   return 0;
1356: }

1358: /* Remap geometry to cylinder
1359:    TODO: This only works for a single refinement, then it is broken

1361:      Interior square: Linear interpolation is correct
1362:      The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1363:      such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance

1365:        phi     = arctan(y/x)
1366:        d_close = sqrt(1/8 + 1/4 sin^2(phi))
1367:        d_far   = sqrt(1/2 + sin^2(phi))

1369:      so we remap them using

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

1374:      If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1375: */
1376: static void snapToCylinder(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1377:                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1378:                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1379:                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1380: {
1381:   const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1382:   const PetscReal ds2 = 0.5*dis;

1384:   if ((PetscAbsScalar(u[0]) <= ds2) && (PetscAbsScalar(u[1]) <= ds2)) {
1385:     f0[0] = u[0];
1386:     f0[1] = u[1];
1387:   } else {
1388:     PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;

1390:     x    = PetscRealPart(u[0]);
1391:     y    = PetscRealPart(u[1]);
1392:     phi  = PetscAtan2Real(y, x);
1393:     sinp = PetscSinReal(phi);
1394:     cosp = PetscCosReal(phi);
1395:     if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
1396:       dc = PetscAbsReal(ds2/sinp);
1397:       df = PetscAbsReal(dis/sinp);
1398:       xc = ds2*x/PetscAbsReal(y);
1399:       yc = ds2*PetscSignReal(y);
1400:     } else {
1401:       dc = PetscAbsReal(ds2/cosp);
1402:       df = PetscAbsReal(dis/cosp);
1403:       xc = ds2*PetscSignReal(x);
1404:       yc = ds2*y/PetscAbsReal(x);
1405:     }
1406:     f0[0] = xc + (u[0] - xc)*(1.0 - dc)/(df - dc);
1407:     f0[1] = yc + (u[1] - yc)*(1.0 - dc)/(df - dc);
1408:   }
1409:   f0[2] = u[2];
1410: }

1412: static PetscErrorCode DMPlexCreateHexCylinderMesh_Internal(DM dm, DMBoundaryType periodicZ)
1413: {
1414:   const PetscInt dim = 3;
1415:   PetscInt       numCells, numVertices;
1416:   PetscMPIInt    rank;

1418:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1419:   DMSetDimension(dm, dim);
1420:   /* Create topology */
1421:   {
1422:     PetscInt cone[8], c;

1424:     numCells    = rank == 0 ?  5 : 0;
1425:     numVertices = rank == 0 ? 16 : 0;
1426:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1427:       numCells   *= 3;
1428:       numVertices = rank == 0 ? 24 : 0;
1429:     }
1430:     DMPlexSetChart(dm, 0, numCells+numVertices);
1431:     for (c = 0; c < numCells; c++) DMPlexSetConeSize(dm, c, 8);
1432:     DMSetUp(dm);
1433:     if (rank == 0) {
1434:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1435:         cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1436:         cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1437:         DMPlexSetCone(dm, 0, cone);
1438:         cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1439:         cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1440:         DMPlexSetCone(dm, 1, cone);
1441:         cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1442:         cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1443:         DMPlexSetCone(dm, 2, cone);
1444:         cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1445:         cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1446:         DMPlexSetCone(dm, 3, cone);
1447:         cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1448:         cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1449:         DMPlexSetCone(dm, 4, cone);

1451:         cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1452:         cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1453:         DMPlexSetCone(dm, 5, cone);
1454:         cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1455:         cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1456:         DMPlexSetCone(dm, 6, cone);
1457:         cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1458:         cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1459:         DMPlexSetCone(dm, 7, cone);
1460:         cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1461:         cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1462:         DMPlexSetCone(dm, 8, cone);
1463:         cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1464:         cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1465:         DMPlexSetCone(dm, 9, cone);

1467:         cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1468:         cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1469:         DMPlexSetCone(dm, 10, cone);
1470:         cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1471:         cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1472:         DMPlexSetCone(dm, 11, cone);
1473:         cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1474:         cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1475:         DMPlexSetCone(dm, 12, cone);
1476:         cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1477:         cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1478:         DMPlexSetCone(dm, 13, cone);
1479:         cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1480:         cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1481:         DMPlexSetCone(dm, 14, cone);
1482:       } else {
1483:         cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
1484:         cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
1485:         DMPlexSetCone(dm, 0, cone);
1486:         cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
1487:         cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
1488:         DMPlexSetCone(dm, 1, cone);
1489:         cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
1490:         cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
1491:         DMPlexSetCone(dm, 2, cone);
1492:         cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
1493:         cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
1494:         DMPlexSetCone(dm, 3, cone);
1495:         cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
1496:         cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
1497:         DMPlexSetCone(dm, 4, cone);
1498:       }
1499:     }
1500:     DMPlexSymmetrize(dm);
1501:     DMPlexStratify(dm);
1502:   }
1503:   /* Create cube geometry */
1504:   {
1505:     Vec             coordinates;
1506:     PetscSection    coordSection;
1507:     PetscScalar    *coords;
1508:     PetscInt        coordSize, v;
1509:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1510:     const PetscReal ds2 = dis/2.0;

1512:     /* Build coordinates */
1513:     DMGetCoordinateSection(dm, &coordSection);
1514:     PetscSectionSetNumFields(coordSection, 1);
1515:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1516:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1517:     for (v = numCells; v < numCells+numVertices; ++v) {
1518:       PetscSectionSetDof(coordSection, v, dim);
1519:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1520:     }
1521:     PetscSectionSetUp(coordSection);
1522:     PetscSectionGetStorageSize(coordSection, &coordSize);
1523:     VecCreate(PETSC_COMM_SELF, &coordinates);
1524:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1525:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1526:     VecSetBlockSize(coordinates, dim);
1527:     VecSetType(coordinates,VECSTANDARD);
1528:     VecGetArray(coordinates, &coords);
1529:     if (rank == 0) {
1530:       coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
1531:       coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
1532:       coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
1533:       coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
1534:       coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
1535:       coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
1536:       coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
1537:       coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
1538:       coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
1539:       coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
1540:       coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
1541:       coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
1542:       coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
1543:       coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
1544:       coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
1545:       coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1546:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1547:         /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1548:         /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1549:         /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1550:         /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1551:         /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1552:         /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1553:         /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1554:         /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1555:       }
1556:     }
1557:     VecRestoreArray(coordinates, &coords);
1558:     DMSetCoordinatesLocal(dm, coordinates);
1559:     VecDestroy(&coordinates);
1560:   }
1561:   /* Create periodicity */
1562:   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1563:     PetscReal      L[3];
1564:     PetscReal      maxCell[3];
1565:     DMBoundaryType bdType[3];
1566:     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1567:     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1568:     PetscInt       i, numZCells = 3;

1570:     bdType[0] = DM_BOUNDARY_NONE;
1571:     bdType[1] = DM_BOUNDARY_NONE;
1572:     bdType[2] = periodicZ;
1573:     for (i = 0; i < dim; i++) {
1574:       L[i]       = upper[i] - lower[i];
1575:       maxCell[i] = 1.1 * (L[i] / numZCells);
1576:     }
1577:     DMSetPeriodicity(dm, PETSC_TRUE, maxCell, L, bdType);
1578:   }
1579:   {
1580:     DM          cdm;
1581:     PetscDS     cds;
1582:     PetscScalar c[2] = {1.0, 1.0};

1584:     DMPlexCreateCoordinateSpace(dm, 1, snapToCylinder);
1585:     DMGetCoordinateDM(dm, &cdm);
1586:     DMGetDS(cdm, &cds);
1587:     PetscDSSetConstants(cds, 2, c);
1588:   }
1589:   /* Wait for coordinate creation before doing in-place modification */
1590:   DMPlexInterpolateInPlace_Internal(dm);
1591:   return 0;
1592: }

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

1597:   Collective

1599:   Input Parameters:
1600: + comm      - The communicator for the DM object
1601: - periodicZ - The boundary type for the Z direction

1603:   Output Parameter:
1604: . dm  - The DM object

1606:   Note:
1607:   Here is the output numbering looking from the bottom of the cylinder:
1608: $       17-----14
1609: $        |     |
1610: $        |  2  |
1611: $        |     |
1612: $ 17-----8-----7-----14
1613: $  |     |     |     |
1614: $  |  3  |  0  |  1  |
1615: $  |     |     |     |
1616: $ 19-----5-----6-----13
1617: $        |     |
1618: $        |  4  |
1619: $        |     |
1620: $       19-----13
1621: $
1622: $ and up through the top
1623: $
1624: $       18-----16
1625: $        |     |
1626: $        |  2  |
1627: $        |     |
1628: $ 18----10----11-----16
1629: $  |     |     |     |
1630: $  |  3  |  0  |  1  |
1631: $  |     |     |     |
1632: $ 20-----9----12-----15
1633: $        |     |
1634: $        |  4  |
1635: $        |     |
1636: $       20-----15

1638:   Level: beginner

1640: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1641: @*/
1642: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, DMBoundaryType periodicZ, DM *dm)
1643: {
1645:   DMCreate(comm, dm);
1646:   DMSetType(*dm, DMPLEX);
1647:   DMPlexCreateHexCylinderMesh_Internal(*dm, periodicZ);
1648:   return 0;
1649: }

1651: static PetscErrorCode DMPlexCreateWedgeCylinderMesh_Internal(DM dm, PetscInt n, PetscBool interpolate)
1652: {
1653:   const PetscInt dim = 3;
1654:   PetscInt       numCells, numVertices, v;
1655:   PetscMPIInt    rank;

1658:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1659:   DMSetDimension(dm, dim);
1660:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1661:   DMCreateLabel(dm, "celltype");
1662:   /* Create topology */
1663:   {
1664:     PetscInt cone[6], c;

1666:     numCells    = rank == 0 ?        n : 0;
1667:     numVertices = rank == 0 ?  2*(n+1) : 0;
1668:     DMPlexSetChart(dm, 0, numCells+numVertices);
1669:     for (c = 0; c < numCells; c++) DMPlexSetConeSize(dm, c, 6);
1670:     DMSetUp(dm);
1671:     for (c = 0; c < numCells; c++) {
1672:       cone[0] =  c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n;
1673:       cone[3] =  c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n;
1674:       DMPlexSetCone(dm, c, cone);
1675:       DMPlexSetCellType(dm, c, DM_POLYTOPE_TRI_PRISM_TENSOR);
1676:     }
1677:     DMPlexSymmetrize(dm);
1678:     DMPlexStratify(dm);
1679:   }
1680:   for (v = numCells; v < numCells+numVertices; ++v) {
1681:     DMPlexSetCellType(dm, v, DM_POLYTOPE_POINT);
1682:   }
1683:   /* Create cylinder geometry */
1684:   {
1685:     Vec          coordinates;
1686:     PetscSection coordSection;
1687:     PetscScalar *coords;
1688:     PetscInt     coordSize, c;

1690:     /* Build coordinates */
1691:     DMGetCoordinateSection(dm, &coordSection);
1692:     PetscSectionSetNumFields(coordSection, 1);
1693:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1694:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1695:     for (v = numCells; v < numCells+numVertices; ++v) {
1696:       PetscSectionSetDof(coordSection, v, dim);
1697:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1698:     }
1699:     PetscSectionSetUp(coordSection);
1700:     PetscSectionGetStorageSize(coordSection, &coordSize);
1701:     VecCreate(PETSC_COMM_SELF, &coordinates);
1702:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1703:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1704:     VecSetBlockSize(coordinates, dim);
1705:     VecSetType(coordinates,VECSTANDARD);
1706:     VecGetArray(coordinates, &coords);
1707:     for (c = 0; c < numCells; c++) {
1708:       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;
1709:       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;
1710:     }
1711:     if (rank == 0) {
1712:       coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0;
1713:       coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0;
1714:     }
1715:     VecRestoreArray(coordinates, &coords);
1716:     DMSetCoordinatesLocal(dm, coordinates);
1717:     VecDestroy(&coordinates);
1718:   }
1719:   /* Interpolate */
1720:   if (interpolate) DMPlexInterpolateInPlace_Internal(dm);
1721:   return 0;
1722: }

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

1727:   Collective

1729:   Input Parameters:
1730: + comm - The communicator for the DM object
1731: . n    - The number of wedges around the origin
1732: - interpolate - Create edges and faces

1734:   Output Parameter:
1735: . dm  - The DM object

1737:   Level: beginner

1739: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1740: @*/
1741: PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
1742: {
1744:   DMCreate(comm, dm);
1745:   DMSetType(*dm, DMPLEX);
1746:   DMPlexCreateWedgeCylinderMesh_Internal(*dm, n, interpolate);
1747:   return 0;
1748: }

1750: static inline PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1751: {
1752:   PetscReal prod = 0.0;
1753:   PetscInt  i;
1754:   for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
1755:   return PetscSqrtReal(prod);
1756: }
1757: static inline PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1758: {
1759:   PetscReal prod = 0.0;
1760:   PetscInt  i;
1761:   for (i = 0; i < dim; ++i) prod += x[i]*y[i];
1762:   return prod;
1763: }

1765: /* The first constant is the sphere radius */
1766: static void snapToSphere(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1767:                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1768:                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1769:                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1770: {
1771:   PetscReal r = PetscRealPart(constants[0]);
1772:   PetscReal norm2 = 0.0, fac;
1773:   PetscInt  n = uOff[1] - uOff[0], d;

1775:   for (d = 0; d < n; ++d) norm2 += PetscSqr(PetscRealPart(u[d]));
1776:   fac = r/PetscSqrtReal(norm2);
1777:   for (d = 0; d < n; ++d) f0[d] = u[d]*fac;
1778: }

1780: static PetscErrorCode DMPlexCreateSphereMesh_Internal(DM dm, PetscInt dim, PetscBool simplex, PetscReal R)
1781: {
1782:   const PetscInt  embedDim = dim+1;
1783:   PetscSection    coordSection;
1784:   Vec             coordinates;
1785:   PetscScalar    *coords;
1786:   PetscReal      *coordsIn;
1787:   PetscInt        numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e;
1788:   PetscMPIInt     rank;

1791:   DMSetDimension(dm, dim);
1792:   DMSetCoordinateDim(dm, dim+1);
1793:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1794:   switch (dim) {
1795:   case 2:
1796:     if (simplex) {
1797:       const PetscReal radius    = PetscSqrtReal(1 + PETSC_PHI*PETSC_PHI)/(1.0 + PETSC_PHI);
1798:       const PetscReal edgeLen   = 2.0/(1.0 + PETSC_PHI) * (R/radius);
1799:       const PetscInt  degree    = 5;
1800:       PetscReal       vertex[3] = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)};
1801:       PetscInt        s[3]      = {1, 1, 1};
1802:       PetscInt        cone[3];
1803:       PetscInt       *graph, p, i, j, k;

1805:       vertex[0] *= R/radius; vertex[1] *= R/radius; vertex[2] *= R/radius;
1806:       numCells    = rank == 0 ? 20 : 0;
1807:       numVerts    = rank == 0 ? 12 : 0;
1808:       firstVertex = numCells;
1809:       /* Use icosahedron, which for a R-sphere has coordinates which are all cyclic permutations of

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

1813:          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
1814:          length is then given by 2/(1+\phi) = 2 * 0.38197 = 0.76393.
1815:       */
1816:       /* Construct vertices */
1817:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1818:       if (rank == 0) {
1819:         for (p = 0, i = 0; p < embedDim; ++p) {
1820:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
1821:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
1822:               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim];
1823:               ++i;
1824:             }
1825:           }
1826:         }
1827:       }
1828:       /* Construct graph */
1829:       PetscCalloc1(numVerts * numVerts, &graph);
1830:       for (i = 0; i < numVerts; ++i) {
1831:         for (j = 0, k = 0; j < numVerts; ++j) {
1832:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
1833:         }
1835:       }
1836:       /* Build Topology */
1837:       DMPlexSetChart(dm, 0, numCells+numVerts);
1838:       for (c = 0; c < numCells; c++) {
1839:         DMPlexSetConeSize(dm, c, embedDim);
1840:       }
1841:       DMSetUp(dm); /* Allocate space for cones */
1842:       /* Cells */
1843:       for (i = 0, c = 0; i < numVerts; ++i) {
1844:         for (j = 0; j < i; ++j) {
1845:           for (k = 0; k < j; ++k) {
1846:             if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) {
1847:               cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k;
1848:               /* Check orientation */
1849:               {
1850:                 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}}};
1851:                 PetscReal normal[3];
1852:                 PetscInt  e, f;

1854:                 for (d = 0; d < embedDim; ++d) {
1855:                   normal[d] = 0.0;
1856:                   for (e = 0; e < embedDim; ++e) {
1857:                     for (f = 0; f < embedDim; ++f) {
1858:                       normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]);
1859:                     }
1860:                   }
1861:                 }
1862:                 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
1863:               }
1864:               DMPlexSetCone(dm, c++, cone);
1865:             }
1866:           }
1867:         }
1868:       }
1869:       DMPlexSymmetrize(dm);
1870:       DMPlexStratify(dm);
1871:       PetscFree(graph);
1872:     } else {
1873:       /*
1874:         12-21--13
1875:          |     |
1876:         25  4  24
1877:          |     |
1878:   12-25--9-16--8-24--13
1879:    |     |     |     |
1880:   23  5 17  0 15  3  22
1881:    |     |     |     |
1882:   10-20--6-14--7-19--11
1883:          |     |
1884:         20  1  19
1885:          |     |
1886:         10-18--11
1887:          |     |
1888:         23  2  22
1889:          |     |
1890:         12-21--13
1891:        */
1892:       PetscInt cone[4], ornt[4];

1894:       numCells    = rank == 0 ?  6 : 0;
1895:       numEdges    = rank == 0 ? 12 : 0;
1896:       numVerts    = rank == 0 ?  8 : 0;
1897:       firstVertex = numCells;
1898:       firstEdge   = numCells + numVerts;
1899:       /* Build Topology */
1900:       DMPlexSetChart(dm, 0, numCells+numEdges+numVerts);
1901:       for (c = 0; c < numCells; c++) {
1902:         DMPlexSetConeSize(dm, c, 4);
1903:       }
1904:       for (e = firstEdge; e < firstEdge+numEdges; ++e) {
1905:         DMPlexSetConeSize(dm, e, 2);
1906:       }
1907:       DMSetUp(dm); /* Allocate space for cones */
1908:       if (rank == 0) {
1909:         /* Cell 0 */
1910:         cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17;
1911:         DMPlexSetCone(dm, 0, cone);
1912:         ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0;
1913:         DMPlexSetConeOrientation(dm, 0, ornt);
1914:         /* Cell 1 */
1915:         cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20;
1916:         DMPlexSetCone(dm, 1, cone);
1917:         ornt[0] = 0; ornt[1] = 0; ornt[2] = -1; ornt[3] = 0;
1918:         DMPlexSetConeOrientation(dm, 1, ornt);
1919:         /* Cell 2 */
1920:         cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23;
1921:         DMPlexSetCone(dm, 2, cone);
1922:         ornt[0] = 0; ornt[1] = 0; ornt[2] = -1; ornt[3] = 0;
1923:         DMPlexSetConeOrientation(dm, 2, ornt);
1924:         /* Cell 3 */
1925:         cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15;
1926:         DMPlexSetCone(dm, 3, cone);
1927:         ornt[0] = -1; ornt[1] = -1; ornt[2] = 0; ornt[3] = -1;
1928:         DMPlexSetConeOrientation(dm, 3, ornt);
1929:         /* Cell 4 */
1930:         cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25;
1931:         DMPlexSetCone(dm, 4, cone);
1932:         ornt[0] = -1; ornt[1] = -1; ornt[2] = -1; ornt[3] = 0;
1933:         DMPlexSetConeOrientation(dm, 4, ornt);
1934:         /* Cell 5 */
1935:         cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23;
1936:         DMPlexSetCone(dm, 5, cone);
1937:         ornt[0] = -1; ornt[1] = -1; ornt[2] = -1; ornt[3] = -1;
1938:         DMPlexSetConeOrientation(dm, 5, ornt);
1939:         /* Edges */
1940:         cone[0] =  6; cone[1] =  7;
1941:         DMPlexSetCone(dm, 14, cone);
1942:         cone[0] =  7; cone[1] =  8;
1943:         DMPlexSetCone(dm, 15, cone);
1944:         cone[0] =  8; cone[1] =  9;
1945:         DMPlexSetCone(dm, 16, cone);
1946:         cone[0] =  9; cone[1] =  6;
1947:         DMPlexSetCone(dm, 17, cone);
1948:         cone[0] = 10; cone[1] = 11;
1949:         DMPlexSetCone(dm, 18, cone);
1950:         cone[0] = 11; cone[1] =  7;
1951:         DMPlexSetCone(dm, 19, cone);
1952:         cone[0] =  6; cone[1] = 10;
1953:         DMPlexSetCone(dm, 20, cone);
1954:         cone[0] = 12; cone[1] = 13;
1955:         DMPlexSetCone(dm, 21, cone);
1956:         cone[0] = 13; cone[1] = 11;
1957:         DMPlexSetCone(dm, 22, cone);
1958:         cone[0] = 10; cone[1] = 12;
1959:         DMPlexSetCone(dm, 23, cone);
1960:         cone[0] = 13; cone[1] =  8;
1961:         DMPlexSetCone(dm, 24, cone);
1962:         cone[0] = 12; cone[1] =  9;
1963:         DMPlexSetCone(dm, 25, cone);
1964:       }
1965:       DMPlexSymmetrize(dm);
1966:       DMPlexStratify(dm);
1967:       /* Build coordinates */
1968:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1969:       if (rank == 0) {
1970:         coordsIn[0*embedDim+0] = -R; coordsIn[0*embedDim+1] =  R; coordsIn[0*embedDim+2] = -R;
1971:         coordsIn[1*embedDim+0] =  R; coordsIn[1*embedDim+1] =  R; coordsIn[1*embedDim+2] = -R;
1972:         coordsIn[2*embedDim+0] =  R; coordsIn[2*embedDim+1] = -R; coordsIn[2*embedDim+2] = -R;
1973:         coordsIn[3*embedDim+0] = -R; coordsIn[3*embedDim+1] = -R; coordsIn[3*embedDim+2] = -R;
1974:         coordsIn[4*embedDim+0] = -R; coordsIn[4*embedDim+1] =  R; coordsIn[4*embedDim+2] =  R;
1975:         coordsIn[5*embedDim+0] =  R; coordsIn[5*embedDim+1] =  R; coordsIn[5*embedDim+2] =  R;
1976:         coordsIn[6*embedDim+0] = -R; coordsIn[6*embedDim+1] = -R; coordsIn[6*embedDim+2] =  R;
1977:         coordsIn[7*embedDim+0] =  R; coordsIn[7*embedDim+1] = -R; coordsIn[7*embedDim+2] =  R;
1978:       }
1979:     }
1980:     break;
1981:   case 3:
1982:     if (simplex) {
1983:       const PetscReal edgeLen         = 1.0/PETSC_PHI;
1984:       PetscReal       vertexA[4]      = {0.5, 0.5, 0.5, 0.5};
1985:       PetscReal       vertexB[4]      = {1.0, 0.0, 0.0, 0.0};
1986:       PetscReal       vertexC[4]      = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0};
1987:       const PetscInt  degree          = 12;
1988:       PetscInt        s[4]            = {1, 1, 1};
1989:       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},
1990:                                          {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}};
1991:       PetscInt        cone[4];
1992:       PetscInt       *graph, p, i, j, k, l;

1994:       vertexA[0] *= R; vertexA[1] *= R; vertexA[2] *= R; vertexA[3] *= R;
1995:       vertexB[0] *= R; vertexB[1] *= R; vertexB[2] *= R; vertexB[3] *= R;
1996:       vertexC[0] *= R; vertexC[1] *= R; vertexC[2] *= R; vertexC[3] *= R;
1997:       numCells    = rank == 0 ? 600 : 0;
1998:       numVerts    = rank == 0 ? 120 : 0;
1999:       firstVertex = numCells;
2000:       /* Use the 600-cell, which for a unit sphere has coordinates which are

2002:            1/2 (\pm 1, \pm 1,    \pm 1, \pm 1)                          16
2003:                (\pm 1,    0,       0,      0)  all cyclic permutations   8
2004:            1/2 (\pm 1, \pm phi, \pm 1/phi, 0)  all even permutations    96

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

2009:          http://buzzard.pugetsound.edu/sage-practice/ch03s03.html
2010:          http://mathworld.wolfram.com/600-Cell.html
2011:       */
2012:       /* Construct vertices */
2013:       PetscCalloc1(numVerts * embedDim, &coordsIn);
2014:       i    = 0;
2015:       if (rank == 0) {
2016:         for (s[0] = -1; s[0] < 2; s[0] += 2) {
2017:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
2018:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
2019:               for (s[3] = -1; s[3] < 2; s[3] += 2) {
2020:                 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d];
2021:                 ++i;
2022:               }
2023:             }
2024:           }
2025:         }
2026:         for (p = 0; p < embedDim; ++p) {
2027:           s[1] = s[2] = s[3] = 1;
2028:           for (s[0] = -1; s[0] < 2; s[0] += 2) {
2029:             for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim];
2030:             ++i;
2031:           }
2032:         }
2033:         for (p = 0; p < 12; ++p) {
2034:           s[3] = 1;
2035:           for (s[0] = -1; s[0] < 2; s[0] += 2) {
2036:             for (s[1] = -1; s[1] < 2; s[1] += 2) {
2037:               for (s[2] = -1; s[2] < 2; s[2] += 2) {
2038:                 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]];
2039:                 ++i;
2040:               }
2041:             }
2042:           }
2043:         }
2044:       }
2046:       /* Construct graph */
2047:       PetscCalloc1(numVerts * numVerts, &graph);
2048:       for (i = 0; i < numVerts; ++i) {
2049:         for (j = 0, k = 0; j < numVerts; ++j) {
2050:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
2051:         }
2053:       }
2054:       /* Build Topology */
2055:       DMPlexSetChart(dm, 0, numCells+numVerts);
2056:       for (c = 0; c < numCells; c++) {
2057:         DMPlexSetConeSize(dm, c, embedDim);
2058:       }
2059:       DMSetUp(dm); /* Allocate space for cones */
2060:       /* Cells */
2061:       if (rank == 0) {
2062:         for (i = 0, c = 0; i < numVerts; ++i) {
2063:           for (j = 0; j < i; ++j) {
2064:             for (k = 0; k < j; ++k) {
2065:               for (l = 0; l < k; ++l) {
2066:                 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] &&
2067:                     graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) {
2068:                   cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l;
2069:                   /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */
2070:                   {
2071:                     const PetscInt epsilon[4][4][4][4] = {{{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2072:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  1}, { 0,  0, -1, 0}},
2073:                                                            {{0,  0,  0,  0}, { 0, 0,  0, -1}, { 0,  0, 0,  0}, { 0,  1,  0, 0}},
2074:                                                            {{0,  0,  0,  0}, { 0, 0,  1,  0}, { 0, -1, 0,  0}, { 0,  0,  0, 0}}},

2076:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0, -1}, { 0,  0,  1, 0}},
2077:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2078:                                                            {{0,  0,  0,  1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, {-1,  0,  0, 0}},
2079:                                                            {{0,  0, -1,  0}, { 0, 0,  0,  0}, { 1,  0, 0,  0}, { 0,  0,  0, 0}}},

2081:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  1}, { 0,  0, 0,  0}, { 0, -1,  0, 0}},
2082:                                                            {{0,  0,  0, -1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 1,  0,  0, 0}},
2083:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2084:                                                            {{0,  1,  0,  0}, {-1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}},

2086:                                                           {{{0,  0,  0,  0}, { 0, 0, -1,  0}, { 0,  1, 0,  0}, { 0,  0,  0, 0}},
2087:                                                            {{0,  0,  1,  0}, { 0, 0,  0,  0}, {-1,  0, 0,  0}, { 0,  0,  0, 0}},
2088:                                                            {{0, -1,  0,  0}, { 1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2089:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}}};
2090:                     PetscReal normal[4];
2091:                     PetscInt  e, f, g;

2093:                     for (d = 0; d < embedDim; ++d) {
2094:                       normal[d] = 0.0;
2095:                       for (e = 0; e < embedDim; ++e) {
2096:                         for (f = 0; f < embedDim; ++f) {
2097:                           for (g = 0; g < embedDim; ++g) {
2098:                             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]);
2099:                           }
2100:                         }
2101:                       }
2102:                     }
2103:                     if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
2104:                   }
2105:                   DMPlexSetCone(dm, c++, cone);
2106:                 }
2107:               }
2108:             }
2109:           }
2110:         }
2111:       }
2112:       DMPlexSymmetrize(dm);
2113:       DMPlexStratify(dm);
2114:       PetscFree(graph);
2115:       break;
2116:     }
2117:   default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim);
2118:   }
2119:   /* Create coordinates */
2120:   DMGetCoordinateSection(dm, &coordSection);
2121:   PetscSectionSetNumFields(coordSection, 1);
2122:   PetscSectionSetFieldComponents(coordSection, 0, embedDim);
2123:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);
2124:   for (v = firstVertex; v < firstVertex+numVerts; ++v) {
2125:     PetscSectionSetDof(coordSection, v, embedDim);
2126:     PetscSectionSetFieldDof(coordSection, v, 0, embedDim);
2127:   }
2128:   PetscSectionSetUp(coordSection);
2129:   PetscSectionGetStorageSize(coordSection, &coordSize);
2130:   VecCreate(PETSC_COMM_SELF, &coordinates);
2131:   VecSetBlockSize(coordinates, embedDim);
2132:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2133:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2134:   VecSetType(coordinates,VECSTANDARD);
2135:   VecGetArray(coordinates, &coords);
2136:   for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];}
2137:   VecRestoreArray(coordinates, &coords);
2138:   DMSetCoordinatesLocal(dm, coordinates);
2139:   VecDestroy(&coordinates);
2140:   PetscFree(coordsIn);
2141:   {
2142:     DM          cdm;
2143:     PetscDS     cds;
2144:     PetscScalar c = R;

2146:     DMPlexCreateCoordinateSpace(dm, 1, snapToSphere);
2147:     DMGetCoordinateDM(dm, &cdm);
2148:     DMGetDS(cdm, &cds);
2149:     PetscDSSetConstants(cds, 1, &c);
2150:   }
2151:   /* Wait for coordinate creation before doing in-place modification */
2152:   if (simplex) DMPlexInterpolateInPlace_Internal(dm);
2153:   return 0;
2154: }

2156: typedef void (*TPSEvaluateFunc)(const PetscReal[], PetscReal*, PetscReal[], PetscReal(*)[3]);

2158: /*
2159:  The Schwarz P implicit surface is

2161:      f(x) = cos(x0) + cos(x1) + cos(x2) = 0
2162: */
2163: static void TPSEvaluate_SchwarzP(const PetscReal y[3], PetscReal *f, PetscReal grad[], PetscReal (*hess)[3])
2164: {
2165:   PetscReal c[3] = {PetscCosReal(y[0] * PETSC_PI), PetscCosReal(y[1] * PETSC_PI), PetscCosReal(y[2] * PETSC_PI)};
2166:   PetscReal g[3] = {-PetscSinReal(y[0] * PETSC_PI), -PetscSinReal(y[1] * PETSC_PI), -PetscSinReal(y[2] * PETSC_PI)};
2167:   f[0] = c[0] + c[1] + c[2];
2168:   for (PetscInt i=0; i<3; i++) {
2169:     grad[i] = PETSC_PI * g[i];
2170:     for (PetscInt j=0; j<3; j++) {
2171:       hess[i][j] = (i == j) ? -PetscSqr(PETSC_PI) * c[i] : 0.;
2172:     }
2173:   }
2174: }

2176: // u[] is a tentative normal on input. Replace with the implicit function gradient in the same direction
2177: static PetscErrorCode TPSExtrudeNormalFunc_SchwarzP(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) {
2178:   for (PetscInt i=0; i<3; i++) {
2179:     u[i] = -PETSC_PI * PetscSinReal(x[i] * PETSC_PI);
2180:   }
2181:   return 0;
2182: }

2184: /*
2185:  The Gyroid implicit surface is

2187:  f(x,y,z) = sin(pi * x) * cos (pi * (y + 1/2))  + sin(pi * (y + 1/2)) * cos(pi * (z + 1/4)) + sin(pi * (z + 1/4)) * cos(pi * x)

2189: */
2190: static void TPSEvaluate_Gyroid(const PetscReal y[3], PetscReal *f, PetscReal grad[], PetscReal (*hess)[3])
2191: {
2192:   PetscReal s[3] = {PetscSinReal(PETSC_PI * y[0]), PetscSinReal(PETSC_PI * (y[1] + .5)), PetscSinReal(PETSC_PI * (y[2] + .25))};
2193:   PetscReal c[3] = {PetscCosReal(PETSC_PI * y[0]), PetscCosReal(PETSC_PI * (y[1] + .5)), PetscCosReal(PETSC_PI * (y[2] + .25))};
2194:   f[0] = s[0] * c[1] + s[1] * c[2] + s[2] * c[0];
2195:   grad[0] = PETSC_PI * (c[0] * c[1] - s[2] * s[0]);
2196:   grad[1] = PETSC_PI * (c[1] * c[2] - s[0] * s[1]);
2197:   grad[2] = PETSC_PI * (c[2] * c[0] - s[1] * s[2]);
2198:   hess[0][0] = -PetscSqr(PETSC_PI) * (s[0] * c[1] + s[2] * c[0]);
2199:   hess[0][1] = -PetscSqr(PETSC_PI) * (c[0] * s[1]);
2200:   hess[0][2] = -PetscSqr(PETSC_PI) * (c[2] * s[0]);
2201:   hess[1][0] = -PetscSqr(PETSC_PI) * (s[1] * c[2] + s[0] * c[1]);
2202:   hess[1][1] = -PetscSqr(PETSC_PI) * (c[1] * s[2]);
2203:   hess[2][2] = -PetscSqr(PETSC_PI) * (c[0] * s[1]);
2204:   hess[2][0] = -PetscSqr(PETSC_PI) * (s[2] * c[0] + s[1] * c[2]);
2205:   hess[2][1] = -PetscSqr(PETSC_PI) * (c[2] * s[0]);
2206:   hess[2][2] = -PetscSqr(PETSC_PI) * (c[1] * s[2]);
2207: }

2209: // u[] is a tentative normal on input. Replace with the implicit function gradient in the same direction
2210: static PetscErrorCode TPSExtrudeNormalFunc_Gyroid(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) {
2211:   PetscReal s[3] = {PetscSinReal(PETSC_PI * x[0]), PetscSinReal(PETSC_PI * (x[1] + .5)), PetscSinReal(PETSC_PI * (x[2] + .25))};
2212:   PetscReal c[3] = {PetscCosReal(PETSC_PI * x[0]), PetscCosReal(PETSC_PI * (x[1] + .5)), PetscCosReal(PETSC_PI * (x[2] + .25))};
2213:   u[0] = PETSC_PI * (c[0] * c[1] - s[2] * s[0]);
2214:   u[1] = PETSC_PI * (c[1] * c[2] - s[0] * s[1]);
2215:   u[2] = PETSC_PI * (c[2] * c[0] - s[1] * s[2]);
2216:   return 0;
2217: }

2219: /*
2220:    We wish to solve

2222:          min_y || y - x ||^2  subject to f(y) = 0

2224:    Let g(y) = grad(f).  The minimization problem is equivalent to asking to satisfy
2225:    f(y) = 0 and (y-x) is parallel to g(y).  We do this by using Householder QR to obtain a basis for the
2226:    tangent space and ask for both components in the tangent space to be zero.

2228:    Take g to be a column vector and compute the "full QR" factorization Q R = g,
2229:    where Q = I - 2 n n^T is a symmetric orthogonal matrix.
2230:    The first column of Q is parallel to g so the remaining two columns span the null space.
2231:    Let Qn = Q[:,1:] be those remaining columns.  Then Qn Qn^T is an orthogonal projector into the tangent space.
2232:    Since Q is symmetric, this is equivalent to multipyling by Q and taking the last two entries.
2233:    In total, we have a system of 3 equations in 3 unknowns:

2235:      f(y) = 0                       1 equation
2236:      Qn^T (y - x) = 0               2 equations

2238:    Here, we compute the residual and Jacobian of this system.
2239: */
2240: static void TPSNearestPointResJac(TPSEvaluateFunc feval, const PetscScalar x[], const PetscScalar y[], PetscScalar res[], PetscScalar J[])
2241: {
2242:   PetscReal yreal[3] = {PetscRealPart(y[0]), PetscRealPart(y[1]), PetscRealPart(y[2])};
2243:   PetscReal d[3] = {PetscRealPart(y[0] - x[0]), PetscRealPart(y[1] - x[1]), PetscRealPart(y[2] - x[2])};
2244:   PetscReal f, grad[3], n[3], norm, norm_y[3], nd, nd_y[3], sign;
2245:   PetscReal n_y[3][3] = {{0,0,0},{0,0,0},{0,0,0}};

2247:   feval(yreal, &f, grad, n_y);

2249:   for (PetscInt i=0; i<3; i++) n[i] = grad[i];
2250:   norm = PetscSqrtReal(PetscSqr(n[0]) + PetscSqr(n[1]) + PetscSqr(n[2]));
2251:   for (PetscInt i=0; i<3; i++) {
2252:     norm_y[i] = 1. / norm * n[i] * n_y[i][i];
2253:   }

2255:   // Define the Householder reflector
2256:   sign = n[0] >= 0 ? 1. : -1.;
2257:   n[0] += norm * sign;
2258:   for (PetscInt i=0; i<3; i++) n_y[0][i] += norm_y[i] * sign;

2260:   norm = PetscSqrtReal(PetscSqr(n[0]) + PetscSqr(n[1]) + PetscSqr(n[2]));
2261:   norm_y[0] = 1. / norm * (n[0] * n_y[0][0]);
2262:   norm_y[1] = 1. / norm * (n[0] * n_y[0][1] + n[1] * n_y[1][1]);
2263:   norm_y[2] = 1. / norm * (n[0] * n_y[0][2] + n[2] * n_y[2][2]);

2265:   for (PetscInt i=0; i<3; i++) {
2266:     n[i] /= norm;
2267:     for (PetscInt j=0; j<3; j++) {
2268:       // note that n[i] is n_old[i]/norm when executing the code below
2269:       n_y[i][j] = n_y[i][j] / norm - n[i] / norm * norm_y[j];
2270:     }
2271:   }

2273:   nd = n[0] * d[0] + n[1] * d[1] + n[2] * d[2];
2274:   for (PetscInt i=0; i<3; i++) nd_y[i] = n[i] + n_y[0][i] * d[0] + n_y[1][i] * d[1] + n_y[2][i] * d[2];

2276:   res[0] = f;
2277:   res[1] = d[1] - 2 * n[1] * nd;
2278:   res[2] = d[2] - 2 * n[2] * nd;
2279:   // J[j][i] is J_{ij} (column major)
2280:   for (PetscInt j=0; j<3; j++) {
2281:     J[0 + j*3] = grad[j];
2282:     J[1 + j*3] = (j == 1)*1. - 2 * (n_y[1][j] * nd + n[1] * nd_y[j]);
2283:     J[2 + j*3] = (j == 2)*1. - 2 * (n_y[2][j] * nd + n[2] * nd_y[j]);
2284:   }
2285: }

2287: /*
2288:    Project x to the nearest point on the implicit surface using Newton's method.
2289: */
2290: static PetscErrorCode TPSNearestPoint(TPSEvaluateFunc feval, PetscScalar x[])
2291: {
2292:   PetscScalar y[3] = {x[0], x[1], x[2]}; // Initial guess

2294:   for (PetscInt iter=0; iter<10; iter++) {
2295:     PetscScalar res[3], J[9];
2296:     PetscReal resnorm;
2297:     TPSNearestPointResJac(feval, x, y, res, J);
2298:     resnorm = PetscSqrtReal(PetscSqr(PetscRealPart(res[0])) + PetscSqr(PetscRealPart(res[1])) + PetscSqr(PetscRealPart(res[2])));
2299:     if (0) { // Turn on this monitor if you need to confirm quadratic convergence
2300:       PetscPrintf(PETSC_COMM_SELF, "[%D] res [%g %g %g]\n", iter, PetscRealPart(res[0]), PetscRealPart(res[1]), PetscRealPart(res[2]));
2301:     }
2302:     if (resnorm < PETSC_SMALL) break;

2304:     // Take the Newton step
2305:     PetscKernel_A_gets_inverse_A_3(J, 0., PETSC_FALSE, NULL);
2306:     PetscKernel_v_gets_v_minus_A_times_w_3(y, J, res);
2307:   }
2308:   for (PetscInt i=0; i<3; i++) x[i] = y[i];
2309:   return 0;
2310: }

2312: const char *const DMPlexTPSTypes[] = {"SCHWARZ_P", "GYROID", "DMPlexTPSType", "DMPLEX_TPS_", NULL};

2314: static PetscErrorCode DMPlexCreateTPSMesh_Internal(DM dm, DMPlexTPSType tpstype, const PetscInt extent[], const DMBoundaryType periodic[], PetscBool tps_distribute, PetscInt refinements, PetscInt layers, PetscReal thickness)
2315: {
2316:   PetscMPIInt rank;
2317:   PetscInt topoDim = 2, spaceDim = 3, numFaces = 0, numVertices = 0, numEdges = 0;
2318:   PetscInt (*edges)[2] = NULL, *edgeSets = NULL;
2319:   PetscInt *cells_flat = NULL;
2320:   PetscReal *vtxCoords = NULL;
2321:   TPSEvaluateFunc evalFunc = NULL;
2322:   PetscSimplePointFunc normalFunc = NULL;
2323:   DMLabel label;

2325:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
2327:   switch (tpstype) {
2328:   case DMPLEX_TPS_SCHWARZ_P:
2330:     if (!rank) {
2331:       PetscInt (*cells)[6][4][4] = NULL; // [junction, junction-face, cell, conn]
2332:       PetscInt Njunctions = 0, Ncuts = 0, Npipes[3], vcount;
2333:       PetscReal L = 1;

2335:       Npipes[0] = (extent[0] + 1) * extent[1] * extent[2];
2336:       Npipes[1] = extent[0] * (extent[1] + 1) * extent[2];
2337:       Npipes[2] = extent[0] * extent[1] * (extent[2] + 1);
2338:       Njunctions = extent[0] * extent[1] * extent[2];
2339:       Ncuts = 2 * (extent[0] * extent[1] + extent[1] * extent[2] + extent[2] * extent[0]);
2340:       numVertices = 4 * (Npipes[0] + Npipes[1] + Npipes[2]) + 8 * Njunctions;
2341:       PetscMalloc1(3*numVertices, &vtxCoords);
2342:       PetscMalloc1(Njunctions, &cells);
2343:       PetscMalloc1(Ncuts*4, &edges);
2344:       PetscMalloc1(Ncuts*4, &edgeSets);
2345:       // x-normal pipes
2346:       vcount = 0;
2347:       for (PetscInt i=0; i<extent[0]+1; i++) {
2348:         for (PetscInt j=0; j<extent[1]; j++) {
2349:           for (PetscInt k=0; k<extent[2]; k++) {
2350:             for (PetscInt l=0; l<4; l++) {
2351:               vtxCoords[vcount++] = (2*i - 1) * L;
2352:               vtxCoords[vcount++] = 2 * j * L + PetscCosReal((2*l + 1) * PETSC_PI / 4) * L / 2;
2353:               vtxCoords[vcount++] = 2 * k * L + PetscSinReal((2*l + 1) * PETSC_PI / 4) * L / 2;
2354:             }
2355:           }
2356:         }
2357:       }
2358:       // y-normal pipes
2359:       for (PetscInt i=0; i<extent[0]; i++) {
2360:         for (PetscInt j=0; j<extent[1]+1; j++) {
2361:           for (PetscInt k=0; k<extent[2]; k++) {
2362:             for (PetscInt l=0; l<4; l++) {
2363:               vtxCoords[vcount++] = 2 * i * L + PetscSinReal((2*l + 1) * PETSC_PI / 4) * L / 2;
2364:               vtxCoords[vcount++] = (2*j - 1) * L;
2365:               vtxCoords[vcount++] = 2 * k * L + PetscCosReal((2*l + 1) * PETSC_PI / 4) * L / 2;
2366:             }
2367:           }
2368:         }
2369:       }
2370:       // z-normal pipes
2371:       for (PetscInt i=0; i<extent[0]; i++) {
2372:         for (PetscInt j=0; j<extent[1]; j++) {
2373:           for (PetscInt k=0; k<extent[2]+1; k++) {
2374:             for (PetscInt l=0; l<4; l++) {
2375:               vtxCoords[vcount++] = 2 * i * L + PetscCosReal((2*l + 1) * PETSC_PI / 4) * L / 2;
2376:               vtxCoords[vcount++] = 2 * j * L + PetscSinReal((2*l + 1) * PETSC_PI / 4) * L / 2;
2377:               vtxCoords[vcount++] = (2*k - 1) * L;
2378:             }
2379:           }
2380:         }
2381:       }
2382:       // junctions
2383:       for (PetscInt i=0; i<extent[0]; i++) {
2384:         for (PetscInt j=0; j<extent[1]; j++) {
2385:           for (PetscInt k=0; k<extent[2]; k++) {
2386:             const PetscInt J = (i*extent[1] + j)*extent[2] + k, Jvoff = (Npipes[0] + Npipes[1] + Npipes[2])*4 + J*8;
2388:             for (PetscInt ii=0; ii<2; ii++) {
2389:               for (PetscInt jj=0; jj<2; jj++) {
2390:                 for (PetscInt kk=0; kk<2; kk++) {
2391:                   double Ls = (1 - sqrt(2) / 4) * L;
2392:                   vtxCoords[vcount++] = 2*i*L + (2*ii-1) * Ls;
2393:                   vtxCoords[vcount++] = 2*j*L + (2*jj-1) * Ls;
2394:                   vtxCoords[vcount++] = 2*k*L + (2*kk-1) * Ls;
2395:                 }
2396:               }
2397:             }
2398:             const PetscInt jfaces[3][2][4] = {
2399:               {{3,1,0,2}, {7,5,4,6}}, // x-aligned
2400:               {{5,4,0,1}, {7,6,2,3}}, // y-aligned
2401:               {{6,2,0,4}, {7,3,1,5}}  // z-aligned
2402:             };
2403:             const PetscInt pipe_lo[3] = { // vertex numbers of pipes
2404:               ((i * extent[1] + j) * extent[2] + k)*4,
2405:               ((i * (extent[1] + 1) + j) * extent[2] + k + Npipes[0])*4,
2406:               ((i * extent[1] + j) * (extent[2]+1) + k + Npipes[0] + Npipes[1])*4
2407:             };
2408:             const PetscInt pipe_hi[3] = { // vertex numbers of pipes
2409:               (((i + 1) * extent[1] + j) * extent[2] + k)*4,
2410:               ((i * (extent[1] + 1) + j + 1) * extent[2] + k + Npipes[0])*4,
2411:               ((i * extent[1] + j) * (extent[2]+1) + k + 1 + Npipes[0] + Npipes[1])*4
2412:             };
2413:             for (PetscInt dir=0; dir<3; dir++) { // x,y,z
2414:               const PetscInt ijk[3] = {i, j, k};
2415:               for (PetscInt l=0; l<4; l++) { // rotations
2416:                 cells[J][dir*2+0][l][0] = pipe_lo[dir] + l;
2417:                 cells[J][dir*2+0][l][1] = Jvoff + jfaces[dir][0][l];
2418:                 cells[J][dir*2+0][l][2] = Jvoff + jfaces[dir][0][(l-1+4)%4];
2419:                 cells[J][dir*2+0][l][3] = pipe_lo[dir] + (l-1+4)%4;
2420:                 cells[J][dir*2+1][l][0] = Jvoff + jfaces[dir][1][l];
2421:                 cells[J][dir*2+1][l][1] = pipe_hi[dir] + l;
2422:                 cells[J][dir*2+1][l][2] = pipe_hi[dir] + (l-1+4)%4;
2423:                 cells[J][dir*2+1][l][3] = Jvoff + jfaces[dir][1][(l-1+4)%4];
2424:                 if (ijk[dir] == 0) {
2425:                   edges[numEdges][0] = pipe_lo[dir] + l;
2426:                   edges[numEdges][1] = pipe_lo[dir] + (l+1) % 4;
2427:                   edgeSets[numEdges] = dir*2 + 1;
2428:                   numEdges++;
2429:                 }
2430:                 if (ijk[dir] + 1 == extent[dir]) {
2431:                   edges[numEdges][0] = pipe_hi[dir] + l;
2432:                   edges[numEdges][1] = pipe_hi[dir] + (l+1) % 4;
2433:                   edgeSets[numEdges] = dir*2 + 2;
2434:                   numEdges++;
2435:                 }
2436:               }
2437:             }
2438:           }
2439:         }
2440:       }
2442:       numFaces = 24 * Njunctions;
2443:       cells_flat = cells[0][0][0];
2444:     }
2445:     evalFunc = TPSEvaluate_SchwarzP;
2446:     normalFunc = TPSExtrudeNormalFunc_SchwarzP;
2447:     break;
2448:   case DMPLEX_TPS_GYROID:
2449:     if (!rank) {
2450:       // This is a coarse mesh approximation of the gyroid shifted to being the zero of the level set
2451:       //
2452:       //     sin(pi*x)*cos(pi*(y+1/2)) + sin(pi*(y+1/2))*cos(pi*(z+1/4)) + sin(pi*(z+1/4))*cos(x)
2453:       //
2454:       // on the cell [0,2]^3.
2455:       //
2456:       // Think about dividing that cell into four columns, and focus on the column [0,1]x[0,1]x[0,2].
2457:       // If you looked at the gyroid in that column at different slices of z you would see that it kind of spins
2458:       // like a boomerang:
2459:       //
2460:       //     z = 0          z = 1/4        z = 1/2        z = 3/4     //
2461:       //     -----          -------        -------        -------     //
2462:       //                                                              //
2463:       //     +       +      +       +      +       +      +   \   +   //
2464:       //      \                                   /            \      //
2465:       //       \            `-_   _-'            /              }     //
2466:       //        *-_            `-'            _-'              /      //
2467:       //     +     `-+      +       +      +-'     +      +   /   +   //
2468:       //                                                              //
2469:       //                                                              //
2470:       //     z = 1          z = 5/4        z = 3/2        z = 7/4     //
2471:       //     -----          -------        -------        -------     //
2472:       //                                                              //
2473:       //     +-_     +      +       +      +     _-+      +   /   +   //
2474:       //        `-_            _-_            _-`            /        //
2475:       //           \        _-'   `-_        /              {         //
2476:       //            \                       /                \        //
2477:       //     +       +      +       +      +       +      +   \   +   //
2478:       //
2479:       //
2480:       // This course mesh approximates each of these slices by two line segments,
2481:       // and then connects the segments in consecutive layers with quadrilateral faces.
2482:       // All of the end points of the segments are multiples of 1/4 except for the
2483:       // point * in the picture for z = 0 above and the similar points in other layers.
2484:       // That point is at (gamma, gamma, 0), where gamma is calculated below.
2485:       //
2486:       // The column  [1,2]x[1,2]x[0,2] looks the same as this column;
2487:       // The columns [1,2]x[0,1]x[0,2] and [0,1]x[1,2]x[0,2] are mirror images.
2488:       //
2489:       // As for how this method turned into the names given to the vertices:
2490:       // that was not systematic, it was just the way it worked out in my handwritten notes.

2492:       PetscInt facesPerBlock = 64;
2493:       PetscInt vertsPerBlock = 56;
2494:       PetscInt extentPlus[3];
2495:       PetscInt numBlocks, numBlocksPlus;
2496:       const PetscInt A =  0,   B =  1,   C =  2,   D =  3,   E =  4,   F =  5,   G =  6,   H =  7,
2497:         II =  8,   J =  9,   K = 10,   L = 11,   M = 12,   N = 13,   O = 14,   P = 15,
2498:         Q = 16,   R = 17,   S = 18,   T = 19,   U = 20,   V = 21,   W = 22,   X = 23,
2499:         Y = 24,   Z = 25,  Ap = 26,  Bp = 27,  Cp = 28,  Dp = 29,  Ep = 30,  Fp = 31,
2500:         Gp = 32,  Hp = 33,  Ip = 34,  Jp = 35,  Kp = 36,  Lp = 37,  Mp = 38,  Np = 39,
2501:         Op = 40,  Pp = 41,  Qp = 42,  Rp = 43,  Sp = 44,  Tp = 45,  Up = 46,  Vp = 47,
2502:         Wp = 48,  Xp = 49,  Yp = 50,  Zp = 51,  Aq = 52,  Bq = 53,  Cq = 54,  Dq = 55;
2503:       const PetscInt pattern[64][4] =
2504:         { /* face to vertex within the coarse discretization of a single gyroid block */
2505:           /* layer 0 */
2506:           {A,C,K,G},{C,B,II,K},{D,A,H,L},{B+56*1,D,L,J},{E,B+56*1,J,N},{A+56*2,E,N,H+56*2},{F,A+56*2,G+56*2,M},{B,F,M,II},
2507:           /* layer 1 */
2508:           {G,K,Q,O},{K,II,P,Q},{L,H,O+56*1,R},{J,L,R,P},{N,J,P,S},{H+56*2,N,S,O+56*3},{M,G+56*2,O+56*2,T},{II,M,T,P},
2509:           /* layer 2 */
2510:           {O,Q,Y,U},{Q,P,W,Y},{R,O+56*1,U+56*1,Ap},{P,R,Ap,W},{S,P,X,Bp},{O+56*3,S,Bp,V+56*1},{T,O+56*2,V,Z},{P,T,Z,X},
2511:           /* layer 3 */
2512:           {U,Y,Ep,Dp},{Y,W,Cp,Ep},{Ap,U+56*1,Dp+56*1,Gp},{W,Ap,Gp,Cp},{Bp,X,Cp+56*2,Fp},{V+56*1,Bp,Fp,Dp+56*1},{Z,V,Dp,Hp},{X,Z,Hp,Cp+56*2},
2513:           /* layer 4 */
2514:           {Dp,Ep,Mp,Kp},{Ep,Cp,Ip,Mp},{Gp,Dp+56*1,Lp,Np},{Cp,Gp,Np,Jp},{Fp,Cp+56*2,Jp+56*2,Pp},{Dp+56*1,Fp,Pp,Lp},{Hp,Dp,Kp,Op},{Cp+56*2,Hp,Op,Ip+56*2},
2515:           /* layer 5 */
2516:           {Kp,Mp,Sp,Rp},{Mp,Ip,Qp,Sp},{Np,Lp,Rp,Tp},{Jp,Np,Tp,Qp+56*1},{Pp,Jp+56*2,Qp+56*3,Up},{Lp,Pp,Up,Rp},{Op,Kp,Rp,Vp},{Ip+56*2,Op,Vp,Qp+56*2},
2517:           /* layer 6 */
2518:           {Rp,Sp,Aq,Yp},{Sp,Qp,Wp,Aq},{Tp,Rp,Yp,Cq},{Qp+56*1,Tp,Cq,Wp+56*1},{Up,Qp+56*3,Xp+56*1,Dq},{Rp,Up,Dq,Zp},{Vp,Rp,Zp,Bq},{Qp+56*2,Vp,Bq,Xp},
2519:           /* layer 7 (the top is the periodic image of the bottom of layer 0) */
2520:           {Yp,Aq,C+56*4,A+56*4},{Aq,Wp,B+56*4,C+56*4},{Cq,Yp,A+56*4,D+56*4},{Wp+56*1,Cq,D+56*4,B+56*5},{Dq,Xp+56*1,B+56*5,E+56*4},{Zp,Dq,E+56*4,A+56*6},{Bq,Zp,A+56*6,F+56*4},{Xp,Bq,F+56*4,B+56*4}
2521:         };
2522:       const PetscReal gamma = PetscAcosReal((PetscSqrtReal(3.)-1.) / PetscSqrtReal(2.)) / PETSC_PI;
2523:       const PetscReal patternCoords[56][3] =
2524:         {
2525:           /* A  */ {1.,0.,0.},
2526:           /* B  */ {0.,1.,0.},
2527:           /* C  */ {gamma,gamma,0.},
2528:           /* D  */ {1+gamma,1-gamma,0.},
2529:           /* E  */ {2-gamma,2-gamma,0.},
2530:           /* F  */ {1-gamma,1+gamma,0.},

2532:           /* G  */ {.5,0,.25},
2533:           /* H  */ {1.5,0.,.25},
2534:           /* II */ {.5,1.,.25},
2535:           /* J  */ {1.5,1.,.25},
2536:           /* K  */ {.25,.5,.25},
2537:           /* L  */ {1.25,.5,.25},
2538:           /* M  */ {.75,1.5,.25},
2539:           /* N  */ {1.75,1.5,.25},

2541:           /* O  */ {0.,0.,.5},
2542:           /* P  */ {1.,1.,.5},
2543:           /* Q  */ {gamma,1-gamma,.5},
2544:           /* R  */ {1+gamma,gamma,.5},
2545:           /* S  */ {2-gamma,1+gamma,.5},
2546:           /* T  */ {1-gamma,2-gamma,.5},

2548:           /* U  */ {0.,.5,.75},
2549:           /* V  */ {0.,1.5,.75},
2550:           /* W  */ {1.,.5,.75},
2551:           /* X  */ {1.,1.5,.75},
2552:           /* Y  */ {.5,.75,.75},
2553:           /* Z  */ {.5,1.75,.75},
2554:           /* Ap */ {1.5,.25,.75},
2555:           /* Bp */ {1.5,1.25,.75},

2557:           /* Cp */ {1.,0.,1.},
2558:           /* Dp */ {0.,1.,1.},
2559:           /* Ep */ {1-gamma,1-gamma,1.},
2560:           /* Fp */ {1+gamma,1+gamma,1.},
2561:           /* Gp */ {2-gamma,gamma,1.},
2562:           /* Hp */ {gamma,2-gamma,1.},

2564:           /* Ip */ {.5,0.,1.25},
2565:           /* Jp */ {1.5,0.,1.25},
2566:           /* Kp */ {.5,1.,1.25},
2567:           /* Lp */ {1.5,1.,1.25},
2568:           /* Mp */ {.75,.5,1.25},
2569:           /* Np */ {1.75,.5,1.25},
2570:           /* Op */ {.25,1.5,1.25},
2571:           /* Pp */ {1.25,1.5,1.25},

2573:           /* Qp */ {0.,0.,1.5},
2574:           /* Rp */ {1.,1.,1.5},
2575:           /* Sp */ {1-gamma,gamma,1.5},
2576:           /* Tp */ {2-gamma,1-gamma,1.5},
2577:           /* Up */ {1+gamma,2-gamma,1.5},
2578:           /* Vp */ {gamma,1+gamma,1.5},

2580:           /* Wp */ {0.,.5,1.75},
2581:           /* Xp */ {0.,1.5,1.75},
2582:           /* Yp */ {1.,.5,1.75},
2583:           /* Zp */ {1.,1.5,1.75},
2584:           /* Aq */ {.5,.25,1.75},
2585:           /* Bq */ {.5,1.25,1.75},
2586:           /* Cq */ {1.5,.75,1.75},
2587:           /* Dq */ {1.5,1.75,1.75},
2588:         };
2589:       PetscInt  (*cells)[64][4] = NULL;
2590:       PetscBool *seen;
2591:       PetscInt  *vertToTrueVert;
2592:       PetscInt  count;

2594:       for (PetscInt i = 0; i < 3; i++) extentPlus[i]  = extent[i] + 1;
2595:       numBlocks = 1;
2596:       for (PetscInt i = 0; i < 3; i++)     numBlocks *= extent[i];
2597:       numBlocksPlus = 1;
2598:       for (PetscInt i = 0; i < 3; i++) numBlocksPlus *= extentPlus[i];
2599:       numFaces = numBlocks * facesPerBlock;
2600:       PetscMalloc1(numBlocks, &cells);
2601:       PetscCalloc1(numBlocksPlus * vertsPerBlock,&seen);
2602:       for (PetscInt k = 0; k < extent[2]; k++) {
2603:         for (PetscInt j = 0; j < extent[1]; j++) {
2604:           for (PetscInt i = 0; i < extent[0]; i++) {
2605:             for (PetscInt f = 0; f < facesPerBlock; f++) {
2606:               for (PetscInt v = 0; v < 4; v++) {
2607:                 PetscInt vertRaw = pattern[f][v];
2608:                 PetscInt blockidx = vertRaw / 56;
2609:                 PetscInt patternvert = vertRaw % 56;
2610:                 PetscInt xplus = (blockidx & 1);
2611:                 PetscInt yplus = (blockidx & 2) >> 1;
2612:                 PetscInt zplus = (blockidx & 4) >> 2;
2613:                 PetscInt zcoord = (periodic && periodic[2] == DM_BOUNDARY_PERIODIC) ? ((k + zplus) % extent[2]) : (k + zplus);
2614:                 PetscInt ycoord = (periodic && periodic[1] == DM_BOUNDARY_PERIODIC) ? ((j + yplus) % extent[1]) : (j + yplus);
2615:                 PetscInt xcoord = (periodic && periodic[0] == DM_BOUNDARY_PERIODIC) ? ((i + xplus) % extent[0]) : (i + xplus);
2616:                 PetscInt vert = ((zcoord * extentPlus[1] + ycoord) * extentPlus[0] + xcoord) * 56 + patternvert;

2618:                 cells[(k * extent[1] + j) * extent[0] + i][f][v] = vert;
2619:                 seen[vert] = PETSC_TRUE;
2620:               }
2621:             }
2622:           }
2623:         }
2624:       }
2625:       for (PetscInt i = 0; i < numBlocksPlus * vertsPerBlock; i++) if (seen[i]) numVertices++;
2626:       count = 0;
2627:       PetscMalloc1(numBlocksPlus * vertsPerBlock, &vertToTrueVert);
2628:       PetscMalloc1(numVertices * 3, &vtxCoords);
2629:       for (PetscInt i = 0; i < numBlocksPlus * vertsPerBlock; i++) vertToTrueVert[i] = -1;
2630:       for (PetscInt k = 0; k < extentPlus[2]; k++) {
2631:         for (PetscInt j = 0; j < extentPlus[1]; j++) {
2632:           for (PetscInt i = 0; i < extentPlus[0]; i++) {
2633:             for (PetscInt v = 0; v < vertsPerBlock; v++) {
2634:               PetscInt vIdx = ((k * extentPlus[1] + j) * extentPlus[0] + i) * vertsPerBlock + v;

2636:               if (seen[vIdx]) {
2637:                 PetscInt thisVert;

2639:                 vertToTrueVert[vIdx] = thisVert = count++;

2641:                 for (PetscInt d = 0; d < 3; d++) vtxCoords[3 * thisVert + d] = patternCoords[v][d];
2642:                 vtxCoords[3 * thisVert + 0] += i * 2;
2643:                 vtxCoords[3 * thisVert + 1] += j * 2;
2644:                 vtxCoords[3 * thisVert + 2] += k * 2;
2645:               }
2646:             }
2647:           }
2648:         }
2649:       }
2650:       for (PetscInt i = 0; i < numBlocks; i++) {
2651:         for (PetscInt f = 0; f < facesPerBlock; f++) {
2652:           for (PetscInt v = 0; v < 4; v++) {
2653:             cells[i][f][v] = vertToTrueVert[cells[i][f][v]];
2654:           }
2655:         }
2656:       }
2657:       PetscFree(vertToTrueVert);
2658:       PetscFree(seen);
2659:       cells_flat = cells[0][0];
2660:       numEdges = 0;
2661:       for (PetscInt i = 0; i < numFaces; i++) {
2662:         for (PetscInt e = 0; e < 4; e++) {
2663:           PetscInt ev[] = {cells_flat[i*4 + e], cells_flat[i*4 + ((e+1)%4)]};
2664:           const PetscReal *evCoords[] = {&vtxCoords[3*ev[0]], &vtxCoords[3*ev[1]]};

2666:           for (PetscInt d = 0; d < 3; d++) {
2667:             if (!periodic || periodic[0] != DM_BOUNDARY_PERIODIC) {
2668:               if (evCoords[0][d] == 0. && evCoords[1][d] == 0.) numEdges++;
2669:               if (evCoords[0][d] == 2.*extent[d] && evCoords[1][d] == 2.*extent[d]) numEdges++;
2670:             }
2671:           }
2672:         }
2673:       }
2674:       PetscMalloc1(numEdges, &edges);
2675:       PetscMalloc1(numEdges, &edgeSets);
2676:       for (PetscInt edge = 0, i = 0; i < numFaces; i++) {
2677:         for (PetscInt e = 0; e < 4; e++) {
2678:           PetscInt ev[] = {cells_flat[i*4 + e], cells_flat[i*4 + ((e+1)%4)]};
2679:           const PetscReal *evCoords[] = {&vtxCoords[3*ev[0]], &vtxCoords[3*ev[1]]};

2681:           for (PetscInt d = 0; d < 3; d++) {
2682:             if (!periodic || periodic[d] != DM_BOUNDARY_PERIODIC) {
2683:               if (evCoords[0][d] == 0. && evCoords[1][d] == 0.) {
2684:                 edges[edge][0] = ev[0];
2685:                 edges[edge][1] = ev[1];
2686:                 edgeSets[edge++] = 2 * d;
2687:               }
2688:               if (evCoords[0][d] == 2.*extent[d] && evCoords[1][d] == 2.*extent[d]) {
2689:                 edges[edge][0] = ev[0];
2690:                 edges[edge][1] = ev[1];
2691:                 edgeSets[edge++] = 2 * d + 1;
2692:               }
2693:             }
2694:           }
2695:         }
2696:       }
2697:     }
2698:     evalFunc = TPSEvaluate_Gyroid;
2699:     normalFunc = TPSExtrudeNormalFunc_Gyroid;
2700:     break;
2701:   }

2703:   DMSetDimension(dm, topoDim);
2704:   if (!rank) DMPlexBuildFromCellList(dm, numFaces, numVertices, 4, cells_flat);
2705:   else       DMPlexBuildFromCellList(dm, 0, 0, 0, NULL);
2706:   PetscFree(cells_flat);
2707:   {
2708:     DM idm;
2709:     DMPlexInterpolate(dm, &idm);
2710:     DMPlexReplace_Static(dm, &idm);
2711:   }
2712:   if (!rank) DMPlexBuildCoordinatesFromCellList(dm, spaceDim, vtxCoords);
2713:   else       DMPlexBuildCoordinatesFromCellList(dm, spaceDim, NULL);
2714:   PetscFree(vtxCoords);

2716:   DMCreateLabel(dm, "Face Sets");
2717:   DMGetLabel(dm, "Face Sets", &label);
2718:   for (PetscInt e=0; e<numEdges; e++) {
2719:     PetscInt njoin;
2720:     const PetscInt *join, verts[] = {numFaces + edges[e][0], numFaces + edges[e][1]};
2721:     DMPlexGetJoin(dm, 2, verts, &njoin, &join);
2723:     DMLabelSetValue(label, join[0], edgeSets[e]);
2724:     DMPlexRestoreJoin(dm, 2, verts, &njoin, &join);
2725:   }
2726:   PetscFree(edges);
2727:   PetscFree(edgeSets);
2728:   if (tps_distribute) {
2729:     DM               pdm = NULL;
2730:     PetscPartitioner part;

2732:     DMPlexGetPartitioner(dm, &part);
2733:     PetscPartitionerSetFromOptions(part);
2734:     DMPlexDistribute(dm, 0, NULL, &pdm);
2735:     if (pdm) {
2736:       DMPlexReplace_Static(dm, &pdm);
2737:     }
2738:     // Do not auto-distribute again
2739:     DMPlexDistributeSetDefault(dm, PETSC_FALSE);
2740:   }

2742:   DMPlexSetRefinementUniform(dm, PETSC_TRUE);
2743:   for (PetscInt refine=0; refine<refinements; refine++) {
2744:     PetscInt m;
2745:     DM dmf;
2746:     Vec X;
2747:     PetscScalar *x;
2748:     DMRefine(dm, MPI_COMM_NULL, &dmf);
2749:     DMPlexReplace_Static(dm, &dmf);

2751:     DMGetCoordinatesLocal(dm, &X);
2752:     VecGetLocalSize(X, &m);
2753:     VecGetArray(X, &x);
2754:     for (PetscInt i=0; i<m; i+=3) {
2755:       TPSNearestPoint(evalFunc, &x[i]);
2756:     }
2757:     VecRestoreArray(X, &x);
2758:   }

2760:   // Face Sets has already been propagated to new vertices during refinement; this propagates to the initial vertices.
2761:   DMGetLabel(dm, "Face Sets", &label);
2762:   DMPlexLabelComplete(dm, label);

2764:   if (thickness > 0) {
2765:     DM edm,cdm,ecdm;
2766:     DMPlexTransform tr;
2767:     const char *prefix;
2768:     PetscOptions options;
2769:     // Code from DMPlexExtrude
2770:     DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr);
2771:     DMPlexTransformSetDM(tr, dm);
2772:     DMPlexTransformSetType(tr, DMPLEXEXTRUDE);
2773:     PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix);
2774:     PetscObjectSetOptionsPrefix((PetscObject) tr,  prefix);
2775:     PetscObjectGetOptions((PetscObject) dm, &options);
2776:     PetscObjectSetOptions((PetscObject) tr, options);
2777:     DMPlexTransformExtrudeSetLayers(tr, layers);
2778:     DMPlexTransformExtrudeSetThickness(tr, thickness);
2779:     DMPlexTransformExtrudeSetTensor(tr, PETSC_FALSE);
2780:     DMPlexTransformExtrudeSetSymmetric(tr, PETSC_TRUE);
2781:     DMPlexTransformExtrudeSetNormalFunction(tr, normalFunc);
2782:     DMPlexTransformSetFromOptions(tr);
2783:     PetscObjectSetOptions((PetscObject) tr, NULL);
2784:     DMPlexTransformSetUp(tr);
2785:     PetscObjectViewFromOptions((PetscObject) tr, NULL, "-dm_plex_tps_transform_view");
2786:     DMPlexTransformApply(tr, dm, &edm);
2787:     DMCopyDisc(dm, edm);
2788:     DMGetCoordinateDM(dm, &cdm);
2789:     DMGetCoordinateDM(edm, &ecdm);
2790:     DMCopyDisc(cdm, ecdm);
2791:     DMPlexTransformCreateDiscLabels(tr, edm);
2792:     DMPlexTransformDestroy(&tr);
2793:     if (edm) {
2794:       ((DM_Plex *)edm->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
2795:       ((DM_Plex *)edm->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
2796:     }
2797:     DMPlexReplace_Static(dm, &edm);
2798:   }
2799:   return 0;
2800: }

2802: /*@
2803:   DMPlexCreateTPSMesh - Create a distributed, interpolated mesh of a triply-periodic surface

2805:   Collective

2807:   Input Parameters:
2808: + comm   - The communicator for the DM object
2809: . tpstype - Type of triply-periodic surface
2810: . extent - Array of length 3 containing number of periods in each direction
2811: . periodic - array of length 3 with periodicity, or NULL for non-periodic
2812: . tps_distribute - Distribute 2D manifold mesh prior to refinement and extrusion (more scalable)
2813: . refinements - Number of factor-of-2 refinements of 2D manifold mesh
2814: . layers - Number of cell layers extruded in normal direction
2815: - thickness - Thickness in normal direction

2817:   Output Parameter:
2818: . dm  - The DM object

2820:   Notes:
2821:   This meshes the surface of the Schwarz P or Gyroid surfaces.  Schwarz P is is the simplest member of the triply-periodic minimal surfaces.
2822:   https://en.wikipedia.org/wiki/Schwarz_minimal_surface#Schwarz_P_(%22Primitive%22) and can be cut with "clean" boundaries.
2823:   The Gyroid (https://en.wikipedia.org/wiki/Gyroid) is another triply-periodic minimal surface with applications in additive manufacturing; it is much more difficult to "cut" since there are no planes of symmetry.
2824:   Our implementation creates a very coarse mesh of the surface and refines (by 4-way splitting) as many times as requested.
2825:   On each refinement, all vertices are projected to their nearest point on the surface.
2826:   This projection could readily be extended to related surfaces.

2828:   The face (edge) sets for the Schwarz P surface are numbered 1(-x), 2(+x), 3(-y), 4(+y), 5(-z), 6(+z).
2829:   When the mesh is refined, "Face Sets" contain the new vertices (created during refinement).  Use DMPlexLabelComplete() to propagate to coarse-level vertices.

2831:   References:
2832: . * - Maskery et al, Insights into the mechanical properties of several triply periodic minimal surface lattice structures made by polymer additive manufacturing, 2017. https://doi.org/10.1016/j.polymer.2017.11.049

2834:   Developer Notes:
2835:   The Gyroid mesh does not currently mark boundary sets.

2837:   Level: beginner

2839: .seealso: DMPlexCreateSphereMesh(), DMSetType(), DMCreate()
2840: @*/
2841: PetscErrorCode DMPlexCreateTPSMesh(MPI_Comm comm, DMPlexTPSType tpstype, const PetscInt extent[], const DMBoundaryType periodic[], PetscBool tps_distribute, PetscInt refinements, PetscInt layers, PetscReal thickness, DM *dm)
2842: {
2843:   DMCreate(comm, dm);
2844:   DMSetType(*dm, DMPLEX);
2845:   DMPlexCreateTPSMesh_Internal(*dm, tpstype, extent, periodic, tps_distribute, refinements, layers, thickness);
2846:   return 0;
2847: }

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

2852:   Collective

2854:   Input Parameters:
2855: + comm    - The communicator for the DM object
2856: . dim     - The dimension
2857: . simplex - Use simplices, or tensor product cells
2858: - R       - The radius

2860:   Output Parameter:
2861: . dm  - The DM object

2863:   Level: beginner

2865: .seealso: DMPlexCreateBallMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
2866: @*/
2867: PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscReal R, DM *dm)
2868: {
2870:   DMCreate(comm, dm);
2871:   DMSetType(*dm, DMPLEX);
2872:   DMPlexCreateSphereMesh_Internal(*dm, dim, simplex, R);
2873:   return 0;
2874: }

2876: static PetscErrorCode DMPlexCreateBallMesh_Internal(DM dm, PetscInt dim, PetscReal R)
2877: {
2878:   DM             sdm, vol;
2879:   DMLabel        bdlabel;

2881:   DMCreate(PetscObjectComm((PetscObject) dm), &sdm);
2882:   DMSetType(sdm, DMPLEX);
2883:   PetscObjectSetOptionsPrefix((PetscObject) sdm, "bd_");
2884:   DMPlexCreateSphereMesh_Internal(sdm, dim-1, PETSC_TRUE, R);
2885:   DMSetFromOptions(sdm);
2886:   DMViewFromOptions(sdm, NULL, "-dm_view");
2887:   DMPlexGenerate(sdm, NULL, PETSC_TRUE, &vol);
2888:   DMDestroy(&sdm);
2889:   DMPlexReplace_Static(dm, &vol);
2890:   DMCreateLabel(dm, "marker");
2891:   DMGetLabel(dm, "marker", &bdlabel);
2892:   DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, bdlabel);
2893:   DMPlexLabelComplete(dm, bdlabel);
2894:   return 0;
2895: }

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

2900:   Collective

2902:   Input Parameters:
2903: + comm  - The communicator for the DM object
2904: . dim   - The dimension
2905: - R     - The radius

2907:   Output Parameter:
2908: . dm  - The DM object

2910:   Options Database Keys:
2911: - bd_dm_refine - This will refine the surface mesh preserving the sphere geometry

2913:   Level: beginner

2915: .seealso: DMPlexCreateSphereMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
2916: @*/
2917: PetscErrorCode DMPlexCreateBallMesh(MPI_Comm comm, PetscInt dim, PetscReal R, DM *dm)
2918: {
2919:   DMCreate(comm, dm);
2920:   DMSetType(*dm, DMPLEX);
2921:   DMPlexCreateBallMesh_Internal(*dm, dim, R);
2922:   return 0;
2923: }

2925: static PetscErrorCode DMPlexCreateReferenceCell_Internal(DM rdm, DMPolytopeType ct)
2926: {
2927:   switch (ct) {
2928:     case DM_POLYTOPE_POINT:
2929:     {
2930:       PetscInt    numPoints[1]        = {1};
2931:       PetscInt    coneSize[1]         = {0};
2932:       PetscInt    cones[1]            = {0};
2933:       PetscInt    coneOrientations[1] = {0};
2934:       PetscScalar vertexCoords[1]     = {0.0};

2936:       DMSetDimension(rdm, 0);
2937:       DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2938:     }
2939:     break;
2940:     case DM_POLYTOPE_SEGMENT:
2941:     {
2942:       PetscInt    numPoints[2]        = {2, 1};
2943:       PetscInt    coneSize[3]         = {2, 0, 0};
2944:       PetscInt    cones[2]            = {1, 2};
2945:       PetscInt    coneOrientations[2] = {0, 0};
2946:       PetscScalar vertexCoords[2]     = {-1.0,  1.0};

2948:       DMSetDimension(rdm, 1);
2949:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2950:     }
2951:     break;
2952:     case DM_POLYTOPE_POINT_PRISM_TENSOR:
2953:     {
2954:       PetscInt    numPoints[2]        = {2, 1};
2955:       PetscInt    coneSize[3]         = {2, 0, 0};
2956:       PetscInt    cones[2]            = {1, 2};
2957:       PetscInt    coneOrientations[2] = {0, 0};
2958:       PetscScalar vertexCoords[2]     = {-1.0,  1.0};

2960:       DMSetDimension(rdm, 1);
2961:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2962:     }
2963:     break;
2964:     case DM_POLYTOPE_TRIANGLE:
2965:     {
2966:       PetscInt    numPoints[2]        = {3, 1};
2967:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
2968:       PetscInt    cones[3]            = {1, 2, 3};
2969:       PetscInt    coneOrientations[3] = {0, 0, 0};
2970:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

2972:       DMSetDimension(rdm, 2);
2973:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2974:     }
2975:     break;
2976:     case DM_POLYTOPE_QUADRILATERAL:
2977:     {
2978:       PetscInt    numPoints[2]        = {4, 1};
2979:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2980:       PetscInt    cones[4]            = {1, 2, 3, 4};
2981:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2982:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

2984:       DMSetDimension(rdm, 2);
2985:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2986:     }
2987:     break;
2988:     case DM_POLYTOPE_SEG_PRISM_TENSOR:
2989:     {
2990:       PetscInt    numPoints[2]        = {4, 1};
2991:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2992:       PetscInt    cones[4]            = {1, 2, 3, 4};
2993:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2994:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0,  1.0, 1.0};

2996:       DMSetDimension(rdm, 2);
2997:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2998:     }
2999:     break;
3000:     case DM_POLYTOPE_TETRAHEDRON:
3001:     {
3002:       PetscInt    numPoints[2]        = {4, 1};
3003:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3004:       PetscInt    cones[4]            = {1, 2, 3, 4};
3005:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3006:       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};

3008:       DMSetDimension(rdm, 3);
3009:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3010:     }
3011:     break;
3012:     case DM_POLYTOPE_HEXAHEDRON:
3013:     {
3014:       PetscInt    numPoints[2]        = {8, 1};
3015:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3016:       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
3017:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3018:       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,
3019:                                          -1.0, -1.0,  1.0,   1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0,  1.0,  1.0};

3021:       DMSetDimension(rdm, 3);
3022:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3023:     }
3024:     break;
3025:     case DM_POLYTOPE_TRI_PRISM:
3026:     {
3027:       PetscInt    numPoints[2]        = {6, 1};
3028:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3029:       PetscInt    cones[6]            = {1, 2, 3, 4, 5, 6};
3030:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3031:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0, -1.0,  1.0, -1.0,   1.0, -1.0, -1.0,
3032:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0,  1.0,  1.0};

3034:       DMSetDimension(rdm, 3);
3035:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3036:     }
3037:     break;
3038:     case DM_POLYTOPE_TRI_PRISM_TENSOR:
3039:     {
3040:       PetscInt    numPoints[2]        = {6, 1};
3041:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3042:       PetscInt    cones[6]            = {1, 2, 3, 4, 5, 6};
3043:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3044:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3045:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3047:       DMSetDimension(rdm, 3);
3048:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3049:     }
3050:     break;
3051:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
3052:     {
3053:       PetscInt    numPoints[2]        = {8, 1};
3054:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3055:       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
3056:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3057:       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,
3058:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3060:       DMSetDimension(rdm, 3);
3061:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3062:     }
3063:     break;
3064:     case DM_POLYTOPE_PYRAMID:
3065:     {
3066:       PetscInt    numPoints[2]        = {5, 1};
3067:       PetscInt    coneSize[6]         = {5, 0, 0, 0, 0, 0};
3068:       PetscInt    cones[5]            = {1, 2, 3, 4, 5};
3069:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3070:       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,
3071:                                           0.0,  0.0,  1.0};

3073:       DMSetDimension(rdm, 3);
3074:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3075:     }
3076:     break;
3077:     default: SETERRQ(PetscObjectComm((PetscObject) rdm), PETSC_ERR_ARG_WRONG, "Cannot create reference cell for cell type %s", DMPolytopeTypes[ct]);
3078:   }
3079:   {
3080:     PetscInt Nv, v;

3082:     /* Must create the celltype label here so that we do not automatically try to compute the types */
3083:     DMCreateLabel(rdm, "celltype");
3084:     DMPlexSetCellType(rdm, 0, ct);
3085:     DMPlexGetChart(rdm, NULL, &Nv);
3086:     for (v = 1; v < Nv; ++v) DMPlexSetCellType(rdm, v, DM_POLYTOPE_POINT);
3087:   }
3088:   DMPlexInterpolateInPlace_Internal(rdm);
3089:   PetscObjectSetName((PetscObject) rdm, DMPolytopeTypes[ct]);
3090:   return 0;
3091: }

3093: /*@
3094:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

3096:   Collective

3098:   Input Parameters:
3099: + comm - The communicator
3100: - ct   - The cell type of the reference cell

3102:   Output Parameter:
3103: . refdm - The reference cell

3105:   Level: intermediate

3107: .seealso: DMPlexCreateReferenceCell(), DMPlexCreateBoxMesh()
3108: @*/
3109: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, DMPolytopeType ct, DM *refdm)
3110: {
3111:   DMCreate(comm, refdm);
3112:   DMSetType(*refdm, DMPLEX);
3113:   DMPlexCreateReferenceCell_Internal(*refdm, ct);
3114:   return 0;
3115: }

3117: static PetscErrorCode DMPlexCreateBoundaryLabel_Private(DM dm, const char name[])
3118: {
3119:   DM             plex;
3120:   DMLabel        label;
3121:   PetscBool      hasLabel;

3124:   DMHasLabel(dm, name, &hasLabel);
3125:   if (hasLabel) return 0;
3126:   DMCreateLabel(dm, name);
3127:   DMGetLabel(dm, name, &label);
3128:   DMConvert(dm, DMPLEX, &plex);
3129:   DMPlexMarkBoundaryFaces(plex, 1, label);
3130:   DMDestroy(&plex);
3131:   return 0;
3132: }

3134: const char * const DMPlexShapes[] = {"box", "box_surface", "ball", "sphere", "cylinder", "schwarz_p", "gyroid", "doublet", "unknown", "DMPlexShape", "DM_SHAPE_", NULL};

3136: static PetscErrorCode DMPlexCreateFromOptions_Internal(PetscOptionItems *PetscOptionsObject, PetscBool *useCoordSpace, DM dm)
3137: {
3138:   DMPlexShape    shape = DM_SHAPE_BOX;
3139:   DMPolytopeType cell  = DM_POLYTOPE_TRIANGLE;
3140:   PetscInt       dim   = 2;
3141:   PetscBool      simplex = PETSC_TRUE, interpolate = PETSC_TRUE, adjCone = PETSC_FALSE, adjClosure = PETSC_TRUE, refDomain = PETSC_FALSE;
3142:   PetscBool      flg, flg2, fflg, bdfflg, nameflg;
3143:   MPI_Comm       comm;
3144:   char           filename[PETSC_MAX_PATH_LEN]   = "<unspecified>";
3145:   char           bdFilename[PETSC_MAX_PATH_LEN] = "<unspecified>";
3146:   char           plexname[PETSC_MAX_PATH_LEN]   = "";

3148:   PetscObjectGetComm((PetscObject) dm, &comm);
3149:   /* TODO Turn this into a registration interface */
3150:   PetscOptionsString("-dm_plex_filename", "File containing a mesh", "DMPlexCreateFromFile", filename, filename, sizeof(filename), &fflg);
3151:   PetscOptionsString("-dm_plex_boundary_filename", "File containing a mesh boundary", "DMPlexCreateFromFile", bdFilename, bdFilename, sizeof(bdFilename), &bdfflg);
3152:   PetscOptionsString("-dm_plex_name", "Name of the mesh in the file", "DMPlexCreateFromFile", plexname, plexname, sizeof(plexname), &nameflg);
3153:   PetscOptionsEnum("-dm_plex_cell", "Cell shape", "", DMPolytopeTypes, (PetscEnum) cell, (PetscEnum *) &cell, NULL);
3154:   PetscOptionsBool("-dm_plex_reference_cell_domain", "Use a reference cell domain", "", refDomain, &refDomain, NULL);
3155:   PetscOptionsEnum("-dm_plex_shape", "Shape for built-in mesh", "", DMPlexShapes, (PetscEnum) shape, (PetscEnum *) &shape, &flg);
3156:   PetscOptionsBoundedInt("-dm_plex_dim", "Topological dimension of the mesh", "DMGetDimension", dim, &dim, &flg, 0);
3158:   PetscOptionsBool("-dm_plex_simplex", "Mesh cell shape", "", simplex,  &simplex, &flg);
3159:   PetscOptionsBool("-dm_plex_interpolate", "Flag to create edges and faces automatically", "", interpolate, &interpolate, &flg);
3160:   PetscOptionsBool("-dm_plex_adj_cone", "Set adjacency direction", "DMSetBasicAdjacency", adjCone,  &adjCone, &flg);
3161:   PetscOptionsBool("-dm_plex_adj_closure", "Set adjacency size", "DMSetBasicAdjacency", adjClosure,  &adjClosure, &flg2);
3162:   if (flg || flg2) DMSetBasicAdjacency(dm, adjCone, adjClosure);

3164:   switch (cell) {
3165:     case DM_POLYTOPE_POINT:
3166:     case DM_POLYTOPE_SEGMENT:
3167:     case DM_POLYTOPE_POINT_PRISM_TENSOR:
3168:     case DM_POLYTOPE_TRIANGLE:
3169:     case DM_POLYTOPE_QUADRILATERAL:
3170:     case DM_POLYTOPE_TETRAHEDRON:
3171:     case DM_POLYTOPE_HEXAHEDRON:
3172:       *useCoordSpace = PETSC_TRUE;break;
3173:     default: *useCoordSpace = PETSC_FALSE;break;
3174:   }

3176:   if (fflg) {
3177:     DM dmnew;

3179:     DMPlexCreateFromFile(PetscObjectComm((PetscObject) dm), filename, plexname, interpolate, &dmnew);
3180:     DMPlexCopy_Internal(dm, PETSC_FALSE, dmnew);
3181:     DMPlexReplace_Static(dm, &dmnew);
3182:   } else if (refDomain) {
3183:     DMPlexCreateReferenceCell_Internal(dm, cell);
3184:   } else if (bdfflg) {
3185:     DM bdm, dmnew;

3187:     DMPlexCreateFromFile(PetscObjectComm((PetscObject) dm), bdFilename, plexname, interpolate, &bdm);
3188:     PetscObjectSetOptionsPrefix((PetscObject) bdm, "bd_");
3189:     DMSetFromOptions(bdm);
3190:     DMPlexGenerate(bdm, NULL, interpolate, &dmnew);
3191:     DMDestroy(&bdm);
3192:     DMPlexCopy_Internal(dm, PETSC_FALSE, dmnew);
3193:     DMPlexReplace_Static(dm, &dmnew);
3194:   } else {
3195:     PetscObjectSetName((PetscObject) dm, DMPlexShapes[shape]);
3196:     switch (shape) {
3197:       case DM_SHAPE_BOX:
3198:       {
3199:         PetscInt       faces[3] = {0, 0, 0};
3200:         PetscReal      lower[3] = {0, 0, 0};
3201:         PetscReal      upper[3] = {1, 1, 1};
3202:         DMBoundaryType bdt[3]   = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
3203:         PetscInt       i, n;

3205:         n    = dim;
3206:         for (i = 0; i < dim; ++i) faces[i] = (dim == 1 ? 1 : 4-dim);
3207:         PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", faces, &n, &flg);
3208:         n    = 3;
3209:         PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg);
3211:         n    = 3;
3212:         PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg);
3214:         n    = 3;
3215:         PetscOptionsEnumArray("-dm_plex_box_bd", "Boundary type for each dimension", "", DMBoundaryTypes, (PetscEnum *) bdt, &n, &flg);
3217:         switch (cell) {
3218:           case DM_POLYTOPE_TRI_PRISM_TENSOR:
3219:             DMPlexCreateWedgeBoxMesh_Internal(dm, faces, lower, upper, bdt);
3220:             if (!interpolate) {
3221:               DM udm;

3223:               DMPlexUninterpolate(dm, &udm);
3224:               DMPlexReplace_Static(dm, &udm);
3225:             }
3226:             break;
3227:           default:
3228:             DMPlexCreateBoxMesh_Internal(dm, dim, simplex, faces, lower, upper, bdt, interpolate);
3229:             break;
3230:         }
3231:       }
3232:       break;
3233:       case DM_SHAPE_BOX_SURFACE:
3234:       {
3235:         PetscInt  faces[3] = {0, 0, 0};
3236:         PetscReal lower[3] = {0, 0, 0};
3237:         PetscReal upper[3] = {1, 1, 1};
3238:         PetscInt  i, n;

3240:         n    = dim+1;
3241:         for (i = 0; i < dim+1; ++i) faces[i] = (dim+1 == 1 ? 1 : 4-(dim+1));
3242:         PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", faces, &n, &flg);
3243:         n    = 3;
3244:         PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg);
3246:         n    = 3;
3247:         PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg);
3249:         DMPlexCreateBoxSurfaceMesh_Internal(dm, dim+1, faces, lower, upper, interpolate);
3250:       }
3251:       break;
3252:       case DM_SHAPE_SPHERE:
3253:       {
3254:         PetscReal R = 1.0;

3256:         PetscOptionsReal("-dm_plex_sphere_radius", "Radius of the sphere", "", R,  &R, &flg);
3257:         DMPlexCreateSphereMesh_Internal(dm, dim, simplex, R);
3258:       }
3259:       break;
3260:       case DM_SHAPE_BALL:
3261:       {
3262:         PetscReal R = 1.0;

3264:         PetscOptionsReal("-dm_plex_ball_radius", "Radius of the ball", "", R,  &R, &flg);
3265:         DMPlexCreateBallMesh_Internal(dm, dim, R);
3266:       }
3267:       break;
3268:       case DM_SHAPE_CYLINDER:
3269:       {
3270:         DMBoundaryType bdt = DM_BOUNDARY_NONE;
3271:         PetscInt       Nw  = 6;

3273:         PetscOptionsEnum("-dm_plex_cylinder_bd", "Boundary type in the z direction", "", DMBoundaryTypes, (PetscEnum) bdt, (PetscEnum *) &bdt, NULL);
3274:         PetscOptionsInt("-dm_plex_cylinder_num_wedges", "Number of wedges around the cylinder", "", Nw, &Nw, NULL);
3275:         switch (cell) {
3276:           case DM_POLYTOPE_TRI_PRISM_TENSOR:
3277:             DMPlexCreateWedgeCylinderMesh_Internal(dm, Nw, interpolate);
3278:             break;
3279:           default:
3280:             DMPlexCreateHexCylinderMesh_Internal(dm, bdt);
3281:             break;
3282:         }
3283:       }
3284:       break;
3285:       case DM_SHAPE_SCHWARZ_P: // fallthrough
3286:       case DM_SHAPE_GYROID:
3287:       {
3288:         PetscInt       extent[3] = {1,1,1}, refine = 0, layers = 0, three;
3289:         PetscReal      thickness = 0.;
3290:         DMBoundaryType periodic[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
3291:         DMPlexTPSType  tps_type = shape == DM_SHAPE_SCHWARZ_P ? DMPLEX_TPS_SCHWARZ_P : DMPLEX_TPS_GYROID;
3292:         PetscBool      tps_distribute;
3293:         PetscOptionsIntArray("-dm_plex_tps_extent", "Number of replicas for each of three dimensions", NULL, extent, (three=3, &three), NULL);
3294:         PetscOptionsInt("-dm_plex_tps_refine", "Number of refinements", NULL, refine, &refine, NULL);
3295:         PetscOptionsEnumArray("-dm_plex_tps_periodic", "Periodicity in each of three dimensions", NULL, DMBoundaryTypes, (PetscEnum*)periodic, (three=3, &three), NULL);
3296:         PetscOptionsInt("-dm_plex_tps_layers", "Number of layers in volumetric extrusion (or zero to not extrude)", NULL, layers, &layers, NULL);
3297:         PetscOptionsReal("-dm_plex_tps_thickness", "Thickness of volumetric extrusion", NULL, thickness, &thickness, NULL);
3298:         DMPlexDistributeGetDefault(dm, &tps_distribute);
3299:         PetscOptionsBool("-dm_plex_tps_distribute", "Distribute the 2D mesh prior to refinement and extrusion", NULL, tps_distribute, &tps_distribute, NULL);
3300:         DMPlexCreateTPSMesh_Internal(dm, tps_type, extent, periodic, tps_distribute, refine, layers, thickness);
3301:       }
3302:       break;
3303:       case DM_SHAPE_DOUBLET:
3304:       {
3305:         DM        dmnew;
3306:         PetscReal rl = 0.0;

3308:         PetscOptionsReal("-dm_plex_doublet_refinementlimit", "Refinement limit", NULL, rl, &rl, NULL);
3309:         DMPlexCreateDoublet(PetscObjectComm((PetscObject)dm), dim, simplex, interpolate, rl, &dmnew);
3310:         DMPlexCopy_Internal(dm, PETSC_FALSE, dmnew);
3311:         DMPlexReplace_Static(dm, &dmnew);
3312:       }
3313:       break;
3314:       default: SETERRQ(comm, PETSC_ERR_SUP, "Domain shape %s is unsupported", DMPlexShapes[shape]);
3315:     }
3316:   }
3317:   DMPlexSetRefinementUniform(dm, PETSC_TRUE);
3318:   if (!((PetscObject)dm)->name && nameflg) {
3319:     PetscObjectSetName((PetscObject)dm, plexname);
3320:   }
3321:   return 0;
3322: }

3324: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject, DM dm)
3325: {
3326:   DM_Plex       *mesh = (DM_Plex*) dm->data;
3327:   PetscBool      flg;
3328:   char           bdLabel[PETSC_MAX_PATH_LEN];

3330:   /* Handle viewing */
3331:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL);
3332:   PetscOptionsBoundedInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL,0);
3333:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL);
3334:   PetscOptionsBoundedInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL,0);
3335:   DMMonitorSetFromOptions(dm, "-dm_plex_monitor_throughput", "Monitor the simulation throughput", "DMPlexMonitorThroughput", DMPlexMonitorThroughput, NULL, &flg);
3336:   if (flg) PetscLogDefaultBegin();
3337:   /* Labeling */
3338:   PetscOptionsString("-dm_plex_boundary_label", "Label to mark the mesh boundary", "", bdLabel, bdLabel, sizeof(bdLabel), &flg);
3339:   if (flg) DMPlexCreateBoundaryLabel_Private(dm, bdLabel);
3340:   /* Point Location */
3341:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL);
3342:   /* Partitioning and distribution */
3343:   PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL);
3344:   /* Generation and remeshing */
3345:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL);
3346:   /* Projection behavior */
3347:   PetscOptionsBoundedInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL,0);
3348:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);
3349:   /* Checking structure */
3350:   {
3351:     PetscBool   flg = PETSC_FALSE, flg2 = PETSC_FALSE, all = PETSC_FALSE;

3353:     PetscOptionsBool("-dm_plex_check_all", "Perform all checks", NULL, PETSC_FALSE, &all, &flg2);
3354:     PetscOptionsBool("-dm_plex_check_symmetry", "Check that the adjacency information in the mesh is symmetric", "DMPlexCheckSymmetry", PETSC_FALSE, &flg, &flg2);
3355:     if (all || (flg && flg2)) DMPlexCheckSymmetry(dm);
3356:     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);
3357:     if (all || (flg && flg2)) DMPlexCheckSkeleton(dm, 0);
3358:     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);
3359:     if (all || (flg && flg2)) DMPlexCheckFaces(dm, 0);
3360:     PetscOptionsBool("-dm_plex_check_geometry", "Check that cells have positive volume", "DMPlexCheckGeometry", PETSC_FALSE, &flg, &flg2);
3361:     if (all || (flg && flg2)) DMPlexCheckGeometry(dm);
3362:     PetscOptionsBool("-dm_plex_check_pointsf", "Check some necessary conditions for PointSF", "DMPlexCheckPointSF", PETSC_FALSE, &flg, &flg2);
3363:     if (all || (flg && flg2)) DMPlexCheckPointSF(dm);
3364:     PetscOptionsBool("-dm_plex_check_interface_cones", "Check points on inter-partition interfaces have conforming order of cone points", "DMPlexCheckInterfaceCones", PETSC_FALSE, &flg, &flg2);
3365:     if (all || (flg && flg2)) DMPlexCheckInterfaceCones(dm);
3366:     PetscOptionsBool("-dm_plex_check_cell_shape", "Check cell shape", "DMPlexCheckCellShape", PETSC_FALSE, &flg, &flg2);
3367:     if (flg && flg2) DMPlexCheckCellShape(dm, PETSC_TRUE, PETSC_DETERMINE);
3368:   }
3369:   {
3370:     PetscReal scale = 1.0;

3372:     PetscOptionsReal("-dm_plex_scale", "Scale factor for mesh coordinates", "DMPlexScale", scale, &scale, &flg);
3373:     if (flg) {
3374:       Vec coordinates, coordinatesLocal;

3376:       DMGetCoordinates(dm, &coordinates);
3377:       DMGetCoordinatesLocal(dm, &coordinatesLocal);
3378:       VecScale(coordinates, scale);
3379:       VecScale(coordinatesLocal, scale);
3380:     }
3381:   }
3382:   PetscPartitionerSetFromOptions(mesh->partitioner);
3383:   return 0;
3384: }

3386: static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
3387: {
3388:   PetscFunctionList ordlist;
3389:   char              oname[256];
3390:   PetscReal         volume = -1.0;
3391:   PetscInt          prerefine = 0, refine = 0, r, coarsen = 0, overlap = 0, extLayers = 0, dim;
3392:   PetscBool         uniformOrig, created = PETSC_FALSE, uniform = PETSC_TRUE, distribute, interpolate = PETSC_TRUE, coordSpace = PETSC_TRUE, remap = PETSC_TRUE, ghostCells = PETSC_FALSE, isHierarchy, ignoreModel = PETSC_FALSE, flg;

3395:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
3396:   /* Handle automatic creation */
3397:   DMGetDimension(dm, &dim);
3398:   if (dim < 0) {DMPlexCreateFromOptions_Internal(PetscOptionsObject, &coordSpace, dm);created = PETSC_TRUE;}
3399:   /* Handle interpolation before distribution */
3400:   PetscOptionsBool("-dm_plex_interpolate_pre", "Flag to interpolate mesh before distribution", "", interpolate, &interpolate, &flg);
3401:   if (flg) {
3402:     DMPlexInterpolatedFlag interpolated;

3404:     DMPlexIsInterpolated(dm, &interpolated);
3405:     if (interpolated == DMPLEX_INTERPOLATED_FULL && !interpolate) {
3406:       DM udm;

3408:       DMPlexUninterpolate(dm, &udm);
3409:       DMPlexReplace_Static(dm, &udm);
3410:     } else if (interpolated != DMPLEX_INTERPOLATED_FULL && interpolate) {
3411:       DM idm;

3413:       DMPlexInterpolate(dm, &idm);
3414:       DMPlexReplace_Static(dm, &idm);
3415:     }
3416:   }
3417:   /* Handle DMPlex refinement before distribution */
3418:   PetscOptionsBool("-dm_refine_ignore_model", "Flag to ignore the geometry model when refining", "DMCreate", ignoreModel, &ignoreModel, &flg);
3419:   if (flg) {((DM_Plex *) dm->data)->ignoreModel = ignoreModel;}
3420:   DMPlexGetRefinementUniform(dm, &uniformOrig);
3421:   PetscOptionsBoundedInt("-dm_refine_pre", "The number of refinements before distribution", "DMCreate", prerefine, &prerefine, NULL,0);
3422:   PetscOptionsBool("-dm_refine_remap_pre", "Flag to control coordinate remapping", "DMCreate", remap, &remap, NULL);
3423:   PetscOptionsBool("-dm_refine_uniform_pre", "Flag for uniform refinement before distribution", "DMCreate", uniform, &uniform, &flg);
3424:   if (flg) DMPlexSetRefinementUniform(dm, uniform);
3425:   PetscOptionsReal("-dm_refine_volume_limit_pre", "The maximum cell volume after refinement before distribution", "DMCreate", volume, &volume, &flg);
3426:   if (flg) {
3427:     DMPlexSetRefinementUniform(dm, PETSC_FALSE);
3428:     DMPlexSetRefinementLimit(dm, volume);
3429:     prerefine = PetscMax(prerefine, 1);
3430:   }
3431:   for (r = 0; r < prerefine; ++r) {
3432:     DM             rdm;
3433:     PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

3435:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3436:     DMRefine(dm, PetscObjectComm((PetscObject) dm), &rdm);
3437:     DMPlexReplace_Static(dm, &rdm);
3438:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3439:     if (coordFunc && remap) {
3440:       DMPlexRemapGeometry(dm, 0.0, coordFunc);
3441:       ((DM_Plex*) dm->data)->coordFunc = coordFunc;
3442:     }
3443:   }
3444:   DMPlexSetRefinementUniform(dm, uniformOrig);
3445:   /* Handle DMPlex extrusion before distribution */
3446:   PetscOptionsBoundedInt("-dm_extrude", "The number of layers to extrude", "", extLayers, &extLayers, NULL, 0);
3447:   if (extLayers) {
3448:     DM edm;

3450:     DMExtrude(dm, extLayers, &edm);
3451:     DMPlexReplace_Static(dm, &edm);
3452:     ((DM_Plex *) dm->data)->coordFunc = NULL;
3453:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3454:     extLayers = 0;
3455:   }
3456:   /* Handle DMPlex reordering before distribution */
3457:   MatGetOrderingList(&ordlist);
3458:   PetscOptionsFList("-dm_plex_reorder", "Set mesh reordering type", "DMPlexGetOrdering", ordlist, MATORDERINGNATURAL, oname, sizeof(oname), &flg);
3459:   if (flg) {
3460:     DM pdm;
3461:     IS perm;

3463:     DMPlexGetOrdering(dm, oname, NULL, &perm);
3464:     DMPlexPermute(dm, perm, &pdm);
3465:     ISDestroy(&perm);
3466:     DMPlexReplace_Static(dm, &pdm);
3467:     DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3468:   }
3469:   /* Handle DMPlex distribution */
3470:   DMPlexDistributeGetDefault(dm, &distribute);
3471:   PetscOptionsBool("-dm_distribute", "Flag to redistribute a mesh among processes", "DMCreate", distribute, &distribute, NULL);
3472:   PetscOptionsBoundedInt("-dm_distribute_overlap", "The size of the overlap halo", "DMCreate", overlap, &overlap, NULL, 0);
3473:   if (distribute) {
3474:     DM               pdm = NULL;
3475:     PetscPartitioner part;

3477:     DMPlexGetPartitioner(dm, &part);
3478:     PetscPartitionerSetFromOptions(part);
3479:     DMPlexDistribute(dm, overlap, NULL, &pdm);
3480:     if (pdm) {
3481:       DMPlexReplace_Static(dm, &pdm);
3482:     }
3483:   }
3484:   /* Create coordinate space */
3485:   if (created) {
3486:     DM_Plex  *mesh = (DM_Plex *) dm->data;
3487:     PetscInt  degree = 1;
3488:     PetscBool periodic, flg;

3490:     PetscOptionsBool("-dm_coord_space", "Use an FEM space for coordinates", "", coordSpace, &coordSpace, &flg);
3491:     PetscOptionsInt("-dm_coord_petscspace_degree", "FEM degree for coordinate space", "", degree, &degree, NULL);
3492:     if (coordSpace) DMPlexCreateCoordinateSpace(dm, degree, mesh->coordFunc);
3493:     if (flg && !coordSpace) {
3494:       DM           cdm;
3495:       PetscDS      cds;
3496:       PetscObject  obj;
3497:       PetscClassId id;

3499:       DMGetCoordinateDM(dm, &cdm);
3500:       DMGetDS(cdm, &cds);
3501:       PetscDSGetDiscretization(cds, 0, &obj);
3502:       PetscObjectGetClassId(obj, &id);
3503:       if (id == PETSCFE_CLASSID) {
3504:         PetscContainer dummy;

3506:         PetscContainerCreate(PETSC_COMM_SELF, &dummy);
3507:         PetscObjectSetName((PetscObject) dummy, "coordinates");
3508:         DMSetField(cdm, 0, NULL, (PetscObject) dummy);
3509:         PetscContainerDestroy(&dummy);
3510:         DMClearDS(cdm);
3511:       }
3512:       mesh->coordFunc = NULL;
3513:     }
3514:     DMLocalizeCoordinates(dm);
3515:     DMGetPeriodicity(dm, &periodic, NULL, NULL, NULL);
3516:     if (periodic) DMSetPeriodicity(dm, PETSC_TRUE, NULL, NULL, NULL);
3517:   }
3518:   /* Handle DMPlex refinement */
3519:   remap = PETSC_TRUE;
3520:   PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL,0);
3521:   PetscOptionsBool("-dm_refine_remap", "Flag to control coordinate remapping", "DMCreate", remap, &remap, NULL);
3522:   PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy,0);
3523:   if (refine) DMPlexSetRefinementUniform(dm, PETSC_TRUE);
3524:   if (refine && isHierarchy) {
3525:     DM *dms, coarseDM;

3527:     DMGetCoarseDM(dm, &coarseDM);
3528:     PetscObjectReference((PetscObject)coarseDM);
3529:     PetscMalloc1(refine,&dms);
3530:     DMRefineHierarchy(dm, refine, dms);
3531:     /* Total hack since we do not pass in a pointer */
3532:     DMPlexSwap_Static(dm, dms[refine-1]);
3533:     if (refine == 1) {
3534:       DMSetCoarseDM(dm, dms[0]);
3535:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
3536:     } else {
3537:       DMSetCoarseDM(dm, dms[refine-2]);
3538:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
3539:       DMSetCoarseDM(dms[0], dms[refine-1]);
3540:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
3541:     }
3542:     DMSetCoarseDM(dms[refine-1], coarseDM);
3543:     PetscObjectDereference((PetscObject)coarseDM);
3544:     /* Free DMs */
3545:     for (r = 0; r < refine; ++r) {
3546:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
3547:       DMDestroy(&dms[r]);
3548:     }
3549:     PetscFree(dms);
3550:   } else {
3551:     for (r = 0; r < refine; ++r) {
3552:       DM             rdm;
3553:       PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

3555:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3556:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &rdm);
3557:       /* Total hack since we do not pass in a pointer */
3558:       DMPlexReplace_Static(dm, &rdm);
3559:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3560:       if (coordFunc && remap) {
3561:         DMPlexRemapGeometry(dm, 0.0, coordFunc);
3562:         ((DM_Plex*) dm->data)->coordFunc = coordFunc;
3563:       }
3564:     }
3565:   }
3566:   /* Handle DMPlex coarsening */
3567:   PetscOptionsBoundedInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL,0);
3568:   PetscOptionsBoundedInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy,0);
3569:   if (coarsen && isHierarchy) {
3570:     DM *dms;

3572:     PetscMalloc1(coarsen, &dms);
3573:     DMCoarsenHierarchy(dm, coarsen, dms);
3574:     /* Free DMs */
3575:     for (r = 0; r < coarsen; ++r) {
3576:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
3577:       DMDestroy(&dms[r]);
3578:     }
3579:     PetscFree(dms);
3580:   } else {
3581:     for (r = 0; r < coarsen; ++r) {
3582:       DM             cdm;
3583:       PetscPointFunc coordFunc = ((DM_Plex*) dm->data)->coordFunc;

3585:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3586:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &cdm);
3587:       /* Total hack since we do not pass in a pointer */
3588:       DMPlexReplace_Static(dm, &cdm);
3589:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3590:       if (coordFunc) {
3591:         DMPlexRemapGeometry(dm, 0.0, coordFunc);
3592:         ((DM_Plex*) dm->data)->coordFunc = coordFunc;
3593:       }
3594:     }
3595:   }
3596:   /* Handle ghost cells */
3597:   PetscOptionsBool("-dm_plex_create_fv_ghost_cells", "Flag to create finite volume ghost cells on the boundary", "DMCreate", ghostCells, &ghostCells, NULL);
3598:   if (ghostCells) {
3599:     DM   gdm;
3600:     char lname[PETSC_MAX_PATH_LEN];

3602:     lname[0] = '\0';
3603:     PetscOptionsString("-dm_plex_fv_ghost_cells_label", "Label name for ghost cells boundary", "DMCreate", lname, lname, sizeof(lname), &flg);
3604:     DMPlexConstructGhostCells(dm, flg ? lname : NULL, NULL, &gdm);
3605:     DMPlexReplace_Static(dm, &gdm);
3606:   }
3607:   /* Handle 1D order */
3608:   {
3609:     DM           cdm, rdm;
3610:     PetscDS      cds;
3611:     PetscObject  obj;
3612:     PetscClassId id = PETSC_OBJECT_CLASSID;
3613:     IS           perm;
3614:     PetscInt     dim, Nf;
3615:     PetscBool    distributed;

3617:     DMGetDimension(dm, &dim);
3618:     DMPlexIsDistributed(dm, &distributed);
3619:     DMGetCoordinateDM(dm, &cdm);
3620:     DMGetDS(cdm, &cds);
3621:     PetscDSGetNumFields(cds, &Nf);
3622:     if (Nf) {
3623:       PetscDSGetDiscretization(cds, 0, &obj);
3624:       PetscObjectGetClassId(obj, &id);
3625:     }
3626:     if (dim == 1 && !distributed && id != PETSCFE_CLASSID) {
3627:       DMPlexGetOrdering1D(dm, &perm);
3628:       DMPlexPermute(dm, perm, &rdm);
3629:       DMPlexReplace_Static(dm, &rdm);
3630:       ISDestroy(&perm);
3631:     }
3632:   }
3633:   /* Handle */
3634:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
3635:   PetscOptionsTail();
3636:   return 0;
3637: }

3639: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
3640: {
3641:   DMCreateGlobalVector_Section_Private(dm,vec);
3642:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
3643:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
3644:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
3645:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
3646:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
3647:   return 0;
3648: }

3650: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
3651: {
3652:   DMCreateLocalVector_Section_Private(dm,vec);
3653:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
3654:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
3655:   return 0;
3656: }

3658: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
3659: {
3660:   PetscInt       depth, d;

3662:   DMPlexGetDepth(dm, &depth);
3663:   if (depth == 1) {
3664:     DMGetDimension(dm, &d);
3665:     if (dim == 0)      DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
3666:     else if (dim == d) DMPlexGetDepthStratum(dm, 1, pStart, pEnd);
3667:     else               {*pStart = 0; *pEnd = 0;}
3668:   } else {
3669:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
3670:   }
3671:   return 0;
3672: }

3674: static PetscErrorCode DMGetNeighbors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
3675: {
3676:   PetscSF           sf;
3677:   PetscInt          niranks, njranks, n;
3678:   const PetscMPIInt *iranks, *jranks;
3679:   DM_Plex           *data = (DM_Plex*) dm->data;

3681:   DMGetPointSF(dm, &sf);
3682:   if (!data->neighbors) {
3683:     PetscSFSetUp(sf);
3684:     PetscSFGetRootRanks(sf, &njranks, &jranks, NULL, NULL, NULL);
3685:     PetscSFGetLeafRanks(sf, &niranks, &iranks, NULL, NULL);
3686:     PetscMalloc1(njranks + niranks + 1, &data->neighbors);
3687:     PetscArraycpy(data->neighbors + 1, jranks, njranks);
3688:     PetscArraycpy(data->neighbors + njranks + 1, iranks, niranks);
3689:     n = njranks + niranks;
3690:     PetscSortRemoveDupsMPIInt(&n, data->neighbors + 1);
3691:     /* The following cast should never fail: can't have more neighbors than PETSC_MPI_INT_MAX */
3692:     PetscMPIIntCast(n, data->neighbors);
3693:   }
3694:   if (nranks) *nranks = data->neighbors[0];
3695:   if (ranks) {
3696:     if (data->neighbors[0]) *ranks = data->neighbors + 1;
3697:     else                    *ranks = NULL;
3698:   }
3699:   return 0;
3700: }

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

3704: static PetscErrorCode DMInitialize_Plex(DM dm)
3705: {
3706:   dm->ops->view                            = DMView_Plex;
3707:   dm->ops->load                            = DMLoad_Plex;
3708:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
3709:   dm->ops->clone                           = DMClone_Plex;
3710:   dm->ops->setup                           = DMSetUp_Plex;
3711:   dm->ops->createlocalsection              = DMCreateLocalSection_Plex;
3712:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
3713:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
3714:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
3715:   dm->ops->getlocaltoglobalmapping         = NULL;
3716:   dm->ops->createfieldis                   = NULL;
3717:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
3718:   dm->ops->createcoordinatefield           = DMCreateCoordinateField_Plex;
3719:   dm->ops->getcoloring                     = NULL;
3720:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
3721:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
3722:   dm->ops->createmassmatrix                = DMCreateMassMatrix_Plex;
3723:   dm->ops->createmassmatrixlumped          = DMCreateMassMatrixLumped_Plex;
3724:   dm->ops->createinjection                 = DMCreateInjection_Plex;
3725:   dm->ops->refine                          = DMRefine_Plex;
3726:   dm->ops->coarsen                         = DMCoarsen_Plex;
3727:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
3728:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
3729:   dm->ops->extrude                         = DMExtrude_Plex;
3730:   dm->ops->globaltolocalbegin              = NULL;
3731:   dm->ops->globaltolocalend                = NULL;
3732:   dm->ops->localtoglobalbegin              = NULL;
3733:   dm->ops->localtoglobalend                = NULL;
3734:   dm->ops->destroy                         = DMDestroy_Plex;
3735:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
3736:   dm->ops->createsuperdm                   = DMCreateSuperDM_Plex;
3737:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
3738:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
3739:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
3740:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
3741:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
3742:   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
3743:   dm->ops->projectbdfieldlabellocal        = DMProjectBdFieldLabelLocal_Plex;
3744:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
3745:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
3746:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
3747:   dm->ops->getneighbors                    = DMGetNeighbors_Plex;
3748:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);
3749:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertTimeDerviativeBoundaryValues_C",DMPlexInsertTimeDerivativeBoundaryValues_Plex);
3750:   PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);
3751:   PetscObjectComposeFunction((PetscObject)dm,"DMCreateNeumannOverlap_C",DMCreateNeumannOverlap_Plex);
3752:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexGetOverlap_C",DMPlexGetOverlap_Plex);
3753:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexDistributeGetDefault_C",DMPlexDistributeGetDefault_Plex);
3754:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexDistributeSetDefault_C",DMPlexDistributeSetDefault_Plex);
3755:   PetscObjectComposeFunction((PetscObject)dm,"DMInterpolateSolution_C",DMInterpolateSolution_Plex);
3756:   return 0;
3757: }

3759: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
3760: {
3761:   DM_Plex        *mesh = (DM_Plex *) dm->data;

3763:   mesh->refct++;
3764:   (*newdm)->data = mesh;
3765:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
3766:   DMInitialize_Plex(*newdm);
3767:   return 0;
3768: }

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

3776:   Options Database Keys:
3777: + -dm_refine_pre                     - Refine mesh before distribution
3778: + -dm_refine_uniform_pre             - Choose uniform or generator-based refinement
3779: + -dm_refine_volume_limit_pre        - Cell volume limit after pre-refinement using generator
3780: . -dm_distribute                     - Distribute mesh across processes
3781: . -dm_distribute_overlap             - Number of cells to overlap for distribution
3782: . -dm_refine                         - Refine mesh after distribution
3783: . -dm_plex_hash_location             - Use grid hashing for point location
3784: . -dm_plex_hash_box_faces <n,m,p>    - The number of divisions in each direction of the grid hash
3785: . -dm_plex_partition_balance         - Attempt to evenly divide points on partition boundary between processes
3786: . -dm_plex_remesh_bd                 - Allow changes to the boundary on remeshing
3787: . -dm_plex_max_projection_height     - Maxmimum mesh point height used to project locally
3788: . -dm_plex_regular_refinement        - Use special nested projection algorithm for regular refinement
3789: . -dm_plex_check_all                 - Perform all shecks below
3790: . -dm_plex_check_symmetry            - Check that the adjacency information in the mesh is symmetric
3791: . -dm_plex_check_skeleton <celltype> - Check that each cell has the correct number of vertices
3792: . -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
3793: . -dm_plex_check_geometry            - Check that cells have positive volume
3794: . -dm_view :mesh.tex:ascii_latex     - View the mesh in LaTeX/TikZ
3795: . -dm_plex_view_scale <num>          - Scale the TikZ
3796: - -dm_plex_print_fem <num>           - View FEM assembly information, such as element vectors and matrices

3798:   Level: intermediate

3800: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
3801: M*/

3803: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
3804: {
3805:   DM_Plex       *mesh;
3806:   PetscInt       unit;

3809:   PetscNewLog(dm,&mesh);
3810:   dm->data = mesh;

3812:   mesh->refct             = 1;
3813:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
3814:   mesh->cones             = NULL;
3815:   mesh->coneOrientations  = NULL;
3816:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
3817:   mesh->supports          = NULL;
3818:   mesh->refinementUniform = PETSC_TRUE;
3819:   mesh->refinementLimit   = -1.0;
3820:   mesh->distDefault       = PETSC_TRUE;
3821:   mesh->interpolated      = DMPLEX_INTERPOLATED_INVALID;
3822:   mesh->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID;

3824:   mesh->facesTmp = NULL;

3826:   mesh->tetgenOpts   = NULL;
3827:   mesh->triangleOpts = NULL;
3828:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
3829:   mesh->remeshBd     = PETSC_FALSE;

3831:   mesh->subpointMap = NULL;

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

3835:   mesh->regularRefinement   = PETSC_FALSE;
3836:   mesh->depthState          = -1;
3837:   mesh->celltypeState       = -1;
3838:   mesh->globalVertexNumbers = NULL;
3839:   mesh->globalCellNumbers   = NULL;
3840:   mesh->anchorSection       = NULL;
3841:   mesh->anchorIS            = NULL;
3842:   mesh->createanchors       = NULL;
3843:   mesh->computeanchormatrix = NULL;
3844:   mesh->parentSection       = NULL;
3845:   mesh->parents             = NULL;
3846:   mesh->childIDs            = NULL;
3847:   mesh->childSection        = NULL;
3848:   mesh->children            = NULL;
3849:   mesh->referenceTree       = NULL;
3850:   mesh->getchildsymmetry    = NULL;
3851:   mesh->vtkCellHeight       = 0;
3852:   mesh->useAnchors          = PETSC_FALSE;

3854:   mesh->maxProjectionHeight = 0;

3856:   mesh->neighbors           = NULL;

3858:   mesh->printSetValues = PETSC_FALSE;
3859:   mesh->printFEM       = 0;
3860:   mesh->printTol       = 1.0e-10;

3862:   DMInitialize_Plex(dm);
3863:   return 0;
3864: }

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

3869:   Collective

3871:   Input Parameter:
3872: . comm - The communicator for the DMPlex object

3874:   Output Parameter:
3875: . mesh  - The DMPlex object

3877:   Level: beginner

3879: @*/
3880: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
3881: {
3883:   DMCreate(comm, mesh);
3884:   DMSetType(*mesh, DMPLEX);
3885:   return 0;
3886: }

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

3891:   Input Parameters:
3892: + dm - The DM
3893: . numCells - The number of cells owned by this process
3894: . numVertices - The number of vertices to be owned by this process, or PETSC_DECIDE
3895: . NVertices - The global number of vertices, or PETSC_DETERMINE
3896: . numCorners - The number of vertices for each cell
3897: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell

3899:   Output Parameters:
3900: + vertexSF - (Optional) SF describing complete vertex ownership
3901: - verticesAdjSaved - (Optional) vertex adjacency array

3903:   Notes:
3904:   Two triangles sharing a face
3905: $
3906: $        2
3907: $      / | \
3908: $     /  |  \
3909: $    /   |   \
3910: $   0  0 | 1  3
3911: $    \   |   /
3912: $     \  |  /
3913: $      \ | /
3914: $        1
3915: would have input
3916: $  numCells = 2, numVertices = 4
3917: $  cells = [0 1 2  1 3 2]
3918: $
3919: which would result in the DMPlex
3920: $
3921: $        4
3922: $      / | \
3923: $     /  |  \
3924: $    /   |   \
3925: $   2  0 | 1  5
3926: $    \   |   /
3927: $     \  |  /
3928: $      \ | /
3929: $        3

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

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

3939:   Not currently supported in Fortran.

3941:   Level: advanced

3943: .seealso: DMPlexBuildFromCellList(), DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildCoordinatesFromCellListParallel()
3944: @*/
3945: PetscErrorCode DMPlexBuildFromCellListParallel(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, const PetscInt cells[], PetscSF *vertexSF, PetscInt **verticesAdjSaved)
3946: {
3947:   PetscSF         sfPoint;
3948:   PetscLayout     layout;
3949:   PetscInt        numVerticesAdj, *verticesAdj, *cones, c, p;

3952:   PetscLogEventBegin(DMPLEX_BuildFromCellList,dm,0,0,0);
3953:   /* Get/check global number of vertices */
3954:   {
3955:     PetscInt NVerticesInCells, i;
3956:     const PetscInt len = numCells * numCorners;

3958:     /* NVerticesInCells = max(cells) + 1 */
3959:     NVerticesInCells = PETSC_MIN_INT;
3960:     for (i=0; i<len; i++) if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
3961:     ++NVerticesInCells;
3962:     MPI_Allreduce(MPI_IN_PLACE, &NVerticesInCells, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));

3964:     if (numVertices == PETSC_DECIDE && NVertices == PETSC_DECIDE) NVertices = NVerticesInCells;
3966:   }
3967:   /* Count locally unique vertices */
3968:   {
3969:     PetscHSetI vhash;
3970:     PetscInt off = 0;

3972:     PetscHSetICreate(&vhash);
3973:     for (c = 0; c < numCells; ++c) {
3974:       for (p = 0; p < numCorners; ++p) {
3975:         PetscHSetIAdd(vhash, cells[c*numCorners+p]);
3976:       }
3977:     }
3978:     PetscHSetIGetSize(vhash, &numVerticesAdj);
3979:     if (!verticesAdjSaved) PetscMalloc1(numVerticesAdj, &verticesAdj);
3980:     else { verticesAdj = *verticesAdjSaved; }
3981:     PetscHSetIGetElems(vhash, &off, verticesAdj);
3982:     PetscHSetIDestroy(&vhash);
3984:   }
3985:   PetscSortInt(numVerticesAdj, verticesAdj);
3986:   /* Create cones */
3987:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
3988:   for (c = 0; c < numCells; ++c) DMPlexSetConeSize(dm, c, numCorners);
3989:   DMSetUp(dm);
3990:   DMPlexGetCones(dm,&cones);
3991:   for (c = 0; c < numCells; ++c) {
3992:     for (p = 0; p < numCorners; ++p) {
3993:       const PetscInt gv = cells[c*numCorners+p];
3994:       PetscInt       lv;

3996:       /* Positions within verticesAdj form 0-based local vertex numbering;
3997:          we need to shift it by numCells to get correct DAG points (cells go first) */
3998:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
4000:       cones[c*numCorners+p] = lv+numCells;
4001:     }
4002:   }
4003:   /* Build point sf */
4004:   PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &layout);
4005:   PetscLayoutSetSize(layout, NVertices);
4006:   PetscLayoutSetLocalSize(layout, numVertices);
4007:   PetscLayoutSetBlockSize(layout, 1);
4008:   PetscSFCreateByMatchingIndices(layout, numVerticesAdj, verticesAdj, NULL, numCells, numVerticesAdj, verticesAdj, NULL, numCells, vertexSF, &sfPoint);
4009:   PetscLayoutDestroy(&layout);
4010:   if (!verticesAdjSaved) PetscFree(verticesAdj);
4011:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
4012:   if (dm->sf) {
4013:     const char *prefix;

4015:     PetscObjectGetOptionsPrefix((PetscObject)dm->sf, &prefix);
4016:     PetscObjectSetOptionsPrefix((PetscObject)sfPoint, prefix);
4017:   }
4018:   DMSetPointSF(dm, sfPoint);
4019:   PetscSFDestroy(&sfPoint);
4020:   if (vertexSF) PetscObjectSetName((PetscObject)(*vertexSF), "Vertex Ownership SF");
4021:   /* Fill in the rest of the topology structure */
4022:   DMPlexSymmetrize(dm);
4023:   DMPlexStratify(dm);
4024:   PetscLogEventEnd(DMPLEX_BuildFromCellList,dm,0,0,0);
4025:   return 0;
4026: }

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

4031:   Input Parameters:
4032: + dm - The DM
4033: . spaceDim - The spatial dimension used for coordinates
4034: . sfVert - SF describing complete vertex ownership
4035: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

4037:   Level: advanced

4039:   Notes:
4040:   Not currently supported in Fortran.

4042: .seealso: DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildFromCellListParallel()
4043: @*/
4044: PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM dm, PetscInt spaceDim, PetscSF sfVert, const PetscReal vertexCoords[])
4045: {
4046:   PetscSection   coordSection;
4047:   Vec            coordinates;
4048:   PetscScalar   *coords;
4049:   PetscInt       numVertices, numVerticesAdj, coordSize, v, vStart, vEnd;

4051:   PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
4052:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
4054:   DMSetCoordinateDim(dm, spaceDim);
4055:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
4057:   DMGetCoordinateSection(dm, &coordSection);
4058:   PetscSectionSetNumFields(coordSection, 1);
4059:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
4060:   PetscSectionSetChart(coordSection, vStart, vEnd);
4061:   for (v = vStart; v < vEnd; ++v) {
4062:     PetscSectionSetDof(coordSection, v, spaceDim);
4063:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
4064:   }
4065:   PetscSectionSetUp(coordSection);
4066:   PetscSectionGetStorageSize(coordSection, &coordSize);
4067:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
4068:   VecSetBlockSize(coordinates, spaceDim);
4069:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
4070:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
4071:   VecSetType(coordinates,VECSTANDARD);
4072:   VecGetArray(coordinates, &coords);
4073:   {
4074:     MPI_Datatype coordtype;

4076:     /* Need a temp buffer for coords if we have complex/single */
4077:     MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);
4078:     MPI_Type_commit(&coordtype);
4079: #if defined(PETSC_USE_COMPLEX)
4080:     {
4081:     PetscScalar *svertexCoords;
4082:     PetscInt    i;
4083:     PetscMalloc1(numVertices*spaceDim,&svertexCoords);
4084:     for (i=0; i<numVertices*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
4085:     PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords,MPI_REPLACE);
4086:     PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords,MPI_REPLACE);
4087:     PetscFree(svertexCoords);
4088:     }
4089: #else
4090:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords,MPI_REPLACE);
4091:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords,MPI_REPLACE);
4092: #endif
4093:     MPI_Type_free(&coordtype);
4094:   }
4095:   VecRestoreArray(coordinates, &coords);
4096:   DMSetCoordinatesLocal(dm, coordinates);
4097:   VecDestroy(&coordinates);
4098:   PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
4099:   return 0;
4100: }

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

4105:   Input Parameters:
4106: + comm - The communicator
4107: . dim - The topological dimension of the mesh
4108: . numCells - The number of cells owned by this process
4109: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE
4110: . NVertices - The global number of vertices, or PETSC_DECIDE
4111: . numCorners - The number of vertices for each cell
4112: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
4113: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
4114: . spaceDim - The spatial dimension used for coordinates
4115: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

4117:   Output Parameters:
4118: + dm - The DM
4119: . vertexSF - (Optional) SF describing complete vertex ownership
4120: - verticesAdjSaved - (Optional) vertex adjacency array

4122:   Notes:
4123:   This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(),
4124:   DMPlexBuildFromCellListParallel(), DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellListParallel()

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

4129:   Level: intermediate

4131: .seealso: DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromDAG(), DMPlexCreate()
4132: @*/
4133: 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, PetscInt **verticesAdj, DM *dm)
4134: {
4135:   PetscSF        sfVert;

4137:   DMCreate(comm, dm);
4138:   DMSetType(*dm, DMPLEX);
4141:   DMSetDimension(*dm, dim);
4142:   DMPlexBuildFromCellListParallel(*dm, numCells, numVertices, NVertices, numCorners, cells, &sfVert, verticesAdj);
4143:   if (interpolate) {
4144:     DM idm;

4146:     DMPlexInterpolate(*dm, &idm);
4147:     DMDestroy(dm);
4148:     *dm  = idm;
4149:   }
4150:   DMPlexBuildCoordinatesFromCellListParallel(*dm, spaceDim, sfVert, vertexCoords);
4151:   if (vertexSF) *vertexSF = sfVert;
4152:   else PetscSFDestroy(&sfVert);
4153:   return 0;
4154: }

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

4159:   Input Parameters:
4160: + dm - The DM
4161: . numCells - The number of cells owned by this process
4162: . numVertices - The number of vertices owned by this process, or PETSC_DETERMINE
4163: . numCorners - The number of vertices for each cell
4164: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell

4166:   Level: advanced

4168:   Notes:
4169:   Two triangles sharing a face
4170: $
4171: $        2
4172: $      / | \
4173: $     /  |  \
4174: $    /   |   \
4175: $   0  0 | 1  3
4176: $    \   |   /
4177: $     \  |  /
4178: $      \ | /
4179: $        1
4180: would have input
4181: $  numCells = 2, numVertices = 4
4182: $  cells = [0 1 2  1 3 2]
4183: $
4184: which would result in the DMPlex
4185: $
4186: $        4
4187: $      / | \
4188: $     /  |  \
4189: $    /   |   \
4190: $   2  0 | 1  5
4191: $    \   |   /
4192: $     \  |  /
4193: $      \ | /
4194: $        3

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

4198:   Not currently supported in Fortran.

4200: .seealso: DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromCellListPetsc()
4201: @*/
4202: PetscErrorCode DMPlexBuildFromCellList(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const PetscInt cells[])
4203: {
4204:   PetscInt      *cones, c, p, dim;

4206:   PetscLogEventBegin(DMPLEX_BuildFromCellList,dm,0,0,0);
4207:   DMGetDimension(dm, &dim);
4208:   /* Get/check global number of vertices */
4209:   {
4210:     PetscInt NVerticesInCells, i;
4211:     const PetscInt len = numCells * numCorners;

4213:     /* NVerticesInCells = max(cells) + 1 */
4214:     NVerticesInCells = PETSC_MIN_INT;
4215:     for (i=0; i<len; i++) if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
4216:     ++NVerticesInCells;

4218:     if (numVertices == PETSC_DECIDE) numVertices = NVerticesInCells;
4220:   }
4221:   DMPlexSetChart(dm, 0, numCells+numVertices);
4222:   for (c = 0; c < numCells; ++c) {
4223:     DMPlexSetConeSize(dm, c, numCorners);
4224:   }
4225:   DMSetUp(dm);
4226:   DMPlexGetCones(dm,&cones);
4227:   for (c = 0; c < numCells; ++c) {
4228:     for (p = 0; p < numCorners; ++p) {
4229:       cones[c*numCorners+p] = cells[c*numCorners+p]+numCells;
4230:     }
4231:   }
4232:   DMPlexSymmetrize(dm);
4233:   DMPlexStratify(dm);
4234:   PetscLogEventEnd(DMPLEX_BuildFromCellList,dm,0,0,0);
4235:   return 0;
4236: }

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

4241:   Input Parameters:
4242: + dm - The DM
4243: . spaceDim - The spatial dimension used for coordinates
4244: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

4246:   Level: advanced

4248:   Notes:
4249:   Not currently supported in Fortran.

4251: .seealso: DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellList()
4252: @*/
4253: PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM dm, PetscInt spaceDim, const PetscReal vertexCoords[])
4254: {
4255:   PetscSection   coordSection;
4256:   Vec            coordinates;
4257:   DM             cdm;
4258:   PetscScalar   *coords;
4259:   PetscInt       v, vStart, vEnd, d;

4261:   PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
4262:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
4264:   DMSetCoordinateDim(dm, spaceDim);
4265:   DMGetCoordinateSection(dm, &coordSection);
4266:   PetscSectionSetNumFields(coordSection, 1);
4267:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
4268:   PetscSectionSetChart(coordSection, vStart, vEnd);
4269:   for (v = vStart; v < vEnd; ++v) {
4270:     PetscSectionSetDof(coordSection, v, spaceDim);
4271:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
4272:   }
4273:   PetscSectionSetUp(coordSection);

4275:   DMGetCoordinateDM(dm, &cdm);
4276:   DMCreateLocalVector(cdm, &coordinates);
4277:   VecSetBlockSize(coordinates, spaceDim);
4278:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
4279:   VecGetArrayWrite(coordinates, &coords);
4280:   for (v = 0; v < vEnd-vStart; ++v) {
4281:     for (d = 0; d < spaceDim; ++d) {
4282:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
4283:     }
4284:   }
4285:   VecRestoreArrayWrite(coordinates, &coords);
4286:   DMSetCoordinatesLocal(dm, coordinates);
4287:   VecDestroy(&coordinates);
4288:   PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList,dm,0,0,0);
4289:   return 0;
4290: }

4292: /*@
4293:   DMPlexCreateFromCellListPetsc - Create DMPLEX from a list of vertices for each cell (common mesh generator output), but only process 0 takes in the input

4295:   Collective on comm

4297:   Input Parameters:
4298: + comm - The communicator
4299: . dim - The topological dimension of the mesh
4300: . numCells - The number of cells, only on process 0
4301: . numVertices - The number of vertices owned by this process, or PETSC_DECIDE, only on process 0
4302: . numCorners - The number of vertices for each cell, only on process 0
4303: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
4304: . cells - An array of numCells*numCorners numbers, the vertices for each cell, only on process 0
4305: . spaceDim - The spatial dimension used for coordinates
4306: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex, only on process 0

4308:   Output Parameter:
4309: . dm - The DM, which only has points on process 0

4311:   Notes:
4312:   This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(), DMPlexBuildFromCellList(),
4313:   DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellList()

4315:   See DMPlexBuildFromCellList() for an example and details about the topology-related parameters.
4316:   See DMPlexBuildCoordinatesFromCellList() for details about the geometry-related parameters.
4317:   See DMPlexCreateFromCellListParallelPetsc() for parallel input

4319:   Level: intermediate

4321: .seealso: DMPlexCreateFromCellListParallelPetsc(), DMPlexBuildFromCellList(), DMPlexBuildCoordinatesFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
4322: @*/
4323: 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)
4324: {
4325:   PetscMPIInt    rank;

4328:   MPI_Comm_rank(comm, &rank);
4329:   DMCreate(comm, dm);
4330:   DMSetType(*dm, DMPLEX);
4331:   DMSetDimension(*dm, dim);
4332:   if (!rank) DMPlexBuildFromCellList(*dm, numCells, numVertices, numCorners, cells);
4333:   else       DMPlexBuildFromCellList(*dm, 0, 0, 0, NULL);
4334:   if (interpolate) {
4335:     DM idm;

4337:     DMPlexInterpolate(*dm, &idm);
4338:     DMDestroy(dm);
4339:     *dm  = idm;
4340:   }
4341:   if (!rank) DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, vertexCoords);
4342:   else       DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, NULL);
4343:   return 0;
4344: }

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

4349:   Input Parameters:
4350: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
4351: . depth - The depth of the DAG
4352: . numPoints - Array of size depth + 1 containing the number of points at each depth
4353: . coneSize - The cone size of each point
4354: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
4355: . coneOrientations - The orientation of each cone point
4356: - vertexCoords - An array of numPoints[0]*spacedim numbers representing the coordinates of each vertex, with spacedim the value set via DMSetCoordinateDim()

4358:   Output Parameter:
4359: . dm - The DM

4361:   Note: Two triangles sharing a face would have input
4362: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
4363: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
4364: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
4365: $
4366: which would result in the DMPlex
4367: $
4368: $        4
4369: $      / | \
4370: $     /  |  \
4371: $    /   |   \
4372: $   2  0 | 1  5
4373: $    \   |   /
4374: $     \  |  /
4375: $      \ | /
4376: $        3
4377: $
4378: $ Notice that all points are numbered consecutively, unlike DMPlexCreateFromCellListPetsc()

4380:   Level: advanced

4382: .seealso: DMPlexCreateFromCellListPetsc(), DMPlexCreate()
4383: @*/
4384: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
4385: {
4386:   Vec            coordinates;
4387:   PetscSection   coordSection;
4388:   PetscScalar    *coords;
4389:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

4391:   DMGetDimension(dm, &dim);
4392:   DMGetCoordinateDim(dm, &dimEmbed);
4394:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
4395:   DMPlexSetChart(dm, pStart, pEnd);
4396:   for (p = pStart; p < pEnd; ++p) {
4397:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
4398:     if (firstVertex < 0 && !coneSize[p - pStart]) {
4399:       firstVertex = p - pStart;
4400:     }
4401:   }
4403:   DMSetUp(dm); /* Allocate space for cones */
4404:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
4405:     DMPlexSetCone(dm, p, &cones[off]);
4406:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
4407:   }
4408:   DMPlexSymmetrize(dm);
4409:   DMPlexStratify(dm);
4410:   /* Build coordinates */
4411:   DMGetCoordinateSection(dm, &coordSection);
4412:   PetscSectionSetNumFields(coordSection, 1);
4413:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
4414:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
4415:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
4416:     PetscSectionSetDof(coordSection, v, dimEmbed);
4417:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
4418:   }
4419:   PetscSectionSetUp(coordSection);
4420:   PetscSectionGetStorageSize(coordSection, &coordSize);
4421:   VecCreate(PETSC_COMM_SELF, &coordinates);
4422:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
4423:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
4424:   VecSetBlockSize(coordinates, dimEmbed);
4425:   VecSetType(coordinates,VECSTANDARD);
4426:   if (vertexCoords) {
4427:     VecGetArray(coordinates, &coords);
4428:     for (v = 0; v < numPoints[0]; ++v) {
4429:       PetscInt off;

4431:       PetscSectionGetOffset(coordSection, v+firstVertex, &off);
4432:       for (d = 0; d < dimEmbed; ++d) {
4433:         coords[off+d] = vertexCoords[v*dimEmbed+d];
4434:       }
4435:     }
4436:   }
4437:   VecRestoreArray(coordinates, &coords);
4438:   DMSetCoordinatesLocal(dm, coordinates);
4439:   VecDestroy(&coordinates);
4440:   return 0;
4441: }

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

4446: + comm        - The MPI communicator
4447: . filename    - Name of the .dat file
4448: - interpolate - Create faces and edges in the mesh

4450:   Output Parameter:
4451: . dm  - The DM object representing the mesh

4453:   Note: The format is the simplest possible:
4454: $ Ne
4455: $ v0 v1 ... vk
4456: $ Nv
4457: $ x y z marker

4459:   Level: beginner

4461: .seealso: DMPlexCreateFromFile(), DMPlexCreateMedFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
4462: @*/
4463: PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
4464: {
4465:   DMLabel         marker;
4466:   PetscViewer     viewer;
4467:   Vec             coordinates;
4468:   PetscSection    coordSection;
4469:   PetscScalar    *coords;
4470:   char            line[PETSC_MAX_PATH_LEN];
4471:   PetscInt        dim = 3, cdim = 3, coordSize, v, c, d;
4472:   PetscMPIInt     rank;
4473:   int             snum, Nv, Nc, Ncn, Nl;

4475:   MPI_Comm_rank(comm, &rank);
4476:   PetscViewerCreate(comm, &viewer);
4477:   PetscViewerSetType(viewer, PETSCVIEWERASCII);
4478:   PetscViewerFileSetMode(viewer, FILE_MODE_READ);
4479:   PetscViewerFileSetName(viewer, filename);
4480:   if (rank == 0) {
4481:     PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING);
4482:     snum = sscanf(line, "%d %d %d %d", &Nc, &Nv, &Ncn, &Nl);
4484:   } else {
4485:     Nc = Nv = Ncn = Nl = 0;
4486:   }
4487:   DMCreate(comm, dm);
4488:   DMSetType(*dm, DMPLEX);
4489:   DMPlexSetChart(*dm, 0, Nc+Nv);
4490:   DMSetDimension(*dm, dim);
4491:   DMSetCoordinateDim(*dm, cdim);
4492:   /* Read topology */
4493:   if (rank == 0) {
4494:     char     format[PETSC_MAX_PATH_LEN];
4495:     PetscInt cone[8];
4496:     int      vbuf[8], v;

4498:     for (c = 0; c < Ncn; ++c) {format[c*3+0] = '%'; format[c*3+1] = 'd'; format[c*3+2] = ' ';}
4499:     format[Ncn*3-1] = '\0';
4500:     for (c = 0; c < Nc; ++c) DMPlexSetConeSize(*dm, c, Ncn);
4501:     DMSetUp(*dm);
4502:     for (c = 0; c < Nc; ++c) {
4503:       PetscViewerRead(viewer, line, Ncn, NULL, PETSC_STRING);
4504:       switch (Ncn) {
4505:         case 2: snum = sscanf(line, format, &vbuf[0], &vbuf[1]);break;
4506:         case 3: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2]);break;
4507:         case 4: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3]);break;
4508:         case 6: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5]);break;
4509:         case 8: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5], &vbuf[6], &vbuf[7]);break;
4510:         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No cell shape with %D vertices", Ncn);
4511:       }
4513:       for (v = 0; v < Ncn; ++v) cone[v] = vbuf[v] + Nc;
4514:       /* Hexahedra are inverted */
4515:       if (Ncn == 8) {
4516:         PetscInt tmp = cone[1];
4517:         cone[1] = cone[3];
4518:         cone[3] = tmp;
4519:       }
4520:       DMPlexSetCone(*dm, c, cone);
4521:     }
4522:   }
4523:   DMPlexSymmetrize(*dm);
4524:   DMPlexStratify(*dm);
4525:   /* Read coordinates */
4526:   DMGetCoordinateSection(*dm, &coordSection);
4527:   PetscSectionSetNumFields(coordSection, 1);
4528:   PetscSectionSetFieldComponents(coordSection, 0, cdim);
4529:   PetscSectionSetChart(coordSection, Nc, Nc + Nv);
4530:   for (v = Nc; v < Nc+Nv; ++v) {
4531:     PetscSectionSetDof(coordSection, v, cdim);
4532:     PetscSectionSetFieldDof(coordSection, v, 0, cdim);
4533:   }
4534:   PetscSectionSetUp(coordSection);
4535:   PetscSectionGetStorageSize(coordSection, &coordSize);
4536:   VecCreate(PETSC_COMM_SELF, &coordinates);
4537:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
4538:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
4539:   VecSetBlockSize(coordinates, cdim);
4540:   VecSetType(coordinates, VECSTANDARD);
4541:   VecGetArray(coordinates, &coords);
4542:   if (rank == 0) {
4543:     char   format[PETSC_MAX_PATH_LEN];
4544:     double x[3];
4545:     int    l, val[3];

4547:     if (Nl) {
4548:       for (l = 0; l < Nl; ++l) {format[l*3+0] = '%'; format[l*3+1] = 'd'; format[l*3+2] = ' ';}
4549:       format[Nl*3-1] = '\0';
4550:       DMCreateLabel(*dm, "marker");
4551:       DMGetLabel(*dm, "marker", &marker);
4552:     }
4553:     for (v = 0; v < Nv; ++v) {
4554:       PetscViewerRead(viewer, line, 3+Nl, NULL, PETSC_STRING);
4555:       snum = sscanf(line, "%lg %lg %lg", &x[0], &x[1], &x[2]);
4557:       switch (Nl) {
4558:         case 0: snum = 0;break;
4559:         case 1: snum = sscanf(line, format, &val[0]);break;
4560:         case 2: snum = sscanf(line, format, &val[0], &val[1]);break;
4561:         case 3: snum = sscanf(line, format, &val[0], &val[1], &val[2]);break;
4562:         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Request support for %D labels", Nl);
4563:       }
4565:       for (d = 0; d < cdim; ++d) coords[v*cdim+d] = x[d];
4566:       for (l = 0; l < Nl; ++l) DMLabelSetValue(marker, v+Nc, val[l]);
4567:     }
4568:   }
4569:   VecRestoreArray(coordinates, &coords);
4570:   DMSetCoordinatesLocal(*dm, coordinates);
4571:   VecDestroy(&coordinates);
4572:   PetscViewerDestroy(&viewer);
4573:   if (interpolate) {
4574:     DM      idm;
4575:     DMLabel bdlabel;

4577:     DMPlexInterpolate(*dm, &idm);
4578:     DMDestroy(dm);
4579:     *dm  = idm;

4581:     if (!Nl) {
4582:       DMCreateLabel(*dm, "marker");
4583:       DMGetLabel(*dm, "marker", &bdlabel);
4584:       DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
4585:       DMPlexLabelComplete(*dm, bdlabel);
4586:     }
4587:   }
4588:   return 0;
4589: }

4591: /*@C
4592:   DMPlexCreateFromFile - This takes a filename and produces a DM

4594:   Input Parameters:
4595: + comm - The communicator
4596: . filename - A file name
4597: . plexname - The object name of the resulting DM, also used for intra-datafile lookup by some formats
4598: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

4600:   Output Parameter:
4601: . dm - The DM

4603:   Options Database Keys:
4604: . -dm_plex_create_from_hdf5_xdmf - use the PETSC_VIEWER_HDF5_XDMF format for reading HDF5

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

4609:   Notes:
4610:   Using PETSCVIEWERHDF5 type with PETSC_VIEWER_HDF5_PETSC format, one can save multiple DMPlex
4611:   meshes in a single HDF5 file. This in turn requires one to name the DMPlex object with PetscObjectSetName()
4612:   before saving it with DMView() and before loading it with DMLoad() for identification of the mesh object.
4613:   The input parameter name is thus used to name the DMPlex object when DMPlexCreateFromFile() internally
4614:   calls DMLoad(). Currently, name is ignored for other viewer types and/or formats.

4616:   Level: beginner

4618: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellListPetsc(), DMPlexCreate(), PetscObjectSetName(), DMView(), DMLoad()
4619: @*/
4620: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], const char plexname[], PetscBool interpolate, DM *dm)
4621: {
4622:   const char  extGmsh[]      = ".msh";
4623:   const char  extGmsh2[]     = ".msh2";
4624:   const char  extGmsh4[]     = ".msh4";
4625:   const char  extCGNS[]      = ".cgns";
4626:   const char  extExodus[]    = ".exo";
4627:   const char  extExodus_e[]  = ".e";
4628:   const char  extGenesis[]   = ".gen";
4629:   const char  extFluent[]    = ".cas";
4630:   const char  extHDF5[]      = ".h5";
4631:   const char  extMed[]       = ".med";
4632:   const char  extPLY[]       = ".ply";
4633:   const char  extEGADSLite[] = ".egadslite";
4634:   const char  extEGADS[]     = ".egads";
4635:   const char  extIGES[]      = ".igs";
4636:   const char  extSTEP[]      = ".stp";
4637:   const char  extCV[]        = ".dat";
4638:   size_t      len;
4639:   PetscBool   isGmsh, isGmsh2, isGmsh4, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isEGADSLite, isEGADS, isIGES, isSTEP, isCV;
4640:   PetscMPIInt rank;

4645:   DMInitializePackage();
4646:   PetscLogEventBegin(DMPLEX_CreateFromFile,0,0,0,0);
4647:   MPI_Comm_rank(comm, &rank);
4648:   PetscStrlen(filename, &len);

4651: #define CheckExtension(extension__,is_extension__) do {                                        \
4652:     PetscAssert(sizeof(extension__), comm, PETSC_ERR_PLIB, "Zero-size extension: %s", extension__); \
4653:     /* don't count the null-terminator at the end */                                           \
4654:     const size_t ext_len = sizeof(extension__)-1;                                              \
4655:     if (len < ext_len) {                                                                       \
4656:       is_extension__ = PETSC_FALSE;                                                            \
4657:     } else {                                                                                   \
4658:       PetscStrncmp(filename+len-ext_len, extension__, ext_len, &is_extension__);    \
4659:     }                                                                                          \
4660:   } while (0)

4662:   CheckExtension(extGmsh, isGmsh);
4663:   CheckExtension(extGmsh2, isGmsh2);
4664:   CheckExtension(extGmsh4, isGmsh4);
4665:   CheckExtension(extCGNS, isCGNS);
4666:   CheckExtension(extExodus, isExodus);
4667:   if (!isExodus) CheckExtension(extExodus_e, isExodus);
4668:   CheckExtension(extGenesis, isGenesis);
4669:   CheckExtension(extFluent, isFluent);
4670:   CheckExtension(extHDF5, isHDF5);
4671:   CheckExtension(extMed, isMed);
4672:   CheckExtension(extPLY, isPLY);
4673:   CheckExtension(extEGADSLite, isEGADSLite);
4674:   CheckExtension(extEGADS, isEGADS);
4675:   CheckExtension(extIGES, isIGES);
4676:   CheckExtension(extSTEP, isSTEP);
4677:   CheckExtension(extCV, isCV);

4679: #undef CheckExtension

4681:   if (isGmsh || isGmsh2 || isGmsh4) {
4682:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
4683:   } else if (isCGNS) {
4684:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
4685:   } else if (isExodus || isGenesis) {
4686:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
4687:   } else if (isFluent) {
4688:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
4689:   } else if (isHDF5) {
4690:     PetscBool      load_hdf5_xdmf = PETSC_FALSE;
4691:     PetscViewer viewer;

4693:     /* PETSC_VIEWER_HDF5_XDMF is used if the filename ends with .xdmf.h5, or if -dm_plex_create_from_hdf5_xdmf option is present */
4694:     PetscStrncmp(&filename[PetscMax(0,len-8)], ".xdmf",  5, &load_hdf5_xdmf);
4695:     PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &load_hdf5_xdmf, NULL);
4696:     PetscViewerCreate(comm, &viewer);
4697:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
4698:     PetscViewerSetOptionsPrefix(viewer, "dm_plex_create_");
4699:     PetscViewerSetFromOptions(viewer);
4700:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
4701:     PetscViewerFileSetName(viewer, filename);

4703:     DMCreate(comm, dm);
4704:     PetscObjectSetName((PetscObject)(*dm), plexname);
4705:     DMSetType(*dm, DMPLEX);
4706:     if (load_hdf5_xdmf) PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF);
4707:     DMLoad(*dm, viewer);
4708:     if (load_hdf5_xdmf) PetscViewerPopFormat(viewer);
4709:     PetscViewerDestroy(&viewer);

4711:     if (interpolate) {
4712:       DM idm;

4714:       DMPlexInterpolate(*dm, &idm);
4715:       DMDestroy(dm);
4716:       *dm  = idm;
4717:     }
4718:   } else if (isMed) {
4719:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
4720:   } else if (isPLY) {
4721:     DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);
4722:   } else if (isEGADSLite || isEGADS || isIGES || isSTEP) {
4723:     if (isEGADSLite) DMPlexCreateEGADSLiteFromFile(comm, filename, dm);
4724:     else             DMPlexCreateEGADSFromFile(comm, filename, dm);
4725:     if (!interpolate) {
4726:       DM udm;

4728:       DMPlexUninterpolate(*dm, &udm);
4729:       DMDestroy(dm);
4730:       *dm  = udm;
4731:     }
4732:   } else if (isCV) {
4733:     DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm);
4734:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
4735:   PetscStrlen(plexname, &len);
4736:   if (len) PetscObjectSetName((PetscObject)(*dm), plexname);
4737:   PetscLogEventEnd(DMPLEX_CreateFromFile,0,0,0,0);
4738:   return 0;
4739: }