Actual source code: plexcreate.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: #define PETSCDM_DLL
  2:  #include <petsc/private/dmpleximpl.h>
  3:  #include <petsc/private/hashseti.h>
  4:  #include <petscsf.h>

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

  9:   Collective on MPI_Comm

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

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

 22:   Level: beginner

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

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

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

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

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

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

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

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

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

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

149:   Collective on MPI_Comm

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

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

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

171:   Level: beginner

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

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

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

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

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

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

296:   Collective on MPI_Comm

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

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

307:   Level: beginner

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

326:   if ((faces[0] < 1) || (faces[1] < 1) || (faces[2] < 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side");
327:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
328:   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
329:   numVertices = vertices[0]*vertices[1]*vertices[2];
330:   if (!rank) {
331:     PetscInt f;

333:     DMPlexSetChart(dm, 0, numFaces+numVertices);
334:     for (f = 0; f < numFaces; ++f) {
335:       DMPlexSetConeSize(dm, f, 4);
336:     }
337:     DMSetUp(dm); /* Allocate space for cones */
338:     for (v = 0; v < numFaces+numVertices; ++v) {
339:       DMSetLabelValue(dm, "marker", v, 1);
340:     }

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

352:     /* Side 1 (Bottom) */
353:     for (vy = 0; vy < faces[1]; vy++) {
354:       for (vx = 0; vx < faces[0]; vx++) {
355:         voffset = numFaces + vy*(faces[0]+1) + vx;
356:         cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1;
357:         DMPlexSetCone(dm, iface, cone);
358:         iface++;
359:       }
360:     }

362:     /* Side 2 (Front) */
363:     for (vz = 0; vz < faces[2]; vz++) {
364:       for (vx = 0; vx < faces[0]; vx++) {
365:         voffset = numFaces + vz*vertices[0]*vertices[1] + vx;
366:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1];
367:         DMPlexSetCone(dm, iface, cone);
368:         iface++;
369:       }
370:     }

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

383:     /* Side 4 (Left) */
384:     for (vz = 0; vz < faces[2]; vz++) {
385:       for (vy = 0; vy < faces[1]; vy++) {
386:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
387:         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
388:         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
389:         DMPlexSetCone(dm, iface, cone);
390:         iface++;
391:       }
392:     }

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

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


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

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

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

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

490: static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
491: {
492:   DM             boundary;
493:   PetscInt       i;

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

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

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

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

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

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

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

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

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

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

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

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

921: static PetscErrorCode DMPlexCreateBoxMesh_Tensor_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
922: {
923:   PetscInt       i;

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

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

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

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

963:   Collective on MPI_Comm

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

975:   Output Parameter:
976: . dm  - The DM object

978:   Note: Here is the numbering returned for 2 faces in each direction for tensor cells:
979: $ 10---17---11---18----12
980: $  |         |         |
981: $  |         |         |
982: $ 20    2   22    3    24
983: $  |         |         |
984: $  |         |         |
985: $  7---15----8---16----9
986: $  |         |         |
987: $  |         |         |
988: $ 19    0   21    1   23
989: $  |         |         |
990: $  |         |         |
991: $  4---13----5---14----6

993: and for simplicial cells

995: $ 14----8---15----9----16
996: $  |\     5  |\      7 |
997: $  | \       | \       |
998: $ 13   2    14    3    15
999: $  | 4   \   | 6   \   |
1000: $  |       \ |       \ |
1001: $ 11----6---12----7----13
1002: $  |\        |\        |
1003: $  | \    1  | \     3 |
1004: $ 10   0    11    1    12
1005: $  | 0   \   | 2   \   |
1006: $  |       \ |       \ |
1007: $  8----4----9----5----10

1009:   Level: beginner

1011: .keywords: DM, create
1012: .seealso: DMPlexCreateFromFile(), DMPlexCreateHexCylinderMesh(), DMSetType(), DMCreate()
1013: @*/
1014: 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)
1015: {
1016:   PetscInt       i;
1017:   PetscInt       fac[3] = {0, 0, 0};
1018:   PetscReal      low[3] = {0, 0, 0};
1019:   PetscReal      upp[3] = {1, 1, 1};
1020:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1024:   for (i = 0; i < dim; ++i) fac[i] = faces ? faces[i] : (dim == 1 ? 1 : 4-dim);
1025:   if (lower) for (i = 0; i < dim; ++i) low[i] = lower[i];
1026:   if (upper) for (i = 0; i < dim; ++i) upp[i] = upper[i];
1027:   if (periodicity) for (i = 0; i < dim; ++i) bdt[i] = periodicity[i];

1029:   if (dim == 1)      {DMPlexCreateLineMesh_Internal(comm, fac[0], low[0], upp[0], bdt[0], dm);}
1030:   else if (simplex)  {DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1031:   else               {DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1032:   return(0);
1033: }

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

1038:   Collective on MPI_Comm

1040:   Input Parameters:
1041: + comm        - The communicator for the DM object
1042: . faces       - Number of faces per dimension, or NULL for (1, 1, 1)
1043: . lower       - The lower left corner, or NULL for (0, 0, 0)
1044: . upper       - The upper right corner, or NULL for (1, 1, 1)
1045: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
1046: . ordExt      - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1047: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1049:   Output Parameter:
1050: . dm  - The DM object

1052:   Level: beginner

1054: .keywords: DM, create
1055: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateWedgeCylinderMesh(), DMPlexExtrude(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1056: @*/
1057: PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool ordExt, PetscBool interpolate, DM *dm)
1058: {
1059:   DM             bdm, botdm;
1060:   PetscInt       i;
1061:   PetscInt       fac[3] = {0, 0, 0};
1062:   PetscReal      low[3] = {0, 0, 0};
1063:   PetscReal      upp[3] = {1, 1, 1};
1064:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1068:   for (i = 0; i < 3; ++i) fac[i] = faces ? (faces[i] > 0 ? faces[i] : 1) : 1;
1069:   if (lower) for (i = 0; i < 3; ++i) low[i] = lower[i];
1070:   if (upper) for (i = 0; i < 3; ++i) upp[i] = upper[i];
1071:   if (periodicity) for (i = 0; i < 3; ++i) bdt[i] = periodicity[i];
1072:   for (i = 0; i < 3; ++i) if (bdt[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity not yet supported");

1074:   DMCreate(comm, &bdm);
1075:   DMSetType(bdm, DMPLEX);
1076:   DMSetDimension(bdm, 1);
1077:   DMSetCoordinateDim(bdm, 2);
1078:   DMPlexCreateSquareBoundary(bdm, low, upp, fac);
1079:   DMPlexGenerate(bdm, NULL, PETSC_FALSE, &botdm);
1080:   DMDestroy(&bdm);
1081:   DMPlexExtrude(botdm, fac[2], upp[2] - low[2], ordExt, interpolate, dm);
1082:   if (low[2] != 0.0) {
1083:     Vec         v;
1084:     PetscScalar *x;
1085:     PetscInt    cDim, n;

1087:     DMGetCoordinatesLocal(*dm, &v);
1088:     VecGetBlockSize(v, &cDim);
1089:     VecGetLocalSize(v, &n);
1090:     VecGetArray(v, &x);
1091:     x   += cDim;
1092:     for (i=0; i<n; i+=cDim) x[i] += low[2];
1093:     VecRestoreArray(v,&x);
1094:     DMSetCoordinatesLocal(*dm, v);
1095:   }
1096:   DMDestroy(&botdm);
1097:   return(0);
1098: }

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

1103:   Collective on idm

1105:   Input Parameters:
1106: + idm - The mesh to be extruted
1107: . layers - The number of layers
1108: . height - The height of the extruded layer
1109: . ordExt - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1110: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1112:   Output Parameter:
1113: . dm  - The DM object

1115:   Notes: The object created is an hybrid mesh, the vertex ordering in the cone of the cell is that of the prismatic cells

1117:   Level: advanced

1119: .keywords: DM, create
1120: .seealso: DMPlexCreateWedgeCylinderMesh(), DMPlexCreateWedgeBoxMesh(), DMPlexSetHybridBounds(), DMSetType(), DMCreate()
1121: @*/
1122: PetscErrorCode DMPlexExtrude(DM idm, PetscInt layers, PetscReal height, PetscBool ordExt, PetscBool interpolate, DM* dm)
1123: {
1124:   PetscScalar       *coordsB;
1125:   const PetscScalar *coordsA;
1126:   PetscReal         *normals = NULL;
1127:   Vec               coordinatesA, coordinatesB;
1128:   PetscSection      coordSectionA, coordSectionB;
1129:   PetscInt          dim, cDim, cDimB, c, l, v, coordSize, *newCone;
1130:   PetscInt          cStart, cEnd, vStart, vEnd, cellV, numCells, numVertices;
1131:   PetscErrorCode    ierr;

1138:   DMGetDimension(idm, &dim);
1139:   if (dim < 1 || dim > 3) SETERRQ1(PetscObjectComm((PetscObject)idm), PETSC_ERR_SUP, "Support for dimension %D not coded", dim);

1141:   DMPlexGetHeightStratum(idm, 0, &cStart, &cEnd);
1142:   DMPlexGetDepthStratum(idm, 0, &vStart, &vEnd);
1143:   numCells = (cEnd - cStart)*layers;
1144:   numVertices = (vEnd - vStart)*(layers+1);
1145:   DMCreate(PetscObjectComm((PetscObject)idm), dm);
1146:   DMSetType(*dm, DMPLEX);
1147:   DMSetDimension(*dm, dim+1);
1148:   DMPlexSetChart(*dm, 0, numCells+numVertices);
1149:   for (c = cStart, cellV = 0; c < cEnd; ++c) {
1150:     PetscInt *closure = NULL;
1151:     PetscInt closureSize, numCorners = 0;

1153:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1154:     for (v = 0; v < closureSize*2; v += 2) if ((closure[v] >= vStart) && (closure[v] < vEnd)) numCorners++;
1155:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1156:     for (l = 0; l < layers; ++l) {
1157:       DMPlexSetConeSize(*dm, ordExt ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart, 2*numCorners);
1158:     }
1159:     cellV = PetscMax(numCorners,cellV);
1160:   }
1161:   DMPlexSetHybridBounds(*dm, 0, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);
1162:   DMSetUp(*dm);

1164:   DMGetCoordinateDim(idm, &cDim);
1165:   if (dim != cDim) {
1166:     PetscCalloc1(cDim*(vEnd - vStart), &normals);
1167:   }
1168:   PetscMalloc1(3*cellV,&newCone);
1169:   for (c = cStart; c < cEnd; ++c) {
1170:     PetscInt *closure = NULL;
1171:     PetscInt closureSize, numCorners = 0, l;
1172:     PetscReal normal[3] = {0, 0, 0};

1174:     if (normals) {
1175:       DMPlexComputeCellGeometryFVM(idm, c, NULL, NULL, normal);
1176:     }
1177:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1178:     for (v = 0; v < closureSize*2; v += 2) {
1179:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1180:         PetscInt d;

1182:         newCone[numCorners++] = closure[v] - vStart;
1183:         if (normals) { for (d = 0; d < cDim; ++d) normals[cDim*(closure[v]-vStart)+d] += normal[d]; }
1184:       }
1185:     }
1186:     switch (numCorners) {
1187:     case 4: /* do nothing */
1188:     case 2: /* do nothing */
1189:       break;
1190:     case 3: /* from counter-clockwise to wedge ordering */
1191:       l = newCone[1];
1192:       newCone[1] = newCone[2];
1193:       newCone[2] = l;
1194:       break;
1195:     default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported number of corners: %D", numCorners);
1196:     }
1197:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1198:     for (l = 0; l < layers; ++l) {
1199:       PetscInt i;

1201:       for (i = 0; i < numCorners; ++i) {
1202:         newCone[  numCorners + i] = ordExt ? (layers+1)*newCone[i] + l     + numCells :     l*(vEnd - vStart) + newCone[i] + numCells;
1203:         newCone[2*numCorners + i] = ordExt ? (layers+1)*newCone[i] + l + 1 + numCells : (l+1)*(vEnd - vStart) + newCone[i] + numCells;
1204:       }
1205:       DMPlexSetCone(*dm, ordExt ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart, newCone + numCorners);
1206:     }
1207:   }
1208:   DMPlexSymmetrize(*dm);
1209:   DMPlexStratify(*dm);
1210:   PetscFree(newCone);

1212:   cDimB = cDim == dim ? cDim+1 : cDim;
1213:   DMGetCoordinateSection(*dm, &coordSectionB);
1214:   PetscSectionSetNumFields(coordSectionB, 1);
1215:   PetscSectionSetFieldComponents(coordSectionB, 0, cDimB);
1216:   PetscSectionSetChart(coordSectionB, numCells, numCells+numVertices);
1217:   for (v = numCells; v < numCells+numVertices; ++v) {
1218:     PetscSectionSetDof(coordSectionB, v, cDimB);
1219:     PetscSectionSetFieldDof(coordSectionB, v, 0, cDimB);
1220:   }
1221:   PetscSectionSetUp(coordSectionB);
1222:   PetscSectionGetStorageSize(coordSectionB, &coordSize);
1223:   VecCreate(PETSC_COMM_SELF, &coordinatesB);
1224:   PetscObjectSetName((PetscObject) coordinatesB, "coordinates");
1225:   VecSetSizes(coordinatesB, coordSize, PETSC_DETERMINE);
1226:   VecSetBlockSize(coordinatesB, cDimB);
1227:   VecSetType(coordinatesB,VECSTANDARD);

1229:   DMGetCoordinateSection(idm, &coordSectionA);
1230:   DMGetCoordinatesLocal(idm, &coordinatesA);
1231:   VecGetArray(coordinatesB, &coordsB);
1232:   VecGetArrayRead(coordinatesA, &coordsA);
1233:   for (v = vStart; v < vEnd; ++v) {
1234:     const PetscScalar *cptr;
1235:     PetscReal         ones2[2] = { 0., 1.}, ones3[3] = { 0., 0., 1.};
1236:     PetscReal         *normal, norm, h = height/layers;
1237:     PetscInt          offA, d, cDimA = cDim;

1239:     normal = normals ? normals + cDimB*(v - vStart) : (cDim > 1 ? ones3 : ones2);
1240:     if (normals) {
1241:       for (d = 0, norm = 0.0; d < cDimB; ++d) norm += normal[d]*normal[d];
1242:       for (d = 0; d < cDimB; ++d) normal[d] *= 1./PetscSqrtReal(norm);
1243:     }

1245:     PetscSectionGetOffset(coordSectionA, v, &offA);
1246:     cptr = coordsA + offA;
1247:     for (l = 0; l < layers+1; ++l) {
1248:       PetscInt offB, d, newV;

1250:       newV = ordExt ? (layers+1)*(v -vStart) + l + numCells : (vEnd -vStart)*l + (v -vStart) + numCells;
1251:       PetscSectionGetOffset(coordSectionB, newV, &offB);
1252:       for (d = 0; d < cDimA; ++d) { coordsB[offB+d]  = cptr[d]; }
1253:       for (d = 0; d < cDimB; ++d) { coordsB[offB+d] += l ? normal[d]*h : 0.0; }
1254:       cptr    = coordsB + offB;
1255:       cDimA   = cDimB;
1256:     }
1257:   }
1258:   VecRestoreArrayRead(coordinatesA, &coordsA);
1259:   VecRestoreArray(coordinatesB, &coordsB);
1260:   DMSetCoordinatesLocal(*dm, coordinatesB);
1261:   VecDestroy(&coordinatesB);
1262:   PetscFree(normals);
1263:   if (interpolate) {
1264:     DM idm;

1266:     DMPlexInterpolate(*dm, &idm);
1267:     DMPlexCopyCoordinates(*dm, idm);
1268:     DMDestroy(dm);
1269:     *dm  = idm;
1270:   }
1271:   return(0);
1272: }

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

1277:   Logically Collective on DM

1279:   Input Parameters:
1280: + dm - the DM context
1281: - prefix - the prefix to prepend to all option names

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

1287:   Level: advanced

1289: .seealso: SNESSetFromOptions()
1290: @*/
1291: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1292: {
1293:   DM_Plex       *mesh = (DM_Plex *) dm->data;

1298:   PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);
1299:   PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);
1300:   return(0);
1301: }

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

1306:   Collective on MPI_Comm

1308:   Input Parameters:
1309: + comm      - The communicator for the DM object
1310: . numRefine - The number of regular refinements to the basic 5 cell structure
1311: - periodicZ - The boundary type for the Z direction

1313:   Output Parameter:
1314: . dm  - The DM object

1316:   Note: Here is the output numbering looking from the bottom of the cylinder:
1317: $       17-----14
1318: $        |     |
1319: $        |  2  |
1320: $        |     |
1321: $ 17-----8-----7-----14
1322: $  |     |     |     |
1323: $  |  3  |  0  |  1  |
1324: $  |     |     |     |
1325: $ 19-----5-----6-----13
1326: $        |     |
1327: $        |  4  |
1328: $        |     |
1329: $       19-----13
1330: $
1331: $ and up through the top
1332: $
1333: $       18-----16
1334: $        |     |
1335: $        |  2  |
1336: $        |     |
1337: $ 18----10----11-----16
1338: $  |     |     |     |
1339: $  |  3  |  0  |  1  |
1340: $  |     |     |     |
1341: $ 20-----9----12-----15
1342: $        |     |
1343: $        |  4  |
1344: $        |     |
1345: $       20-----15

1347:   Level: beginner

1349: .keywords: DM, create
1350: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1351: @*/
1352: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
1353: {
1354:   const PetscInt dim = 3;
1355:   PetscInt       numCells, numVertices, r;
1356:   PetscMPIInt    rank;

1361:   MPI_Comm_rank(comm, &rank);
1362:   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
1363:   DMCreate(comm, dm);
1364:   DMSetType(*dm, DMPLEX);
1365:   DMSetDimension(*dm, dim);
1366:   /* Create topology */
1367:   {
1368:     PetscInt cone[8], c;

1370:     numCells    = !rank ?  5 : 0;
1371:     numVertices = !rank ? 16 : 0;
1372:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1373:       numCells   *= 3;
1374:       numVertices = !rank ? 24 : 0;
1375:     }
1376:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1377:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 8);}
1378:     DMSetUp(*dm);
1379:     if (!rank) {
1380:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1381:         cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1382:         cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1383:         DMPlexSetCone(*dm, 0, cone);
1384:         cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1385:         cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1386:         DMPlexSetCone(*dm, 1, cone);
1387:         cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1388:         cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1389:         DMPlexSetCone(*dm, 2, cone);
1390:         cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1391:         cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1392:         DMPlexSetCone(*dm, 3, cone);
1393:         cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1394:         cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1395:         DMPlexSetCone(*dm, 4, cone);

1397:         cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1398:         cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1399:         DMPlexSetCone(*dm, 5, cone);
1400:         cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1401:         cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1402:         DMPlexSetCone(*dm, 6, cone);
1403:         cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1404:         cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1405:         DMPlexSetCone(*dm, 7, cone);
1406:         cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1407:         cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1408:         DMPlexSetCone(*dm, 8, cone);
1409:         cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1410:         cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1411:         DMPlexSetCone(*dm, 9, cone);

1413:         cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1414:         cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1415:         DMPlexSetCone(*dm, 10, cone);
1416:         cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1417:         cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1418:         DMPlexSetCone(*dm, 11, cone);
1419:         cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1420:         cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1421:         DMPlexSetCone(*dm, 12, cone);
1422:         cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1423:         cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1424:         DMPlexSetCone(*dm, 13, cone);
1425:         cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1426:         cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1427:         DMPlexSetCone(*dm, 14, cone);
1428:       } else {
1429:         cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
1430:         cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
1431:         DMPlexSetCone(*dm, 0, cone);
1432:         cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
1433:         cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
1434:         DMPlexSetCone(*dm, 1, cone);
1435:         cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
1436:         cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
1437:         DMPlexSetCone(*dm, 2, cone);
1438:         cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
1439:         cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
1440:         DMPlexSetCone(*dm, 3, cone);
1441:         cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
1442:         cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
1443:         DMPlexSetCone(*dm, 4, cone);
1444:       }
1445:     }
1446:     DMPlexSymmetrize(*dm);
1447:     DMPlexStratify(*dm);
1448:   }
1449:   /* Interpolate */
1450:   {
1451:     DM idm;

1453:     DMPlexInterpolate(*dm, &idm);
1454:     DMDestroy(dm);
1455:     *dm  = idm;
1456:   }
1457:   /* Create cube geometry */
1458:   {
1459:     Vec             coordinates;
1460:     PetscSection    coordSection;
1461:     PetscScalar    *coords;
1462:     PetscInt        coordSize, v;
1463:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1464:     const PetscReal ds2 = dis/2.0;

1466:     /* Build coordinates */
1467:     DMGetCoordinateSection(*dm, &coordSection);
1468:     PetscSectionSetNumFields(coordSection, 1);
1469:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1470:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1471:     for (v = numCells; v < numCells+numVertices; ++v) {
1472:       PetscSectionSetDof(coordSection, v, dim);
1473:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1474:     }
1475:     PetscSectionSetUp(coordSection);
1476:     PetscSectionGetStorageSize(coordSection, &coordSize);
1477:     VecCreate(PETSC_COMM_SELF, &coordinates);
1478:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1479:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1480:     VecSetBlockSize(coordinates, dim);
1481:     VecSetType(coordinates,VECSTANDARD);
1482:     VecGetArray(coordinates, &coords);
1483:     if (!rank) {
1484:       coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
1485:       coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
1486:       coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
1487:       coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
1488:       coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
1489:       coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
1490:       coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
1491:       coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
1492:       coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
1493:       coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
1494:       coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
1495:       coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
1496:       coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
1497:       coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
1498:       coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
1499:       coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1500:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1501:         /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1502:         /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1503:         /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1504:         /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1505:         /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1506:         /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1507:         /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1508:         /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1509:       }
1510:     }
1511:     VecRestoreArray(coordinates, &coords);
1512:     DMSetCoordinatesLocal(*dm, coordinates);
1513:     VecDestroy(&coordinates);
1514:   }
1515:   /* Create periodicity */
1516:   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1517:     PetscReal      L[3];
1518:     PetscReal      maxCell[3];
1519:     DMBoundaryType bdType[3];
1520:     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1521:     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1522:     PetscInt       i, numZCells = 3;

1524:     bdType[0] = DM_BOUNDARY_NONE;
1525:     bdType[1] = DM_BOUNDARY_NONE;
1526:     bdType[2] = periodicZ;
1527:     for (i = 0; i < dim; i++) {
1528:       L[i]       = upper[i] - lower[i];
1529:       maxCell[i] = 1.1 * (L[i] / numZCells);
1530:     }
1531:     DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);
1532:   }
1533:   /* Refine topology */
1534:   for (r = 0; r < numRefine; ++r) {
1535:     DM rdm = NULL;

1537:     DMRefine(*dm, comm, &rdm);
1538:     DMDestroy(dm);
1539:     *dm  = rdm;
1540:   }
1541:   /* Remap geometry to cylinder
1542:        Interior square: Linear interpolation is correct
1543:        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1544:        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance

1546:          phi     = arctan(y/x)
1547:          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1548:          d_far   = sqrt(1/2 + sin^2(phi))

1550:        so we remap them using

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

1555:        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1556:   */
1557:   {
1558:     Vec           coordinates;
1559:     PetscSection  coordSection;
1560:     PetscScalar  *coords;
1561:     PetscInt      vStart, vEnd, v;
1562:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1563:     const PetscReal ds2 = 0.5*dis;

1565:     DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);
1566:     DMGetCoordinateSection(*dm, &coordSection);
1567:     DMGetCoordinatesLocal(*dm, &coordinates);
1568:     VecGetArray(coordinates, &coords);
1569:     for (v = vStart; v < vEnd; ++v) {
1570:       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1571:       PetscInt  off;

1573:       PetscSectionGetOffset(coordSection, v, &off);
1574:       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
1575:       x    = PetscRealPart(coords[off]);
1576:       y    = PetscRealPart(coords[off+1]);
1577:       phi  = PetscAtan2Real(y, x);
1578:       sinp = PetscSinReal(phi);
1579:       cosp = PetscCosReal(phi);
1580:       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
1581:         dc = PetscAbsReal(ds2/sinp);
1582:         df = PetscAbsReal(dis/sinp);
1583:         xc = ds2*x/PetscAbsReal(y);
1584:         yc = ds2*PetscSignReal(y);
1585:       } else {
1586:         dc = PetscAbsReal(ds2/cosp);
1587:         df = PetscAbsReal(dis/cosp);
1588:         xc = ds2*PetscSignReal(x);
1589:         yc = ds2*y/PetscAbsReal(x);
1590:       }
1591:       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1592:       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1593:     }
1594:     VecRestoreArray(coordinates, &coords);
1595:     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1596:       DMLocalizeCoordinates(*dm);
1597:     }
1598:   }
1599:   return(0);
1600: }

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

1605:   Collective on MPI_Comm

1607:   Input Parameters:
1608: + comm - The communicator for the DM object
1609: . n    - The number of wedges around the origin
1610: - interpolate - Create edges and faces

1612:   Output Parameter:
1613: . dm  - The DM object

1615:   Level: beginner

1617: .keywords: DM, create
1618: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1619: @*/
1620: PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
1621: {
1622:   const PetscInt dim = 3;
1623:   PetscInt       numCells, numVertices;
1624:   PetscMPIInt    rank;

1629:   MPI_Comm_rank(comm, &rank);
1630:   if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n);
1631:   DMCreate(comm, dm);
1632:   DMSetType(*dm, DMPLEX);
1633:   DMSetDimension(*dm, dim);
1634:   /* Create topology */
1635:   {
1636:     PetscInt cone[6], c;

1638:     numCells    = !rank ?        n : 0;
1639:     numVertices = !rank ?  2*(n+1) : 0;
1640:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1641:     DMPlexSetHybridBounds(*dm, 0, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);
1642:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 6);}
1643:     DMSetUp(*dm);
1644:     for (c = 0; c < numCells; c++) {
1645:       cone[0] =  c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n;
1646:       cone[3] =  c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n;
1647:       DMPlexSetCone(*dm, c, cone);
1648:     }
1649:     DMPlexSymmetrize(*dm);
1650:     DMPlexStratify(*dm);
1651:   }
1652:   /* Interpolate */
1653:   if (interpolate) {
1654:     DM idm;

1656:     DMPlexInterpolate(*dm, &idm);
1657:     DMDestroy(dm);
1658:     *dm  = idm;
1659:   }
1660:   /* Create cylinder geometry */
1661:   {
1662:     Vec          coordinates;
1663:     PetscSection coordSection;
1664:     PetscScalar *coords;
1665:     PetscInt     coordSize, v, c;

1667:     /* Build coordinates */
1668:     DMGetCoordinateSection(*dm, &coordSection);
1669:     PetscSectionSetNumFields(coordSection, 1);
1670:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1671:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1672:     for (v = numCells; v < numCells+numVertices; ++v) {
1673:       PetscSectionSetDof(coordSection, v, dim);
1674:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1675:     }
1676:     PetscSectionSetUp(coordSection);
1677:     PetscSectionGetStorageSize(coordSection, &coordSize);
1678:     VecCreate(PETSC_COMM_SELF, &coordinates);
1679:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1680:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1681:     VecSetBlockSize(coordinates, dim);
1682:     VecSetType(coordinates,VECSTANDARD);
1683:     VecGetArray(coordinates, &coords);
1684:     for (c = 0; c < numCells; c++) {
1685:       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;
1686:       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;
1687:     }
1688:     if (!rank) {
1689:       coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0;
1690:       coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0;
1691:     }
1692:     VecRestoreArray(coordinates, &coords);
1693:     DMSetCoordinatesLocal(*dm, coordinates);
1694:     VecDestroy(&coordinates);
1695:   }
1696:   return(0);
1697: }

1699: PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1700: {
1701:   PetscReal prod = 0.0;
1702:   PetscInt  i;
1703:   for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
1704:   return PetscSqrtReal(prod);
1705: }
1706: PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1707: {
1708:   PetscReal prod = 0.0;
1709:   PetscInt  i;
1710:   for (i = 0; i < dim; ++i) prod += x[i]*y[i];
1711:   return prod;
1712: }

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

1717:   Collective on MPI_Comm

1719:   Input Parameters:
1720: . comm  - The communicator for the DM object
1721: . dim   - The dimension
1722: - simplex - Use simplices, or tensor product cells

1724:   Output Parameter:
1725: . dm  - The DM object

1727:   Level: beginner

1729: .keywords: DM, create
1730: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1731: @*/
1732: PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *dm)
1733: {
1734:   const PetscInt  embedDim = dim+1;
1735:   PetscSection    coordSection;
1736:   Vec             coordinates;
1737:   PetscScalar    *coords;
1738:   PetscReal      *coordsIn;
1739:   PetscInt        numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e;
1740:   PetscMPIInt     rank;
1741:   PetscErrorCode  ierr;

1745:   DMCreate(comm, dm);
1746:   DMSetType(*dm, DMPLEX);
1747:   DMSetDimension(*dm, dim);
1748:   DMSetCoordinateDim(*dm, dim+1);
1749:   MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank);
1750:   switch (dim) {
1751:   case 2:
1752:     if (simplex) {
1753:       DM              idm;
1754:       const PetscReal edgeLen     = 2.0/(1.0 + PETSC_PHI);
1755:       const PetscReal vertex[3]   = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)};
1756:       const PetscInt  degree      = 5;
1757:       PetscInt        s[3]        = {1, 1, 1};
1758:       PetscInt        cone[3];
1759:       PetscInt       *graph, p, i, j, k;

1761:       numCells    = !rank ? 20 : 0;
1762:       numVerts    = !rank ? 12 : 0;
1763:       firstVertex = numCells;
1764:       /* Use icosahedron, which for a unit sphere has coordinates which are all cyclic permutations of

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

1768:          where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
1769:          length is then given by 2/\phi = 2 * 2.73606 = 5.47214.
1770:        */
1771:       /* Construct vertices */
1772:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1773:       for (p = 0, i = 0; p < embedDim; ++p) {
1774:         for (s[1] = -1; s[1] < 2; s[1] += 2) {
1775:           for (s[2] = -1; s[2] < 2; s[2] += 2) {
1776:             for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim];
1777:             ++i;
1778:           }
1779:         }
1780:       }
1781:       /* Construct graph */
1782:       PetscCalloc1(numVerts * numVerts, &graph);
1783:       for (i = 0; i < numVerts; ++i) {
1784:         for (j = 0, k = 0; j < numVerts; ++j) {
1785:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
1786:         }
1787:         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid icosahedron, vertex %D degree %D != %D", i, k, degree);
1788:       }
1789:       /* Build Topology */
1790:       DMPlexSetChart(*dm, 0, numCells+numVerts);
1791:       for (c = 0; c < numCells; c++) {
1792:         DMPlexSetConeSize(*dm, c, embedDim);
1793:       }
1794:       DMSetUp(*dm); /* Allocate space for cones */
1795:       /* Cells */
1796:       for (i = 0, c = 0; i < numVerts; ++i) {
1797:         for (j = 0; j < i; ++j) {
1798:           for (k = 0; k < j; ++k) {
1799:             if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) {
1800:               cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k;
1801:               /* Check orientation */
1802:               {
1803:                 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}}};
1804:                 PetscReal normal[3];
1805:                 PetscInt  e, f;

1807:                 for (d = 0; d < embedDim; ++d) {
1808:                   normal[d] = 0.0;
1809:                   for (e = 0; e < embedDim; ++e) {
1810:                     for (f = 0; f < embedDim; ++f) {
1811:                       normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]);
1812:                     }
1813:                   }
1814:                 }
1815:                 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
1816:               }
1817:               DMPlexSetCone(*dm, c++, cone);
1818:             }
1819:           }
1820:         }
1821:       }
1822:       DMPlexSymmetrize(*dm);
1823:       DMPlexStratify(*dm);
1824:       PetscFree(graph);
1825:       /* Interpolate mesh */
1826:       DMPlexInterpolate(*dm, &idm);
1827:       DMDestroy(dm);
1828:       *dm  = idm;
1829:     } else {
1830:       /*
1831:         12-21--13
1832:          |     |
1833:         25  4  24
1834:          |     |
1835:   12-25--9-16--8-24--13
1836:    |     |     |     |
1837:   23  5 17  0 15  3  22
1838:    |     |     |     |
1839:   10-20--6-14--7-19--11
1840:          |     |
1841:         20  1  19
1842:          |     |
1843:         10-18--11
1844:          |     |
1845:         23  2  22
1846:          |     |
1847:         12-21--13
1848:        */
1849:       const PetscReal dist = 1.0/PetscSqrtReal(3.0);
1850:       PetscInt        cone[4], ornt[4];

1852:       numCells    = !rank ?  6 : 0;
1853:       numEdges    = !rank ? 12 : 0;
1854:       numVerts    = !rank ?  8 : 0;
1855:       firstVertex = numCells;
1856:       firstEdge   = numCells + numVerts;
1857:       /* Build Topology */
1858:       DMPlexSetChart(*dm, 0, numCells+numEdges+numVerts);
1859:       for (c = 0; c < numCells; c++) {
1860:         DMPlexSetConeSize(*dm, c, 4);
1861:       }
1862:       for (e = firstEdge; e < firstEdge+numEdges; ++e) {
1863:         DMPlexSetConeSize(*dm, e, 2);
1864:       }
1865:       DMSetUp(*dm); /* Allocate space for cones */
1866:       /* Cell 0 */
1867:       cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17;
1868:       DMPlexSetCone(*dm, 0, cone);
1869:       ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0;
1870:       DMPlexSetConeOrientation(*dm, 0, ornt);
1871:       /* Cell 1 */
1872:       cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20;
1873:       DMPlexSetCone(*dm, 1, cone);
1874:       ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
1875:       DMPlexSetConeOrientation(*dm, 1, ornt);
1876:       /* Cell 2 */
1877:       cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23;
1878:       DMPlexSetCone(*dm, 2, cone);
1879:       ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0;
1880:       DMPlexSetConeOrientation(*dm, 2, ornt);
1881:       /* Cell 3 */
1882:       cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15;
1883:       DMPlexSetCone(*dm, 3, cone);
1884:       ornt[0] = -2; ornt[1] = -2; ornt[2] = 0; ornt[3] = -2;
1885:       DMPlexSetConeOrientation(*dm, 3, ornt);
1886:       /* Cell 4 */
1887:       cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25;
1888:       DMPlexSetCone(*dm, 4, cone);
1889:       ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = 0;
1890:       DMPlexSetConeOrientation(*dm, 4, ornt);
1891:       /* Cell 5 */
1892:       cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23;
1893:       DMPlexSetCone(*dm, 5, cone);
1894:       ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = -2;
1895:       DMPlexSetConeOrientation(*dm, 5, ornt);
1896:       /* Edges */
1897:       cone[0] =  6; cone[1] =  7;
1898:       DMPlexSetCone(*dm, 14, cone);
1899:       cone[0] =  7; cone[1] =  8;
1900:       DMPlexSetCone(*dm, 15, cone);
1901:       cone[0] =  8; cone[1] =  9;
1902:       DMPlexSetCone(*dm, 16, cone);
1903:       cone[0] =  9; cone[1] =  6;
1904:       DMPlexSetCone(*dm, 17, cone);
1905:       cone[0] = 10; cone[1] = 11;
1906:       DMPlexSetCone(*dm, 18, cone);
1907:       cone[0] = 11; cone[1] =  7;
1908:       DMPlexSetCone(*dm, 19, cone);
1909:       cone[0] =  6; cone[1] = 10;
1910:       DMPlexSetCone(*dm, 20, cone);
1911:       cone[0] = 12; cone[1] = 13;
1912:       DMPlexSetCone(*dm, 21, cone);
1913:       cone[0] = 13; cone[1] = 11;
1914:       DMPlexSetCone(*dm, 22, cone);
1915:       cone[0] = 10; cone[1] = 12;
1916:       DMPlexSetCone(*dm, 23, cone);
1917:       cone[0] = 13; cone[1] =  8;
1918:       DMPlexSetCone(*dm, 24, cone);
1919:       cone[0] = 12; cone[1] =  9;
1920:       DMPlexSetCone(*dm, 25, cone);
1921:       DMPlexSymmetrize(*dm);
1922:       DMPlexStratify(*dm);
1923:       /* Build coordinates */
1924:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1925:       coordsIn[0*embedDim+0] = -dist; coordsIn[0*embedDim+1] =  dist; coordsIn[0*embedDim+2] = -dist;
1926:       coordsIn[1*embedDim+0] =  dist; coordsIn[1*embedDim+1] =  dist; coordsIn[1*embedDim+2] = -dist;
1927:       coordsIn[2*embedDim+0] =  dist; coordsIn[2*embedDim+1] = -dist; coordsIn[2*embedDim+2] = -dist;
1928:       coordsIn[3*embedDim+0] = -dist; coordsIn[3*embedDim+1] = -dist; coordsIn[3*embedDim+2] = -dist;
1929:       coordsIn[4*embedDim+0] = -dist; coordsIn[4*embedDim+1] =  dist; coordsIn[4*embedDim+2] =  dist;
1930:       coordsIn[5*embedDim+0] =  dist; coordsIn[5*embedDim+1] =  dist; coordsIn[5*embedDim+2] =  dist;
1931:       coordsIn[6*embedDim+0] = -dist; coordsIn[6*embedDim+1] = -dist; coordsIn[6*embedDim+2] =  dist;
1932:       coordsIn[7*embedDim+0] =  dist; coordsIn[7*embedDim+1] = -dist; coordsIn[7*embedDim+2] =  dist;
1933:     }
1934:     break;
1935:   case 3:
1936:     if (simplex) {
1937:       DM              idm;
1938:       const PetscReal edgeLen         = 1.0/PETSC_PHI;
1939:       const PetscReal vertexA[4]      = {0.5, 0.5, 0.5, 0.5};
1940:       const PetscReal vertexB[4]      = {1.0, 0.0, 0.0, 0.0};
1941:       const PetscReal vertexC[4]      = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0};
1942:       const PetscInt  degree          = 12;
1943:       PetscInt        s[4]            = {1, 1, 1};
1944:       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},
1945:                                          {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}};
1946:       PetscInt        cone[4];
1947:       PetscInt       *graph, p, i, j, k, l;

1949:       numCells    = !rank ? 600 : 0;
1950:       numVerts    = !rank ? 120 : 0;
1951:       firstVertex = numCells;
1952:       /* Use the 600-cell, which for a unit sphere has coordinates which are

1954:            1/2 (\pm 1, \pm 1,    \pm 1, \pm 1)                          16
1955:                (\pm 1,    0,       0,      0)  all cyclic permutations   8
1956:            1/2 (\pm 1, \pm phi, \pm 1/phi, 0)  all even permutations    96

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

1961:          http://buzzard.pugetsound.edu/sage-practice/ch03s03.html
1962:          http://mathworld.wolfram.com/600-Cell.html
1963:        */
1964:       /* Construct vertices */
1965:       PetscCalloc1(numVerts * embedDim, &coordsIn);
1966:       i    = 0;
1967:       for (s[0] = -1; s[0] < 2; s[0] += 2) {
1968:         for (s[1] = -1; s[1] < 2; s[1] += 2) {
1969:           for (s[2] = -1; s[2] < 2; s[2] += 2) {
1970:             for (s[3] = -1; s[3] < 2; s[3] += 2) {
1971:               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d];
1972:               ++i;
1973:             }
1974:           }
1975:         }
1976:       }
1977:       for (p = 0; p < embedDim; ++p) {
1978:         s[1] = s[2] = s[3] = 1;
1979:         for (s[0] = -1; s[0] < 2; s[0] += 2) {
1980:           for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim];
1981:           ++i;
1982:         }
1983:       }
1984:       for (p = 0; p < 12; ++p) {
1985:         s[3] = 1;
1986:         for (s[0] = -1; s[0] < 2; s[0] += 2) {
1987:           for (s[1] = -1; s[1] < 2; s[1] += 2) {
1988:             for (s[2] = -1; s[2] < 2; s[2] += 2) {
1989:               for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]];
1990:               ++i;
1991:             }
1992:           }
1993:         }
1994:       }
1995:       if (i != numVerts) SETERRQ2(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertices %D != %D", i, numVerts);
1996:       /* Construct graph */
1997:       PetscCalloc1(numVerts * numVerts, &graph);
1998:       for (i = 0; i < numVerts; ++i) {
1999:         for (j = 0, k = 0; j < numVerts; ++j) {
2000:           if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;}
2001:         }
2002:         if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertex %D degree %D != %D", i, k, degree);
2003:       }
2004:       /* Build Topology */
2005:       DMPlexSetChart(*dm, 0, numCells+numVerts);
2006:       for (c = 0; c < numCells; c++) {
2007:         DMPlexSetConeSize(*dm, c, embedDim);
2008:       }
2009:       DMSetUp(*dm); /* Allocate space for cones */
2010:       /* Cells */
2011:       for (i = 0, c = 0; i < numVerts; ++i) {
2012:         for (j = 0; j < i; ++j) {
2013:           for (k = 0; k < j; ++k) {
2014:             for (l = 0; l < k; ++l) {
2015:               if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] &&
2016:                   graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) {
2017:                 cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l;
2018:                 /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */
2019:                 {
2020:                   const PetscInt epsilon[4][4][4][4] = {{{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2021:                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  1}, { 0,  0, -1, 0}},
2022:                                                          {{0,  0,  0,  0}, { 0, 0,  0, -1}, { 0,  0, 0,  0}, { 0,  1,  0, 0}},
2023:                                                          {{0,  0,  0,  0}, { 0, 0,  1,  0}, { 0, -1, 0,  0}, { 0,  0,  0, 0}}},

2025:                                                         {{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0, -1}, { 0,  0,  1, 0}},
2026:                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2027:                                                          {{0,  0,  0,  1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, {-1,  0,  0, 0}},
2028:                                                          {{0,  0, -1,  0}, { 0, 0,  0,  0}, { 1,  0, 0,  0}, { 0,  0,  0, 0}}},

2030:                                                         {{{0,  0,  0,  0}, { 0, 0,  0,  1}, { 0,  0, 0,  0}, { 0, -1,  0, 0}},
2031:                                                          {{0,  0,  0, -1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 1,  0,  0, 0}},
2032:                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2033:                                                          {{0,  1,  0,  0}, {-1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}},

2035:                                                         {{{0,  0,  0,  0}, { 0, 0, -1,  0}, { 0,  1, 0,  0}, { 0,  0,  0, 0}},
2036:                                                          {{0,  0,  1,  0}, { 0, 0,  0,  0}, {-1,  0, 0,  0}, { 0,  0,  0, 0}},
2037:                                                          {{0, -1,  0,  0}, { 1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2038:                                                          {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}}};
2039:                   PetscReal normal[4];
2040:                   PetscInt  e, f, g;

2042:                   for (d = 0; d < embedDim; ++d) {
2043:                     normal[d] = 0.0;
2044:                     for (e = 0; e < embedDim; ++e) {
2045:                       for (f = 0; f < embedDim; ++f) {
2046:                         for (g = 0; g < embedDim; ++g) {
2047:                           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]);
2048:                         }
2049:                       }
2050:                     }
2051:                   }
2052:                   if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;}
2053:                 }
2054:                 DMPlexSetCone(*dm, c++, cone);
2055:               }
2056:             }
2057:           }
2058:         }
2059:       }
2060:       DMPlexSymmetrize(*dm);
2061:       DMPlexStratify(*dm);
2062:       PetscFree(graph);
2063:       /* Interpolate mesh */
2064:       DMPlexInterpolate(*dm, &idm);
2065:       DMDestroy(dm);
2066:       *dm  = idm;
2067:       break;
2068:     }
2069:   default: SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim);
2070:   }
2071:   /* Create coordinates */
2072:   DMGetCoordinateSection(*dm, &coordSection);
2073:   PetscSectionSetNumFields(coordSection, 1);
2074:   PetscSectionSetFieldComponents(coordSection, 0, embedDim);
2075:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);
2076:   for (v = firstVertex; v < firstVertex+numVerts; ++v) {
2077:     PetscSectionSetDof(coordSection, v, embedDim);
2078:     PetscSectionSetFieldDof(coordSection, v, 0, embedDim);
2079:   }
2080:   PetscSectionSetUp(coordSection);
2081:   PetscSectionGetStorageSize(coordSection, &coordSize);
2082:   VecCreate(PETSC_COMM_SELF, &coordinates);
2083:   VecSetBlockSize(coordinates, embedDim);
2084:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2085:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2086:   VecSetType(coordinates,VECSTANDARD);
2087:   VecGetArray(coordinates, &coords);
2088:   for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];}
2089:   VecRestoreArray(coordinates, &coords);
2090:   DMSetCoordinatesLocal(*dm, coordinates);
2091:   VecDestroy(&coordinates);
2092:   PetscFree(coordsIn);
2093:   return(0);
2094: }

2096: /* External function declarations here */
2097: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
2098: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2099: extern PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2100: extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
2101: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
2102: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
2103: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
2104: extern PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
2105: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
2106: extern PetscErrorCode DMSetUp_Plex(DM dm);
2107: extern PetscErrorCode DMDestroy_Plex(DM dm);
2108: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
2109: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
2110: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
2111: extern PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);
2112: static PetscErrorCode DMInitialize_Plex(DM dm);

2114: /* Replace dm with the contents of dmNew
2115:    - Share the DM_Plex structure
2116:    - Share the coordinates
2117:    - Share the SF
2118: */
2119: static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
2120: {
2121:   PetscSF               sf;
2122:   DM                    coordDM, coarseDM;
2123:   Vec                   coords;
2124:   PetscBool             isper;
2125:   const PetscReal      *maxCell, *L;
2126:   const DMBoundaryType *bd;
2127:   PetscErrorCode        ierr;

2130:   DMGetPointSF(dmNew, &sf);
2131:   DMSetPointSF(dm, sf);
2132:   DMGetCoordinateDM(dmNew, &coordDM);
2133:   DMGetCoordinatesLocal(dmNew, &coords);
2134:   DMSetCoordinateDM(dm, coordDM);
2135:   DMSetCoordinatesLocal(dm, coords);
2136:   DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);
2137:   DMSetPeriodicity(dmNew, isper, maxCell, L, bd);
2138:   DMDestroy_Plex(dm);
2139:   DMInitialize_Plex(dm);
2140:   dm->data = dmNew->data;
2141:   ((DM_Plex *) dmNew->data)->refct++;
2142:   dmNew->labels->refct++;
2143:   if (!--(dm->labels->refct)) {
2144:     DMLabelLink next = dm->labels->next;

2146:     /* destroy the labels */
2147:     while (next) {
2148:       DMLabelLink tmp = next->next;

2150:       DMLabelDestroy(&next->label);
2151:       PetscFree(next);
2152:       next = tmp;
2153:     }
2154:     PetscFree(dm->labels);
2155:   }
2156:   dm->labels = dmNew->labels;
2157:   dm->depthLabel = dmNew->depthLabel;
2158:   DMGetCoarseDM(dmNew,&coarseDM);
2159:   DMSetCoarseDM(dm,coarseDM);
2160:   return(0);
2161: }

2163: /* Swap dm with the contents of dmNew
2164:    - Swap the DM_Plex structure
2165:    - Swap the coordinates
2166:    - Swap the point PetscSF
2167: */
2168: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
2169: {
2170:   DM              coordDMA, coordDMB;
2171:   Vec             coordsA,  coordsB;
2172:   PetscSF         sfA,      sfB;
2173:   void            *tmp;
2174:   DMLabelLinkList listTmp;
2175:   DMLabel         depthTmp;
2176:   PetscErrorCode  ierr;

2179:   DMGetPointSF(dmA, &sfA);
2180:   DMGetPointSF(dmB, &sfB);
2181:   PetscObjectReference((PetscObject) sfA);
2182:   DMSetPointSF(dmA, sfB);
2183:   DMSetPointSF(dmB, sfA);
2184:   PetscObjectDereference((PetscObject) sfA);

2186:   DMGetCoordinateDM(dmA, &coordDMA);
2187:   DMGetCoordinateDM(dmB, &coordDMB);
2188:   PetscObjectReference((PetscObject) coordDMA);
2189:   DMSetCoordinateDM(dmA, coordDMB);
2190:   DMSetCoordinateDM(dmB, coordDMA);
2191:   PetscObjectDereference((PetscObject) coordDMA);

2193:   DMGetCoordinatesLocal(dmA, &coordsA);
2194:   DMGetCoordinatesLocal(dmB, &coordsB);
2195:   PetscObjectReference((PetscObject) coordsA);
2196:   DMSetCoordinatesLocal(dmA, coordsB);
2197:   DMSetCoordinatesLocal(dmB, coordsA);
2198:   PetscObjectDereference((PetscObject) coordsA);

2200:   tmp       = dmA->data;
2201:   dmA->data = dmB->data;
2202:   dmB->data = tmp;
2203:   listTmp   = dmA->labels;
2204:   dmA->labels = dmB->labels;
2205:   dmB->labels = listTmp;
2206:   depthTmp  = dmA->depthLabel;
2207:   dmA->depthLabel = dmB->depthLabel;
2208:   dmB->depthLabel = depthTmp;
2209:   return(0);
2210: }

2212: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2213: {
2214:   DM_Plex       *mesh = (DM_Plex*) dm->data;

2218:   /* Handle viewing */
2219:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL);
2220:   PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL);
2221:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL);
2222:   PetscOptionsInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL);
2223:   /* Point Location */
2224:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL);
2225:   /* Partitioning and distribution */
2226:   PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL);
2227:   /* Generation and remeshing */
2228:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL);
2229:   /* Projection behavior */
2230:   PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);
2231:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);

2233:   PetscPartitionerSetFromOptions(mesh->partitioner);
2234:   return(0);
2235: }

2237: static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2238: {
2239:   PetscInt       refine = 0, coarsen = 0, r;
2240:   PetscBool      isHierarchy;

2245:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
2246:   /* Handle DMPlex refinement */
2247:   PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);
2248:   PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);
2249:   if (refine) {DMPlexSetRefinementUniform(dm, PETSC_TRUE);}
2250:   if (refine && isHierarchy) {
2251:     DM *dms, coarseDM;

2253:     DMGetCoarseDM(dm, &coarseDM);
2254:     PetscObjectReference((PetscObject)coarseDM);
2255:     PetscMalloc1(refine,&dms);
2256:     DMRefineHierarchy(dm, refine, dms);
2257:     /* Total hack since we do not pass in a pointer */
2258:     DMPlexSwap_Static(dm, dms[refine-1]);
2259:     if (refine == 1) {
2260:       DMSetCoarseDM(dm, dms[0]);
2261:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2262:     } else {
2263:       DMSetCoarseDM(dm, dms[refine-2]);
2264:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2265:       DMSetCoarseDM(dms[0], dms[refine-1]);
2266:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
2267:     }
2268:     DMSetCoarseDM(dms[refine-1], coarseDM);
2269:     PetscObjectDereference((PetscObject)coarseDM);
2270:     /* Free DMs */
2271:     for (r = 0; r < refine; ++r) {
2272:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2273:       DMDestroy(&dms[r]);
2274:     }
2275:     PetscFree(dms);
2276:   } else {
2277:     for (r = 0; r < refine; ++r) {
2278:       DM refinedMesh;

2280:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2281:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);
2282:       /* Total hack since we do not pass in a pointer */
2283:       DMPlexReplace_Static(dm, refinedMesh);
2284:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2285:       DMDestroy(&refinedMesh);
2286:     }
2287:   }
2288:   /* Handle DMPlex coarsening */
2289:   PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);
2290:   PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);
2291:   if (coarsen && isHierarchy) {
2292:     DM *dms;

2294:     PetscMalloc1(coarsen, &dms);
2295:     DMCoarsenHierarchy(dm, coarsen, dms);
2296:     /* Free DMs */
2297:     for (r = 0; r < coarsen; ++r) {
2298:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2299:       DMDestroy(&dms[r]);
2300:     }
2301:     PetscFree(dms);
2302:   } else {
2303:     for (r = 0; r < coarsen; ++r) {
2304:       DM coarseMesh;

2306:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2307:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);
2308:       /* Total hack since we do not pass in a pointer */
2309:       DMPlexReplace_Static(dm, coarseMesh);
2310:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2311:       DMDestroy(&coarseMesh);
2312:     }
2313:   }
2314:   /* Handle */
2315:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2316:   PetscOptionsTail();
2317:   return(0);
2318: }

2320: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
2321: {

2325:   DMCreateGlobalVector_Section_Private(dm,vec);
2326:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
2327:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
2328:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
2329:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
2330:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
2331:   return(0);
2332: }

2334: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
2335: {

2339:   DMCreateLocalVector_Section_Private(dm,vec);
2340:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
2341:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
2342:   return(0);
2343: }

2345: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
2346: {
2347:   PetscInt       depth, d;

2351:   DMPlexGetDepth(dm, &depth);
2352:   if (depth == 1) {
2353:     DMGetDimension(dm, &d);
2354:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
2355:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
2356:     else               {*pStart = 0; *pEnd = 0;}
2357:   } else {
2358:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
2359:   }
2360:   return(0);
2361: }

2363: static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
2364: {
2365:   PetscSF        sf;

2369:   DMGetPointSF(dm, &sf);
2370:   PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);
2371:   return(0);
2372: }

2374: static PetscErrorCode DMHasCreateInjection_Plex(DM dm, PetscBool *flg)
2375: {
2379:   *flg = PETSC_TRUE;
2380:   return(0);
2381: }

2383: static PetscErrorCode DMInitialize_Plex(DM dm)
2384: {

2388:   dm->ops->view                            = DMView_Plex;
2389:   dm->ops->load                            = DMLoad_Plex;
2390:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
2391:   dm->ops->clone                           = DMClone_Plex;
2392:   dm->ops->setup                           = DMSetUp_Plex;
2393:   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
2394:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
2395:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
2396:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
2397:   dm->ops->getlocaltoglobalmapping         = NULL;
2398:   dm->ops->createfieldis                   = NULL;
2399:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
2400:   dm->ops->createcoordinatefield           = DMCreateCoordinateField_Plex;
2401:   dm->ops->getcoloring                     = NULL;
2402:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
2403:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
2404:   dm->ops->createmassmatrix                = DMCreateMassMatrix_Plex;
2405:   dm->ops->getaggregates                   = NULL;
2406:   dm->ops->getinjection                    = DMCreateInjection_Plex;
2407:   dm->ops->hascreateinjection              = DMHasCreateInjection_Plex;
2408:   dm->ops->refine                          = DMRefine_Plex;
2409:   dm->ops->coarsen                         = DMCoarsen_Plex;
2410:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
2411:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
2412:   dm->ops->adaptlabel                      = DMAdaptLabel_Plex;
2413:   dm->ops->adaptmetric                     = DMAdaptMetric_Plex;
2414:   dm->ops->globaltolocalbegin              = NULL;
2415:   dm->ops->globaltolocalend                = NULL;
2416:   dm->ops->localtoglobalbegin              = NULL;
2417:   dm->ops->localtoglobalend                = NULL;
2418:   dm->ops->destroy                         = DMDestroy_Plex;
2419:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
2420:   dm->ops->createsuperdm                   = DMCreateSuperDM_Plex;
2421:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
2422:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
2423:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
2424:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
2425:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
2426:   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
2427:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
2428:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
2429:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
2430:   dm->ops->getneighbors                    = DMGetNeighors_Plex;
2431:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);
2432:   PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);
2433:   return(0);
2434: }

2436: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
2437: {
2438:   DM_Plex        *mesh = (DM_Plex *) dm->data;

2442:   mesh->refct++;
2443:   (*newdm)->data = mesh;
2444:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
2445:   DMInitialize_Plex(*newdm);
2446:   return(0);
2447: }

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

2455:   Options Database Keys:
2456: + -dm_view :mesh.tex:ascii_latex - View the mesh in LaTeX/TikZ
2457: . -dm_plex_view_scale <num>      - Scale the TikZ
2458: - -dm_plex_print_fem <num>       - View FEM assembly information, such as element vectors and matrices


2461:   Level: intermediate

2463: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
2464: M*/

2466: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
2467: {
2468:   DM_Plex        *mesh;
2469:   PetscInt       unit, d;

2474:   PetscNewLog(dm,&mesh);
2475:   dm->dim  = 0;
2476:   dm->data = mesh;

2478:   mesh->refct             = 1;
2479:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
2480:   mesh->maxConeSize       = 0;
2481:   mesh->cones             = NULL;
2482:   mesh->coneOrientations  = NULL;
2483:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
2484:   mesh->maxSupportSize    = 0;
2485:   mesh->supports          = NULL;
2486:   mesh->refinementUniform = PETSC_TRUE;
2487:   mesh->refinementLimit   = -1.0;

2489:   mesh->facesTmp = NULL;

2491:   mesh->tetgenOpts   = NULL;
2492:   mesh->triangleOpts = NULL;
2493:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
2494:   mesh->remeshBd     = PETSC_FALSE;

2496:   mesh->subpointMap = NULL;

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

2500:   mesh->regularRefinement   = PETSC_FALSE;
2501:   mesh->depthState          = -1;
2502:   mesh->globalVertexNumbers = NULL;
2503:   mesh->globalCellNumbers   = NULL;
2504:   mesh->anchorSection       = NULL;
2505:   mesh->anchorIS            = NULL;
2506:   mesh->createanchors       = NULL;
2507:   mesh->computeanchormatrix = NULL;
2508:   mesh->parentSection       = NULL;
2509:   mesh->parents             = NULL;
2510:   mesh->childIDs            = NULL;
2511:   mesh->childSection        = NULL;
2512:   mesh->children            = NULL;
2513:   mesh->referenceTree       = NULL;
2514:   mesh->getchildsymmetry    = NULL;
2515:   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
2516:   mesh->vtkCellHeight       = 0;
2517:   mesh->useAnchors          = PETSC_FALSE;

2519:   mesh->maxProjectionHeight = 0;

2521:   mesh->printSetValues = PETSC_FALSE;
2522:   mesh->printFEM       = 0;
2523:   mesh->printTol       = 1.0e-10;

2525:   DMInitialize_Plex(dm);
2526:   return(0);
2527: }

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

2532:   Collective on MPI_Comm

2534:   Input Parameter:
2535: . comm - The communicator for the DMPlex object

2537:   Output Parameter:
2538: . mesh  - The DMPlex object

2540:   Level: beginner

2542: .keywords: DMPlex, create
2543: @*/
2544: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
2545: {

2550:   DMCreate(comm, mesh);
2551:   DMSetType(*mesh, DMPLEX);
2552:   return(0);
2553: }

2555: /*
2556:   This takes as input the common mesh generator output, a list of the vertices for each cell, but vertex numbers are global and an SF is built for them
2557: */
2558: /* TODO: invertCells and spaceDim arguments could be added also to to DMPlexCreateFromCellListParallel(), DMPlexBuildFromCellList_Internal() and DMPlexCreateFromCellList() */
2559: PetscErrorCode DMPlexBuildFromCellList_Parallel_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscBool invertCells, PetscSF *sfVert)
2560: {
2561:   PetscSF         sfPoint;
2562:   PetscLayout     vLayout;
2563:   PetscHSetI      vhash;
2564:   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
2565:   const PetscInt *vrange;
2566:   PetscInt        numVerticesAdj, off = 0, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
2567:   PetscMPIInt     rank, size;
2568:   PetscErrorCode  ierr;

2571:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
2572:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);
2573:   /* Partition vertices */
2574:   PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);
2575:   PetscLayoutSetLocalSize(vLayout, numVertices);
2576:   PetscLayoutSetBlockSize(vLayout, 1);
2577:   PetscLayoutSetUp(vLayout);
2578:   PetscLayoutGetRanges(vLayout, &vrange);
2579:   /* Count vertices and map them to procs */
2580:   PetscHSetICreate(&vhash);
2581:   for (c = 0; c < numCells; ++c) {
2582:     for (p = 0; p < numCorners; ++p) {
2583:       PetscHSetIAdd(vhash, cells[c*numCorners+p]);
2584:     }
2585:   }
2586:   PetscHSetIGetSize(vhash, &numVerticesAdj);
2587:   PetscMalloc1(numVerticesAdj, &verticesAdj);
2588:   PetscHSetIGetElems(vhash, &off, verticesAdj);
2589:   PetscHSetIDestroy(&vhash);
2590:   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
2591:   PetscSortInt(numVerticesAdj, verticesAdj);
2592:   PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);
2593:   for (v = 0; v < numVerticesAdj; ++v) {
2594:     const PetscInt gv = verticesAdj[v];
2595:     PetscInt       vrank;

2597:     PetscFindInt(gv, size+1, vrange, &vrank);
2598:     vrank = vrank < 0 ? -(vrank+2) : vrank;
2599:     remoteVerticesAdj[v].index = gv - vrange[vrank];
2600:     remoteVerticesAdj[v].rank  = vrank;
2601:   }
2602:   /* Create cones */
2603:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
2604:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
2605:   DMSetUp(dm);
2606:   DMGetWorkArray(dm, numCorners, MPIU_INT, &cone);
2607:   for (c = 0; c < numCells; ++c) {
2608:     for (p = 0; p < numCorners; ++p) {
2609:       const PetscInt gv = cells[c*numCorners+p];
2610:       PetscInt       lv;

2612:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
2613:       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
2614:       cone[p] = lv+numCells;
2615:     }
2616:     if (invertCells) { DMPlexInvertCell_Internal(spaceDim, numCorners, cone); }
2617:     DMPlexSetCone(dm, c, cone);
2618:   }
2619:   DMRestoreWorkArray(dm, numCorners, MPIU_INT, &cone);
2620:   /* Create SF for vertices */
2621:   PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);
2622:   PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");
2623:   PetscSFSetFromOptions(*sfVert);
2624:   PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);
2625:   PetscFree(verticesAdj);
2626:   /* Build pointSF */
2627:   PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);
2628:   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
2629:   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
2630:   PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
2631:   PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
2632:   for (v = 0; v < numVertices;    ++v) if (vertexOwner[v].rank < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global vertex %d on rank %d was unclaimed", v + vrange[rank], rank);
2633:   PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);
2634:   PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);
2635:   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
2636:   PetscMalloc1(numVerticesGhost, &localVertex);
2637:   PetscMalloc1(numVerticesGhost, &remoteVertex);
2638:   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
2639:     if (vertexLocal[v].rank != rank) {
2640:       localVertex[g]        = v+numCells;
2641:       remoteVertex[g].index = vertexLocal[v].index;
2642:       remoteVertex[g].rank  = vertexLocal[v].rank;
2643:       ++g;
2644:     }
2645:   }
2646:   PetscFree2(vertexLocal, vertexOwner);
2647:   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
2648:   DMGetPointSF(dm, &sfPoint);
2649:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
2650:   PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);
2651:   PetscLayoutDestroy(&vLayout);
2652:   /* Fill in the rest of the topology structure */
2653:   DMPlexSymmetrize(dm);
2654:   DMPlexStratify(dm);
2655:   return(0);
2656: }

2658: /*
2659:   This takes as input the coordinates for each owned vertex
2660: */
2661: PetscErrorCode DMPlexBuildCoordinates_Parallel_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[])
2662: {
2663:   PetscSection   coordSection;
2664:   Vec            coordinates;
2665:   PetscScalar   *coords;
2666:   PetscInt       numVertices, numVerticesAdj, coordSize, v;

2670:   DMSetCoordinateDim(dm, spaceDim);
2671:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
2672:   DMGetCoordinateSection(dm, &coordSection);
2673:   PetscSectionSetNumFields(coordSection, 1);
2674:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
2675:   PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);
2676:   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
2677:     PetscSectionSetDof(coordSection, v, spaceDim);
2678:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
2679:   }
2680:   PetscSectionSetUp(coordSection);
2681:   PetscSectionGetStorageSize(coordSection, &coordSize);
2682:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
2683:   VecSetBlockSize(coordinates, spaceDim);
2684:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2685:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2686:   VecSetType(coordinates,VECSTANDARD);
2687:   VecGetArray(coordinates, &coords);
2688:   {
2689:     MPI_Datatype coordtype;

2691:     /* Need a temp buffer for coords if we have complex/single */
2692:     MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);
2693:     MPI_Type_commit(&coordtype);
2694: #if defined(PETSC_USE_COMPLEX)
2695:     {
2696:     PetscScalar *svertexCoords;
2697:     PetscInt    i;
2698:     PetscMalloc1(numV*spaceDim,&svertexCoords);
2699:     for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
2700:     PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);
2701:     PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);
2702:     PetscFree(svertexCoords);
2703:     }
2704: #else
2705:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);
2706:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);
2707: #endif
2708:     MPI_Type_free(&coordtype);
2709:   }
2710:   VecRestoreArray(coordinates, &coords);
2711:   DMSetCoordinatesLocal(dm, coordinates);
2712:   VecDestroy(&coordinates);
2713:   return(0);
2714: }

2716: /*@C
2717:   DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM

2719:   Input Parameters:
2720: + comm - The communicator
2721: . dim - The topological dimension of the mesh
2722: . numCells - The number of cells owned by this process
2723: . numVertices - The number of vertices owned by this process
2724: . numCorners - The number of vertices for each cell
2725: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
2726: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
2727: . spaceDim - The spatial dimension used for coordinates
2728: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

2730:   Output Parameter:
2731: + dm - The DM
2732: - vertexSF - Optional, SF describing complete vertex ownership

2734:   Note: Two triangles sharing a face
2735: $
2736: $        2
2737: $      / | \
2738: $     /  |  \
2739: $    /   |   \
2740: $   0  0 | 1  3
2741: $    \   |   /
2742: $     \  |  /
2743: $      \ | /
2744: $        1
2745: would have input
2746: $  numCells = 2, numVertices = 4
2747: $  cells = [0 1 2  1 3 2]
2748: $
2749: which would result in the DMPlex
2750: $
2751: $        4
2752: $      / | \
2753: $     /  |  \
2754: $    /   |   \
2755: $   2  0 | 1  5
2756: $    \   |   /
2757: $     \  |  /
2758: $      \ | /
2759: $        3

2761:   Level: beginner

2763: .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
2764: @*/
2765: PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, DM *dm)
2766: {
2767:   PetscSF        sfVert;

2771:   DMCreate(comm, dm);
2772:   DMSetType(*dm, DMPLEX);
2775:   DMSetDimension(*dm, dim);
2776:   DMPlexBuildFromCellList_Parallel_Internal(*dm, spaceDim, numCells, numVertices, numCorners, cells, PETSC_FALSE, &sfVert);
2777:   if (interpolate) {
2778:     DM idm;

2780:     DMPlexInterpolate(*dm, &idm);
2781:     DMDestroy(dm);
2782:     *dm  = idm;
2783:   }
2784:   DMPlexBuildCoordinates_Parallel_Internal(*dm, spaceDim, numCells, numVertices, sfVert, vertexCoords);
2785:   if (vertexSF) *vertexSF = sfVert;
2786:   else {PetscSFDestroy(&sfVert);}
2787:   return(0);
2788: }

2790: /*
2791:   This takes as input the common mesh generator output, a list of the vertices for each cell
2792: */
2793: PetscErrorCode DMPlexBuildFromCellList_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscBool invertCells)
2794: {
2795:   PetscInt      *cone, c, p;

2799:   DMPlexSetChart(dm, 0, numCells+numVertices);
2800:   for (c = 0; c < numCells; ++c) {
2801:     DMPlexSetConeSize(dm, c, numCorners);
2802:   }
2803:   DMSetUp(dm);
2804:   DMGetWorkArray(dm, numCorners, MPIU_INT, &cone);
2805:   for (c = 0; c < numCells; ++c) {
2806:     for (p = 0; p < numCorners; ++p) {
2807:       cone[p] = cells[c*numCorners+p]+numCells;
2808:     }
2809:     if (invertCells) { DMPlexInvertCell_Internal(spaceDim, numCorners, cone); }
2810:     DMPlexSetCone(dm, c, cone);
2811:   }
2812:   DMRestoreWorkArray(dm, numCorners, MPIU_INT, &cone);
2813:   DMPlexSymmetrize(dm);
2814:   DMPlexStratify(dm);
2815:   return(0);
2816: }

2818: /*
2819:   This takes as input the coordinates for each vertex
2820: */
2821: PetscErrorCode DMPlexBuildCoordinates_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
2822: {
2823:   PetscSection   coordSection;
2824:   Vec            coordinates;
2825:   DM             cdm;
2826:   PetscScalar   *coords;
2827:   PetscInt       v, d;

2831:   DMSetCoordinateDim(dm, spaceDim);
2832:   DMGetCoordinateSection(dm, &coordSection);
2833:   PetscSectionSetNumFields(coordSection, 1);
2834:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
2835:   PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
2836:   for (v = numCells; v < numCells+numVertices; ++v) {
2837:     PetscSectionSetDof(coordSection, v, spaceDim);
2838:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
2839:   }
2840:   PetscSectionSetUp(coordSection);

2842:   DMGetCoordinateDM(dm, &cdm);
2843:   DMCreateLocalVector(cdm, &coordinates);
2844:   VecSetBlockSize(coordinates, spaceDim);
2845:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2846:   VecGetArray(coordinates, &coords);
2847:   for (v = 0; v < numVertices; ++v) {
2848:     for (d = 0; d < spaceDim; ++d) {
2849:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
2850:     }
2851:   }
2852:   VecRestoreArray(coordinates, &coords);
2853:   DMSetCoordinatesLocal(dm, coordinates);
2854:   VecDestroy(&coordinates);
2855:   return(0);
2856: }

2858: /*@C
2859:   DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM

2861:   Input Parameters:
2862: + comm - The communicator
2863: . dim - The topological dimension of the mesh
2864: . numCells - The number of cells
2865: . numVertices - The number of vertices
2866: . numCorners - The number of vertices for each cell
2867: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
2868: . cells - An array of numCells*numCorners numbers, the vertices for each cell
2869: . spaceDim - The spatial dimension used for coordinates
2870: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

2872:   Output Parameter:
2873: . dm - The DM

2875:   Note: Two triangles sharing a face
2876: $
2877: $        2
2878: $      / | \
2879: $     /  |  \
2880: $    /   |   \
2881: $   0  0 | 1  3
2882: $    \   |   /
2883: $     \  |  /
2884: $      \ | /
2885: $        1
2886: would have input
2887: $  numCells = 2, numVertices = 4
2888: $  cells = [0 1 2  1 3 2]
2889: $
2890: which would result in the DMPlex
2891: $
2892: $        4
2893: $      / | \
2894: $     /  |  \
2895: $    /   |   \
2896: $   2  0 | 1  5
2897: $    \   |   /
2898: $     \  |  /
2899: $      \ | /
2900: $        3

2902:   Level: beginner

2904: .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
2905: @*/
2906: PetscErrorCode DMPlexCreateFromCellList(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm)
2907: {

2911:   DMCreate(comm, dm);
2912:   DMSetType(*dm, DMPLEX);
2913:   DMSetDimension(*dm, dim);
2914:   DMPlexBuildFromCellList_Internal(*dm, spaceDim, numCells, numVertices, numCorners, cells, PETSC_FALSE);
2915:   if (interpolate) {
2916:     DM idm;

2918:     DMPlexInterpolate(*dm, &idm);
2919:     DMDestroy(dm);
2920:     *dm  = idm;
2921:   }
2922:   DMPlexBuildCoordinates_Internal(*dm, spaceDim, numCells, numVertices, vertexCoords);
2923:   return(0);
2924: }

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

2929:   Input Parameters:
2930: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
2931: . depth - The depth of the DAG
2932: . numPoints - The number of points at each depth
2933: . coneSize - The cone size of each point
2934: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
2935: . coneOrientations - The orientation of each cone point
2936: - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex

2938:   Output Parameter:
2939: . dm - The DM

2941:   Note: Two triangles sharing a face would have input
2942: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
2943: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
2944: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
2945: $
2946: which would result in the DMPlex
2947: $
2948: $        4
2949: $      / | \
2950: $     /  |  \
2951: $    /   |   \
2952: $   2  0 | 1  5
2953: $    \   |   /
2954: $     \  |  /
2955: $      \ | /
2956: $        3
2957: $
2958: $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()

2960:   Level: advanced

2962: .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
2963: @*/
2964: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
2965: {
2966:   Vec            coordinates;
2967:   PetscSection   coordSection;
2968:   PetscScalar    *coords;
2969:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

2973:   DMGetDimension(dm, &dim);
2974:   DMGetCoordinateDim(dm, &dimEmbed);
2975:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
2976:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
2977:   DMPlexSetChart(dm, pStart, pEnd);
2978:   for (p = pStart; p < pEnd; ++p) {
2979:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
2980:     if (firstVertex < 0 && !coneSize[p - pStart]) {
2981:       firstVertex = p - pStart;
2982:     }
2983:   }
2984:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
2985:   DMSetUp(dm); /* Allocate space for cones */
2986:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
2987:     DMPlexSetCone(dm, p, &cones[off]);
2988:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
2989:   }
2990:   DMPlexSymmetrize(dm);
2991:   DMPlexStratify(dm);
2992:   /* Build coordinates */
2993:   DMGetCoordinateSection(dm, &coordSection);
2994:   PetscSectionSetNumFields(coordSection, 1);
2995:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
2996:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
2997:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
2998:     PetscSectionSetDof(coordSection, v, dimEmbed);
2999:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
3000:   }
3001:   PetscSectionSetUp(coordSection);
3002:   PetscSectionGetStorageSize(coordSection, &coordSize);
3003:   VecCreate(PETSC_COMM_SELF, &coordinates);
3004:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3005:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3006:   VecSetBlockSize(coordinates, dimEmbed);
3007:   VecSetType(coordinates,VECSTANDARD);
3008:   VecGetArray(coordinates, &coords);
3009:   for (v = 0; v < numPoints[0]; ++v) {
3010:     PetscInt off;

3012:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
3013:     for (d = 0; d < dimEmbed; ++d) {
3014:       coords[off+d] = vertexCoords[v*dimEmbed+d];
3015:     }
3016:   }
3017:   VecRestoreArray(coordinates, &coords);
3018:   DMSetCoordinatesLocal(dm, coordinates);
3019:   VecDestroy(&coordinates);
3020:   return(0);
3021: }

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

3026: + comm        - The MPI communicator
3027: . filename    - Name of the .dat file
3028: - interpolate - Create faces and edges in the mesh

3030:   Output Parameter:
3031: . dm  - The DM object representing the mesh

3033:   Note: The format is the simplest possible:
3034: $ Ne
3035: $ v0 v1 ... vk
3036: $ Nv
3037: $ x y z marker

3039:   Level: beginner

3041: .seealso: DMPlexCreateFromFile(), DMPlexCreateMedFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
3042: @*/
3043: PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3044: {
3045:   DMLabel         marker;
3046:   PetscViewer     viewer;
3047:   Vec             coordinates;
3048:   PetscSection    coordSection;
3049:   PetscScalar    *coords;
3050:   char            line[PETSC_MAX_PATH_LEN];
3051:   PetscInt        dim = 3, cdim = 3, coordSize, v, c, d;
3052:   PetscMPIInt     rank;
3053:   int             snum, Nv, Nc;
3054:   PetscErrorCode  ierr;

3057:   MPI_Comm_rank(comm, &rank);
3058:   PetscViewerCreate(comm, &viewer);
3059:   PetscViewerSetType(viewer, PETSCVIEWERASCII);
3060:   PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3061:   PetscViewerFileSetName(viewer, filename);
3062:   if (!rank) {
3063:     PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);
3064:     snum = sscanf(line, "%d %d", &Nc, &Nv);
3065:     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3066:   } else {
3067:     Nc = Nv = 0;
3068:   }
3069:   DMCreate(comm, dm);
3070:   DMSetType(*dm, DMPLEX);
3071:   DMPlexSetChart(*dm, 0, Nc+Nv);
3072:   DMSetDimension(*dm, dim);
3073:   DMSetCoordinateDim(*dm, cdim);
3074:   /* Read topology */
3075:   if (!rank) {
3076:     PetscInt cone[8], corners = 8;
3077:     int      vbuf[8], v;

3079:     for (c = 0; c < Nc; ++c) {DMPlexSetConeSize(*dm, c, corners);}
3080:     DMSetUp(*dm);
3081:     for (c = 0; c < Nc; ++c) {
3082:       PetscViewerRead(viewer, line, corners, NULL, PETSC_STRING);
3083:       snum = sscanf(line, "%d %d %d %d %d %d %d %d", &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5], &vbuf[6], &vbuf[7]);
3084:       if (snum != corners) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3085:       for (v = 0; v < corners; ++v) cone[v] = vbuf[v] + Nc;
3086:       /* Hexahedra are inverted */
3087:       {
3088:         PetscInt tmp = cone[1];
3089:         cone[1] = cone[3];
3090:         cone[3] = tmp;
3091:       }
3092:       DMPlexSetCone(*dm, c, cone);
3093:     }
3094:   }
3095:   DMPlexSymmetrize(*dm);
3096:   DMPlexStratify(*dm);
3097:   /* Read coordinates */
3098:   DMGetCoordinateSection(*dm, &coordSection);
3099:   PetscSectionSetNumFields(coordSection, 1);
3100:   PetscSectionSetFieldComponents(coordSection, 0, cdim);
3101:   PetscSectionSetChart(coordSection, Nc, Nc + Nv);
3102:   for (v = Nc; v < Nc+Nv; ++v) {
3103:     PetscSectionSetDof(coordSection, v, cdim);
3104:     PetscSectionSetFieldDof(coordSection, v, 0, cdim);
3105:   }
3106:   PetscSectionSetUp(coordSection);
3107:   PetscSectionGetStorageSize(coordSection, &coordSize);
3108:   VecCreate(PETSC_COMM_SELF, &coordinates);
3109:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3110:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3111:   VecSetBlockSize(coordinates, cdim);
3112:   VecSetType(coordinates, VECSTANDARD);
3113:   VecGetArray(coordinates, &coords);
3114:   if (!rank) {
3115:     double x[3];
3116:     int    val;

3118:     DMCreateLabel(*dm, "marker");
3119:     DMGetLabel(*dm, "marker", &marker);
3120:     for (v = 0; v < Nv; ++v) {
3121:       PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING);
3122:       snum = sscanf(line, "%lg %lg %lg %d", &x[0], &x[1], &x[2], &val);
3123:       if (snum != 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3124:       for (d = 0; d < cdim; ++d) coords[v*cdim+d] = x[d];
3125:       if (val) {DMLabelSetValue(marker, v+Nc, val);}
3126:     }
3127:   }
3128:   VecRestoreArray(coordinates, &coords);
3129:   DMSetCoordinatesLocal(*dm, coordinates);
3130:   VecDestroy(&coordinates);
3131:   PetscViewerDestroy(&viewer);
3132:   if (interpolate) {
3133:     DM      idm;
3134:     DMLabel bdlabel;

3136:     DMPlexInterpolate(*dm, &idm);
3137:     DMDestroy(dm);
3138:     *dm  = idm;

3140:     DMGetLabel(*dm, "marker", &bdlabel);
3141:     DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
3142:     DMPlexLabelComplete(*dm, bdlabel);
3143:   }
3144:   return(0);
3145: }

3147: /*@C
3148:   DMPlexCreateFromFile - This takes a filename and produces a DM

3150:   Input Parameters:
3151: + comm - The communicator
3152: . filename - A file name
3153: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

3155:   Output Parameter:
3156: . dm - The DM

3158:   Options Database Keys:
3159: . -dm_plex_create_from_hdf5_xdmf - use the PETSC_VIEWER_HDF5_XDMF format for reading HDF5 

3161:   Level: beginner

3163: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
3164: @*/
3165: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3166: {
3167:   const char    *extGmsh    = ".msh";
3168:   const char    *extCGNS    = ".cgns";
3169:   const char    *extExodus  = ".exo";
3170:   const char    *extGenesis = ".gen";
3171:   const char    *extFluent  = ".cas";
3172:   const char    *extHDF5    = ".h5";
3173:   const char    *extMed     = ".med";
3174:   const char    *extPLY     = ".ply";
3175:   const char    *extCV      = ".dat";
3176:   size_t         len;
3177:   PetscBool      isGmsh, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isCV;
3178:   PetscMPIInt    rank;

3184:   MPI_Comm_rank(comm, &rank);
3185:   PetscStrlen(filename, &len);
3186:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
3187:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,    4, &isGmsh);
3188:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,    5, &isCGNS);
3189:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus,  4, &isExodus);
3190:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);
3191:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent,  4, &isFluent);
3192:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,    3, &isHDF5);
3193:   PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,     4, &isMed);
3194:   PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY,     4, &isPLY);
3195:   PetscStrncmp(&filename[PetscMax(0,len-4)], extCV,      4, &isCV);
3196:   if (isGmsh) {
3197:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
3198:   } else if (isCGNS) {
3199:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
3200:   } else if (isExodus || isGenesis) {
3201:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
3202:   } else if (isFluent) {
3203:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
3204:   } else if (isHDF5) {
3205:     PetscBool      load_hdf5_xdmf = PETSC_FALSE;
3206:     PetscViewer viewer;

3208:     PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &load_hdf5_xdmf, NULL);
3209:     PetscViewerCreate(comm, &viewer);
3210:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
3211:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3212:     PetscViewerFileSetName(viewer, filename);
3213:     DMCreate(comm, dm);
3214:     DMSetType(*dm, DMPLEX);
3215:     if (load_hdf5_xdmf) {PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF);}
3216:     DMLoad(*dm, viewer);
3217:     if (load_hdf5_xdmf) {PetscViewerPopFormat(viewer);}
3218:     PetscViewerDestroy(&viewer);

3220:     if (interpolate) {
3221:       DM idm;

3223:       DMPlexInterpolate(*dm, &idm);
3224:       DMDestroy(dm);
3225:       *dm  = idm;
3226:     }
3227:   } else if (isMed) {
3228:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
3229:   } else if (isPLY) {
3230:     DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);
3231:   } else if (isCV) {
3232:     DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm);
3233:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
3234:   return(0);
3235: }

3237: /*@
3238:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

3240:   Collective on comm

3242:   Input Parameters:
3243: + comm    - The communicator
3244: . dim     - The spatial dimension
3245: - simplex - Flag for simplex, otherwise use a tensor-product cell

3247:   Output Parameter:
3248: . refdm - The reference cell

3250:   Level: intermediate

3252: .keywords: reference cell
3253: .seealso:
3254: @*/
3255: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
3256: {
3257:   DM             rdm;
3258:   Vec            coords;

3262:   DMCreate(comm, &rdm);
3263:   DMSetType(rdm, DMPLEX);
3264:   DMSetDimension(rdm, dim);
3265:   switch (dim) {
3266:   case 0:
3267:   {
3268:     PetscInt    numPoints[1]        = {1};
3269:     PetscInt    coneSize[1]         = {0};
3270:     PetscInt    cones[1]            = {0};
3271:     PetscInt    coneOrientations[1] = {0};
3272:     PetscScalar vertexCoords[1]     = {0.0};

3274:     DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3275:   }
3276:   break;
3277:   case 1:
3278:   {
3279:     PetscInt    numPoints[2]        = {2, 1};
3280:     PetscInt    coneSize[3]         = {2, 0, 0};
3281:     PetscInt    cones[2]            = {1, 2};
3282:     PetscInt    coneOrientations[2] = {0, 0};
3283:     PetscScalar vertexCoords[2]     = {-1.0,  1.0};

3285:     DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3286:   }
3287:   break;
3288:   case 2:
3289:     if (simplex) {
3290:       PetscInt    numPoints[2]        = {3, 1};
3291:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
3292:       PetscInt    cones[3]            = {1, 2, 3};
3293:       PetscInt    coneOrientations[3] = {0, 0, 0};
3294:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

3296:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3297:     } else {
3298:       PetscInt    numPoints[2]        = {4, 1};
3299:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3300:       PetscInt    cones[4]            = {1, 2, 3, 4};
3301:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3302:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

3304:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3305:     }
3306:   break;
3307:   case 3:
3308:     if (simplex) {
3309:       PetscInt    numPoints[2]        = {4, 1};
3310:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3311:       PetscInt    cones[4]            = {1, 3, 2, 4};
3312:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3313:       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};

3315:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3316:     } else {
3317:       PetscInt    numPoints[2]        = {8, 1};
3318:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3319:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
3320:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3321:       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,
3322:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3324:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3325:     }
3326:   break;
3327:   default:
3328:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
3329:   }
3330:   DMPlexInterpolate(rdm, refdm);
3331:   if (rdm->coordinateDM) {
3332:     DM           ncdm;
3333:     PetscSection cs;
3334:     PetscInt     pEnd = -1;

3336:     DMGetSection(rdm->coordinateDM, &cs);
3337:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
3338:     if (pEnd >= 0) {
3339:       DMClone(rdm->coordinateDM, &ncdm);
3340:       DMSetSection(ncdm, cs);
3341:       DMSetCoordinateDM(*refdm, ncdm);
3342:       DMDestroy(&ncdm);
3343:     }
3344:   }
3345:   DMGetCoordinatesLocal(rdm, &coords);
3346:   if (coords) {
3347:     DMSetCoordinatesLocal(*refdm, coords);
3348:   } else {
3349:     DMGetCoordinates(rdm, &coords);
3350:     if (coords) {DMSetCoordinates(*refdm, coords);}
3351:   }
3352:   DMDestroy(&rdm);
3353:   return(0);
3354: }