Actual source code: plexcreate.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: #define PETSCDM_DLL
  2:  #include <petsc/private/dmpleximpl.h>
  3:  #include <petsc/private/hashseti.h>
  4:  #include <petscsf.h>

  6: PetscLogEvent DMPLEX_CreateFromFile, DMPLEX_CreateFromCellList, DMPLEX_CreateFromCellList_Coordinates;

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

 11:   Collective

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

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

 24:   Level: beginner

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

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

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

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

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

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

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

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

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

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

150:   Collective

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

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

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

172:   Level: beginner

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

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

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

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

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

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

296:   Collective

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

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

307:   Level: beginner

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

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

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

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

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

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

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

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

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

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


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

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

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

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

516: static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
517: {
518:   DM             boundary;
519:   PetscInt       i;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

989:   Collective

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

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

1004:   Options Database Keys:
1005: + -dm_plex_box_lower <x,y,z> - Specify lower-left-bottom coordinates for the box
1006: . -dm_plex_box_upper <x,y,z> - Specify upper-right-top coordinates for the box
1007: - -dm_plex_box_faces <m,n,p> - Number of faces in each linear direction

1009:   Notes:
1010:   The options database keys above take lists of length d in d dimensions.

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

1027: and for simplicial cells

1029: $ 14----8---15----9----16
1030: $  |\     5  |\      7 |
1031: $  | \       | \       |
1032: $ 13   2    14    3    15
1033: $  | 4   \   | 6   \   |
1034: $  |       \ |       \ |
1035: $ 11----6---12----7----13
1036: $  |\        |\        |
1037: $  | \    1  | \     3 |
1038: $ 10   0    11    1    12
1039: $  | 0   \   | 2   \   |
1040: $  |       \ |       \ |
1041: $  8----4----9----5----10

1043:   Level: beginner

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

1058:   n    = 3;
1059:   PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", fac, &n, &flg);
1060:   for (i = 0; i < dim; ++i) fac[i] = faces ? faces[i] : (flg && i < n ? fac[i] : (dim == 1 ? 1 : 4-dim));
1061:   if (lower) for (i = 0; i < dim; ++i) low[i] = lower[i];
1062:   if (upper) for (i = 0; i < dim; ++i) upp[i] = upper[i];
1063:   if (periodicity) for (i = 0; i < dim; ++i) bdt[i] = periodicity[i];
1064:   /* Allow bounds to be specified from the command line */
1065:   n    = 3;
1066:   PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_lower", low, &n, &flg);
1067:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Lower box point had %D values, should have been %D", n, dim);
1068:   n    = 3;
1069:   PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_upper", upp, &n, &flg);
1070:   if (flg && (n != dim)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Upper box point had %D values, should have been %D", n, dim);

1072:   if (dim == 1)      {DMPlexCreateLineMesh_Internal(comm, fac[0], low[0], upp[0], bdt[0], dm);}
1073:   else if (simplex)  {DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1074:   else               {DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);}
1075:   return(0);
1076: }

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

1081:   Collective

1083:   Input Parameters:
1084: + comm        - The communicator for the DM object
1085: . faces       - Number of faces per dimension, or NULL for (1, 1, 1)
1086: . lower       - The lower left corner, or NULL for (0, 0, 0)
1087: . upper       - The upper right corner, or NULL for (1, 1, 1)
1088: . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE
1089: . ordExt      - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1090: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1092:   Output Parameter:
1093: . dm  - The DM object

1095:   Level: beginner

1097: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateWedgeCylinderMesh(), DMPlexExtrude(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1098: @*/
1099: PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool ordExt, PetscBool interpolate, DM *dm)
1100: {
1101:   DM             bdm, botdm;
1102:   PetscInt       i;
1103:   PetscInt       fac[3] = {0, 0, 0};
1104:   PetscReal      low[3] = {0, 0, 0};
1105:   PetscReal      upp[3] = {1, 1, 1};
1106:   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};

1110:   for (i = 0; i < 3; ++i) fac[i] = faces ? (faces[i] > 0 ? faces[i] : 1) : 1;
1111:   if (lower) for (i = 0; i < 3; ++i) low[i] = lower[i];
1112:   if (upper) for (i = 0; i < 3; ++i) upp[i] = upper[i];
1113:   if (periodicity) for (i = 0; i < 3; ++i) bdt[i] = periodicity[i];
1114:   for (i = 0; i < 3; ++i) if (bdt[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity not yet supported");

1116:   DMCreate(comm, &bdm);
1117:   DMSetType(bdm, DMPLEX);
1118:   DMSetDimension(bdm, 1);
1119:   DMSetCoordinateDim(bdm, 2);
1120:   DMPlexCreateSquareBoundary(bdm, low, upp, fac);
1121:   DMPlexGenerate(bdm, NULL, PETSC_FALSE, &botdm);
1122:   DMDestroy(&bdm);
1123:   DMPlexExtrude(botdm, fac[2], upp[2] - low[2], ordExt, interpolate, dm);
1124:   if (low[2] != 0.0) {
1125:     Vec         v;
1126:     PetscScalar *x;
1127:     PetscInt    cDim, n;

1129:     DMGetCoordinatesLocal(*dm, &v);
1130:     VecGetBlockSize(v, &cDim);
1131:     VecGetLocalSize(v, &n);
1132:     VecGetArray(v, &x);
1133:     x   += cDim;
1134:     for (i=0; i<n; i+=cDim) x[i] += low[2];
1135:     VecRestoreArray(v,&x);
1136:     DMSetCoordinatesLocal(*dm, v);
1137:   }
1138:   DMDestroy(&botdm);
1139:   return(0);
1140: }

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

1145:   Collective on idm

1147:   Input Parameters:
1148: + idm - The mesh to be extruted
1149: . layers - The number of layers
1150: . height - The height of the extruded layer
1151: . ordExt - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1152: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1154:   Output Parameter:
1155: . dm  - The DM object

1157:   Notes: The mesh created has prismatic cells, and the vertex ordering in the cone of the cell is that of the tensor prismatic cells.

1159:   Level: advanced

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

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

1182:   DMPlexGetHeightStratum(idm, 0, &cStart, &cEnd);
1183:   DMPlexGetDepthStratum(idm, 0, &vStart, &vEnd);
1184:   numCells = (cEnd - cStart)*layers;
1185:   numVertices = (vEnd - vStart)*(layers+1);
1186:   DMCreate(PetscObjectComm((PetscObject)idm), dm);
1187:   DMSetType(*dm, DMPLEX);
1188:   DMSetDimension(*dm, dim+1);
1189:   DMPlexSetChart(*dm, 0, numCells+numVertices);
1190:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1191:   DMCreateLabel(*dm, "celltype");
1192:   for (c = cStart, cellV = 0; c < cEnd; ++c) {
1193:     DMPolytopeType ct, nct;
1194:     PetscInt      *closure = NULL;
1195:     PetscInt       closureSize, numCorners = 0;

1197:     DMPlexGetCellType(idm, c, &ct);
1198:     switch (ct) {
1199:       case DM_POLYTOPE_SEGMENT:       nct = DM_POLYTOPE_SEG_PRISM_TENSOR;break;
1200:       case DM_POLYTOPE_TRIANGLE:      nct = DM_POLYTOPE_TRI_PRISM_TENSOR;break;
1201:       case DM_POLYTOPE_QUADRILATERAL: nct = DM_POLYTOPE_QUAD_PRISM_TENSOR;break;
1202:       default: nct = DM_POLYTOPE_UNKNOWN;
1203:     }
1204:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1205:     for (v = 0; v < closureSize*2; v += 2) if ((closure[v] >= vStart) && (closure[v] < vEnd)) numCorners++;
1206:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1207:     for (l = 0; l < layers; ++l) {
1208:       const PetscInt cell = ordExt ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart;

1210:       DMPlexSetConeSize(*dm, cell, 2*numCorners);
1211:       DMPlexSetCellType(*dm, cell, nct);
1212:     }
1213:     cellV = PetscMax(numCorners,cellV);
1214:   }
1215:   DMSetUp(*dm);

1217:   DMGetCoordinateDim(idm, &cDim);
1218:   if (dim != cDim) {
1219:     PetscCalloc1(cDim*(vEnd - vStart), &normals);
1220:   }
1221:   PetscMalloc1(3*cellV,&newCone);
1222:   for (c = cStart; c < cEnd; ++c) {
1223:     PetscInt *closure = NULL;
1224:     PetscInt closureSize, numCorners = 0, l;
1225:     PetscReal normal[3] = {0, 0, 0};

1227:     if (normals) {
1228:       DMPlexComputeCellGeometryFVM(idm, c, NULL, NULL, normal);
1229:     }
1230:     DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1231:     for (v = 0; v < closureSize*2; v += 2) {
1232:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1233:         PetscInt d;

1235:         newCone[numCorners++] = closure[v] - vStart;
1236:         if (normals) { for (d = 0; d < cDim; ++d) normals[cDim*(closure[v]-vStart)+d] += normal[d]; }
1237:       }
1238:     }
1239:     DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &closureSize, &closure);
1240:     for (l = 0; l < layers; ++l) {
1241:       PetscInt i;

1243:       for (i = 0; i < numCorners; ++i) {
1244:         newCone[  numCorners + i] = ordExt ? (layers+1)*newCone[i] + l     + numCells :     l*(vEnd - vStart) + newCone[i] + numCells;
1245:         newCone[2*numCorners + i] = ordExt ? (layers+1)*newCone[i] + l + 1 + numCells : (l+1)*(vEnd - vStart) + newCone[i] + numCells;
1246:       }
1247:       DMPlexSetCone(*dm, ordExt ? layers*(c - cStart) + l : l*(cEnd - cStart) + c - cStart, newCone + numCorners);
1248:     }
1249:   }
1250:   DMPlexSymmetrize(*dm);
1251:   DMPlexStratify(*dm);
1252:   PetscFree(newCone);

1254:   cDimB = cDim == dim ? cDim+1 : cDim;
1255:   DMGetCoordinateSection(*dm, &coordSectionB);
1256:   PetscSectionSetNumFields(coordSectionB, 1);
1257:   PetscSectionSetFieldComponents(coordSectionB, 0, cDimB);
1258:   PetscSectionSetChart(coordSectionB, numCells, numCells+numVertices);
1259:   for (v = numCells; v < numCells+numVertices; ++v) {
1260:     PetscSectionSetDof(coordSectionB, v, cDimB);
1261:     PetscSectionSetFieldDof(coordSectionB, v, 0, cDimB);
1262:     DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);
1263:   }
1264:   PetscSectionSetUp(coordSectionB);
1265:   PetscSectionGetStorageSize(coordSectionB, &coordSize);
1266:   VecCreate(PETSC_COMM_SELF, &coordinatesB);
1267:   PetscObjectSetName((PetscObject) coordinatesB, "coordinates");
1268:   VecSetSizes(coordinatesB, coordSize, PETSC_DETERMINE);
1269:   VecSetBlockSize(coordinatesB, cDimB);
1270:   VecSetType(coordinatesB,VECSTANDARD);

1272:   DMGetCoordinateSection(idm, &coordSectionA);
1273:   DMGetCoordinatesLocal(idm, &coordinatesA);
1274:   VecGetArray(coordinatesB, &coordsB);
1275:   VecGetArrayRead(coordinatesA, &coordsA);
1276:   for (v = vStart; v < vEnd; ++v) {
1277:     const PetscScalar *cptr;
1278:     PetscReal         ones2[2] = { 0., 1.}, ones3[3] = { 0., 0., 1.};
1279:     PetscReal         *normal, norm, h = height/layers;
1280:     PetscInt          offA, d, cDimA = cDim;

1282:     normal = normals ? normals + cDimB*(v - vStart) : (cDim > 1 ? ones3 : ones2);
1283:     if (normals) {
1284:       for (d = 0, norm = 0.0; d < cDimB; ++d) norm += normal[d]*normal[d];
1285:       for (d = 0; d < cDimB; ++d) normal[d] *= 1./PetscSqrtReal(norm);
1286:     }

1288:     PetscSectionGetOffset(coordSectionA, v, &offA);
1289:     cptr = coordsA + offA;
1290:     for (l = 0; l < layers+1; ++l) {
1291:       PetscInt offB, d, newV;

1293:       newV = ordExt ? (layers+1)*(v -vStart) + l + numCells : (vEnd -vStart)*l + (v -vStart) + numCells;
1294:       PetscSectionGetOffset(coordSectionB, newV, &offB);
1295:       for (d = 0; d < cDimA; ++d) { coordsB[offB+d]  = cptr[d]; }
1296:       for (d = 0; d < cDimB; ++d) { coordsB[offB+d] += l ? normal[d]*h : 0.0; }
1297:       cptr    = coordsB + offB;
1298:       cDimA   = cDimB;
1299:     }
1300:   }
1301:   VecRestoreArrayRead(coordinatesA, &coordsA);
1302:   VecRestoreArray(coordinatesB, &coordsB);
1303:   DMSetCoordinatesLocal(*dm, coordinatesB);
1304:   VecDestroy(&coordinatesB);
1305:   PetscFree(normals);
1306:   if (interpolate) {
1307:     DM idm;

1309:     DMPlexInterpolate(*dm, &idm);
1310:     DMPlexCopyCoordinates(*dm, idm);
1311:     DMDestroy(dm);
1312:     *dm  = idm;
1313:   }
1314:   return(0);
1315: }

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

1320:   Logically Collective on dm

1322:   Input Parameters:
1323: + dm - the DM context
1324: - prefix - the prefix to prepend to all option names

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

1330:   Level: advanced

1332: .seealso: SNESSetFromOptions()
1333: @*/
1334: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1335: {
1336:   DM_Plex       *mesh = (DM_Plex *) dm->data;

1341:   PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);
1342:   PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);
1343:   return(0);
1344: }

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

1349:   Collective

1351:   Input Parameters:
1352: + comm      - The communicator for the DM object
1353: . numRefine - The number of regular refinements to the basic 5 cell structure
1354: - periodicZ - The boundary type for the Z direction

1356:   Output Parameter:
1357: . dm  - The DM object

1359:   Note: Here is the output numbering looking from the bottom of the cylinder:
1360: $       17-----14
1361: $        |     |
1362: $        |  2  |
1363: $        |     |
1364: $ 17-----8-----7-----14
1365: $  |     |     |     |
1366: $  |  3  |  0  |  1  |
1367: $  |     |     |     |
1368: $ 19-----5-----6-----13
1369: $        |     |
1370: $        |  4  |
1371: $        |     |
1372: $       19-----13
1373: $
1374: $ and up through the top
1375: $
1376: $       18-----16
1377: $        |     |
1378: $        |  2  |
1379: $        |     |
1380: $ 18----10----11-----16
1381: $  |     |     |     |
1382: $  |  3  |  0  |  1  |
1383: $  |     |     |     |
1384: $ 20-----9----12-----15
1385: $        |     |
1386: $        |  4  |
1387: $        |     |
1388: $       20-----15

1390:   Level: beginner

1392: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1393: @*/
1394: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
1395: {
1396:   const PetscInt dim = 3;
1397:   PetscInt       numCells, numVertices, r;
1398:   PetscMPIInt    rank;

1403:   MPI_Comm_rank(comm, &rank);
1404:   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
1405:   DMCreate(comm, dm);
1406:   DMSetType(*dm, DMPLEX);
1407:   DMSetDimension(*dm, dim);
1408:   /* Create topology */
1409:   {
1410:     PetscInt cone[8], c;

1412:     numCells    = !rank ?  5 : 0;
1413:     numVertices = !rank ? 16 : 0;
1414:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1415:       numCells   *= 3;
1416:       numVertices = !rank ? 24 : 0;
1417:     }
1418:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1419:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 8);}
1420:     DMSetUp(*dm);
1421:     if (!rank) {
1422:       if (periodicZ == DM_BOUNDARY_PERIODIC) {
1423:         cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1424:         cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1425:         DMPlexSetCone(*dm, 0, cone);
1426:         cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1427:         cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1428:         DMPlexSetCone(*dm, 1, cone);
1429:         cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1430:         cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1431:         DMPlexSetCone(*dm, 2, cone);
1432:         cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1433:         cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1434:         DMPlexSetCone(*dm, 3, cone);
1435:         cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1436:         cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1437:         DMPlexSetCone(*dm, 4, cone);

1439:         cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1440:         cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1441:         DMPlexSetCone(*dm, 5, cone);
1442:         cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1443:         cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1444:         DMPlexSetCone(*dm, 6, cone);
1445:         cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1446:         cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1447:         DMPlexSetCone(*dm, 7, cone);
1448:         cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1449:         cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1450:         DMPlexSetCone(*dm, 8, cone);
1451:         cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1452:         cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1453:         DMPlexSetCone(*dm, 9, cone);

1455:         cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1456:         cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1457:         DMPlexSetCone(*dm, 10, cone);
1458:         cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1459:         cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1460:         DMPlexSetCone(*dm, 11, cone);
1461:         cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1462:         cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1463:         DMPlexSetCone(*dm, 12, cone);
1464:         cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1465:         cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1466:         DMPlexSetCone(*dm, 13, cone);
1467:         cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1468:         cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1469:         DMPlexSetCone(*dm, 14, cone);
1470:       } else {
1471:         cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
1472:         cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
1473:         DMPlexSetCone(*dm, 0, cone);
1474:         cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
1475:         cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
1476:         DMPlexSetCone(*dm, 1, cone);
1477:         cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
1478:         cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
1479:         DMPlexSetCone(*dm, 2, cone);
1480:         cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
1481:         cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
1482:         DMPlexSetCone(*dm, 3, cone);
1483:         cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
1484:         cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
1485:         DMPlexSetCone(*dm, 4, cone);
1486:       }
1487:     }
1488:     DMPlexSymmetrize(*dm);
1489:     DMPlexStratify(*dm);
1490:   }
1491:   /* Interpolate */
1492:   {
1493:     DM idm;

1495:     DMPlexInterpolate(*dm, &idm);
1496:     DMDestroy(dm);
1497:     *dm  = idm;
1498:   }
1499:   /* Create cube geometry */
1500:   {
1501:     Vec             coordinates;
1502:     PetscSection    coordSection;
1503:     PetscScalar    *coords;
1504:     PetscInt        coordSize, v;
1505:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1506:     const PetscReal ds2 = dis/2.0;

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

1566:     bdType[0] = DM_BOUNDARY_NONE;
1567:     bdType[1] = DM_BOUNDARY_NONE;
1568:     bdType[2] = periodicZ;
1569:     for (i = 0; i < dim; i++) {
1570:       L[i]       = upper[i] - lower[i];
1571:       maxCell[i] = 1.1 * (L[i] / numZCells);
1572:     }
1573:     DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);
1574:   }
1575:   /* Refine topology */
1576:   for (r = 0; r < numRefine; ++r) {
1577:     DM rdm = NULL;

1579:     DMRefine(*dm, comm, &rdm);
1580:     DMDestroy(dm);
1581:     *dm  = rdm;
1582:   }
1583:   /* Remap geometry to cylinder
1584:        Interior square: Linear interpolation is correct
1585:        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1586:        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance

1588:          phi     = arctan(y/x)
1589:          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1590:          d_far   = sqrt(1/2 + sin^2(phi))

1592:        so we remap them using

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

1597:        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1598:   */
1599:   {
1600:     Vec           coordinates;
1601:     PetscSection  coordSection;
1602:     PetscScalar  *coords;
1603:     PetscInt      vStart, vEnd, v;
1604:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1605:     const PetscReal ds2 = 0.5*dis;

1607:     DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);
1608:     DMGetCoordinateSection(*dm, &coordSection);
1609:     DMGetCoordinatesLocal(*dm, &coordinates);
1610:     VecGetArray(coordinates, &coords);
1611:     for (v = vStart; v < vEnd; ++v) {
1612:       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1613:       PetscInt  off;

1615:       PetscSectionGetOffset(coordSection, v, &off);
1616:       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
1617:       x    = PetscRealPart(coords[off]);
1618:       y    = PetscRealPart(coords[off+1]);
1619:       phi  = PetscAtan2Real(y, x);
1620:       sinp = PetscSinReal(phi);
1621:       cosp = PetscCosReal(phi);
1622:       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
1623:         dc = PetscAbsReal(ds2/sinp);
1624:         df = PetscAbsReal(dis/sinp);
1625:         xc = ds2*x/PetscAbsReal(y);
1626:         yc = ds2*PetscSignReal(y);
1627:       } else {
1628:         dc = PetscAbsReal(ds2/cosp);
1629:         df = PetscAbsReal(dis/cosp);
1630:         xc = ds2*PetscSignReal(x);
1631:         yc = ds2*y/PetscAbsReal(x);
1632:       }
1633:       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1634:       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1635:     }
1636:     VecRestoreArray(coordinates, &coords);
1637:     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1638:       DMLocalizeCoordinates(*dm);
1639:     }
1640:   }
1641:   return(0);
1642: }

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

1647:   Collective

1649:   Input Parameters:
1650: + comm - The communicator for the DM object
1651: . n    - The number of wedges around the origin
1652: - interpolate - Create edges and faces

1654:   Output Parameter:
1655: . dm  - The DM object

1657:   Level: beginner

1659: .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1660: @*/
1661: PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
1662: {
1663:   const PetscInt dim = 3;
1664:   PetscInt       numCells, numVertices, v;
1665:   PetscMPIInt    rank;

1670:   MPI_Comm_rank(comm, &rank);
1671:   if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n);
1672:   DMCreate(comm, dm);
1673:   DMSetType(*dm, DMPLEX);
1674:   DMSetDimension(*dm, dim);
1675:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1676:   DMCreateLabel(*dm, "celltype");
1677:   /* Create topology */
1678:   {
1679:     PetscInt cone[6], c;

1681:     numCells    = !rank ?        n : 0;
1682:     numVertices = !rank ?  2*(n+1) : 0;
1683:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1684:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 6);}
1685:     DMSetUp(*dm);
1686:     for (c = 0; c < numCells; c++) {
1687:       cone[0] =  c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n;
1688:       cone[3] =  c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n;
1689:       DMPlexSetCone(*dm, c, cone);
1690:       DMPlexSetCellType(*dm, c, DM_POLYTOPE_TRI_PRISM_TENSOR);
1691:     }
1692:     DMPlexSymmetrize(*dm);
1693:     DMPlexStratify(*dm);
1694:   }
1695:   for (v = numCells; v < numCells+numVertices; ++v) {
1696:     DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);
1697:   }
1698:   /* Interpolate */
1699:   if (interpolate) {
1700:     DM idm;

1702:     DMPlexInterpolate(*dm, &idm);
1703:     DMDestroy(dm);
1704:     *dm  = idm;
1705:   }
1706:   /* Create cylinder geometry */
1707:   {
1708:     Vec          coordinates;
1709:     PetscSection coordSection;
1710:     PetscScalar *coords;
1711:     PetscInt     coordSize, c;

1713:     /* Build coordinates */
1714:     DMGetCoordinateSection(*dm, &coordSection);
1715:     PetscSectionSetNumFields(coordSection, 1);
1716:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1717:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1718:     for (v = numCells; v < numCells+numVertices; ++v) {
1719:       PetscSectionSetDof(coordSection, v, dim);
1720:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1721:     }
1722:     PetscSectionSetUp(coordSection);
1723:     PetscSectionGetStorageSize(coordSection, &coordSize);
1724:     VecCreate(PETSC_COMM_SELF, &coordinates);
1725:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1726:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1727:     VecSetBlockSize(coordinates, dim);
1728:     VecSetType(coordinates,VECSTANDARD);
1729:     VecGetArray(coordinates, &coords);
1730:     for (c = 0; c < numCells; c++) {
1731:       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;
1732:       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;
1733:     }
1734:     if (!rank) {
1735:       coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0;
1736:       coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0;
1737:     }
1738:     VecRestoreArray(coordinates, &coords);
1739:     DMSetCoordinatesLocal(*dm, coordinates);
1740:     VecDestroy(&coordinates);
1741:   }
1742:   return(0);
1743: }

1745: PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1746: {
1747:   PetscReal prod = 0.0;
1748:   PetscInt  i;
1749:   for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
1750:   return PetscSqrtReal(prod);
1751: }
1752: PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
1753: {
1754:   PetscReal prod = 0.0;
1755:   PetscInt  i;
1756:   for (i = 0; i < dim; ++i) prod += x[i]*y[i];
1757:   return prod;
1758: }

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

1763:   Collective

1765:   Input Parameters:
1766: + comm  - The communicator for the DM object
1767: . dim   - The dimension
1768: - simplex - Use simplices, or tensor product cells

1770:   Output Parameter:
1771: . dm  - The DM object

1773:   Level: beginner

1775: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1776: @*/
1777: PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *dm)
1778: {
1779:   const PetscInt  embedDim = dim+1;
1780:   PetscSection    coordSection;
1781:   Vec             coordinates;
1782:   PetscScalar    *coords;
1783:   PetscReal      *coordsIn;
1784:   PetscInt        numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e;
1785:   PetscMPIInt     rank;
1786:   PetscErrorCode  ierr;

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

1806:       numCells    = !rank ? 20 : 0;
1807:       numVerts    = !rank ? 12 : 0;
1808:       firstVertex = numCells;
1809:       /* Use icosahedron, which for a unit sphere has coordinates which are all cyclic permutations of

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

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

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

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

2000:       numCells    = !rank ? 600 : 0;
2001:       numVerts    = !rank ? 120 : 0;
2002:       firstVertex = numCells;
2003:       /* Use the 600-cell, which for a unit sphere has coordinates which are

2005:            1/2 (\pm 1, \pm 1,    \pm 1, \pm 1)                          16
2006:                (\pm 1,    0,       0,      0)  all cyclic permutations   8
2007:            1/2 (\pm 1, \pm phi, \pm 1/phi, 0)  all even permutations    96

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

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

2079:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0, -1}, { 0,  0,  1, 0}},
2080:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2081:                                                            {{0,  0,  0,  1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, {-1,  0,  0, 0}},
2082:                                                            {{0,  0, -1,  0}, { 0, 0,  0,  0}, { 1,  0, 0,  0}, { 0,  0,  0, 0}}},

2084:                                                           {{{0,  0,  0,  0}, { 0, 0,  0,  1}, { 0,  0, 0,  0}, { 0, -1,  0, 0}},
2085:                                                            {{0,  0,  0, -1}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 1,  0,  0, 0}},
2086:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2087:                                                            {{0,  1,  0,  0}, {-1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}},

2089:                                                           {{{0,  0,  0,  0}, { 0, 0, -1,  0}, { 0,  1, 0,  0}, { 0,  0,  0, 0}},
2090:                                                            {{0,  0,  1,  0}, { 0, 0,  0,  0}, {-1,  0, 0,  0}, { 0,  0,  0, 0}},
2091:                                                            {{0, -1,  0,  0}, { 1, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}},
2092:                                                            {{0,  0,  0,  0}, { 0, 0,  0,  0}, { 0,  0, 0,  0}, { 0,  0,  0, 0}}}};
2093:                     PetscReal normal[4];
2094:                     PetscInt  e, f, g;

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

2151: /* External function declarations here */
2152: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
2153: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2154: extern PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
2155: extern PetscErrorCode DMCreateLocalSection_Plex(DM dm);
2156: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
2157: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
2158: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
2159: extern PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
2160: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
2161: extern PetscErrorCode DMSetUp_Plex(DM dm);
2162: extern PetscErrorCode DMDestroy_Plex(DM dm);
2163: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
2164: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
2165: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
2166: extern PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);
2167: static PetscErrorCode DMInitialize_Plex(DM dm);

2169: /* Replace dm with the contents of dmNew
2170:    - Share the DM_Plex structure
2171:    - Share the coordinates
2172:    - Share the SF
2173: */
2174: static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
2175: {
2176:   PetscSF               sf;
2177:   DM                    coordDM, coarseDM;
2178:   Vec                   coords;
2179:   PetscBool             isper;
2180:   const PetscReal      *maxCell, *L;
2181:   const DMBoundaryType *bd;
2182:   PetscErrorCode        ierr;

2185:   DMGetPointSF(dmNew, &sf);
2186:   DMSetPointSF(dm, sf);
2187:   DMGetCoordinateDM(dmNew, &coordDM);
2188:   DMGetCoordinatesLocal(dmNew, &coords);
2189:   DMSetCoordinateDM(dm, coordDM);
2190:   DMSetCoordinatesLocal(dm, coords);
2191:   DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);
2192:   DMSetPeriodicity(dmNew, isper, maxCell, L, bd);
2193:   DMDestroy_Plex(dm);
2194:   DMInitialize_Plex(dm);
2195:   dm->data = dmNew->data;
2196:   ((DM_Plex *) dmNew->data)->refct++;
2197:   DMDestroyLabelLinkList_Internal(dm);
2198:   DMCopyLabels(dmNew, dm, PETSC_OWN_POINTER, PETSC_TRUE);
2199:   DMGetCoarseDM(dmNew,&coarseDM);
2200:   DMSetCoarseDM(dm,coarseDM);
2201:   return(0);
2202: }

2204: /* Swap dm with the contents of dmNew
2205:    - Swap the DM_Plex structure
2206:    - Swap the coordinates
2207:    - Swap the point PetscSF
2208: */
2209: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
2210: {
2211:   DM              coordDMA, coordDMB;
2212:   Vec             coordsA,  coordsB;
2213:   PetscSF         sfA,      sfB;
2214:   void            *tmp;
2215:   DMLabelLink     listTmp;
2216:   DMLabel         depthTmp;
2217:   PetscInt        tmpI;
2218:   PetscErrorCode  ierr;

2221:   DMGetPointSF(dmA, &sfA);
2222:   DMGetPointSF(dmB, &sfB);
2223:   PetscObjectReference((PetscObject) sfA);
2224:   DMSetPointSF(dmA, sfB);
2225:   DMSetPointSF(dmB, sfA);
2226:   PetscObjectDereference((PetscObject) sfA);

2228:   DMGetCoordinateDM(dmA, &coordDMA);
2229:   DMGetCoordinateDM(dmB, &coordDMB);
2230:   PetscObjectReference((PetscObject) coordDMA);
2231:   DMSetCoordinateDM(dmA, coordDMB);
2232:   DMSetCoordinateDM(dmB, coordDMA);
2233:   PetscObjectDereference((PetscObject) coordDMA);

2235:   DMGetCoordinatesLocal(dmA, &coordsA);
2236:   DMGetCoordinatesLocal(dmB, &coordsB);
2237:   PetscObjectReference((PetscObject) coordsA);
2238:   DMSetCoordinatesLocal(dmA, coordsB);
2239:   DMSetCoordinatesLocal(dmB, coordsA);
2240:   PetscObjectDereference((PetscObject) coordsA);

2242:   tmp       = dmA->data;
2243:   dmA->data = dmB->data;
2244:   dmB->data = tmp;
2245:   listTmp   = dmA->labels;
2246:   dmA->labels = dmB->labels;
2247:   dmB->labels = listTmp;
2248:   depthTmp  = dmA->depthLabel;
2249:   dmA->depthLabel = dmB->depthLabel;
2250:   dmB->depthLabel = depthTmp;
2251:   depthTmp  = dmA->celltypeLabel;
2252:   dmA->celltypeLabel = dmB->celltypeLabel;
2253:   dmB->celltypeLabel = depthTmp;
2254:   tmpI         = dmA->levelup;
2255:   dmA->levelup = dmB->levelup;
2256:   dmB->levelup = tmpI;
2257:   return(0);
2258: }

2260: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2261: {
2262:   DM_Plex       *mesh = (DM_Plex*) dm->data;
2263:   PetscBool      flg;

2267:   /* Handle viewing */
2268:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL);
2269:   PetscOptionsBoundedInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL,0);
2270:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL);
2271:   PetscOptionsBoundedInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL,0);
2272:   DMMonitorSetFromOptions(dm, "-dm_plex_monitor_throughput", "Monitor the simulation throughput", "DMPlexMonitorThroughput", DMPlexMonitorThroughput, NULL, &flg);
2273:   if (flg) {PetscLogDefaultBegin();}
2274:   /* Point Location */
2275:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL);
2276:   /* Partitioning and distribution */
2277:   PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL);
2278:   /* Generation and remeshing */
2279:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL);
2280:   /* Projection behavior */
2281:   PetscOptionsBoundedInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL,0);
2282:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);
2283:   PetscOptionsEnum("-dm_plex_cell_refiner", "Strategy for cell refinment", "ex40.c", DMPlexCellRefinerTypes, (PetscEnum) mesh->cellRefiner, (PetscEnum *) &mesh->cellRefiner, NULL);
2284:   /* Checking structure */
2285:   {
2286:     PetscBool   flg = PETSC_FALSE, flg2 = PETSC_FALSE, all = PETSC_FALSE;

2288:     PetscOptionsBool("-dm_plex_check_all", "Perform all checks", NULL, PETSC_FALSE, &all, &flg2);
2289:     PetscOptionsBool("-dm_plex_check_symmetry", "Check that the adjacency information in the mesh is symmetric", "DMPlexCheckSymmetry", PETSC_FALSE, &flg, &flg2);
2290:     if (all || (flg && flg2)) {DMPlexCheckSymmetry(dm);}
2291:     PetscOptionsBool("-dm_plex_check_skeleton", "Check that each cell has the correct number of vertices (only for homogeneous simplex or tensor meshes)", "DMPlexCheckSkeleton", PETSC_FALSE, &flg, &flg2);
2292:     if (all || (flg && flg2)) {DMPlexCheckSkeleton(dm, 0);}
2293:     PetscOptionsBool("-dm_plex_check_faces", "Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type", "DMPlexCheckFaces", PETSC_FALSE, &flg, &flg2);
2294:     if (all || (flg && flg2)) {DMPlexCheckFaces(dm, 0);}
2295:     PetscOptionsBool("-dm_plex_check_geometry", "Check that cells have positive volume", "DMPlexCheckGeometry", PETSC_FALSE, &flg, &flg2);
2296:     if (all || (flg && flg2)) {DMPlexCheckGeometry(dm);}
2297:     PetscOptionsBool("-dm_plex_check_pointsf", "Check some necessary conditions for PointSF", "DMPlexCheckPointSF", PETSC_FALSE, &flg, &flg2);
2298:     if (all || (flg && flg2)) {DMPlexCheckPointSF(dm);}
2299:     PetscOptionsBool("-dm_plex_check_interface_cones", "Check points on inter-partition interfaces have conforming order of cone points", "DMPlexCheckInterfaceCones", PETSC_FALSE, &flg, &flg2);
2300:     if (all || (flg && flg2)) {DMPlexCheckInterfaceCones(dm);}
2301:     PetscOptionsBool("-dm_plex_check_cell_shape", "Check cell shape", "DMPlexCheckCellShape", PETSC_FALSE, &flg, &flg2);
2302:     if (flg && flg2) {DMPlexCheckCellShape(dm, PETSC_TRUE, PETSC_DETERMINE);}
2303:   }

2305:   PetscPartitionerSetFromOptions(mesh->partitioner);
2306:   return(0);
2307: }

2309: static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
2310: {
2311:   PetscInt       refine = 0, coarsen = 0, r;
2312:   PetscBool      isHierarchy;

2317:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
2318:   /* Handle DMPlex refinement */
2319:   PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL,0);
2320:   PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy,0);
2321:   if (refine) {DMPlexSetRefinementUniform(dm, PETSC_TRUE);}
2322:   if (refine && isHierarchy) {
2323:     DM *dms, coarseDM;

2325:     DMGetCoarseDM(dm, &coarseDM);
2326:     PetscObjectReference((PetscObject)coarseDM);
2327:     PetscMalloc1(refine,&dms);
2328:     DMRefineHierarchy(dm, refine, dms);
2329:     /* Total hack since we do not pass in a pointer */
2330:     DMPlexSwap_Static(dm, dms[refine-1]);
2331:     if (refine == 1) {
2332:       DMSetCoarseDM(dm, dms[0]);
2333:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2334:     } else {
2335:       DMSetCoarseDM(dm, dms[refine-2]);
2336:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
2337:       DMSetCoarseDM(dms[0], dms[refine-1]);
2338:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
2339:     }
2340:     DMSetCoarseDM(dms[refine-1], coarseDM);
2341:     PetscObjectDereference((PetscObject)coarseDM);
2342:     /* Free DMs */
2343:     for (r = 0; r < refine; ++r) {
2344:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2345:       DMDestroy(&dms[r]);
2346:     }
2347:     PetscFree(dms);
2348:   } else {
2349:     for (r = 0; r < refine; ++r) {
2350:       DM refinedMesh;

2352:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2353:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);
2354:       /* Total hack since we do not pass in a pointer */
2355:       DMPlexReplace_Static(dm, refinedMesh);
2356:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2357:       DMDestroy(&refinedMesh);
2358:     }
2359:   }
2360:   /* Handle DMPlex coarsening */
2361:   PetscOptionsBoundedInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL,0);
2362:   PetscOptionsBoundedInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy,0);
2363:   if (coarsen && isHierarchy) {
2364:     DM *dms;

2366:     PetscMalloc1(coarsen, &dms);
2367:     DMCoarsenHierarchy(dm, coarsen, dms);
2368:     /* Free DMs */
2369:     for (r = 0; r < coarsen; ++r) {
2370:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
2371:       DMDestroy(&dms[r]);
2372:     }
2373:     PetscFree(dms);
2374:   } else {
2375:     for (r = 0; r < coarsen; ++r) {
2376:       DM coarseMesh;

2378:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2379:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);
2380:       /* Total hack since we do not pass in a pointer */
2381:       DMPlexReplace_Static(dm, coarseMesh);
2382:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2383:       DMDestroy(&coarseMesh);
2384:     }
2385:   }
2386:   /* Handle */
2387:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
2388:   PetscOptionsTail();
2389:   return(0);
2390: }

2392: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
2393: {

2397:   DMCreateGlobalVector_Section_Private(dm,vec);
2398:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
2399:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
2400:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
2401:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
2402:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
2403:   return(0);
2404: }

2406: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
2407: {

2411:   DMCreateLocalVector_Section_Private(dm,vec);
2412:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
2413:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
2414:   return(0);
2415: }

2417: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
2418: {
2419:   PetscInt       depth, d;

2423:   DMPlexGetDepth(dm, &depth);
2424:   if (depth == 1) {
2425:     DMGetDimension(dm, &d);
2426:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
2427:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
2428:     else               {*pStart = 0; *pEnd = 0;}
2429:   } else {
2430:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
2431:   }
2432:   return(0);
2433: }

2435: static PetscErrorCode DMGetNeighbors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
2436: {
2437:   PetscSF           sf;
2438:   PetscInt          niranks, njranks, n;
2439:   const PetscMPIInt *iranks, *jranks;
2440:   DM_Plex           *data = (DM_Plex*) dm->data;
2441:   PetscErrorCode    ierr;

2444:   DMGetPointSF(dm, &sf);
2445:   if (!data->neighbors) {
2446:     PetscSFGetRootRanks(sf, &njranks, &jranks, NULL, NULL, NULL);
2447:     PetscSFGetLeafRanks(sf, &niranks, &iranks, NULL, NULL);
2448:     PetscMalloc1(njranks + niranks + 1, &data->neighbors);
2449:     PetscArraycpy(data->neighbors + 1, jranks, njranks);
2450:     PetscArraycpy(data->neighbors + njranks + 1, iranks, niranks);
2451:     n = njranks + niranks;
2452:     PetscSortRemoveDupsMPIInt(&n, data->neighbors + 1);
2453:     /* The following cast should never fail: can't have more neighbors than PETSC_MPI_INT_MAX */
2454:     PetscMPIIntCast(n, data->neighbors);
2455:   }
2456:   if (nranks) *nranks = data->neighbors[0];
2457:   if (ranks) {
2458:     if (data->neighbors[0]) *ranks = data->neighbors + 1;
2459:     else                    *ranks = NULL;
2460:   }
2461:   return(0);
2462: }

2464: static PetscErrorCode DMInitialize_Plex(DM dm)
2465: {

2469:   dm->ops->view                            = DMView_Plex;
2470:   dm->ops->load                            = DMLoad_Plex;
2471:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
2472:   dm->ops->clone                           = DMClone_Plex;
2473:   dm->ops->setup                           = DMSetUp_Plex;
2474:   dm->ops->createlocalsection              = DMCreateLocalSection_Plex;
2475:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
2476:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
2477:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
2478:   dm->ops->getlocaltoglobalmapping         = NULL;
2479:   dm->ops->createfieldis                   = NULL;
2480:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
2481:   dm->ops->createcoordinatefield           = DMCreateCoordinateField_Plex;
2482:   dm->ops->getcoloring                     = NULL;
2483:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
2484:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
2485:   dm->ops->createmassmatrix                = DMCreateMassMatrix_Plex;
2486:   dm->ops->createinjection                 = DMCreateInjection_Plex;
2487:   dm->ops->refine                          = DMRefine_Plex;
2488:   dm->ops->coarsen                         = DMCoarsen_Plex;
2489:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
2490:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
2491:   dm->ops->adaptlabel                      = DMAdaptLabel_Plex;
2492:   dm->ops->adaptmetric                     = DMAdaptMetric_Plex;
2493:   dm->ops->globaltolocalbegin              = NULL;
2494:   dm->ops->globaltolocalend                = NULL;
2495:   dm->ops->localtoglobalbegin              = NULL;
2496:   dm->ops->localtoglobalend                = NULL;
2497:   dm->ops->destroy                         = DMDestroy_Plex;
2498:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
2499:   dm->ops->createsuperdm                   = DMCreateSuperDM_Plex;
2500:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
2501:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
2502:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
2503:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
2504:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
2505:   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
2506:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
2507:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
2508:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
2509:   dm->ops->getneighbors                    = DMGetNeighbors_Plex;
2510:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);
2511:   PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);
2512:   PetscObjectComposeFunction((PetscObject)dm,"DMCreateNeumannOverlap_C",DMCreateNeumannOverlap_Plex);
2513:   PetscObjectComposeFunction((PetscObject)dm,"DMPlexGetOverlap_C",DMPlexGetOverlap_Plex);
2514:   return(0);
2515: }

2517: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
2518: {
2519:   DM_Plex        *mesh = (DM_Plex *) dm->data;

2523:   mesh->refct++;
2524:   (*newdm)->data = mesh;
2525:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
2526:   DMInitialize_Plex(*newdm);
2527:   return(0);
2528: }

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

2536:   Options Database Keys:
2537: + -dm_plex_hash_location             - Use grid hashing for point location
2538: . -dm_plex_partition_balance         - Attempt to evenly divide points on partition boundary between processes
2539: . -dm_plex_remesh_bd                 - Allow changes to the boundary on remeshing
2540: . -dm_plex_max_projection_height     - Maxmimum mesh point height used to project locally
2541: . -dm_plex_regular_refinement        - Use special nested projection algorithm for regular refinement
2542: . -dm_plex_check_symmetry            - Check that the adjacency information in the mesh is symmetric
2543: . -dm_plex_check_skeleton <celltype> - Check that each cell has the correct number of vertices
2544: . -dm_plex_check_faces <celltype>    - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type
2545: . -dm_plex_check_geometry            - Check that cells have positive volume
2546: . -dm_view :mesh.tex:ascii_latex     - View the mesh in LaTeX/TikZ
2547: . -dm_plex_view_scale <num>          - Scale the TikZ
2548: - -dm_plex_print_fem <num>           - View FEM assembly information, such as element vectors and matrices


2551:   Level: intermediate

2553: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
2554: M*/

2556: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
2557: {
2558:   DM_Plex       *mesh;
2559:   PetscInt       unit;

2564:   PetscNewLog(dm,&mesh);
2565:   dm->dim  = 0;
2566:   dm->data = mesh;

2568:   mesh->refct             = 1;
2569:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
2570:   mesh->maxConeSize       = 0;
2571:   mesh->cones             = NULL;
2572:   mesh->coneOrientations  = NULL;
2573:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
2574:   mesh->maxSupportSize    = 0;
2575:   mesh->supports          = NULL;
2576:   mesh->refinementUniform = PETSC_TRUE;
2577:   mesh->refinementLimit   = -1.0;
2578:   mesh->interpolated      = DMPLEX_INTERPOLATED_INVALID;
2579:   mesh->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID;

2581:   mesh->facesTmp = NULL;

2583:   mesh->tetgenOpts   = NULL;
2584:   mesh->triangleOpts = NULL;
2585:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
2586:   mesh->remeshBd     = PETSC_FALSE;

2588:   mesh->subpointMap = NULL;

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

2592:   mesh->regularRefinement   = PETSC_FALSE;
2593:   mesh->depthState          = -1;
2594:   mesh->celltypeState       = -1;
2595:   mesh->globalVertexNumbers = NULL;
2596:   mesh->globalCellNumbers   = NULL;
2597:   mesh->anchorSection       = NULL;
2598:   mesh->anchorIS            = NULL;
2599:   mesh->createanchors       = NULL;
2600:   mesh->computeanchormatrix = NULL;
2601:   mesh->parentSection       = NULL;
2602:   mesh->parents             = NULL;
2603:   mesh->childIDs            = NULL;
2604:   mesh->childSection        = NULL;
2605:   mesh->children            = NULL;
2606:   mesh->referenceTree       = NULL;
2607:   mesh->getchildsymmetry    = NULL;
2608:   mesh->vtkCellHeight       = 0;
2609:   mesh->useAnchors          = PETSC_FALSE;

2611:   mesh->maxProjectionHeight = 0;

2613:   mesh->neighbors           = NULL;

2615:   mesh->printSetValues = PETSC_FALSE;
2616:   mesh->printFEM       = 0;
2617:   mesh->printTol       = 1.0e-10;

2619:   DMInitialize_Plex(dm);
2620:   return(0);
2621: }

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

2626:   Collective

2628:   Input Parameter:
2629: . comm - The communicator for the DMPlex object

2631:   Output Parameter:
2632: . mesh  - The DMPlex object

2634:   Level: beginner

2636: @*/
2637: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
2638: {

2643:   DMCreate(comm, mesh);
2644:   DMSetType(*mesh, DMPLEX);
2645:   return(0);
2646: }

2648: static PetscErrorCode DMPlexInvertCell_Internal(PetscInt dim, PetscInt numCorners, PetscInt cone[])
2649: {

2653:   if (dim != 3) return(0);
2654:   switch (numCorners) {
2655:   case 4: DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON,cone); break;
2656:   case 6: DMPlexInvertCell(DM_POLYTOPE_TRI_PRISM,cone); break;
2657:   case 8: DMPlexInvertCell(DM_POLYTOPE_HEXAHEDRON,cone); break;
2658:   default: break;
2659:   }
2660:   return(0);
2661: }

2663: /*
2664:   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
2665: */
2666: /* TODO: invertCells and spaceDim arguments could be added also to to DMPlexCreateFromCellListParallel(), DMPlexBuildFromCellList_Internal() and DMPlexCreateFromCellList() */
2667: PetscErrorCode DMPlexBuildFromCellList_Parallel_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscBool invertCells, PetscSF *sfVert)
2668: {
2669:   PetscSF         sfPoint;
2670:   PetscLayout     vLayout;
2671:   PetscHSetI      vhash;
2672:   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
2673:   const PetscInt *vrange;
2674:   PetscInt        numVerticesAdj, off = 0, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
2675:   PetscMPIInt     rank, size;
2676:   PetscErrorCode  ierr;

2679:   PetscLogEventBegin(DMPLEX_CreateFromCellList,dm,0,0,0);
2680:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
2681:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);
2682:   /* Partition vertices */
2683:   PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);
2684:   PetscLayoutSetLocalSize(vLayout, numVertices);
2685:   PetscLayoutSetBlockSize(vLayout, 1);
2686:   PetscLayoutSetUp(vLayout);
2687:   PetscLayoutGetRanges(vLayout, &vrange);
2688:   /* Count vertices and map them to procs */
2689:   PetscHSetICreate(&vhash);
2690:   for (c = 0; c < numCells; ++c) {
2691:     for (p = 0; p < numCorners; ++p) {
2692:       PetscHSetIAdd(vhash, cells[c*numCorners+p]);
2693:     }
2694:   }
2695:   PetscHSetIGetSize(vhash, &numVerticesAdj);
2696:   PetscMalloc1(numVerticesAdj, &verticesAdj);
2697:   PetscHSetIGetElems(vhash, &off, verticesAdj);
2698:   PetscHSetIDestroy(&vhash);
2699:   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
2700:   PetscSortInt(numVerticesAdj, verticesAdj);
2701:   PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);
2702:   for (v = 0; v < numVerticesAdj; ++v) {
2703:     const PetscInt gv = verticesAdj[v];
2704:     PetscInt       vrank;

2706:     PetscFindInt(gv, size+1, vrange, &vrank);
2707:     vrank = vrank < 0 ? -(vrank+2) : vrank;
2708:     remoteVerticesAdj[v].index = gv - vrange[vrank];
2709:     remoteVerticesAdj[v].rank  = vrank;
2710:   }
2711:   /* Create cones */
2712:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
2713:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
2714:   DMSetUp(dm);
2715:   DMGetWorkArray(dm, numCorners, MPIU_INT, &cone);
2716:   for (c = 0; c < numCells; ++c) {
2717:     for (p = 0; p < numCorners; ++p) {
2718:       const PetscInt gv = cells[c*numCorners+p];
2719:       PetscInt       lv;

2721:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
2722:       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
2723:       cone[p] = lv+numCells;
2724:     }
2725:     if (invertCells) { DMPlexInvertCell_Internal(spaceDim, numCorners, cone); }
2726:     DMPlexSetCone(dm, c, cone);
2727:   }
2728:   DMRestoreWorkArray(dm, numCorners, MPIU_INT, &cone);
2729:   /* Create SF for vertices */
2730:   PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);
2731:   PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");
2732:   PetscSFSetFromOptions(*sfVert);
2733:   PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);
2734:   PetscFree(verticesAdj);
2735:   /* Build pointSF */
2736:   PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);
2737:   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
2738:   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
2739:   PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
2740:   PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
2741:   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);
2742:   PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);
2743:   PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);
2744:   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
2745:   PetscMalloc1(numVerticesGhost, &localVertex);
2746:   PetscMalloc1(numVerticesGhost, &remoteVertex);
2747:   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
2748:     if (vertexLocal[v].rank != rank) {
2749:       localVertex[g]        = v+numCells;
2750:       remoteVertex[g].index = vertexLocal[v].index;
2751:       remoteVertex[g].rank  = vertexLocal[v].rank;
2752:       ++g;
2753:     }
2754:   }
2755:   PetscFree2(vertexLocal, vertexOwner);
2756:   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
2757:   DMGetPointSF(dm, &sfPoint);
2758:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
2759:   PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);
2760:   PetscLayoutDestroy(&vLayout);
2761:   /* Fill in the rest of the topology structure */
2762:   DMPlexSymmetrize(dm);
2763:   DMPlexStratify(dm);
2764:   PetscLogEventEnd(DMPLEX_CreateFromCellList,dm,0,0,0);
2765:   return(0);
2766: }

2768: /*
2769:   This takes as input the coordinates for each owned vertex
2770: */
2771: PetscErrorCode DMPlexBuildCoordinates_Parallel_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[])
2772: {
2773:   PetscSection   coordSection;
2774:   Vec            coordinates;
2775:   PetscScalar   *coords;
2776:   PetscInt       numVertices, numVerticesAdj, coordSize, v;

2780:   PetscLogEventBegin(DMPLEX_CreateFromCellList_Coordinates,dm,0,0,0);
2781:   DMSetCoordinateDim(dm, spaceDim);
2782:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
2783:   DMGetCoordinateSection(dm, &coordSection);
2784:   PetscSectionSetNumFields(coordSection, 1);
2785:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
2786:   PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);
2787:   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
2788:     PetscSectionSetDof(coordSection, v, spaceDim);
2789:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
2790:   }
2791:   PetscSectionSetUp(coordSection);
2792:   PetscSectionGetStorageSize(coordSection, &coordSize);
2793:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
2794:   VecSetBlockSize(coordinates, spaceDim);
2795:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2796:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2797:   VecSetType(coordinates,VECSTANDARD);
2798:   VecGetArray(coordinates, &coords);
2799:   {
2800:     MPI_Datatype coordtype;

2802:     /* Need a temp buffer for coords if we have complex/single */
2803:     MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);
2804:     MPI_Type_commit(&coordtype);
2805: #if defined(PETSC_USE_COMPLEX)
2806:     {
2807:     PetscScalar *svertexCoords;
2808:     PetscInt    i;
2809:     PetscMalloc1(numV*spaceDim,&svertexCoords);
2810:     for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
2811:     PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);
2812:     PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);
2813:     PetscFree(svertexCoords);
2814:     }
2815: #else
2816:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);
2817:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);
2818: #endif
2819:     MPI_Type_free(&coordtype);
2820:   }
2821:   VecRestoreArray(coordinates, &coords);
2822:   DMSetCoordinatesLocal(dm, coordinates);
2823:   VecDestroy(&coordinates);
2824:   PetscLogEventEnd(DMPLEX_CreateFromCellList_Coordinates,dm,0,0,0);
2825:   return(0);
2826: }

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

2831:   Input Parameters:
2832: + comm - The communicator
2833: . dim - The topological dimension of the mesh
2834: . numCells - The number of cells owned by this process
2835: . numVertices - The number of vertices owned by this process
2836: . numCorners - The number of vertices for each cell
2837: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
2838: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
2839: . spaceDim - The spatial dimension used for coordinates
2840: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

2842:   Output Parameter:
2843: + dm - The DM
2844: - vertexSF - Optional, SF describing complete vertex ownership

2846:   Note: Two triangles sharing a face
2847: $
2848: $        2
2849: $      / | \
2850: $     /  |  \
2851: $    /   |   \
2852: $   0  0 | 1  3
2853: $    \   |   /
2854: $     \  |  /
2855: $      \ | /
2856: $        1
2857: would have input
2858: $  numCells = 2, numVertices = 4
2859: $  cells = [0 1 2  1 3 2]
2860: $
2861: which would result in the DMPlex
2862: $
2863: $        4
2864: $      / | \
2865: $     /  |  \
2866: $    /   |   \
2867: $   2  0 | 1  5
2868: $    \   |   /
2869: $     \  |  /
2870: $      \ | /
2871: $        3

2873:   Level: beginner

2875: .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
2876: @*/
2877: 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)
2878: {
2879:   PetscSF        sfVert;

2883:   DMCreate(comm, dm);
2884:   DMSetType(*dm, DMPLEX);
2887:   DMSetDimension(*dm, dim);
2888:   DMPlexBuildFromCellList_Parallel_Internal(*dm, spaceDim, numCells, numVertices, numCorners, cells, PETSC_FALSE, &sfVert);
2889:   if (interpolate) {
2890:     DM idm;

2892:     DMPlexInterpolate(*dm, &idm);
2893:     DMDestroy(dm);
2894:     *dm  = idm;
2895:   }
2896:   DMPlexBuildCoordinates_Parallel_Internal(*dm, spaceDim, numCells, numVertices, sfVert, vertexCoords);
2897:   if (vertexSF) *vertexSF = sfVert;
2898:   else {PetscSFDestroy(&sfVert);}
2899:   return(0);
2900: }

2902: /*
2903:   This takes as input the common mesh generator output, a list of the vertices for each cell
2904: */
2905: PetscErrorCode DMPlexBuildFromCellList_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscBool invertCells)
2906: {
2907:   PetscInt      *cone, c, p;

2911:   PetscLogEventBegin(DMPLEX_CreateFromCellList,dm,0,0,0);
2912:   DMPlexSetChart(dm, 0, numCells+numVertices);
2913:   for (c = 0; c < numCells; ++c) {
2914:     DMPlexSetConeSize(dm, c, numCorners);
2915:   }
2916:   DMSetUp(dm);
2917:   DMGetWorkArray(dm, numCorners, MPIU_INT, &cone);
2918:   for (c = 0; c < numCells; ++c) {
2919:     for (p = 0; p < numCorners; ++p) {
2920:       cone[p] = cells[c*numCorners+p]+numCells;
2921:     }
2922:     if (invertCells) { DMPlexInvertCell_Internal(spaceDim, numCorners, cone); }
2923:     DMPlexSetCone(dm, c, cone);
2924:   }
2925:   DMRestoreWorkArray(dm, numCorners, MPIU_INT, &cone);
2926:   DMPlexSymmetrize(dm);
2927:   DMPlexStratify(dm);
2928:   PetscLogEventEnd(DMPLEX_CreateFromCellList,dm,0,0,0);
2929:   return(0);
2930: }

2932: /*
2933:   This takes as input the coordinates for each vertex
2934: */
2935: PetscErrorCode DMPlexBuildCoordinates_Internal(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
2936: {
2937:   PetscSection   coordSection;
2938:   Vec            coordinates;
2939:   DM             cdm;
2940:   PetscScalar   *coords;
2941:   PetscInt       v, d;

2945:   PetscLogEventBegin(DMPLEX_CreateFromCellList_Coordinates,dm,0,0,0);
2946:   DMSetCoordinateDim(dm, spaceDim);
2947:   DMGetCoordinateSection(dm, &coordSection);
2948:   PetscSectionSetNumFields(coordSection, 1);
2949:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
2950:   PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
2951:   for (v = numCells; v < numCells+numVertices; ++v) {
2952:     PetscSectionSetDof(coordSection, v, spaceDim);
2953:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
2954:   }
2955:   PetscSectionSetUp(coordSection);

2957:   DMGetCoordinateDM(dm, &cdm);
2958:   DMCreateLocalVector(cdm, &coordinates);
2959:   VecSetBlockSize(coordinates, spaceDim);
2960:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2961:   VecGetArray(coordinates, &coords);
2962:   for (v = 0; v < numVertices; ++v) {
2963:     for (d = 0; d < spaceDim; ++d) {
2964:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
2965:     }
2966:   }
2967:   VecRestoreArray(coordinates, &coords);
2968:   DMSetCoordinatesLocal(dm, coordinates);
2969:   VecDestroy(&coordinates);
2970:   PetscLogEventEnd(DMPLEX_CreateFromCellList_Coordinates,dm,0,0,0);
2971:   return(0);
2972: }

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

2977:   Input Parameters:
2978: + comm - The communicator
2979: . dim - The topological dimension of the mesh
2980: . numCells - The number of cells
2981: . numVertices - The number of vertices
2982: . numCorners - The number of vertices for each cell
2983: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
2984: . cells - An array of numCells*numCorners numbers, the vertices for each cell
2985: . spaceDim - The spatial dimension used for coordinates
2986: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

2988:   Output Parameter:
2989: . dm - The DM

2991:   Note: Two triangles sharing a face
2992: $
2993: $        2
2994: $      / | \
2995: $     /  |  \
2996: $    /   |   \
2997: $   0  0 | 1  3
2998: $    \   |   /
2999: $     \  |  /
3000: $      \ | /
3001: $        1
3002: would have input
3003: $  numCells = 2, numVertices = 4
3004: $  cells = [0 1 2  1 3 2]
3005: $
3006: which would result in the DMPlex
3007: $
3008: $        4
3009: $      / | \
3010: $     /  |  \
3011: $    /   |   \
3012: $   2  0 | 1  5
3013: $    \   |   /
3014: $     \  |  /
3015: $      \ | /
3016: $        3

3018:   Level: beginner

3020: .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
3021: @*/
3022: 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)
3023: {

3027:   if (!dim) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "This is not appropriate for 0-dimensional meshes. Consider either creating the DM using DMPlexCreateFromDAG(), by hand, or using DMSwarm.");
3028:   DMCreate(comm, dm);
3029:   DMSetType(*dm, DMPLEX);
3030:   DMSetDimension(*dm, dim);
3031:   DMPlexBuildFromCellList_Internal(*dm, spaceDim, numCells, numVertices, numCorners, cells, PETSC_FALSE);
3032:   if (interpolate) {
3033:     DM idm;

3035:     DMPlexInterpolate(*dm, &idm);
3036:     DMDestroy(dm);
3037:     *dm  = idm;
3038:   }
3039:   DMPlexBuildCoordinates_Internal(*dm, spaceDim, numCells, numVertices, vertexCoords);
3040:   return(0);
3041: }

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

3046:   Input Parameters:
3047: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
3048: . depth - The depth of the DAG
3049: . numPoints - Array of size depth + 1 containing the number of points at each depth
3050: . coneSize - The cone size of each point
3051: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
3052: . coneOrientations - The orientation of each cone point
3053: - vertexCoords - An array of numPoints[0]*spacedim numbers representing the coordinates of each vertex, with spacedim the value set via DMSetCoordinateDim()

3055:   Output Parameter:
3056: . dm - The DM

3058:   Note: Two triangles sharing a face would have input
3059: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
3060: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
3061: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
3062: $
3063: which would result in the DMPlex
3064: $
3065: $        4
3066: $      / | \
3067: $     /  |  \
3068: $    /   |   \
3069: $   2  0 | 1  5
3070: $    \   |   /
3071: $     \  |  /
3072: $      \ | /
3073: $        3
3074: $
3075: $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()

3077:   Level: advanced

3079: .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
3080: @*/
3081: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
3082: {
3083:   Vec            coordinates;
3084:   PetscSection   coordSection;
3085:   PetscScalar    *coords;
3086:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

3090:   DMGetDimension(dm, &dim);
3091:   DMGetCoordinateDim(dm, &dimEmbed);
3092:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
3093:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
3094:   DMPlexSetChart(dm, pStart, pEnd);
3095:   for (p = pStart; p < pEnd; ++p) {
3096:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
3097:     if (firstVertex < 0 && !coneSize[p - pStart]) {
3098:       firstVertex = p - pStart;
3099:     }
3100:   }
3101:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
3102:   DMSetUp(dm); /* Allocate space for cones */
3103:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
3104:     DMPlexSetCone(dm, p, &cones[off]);
3105:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
3106:   }
3107:   DMPlexSymmetrize(dm);
3108:   DMPlexStratify(dm);
3109:   /* Build coordinates */
3110:   DMGetCoordinateSection(dm, &coordSection);
3111:   PetscSectionSetNumFields(coordSection, 1);
3112:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
3113:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
3114:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
3115:     PetscSectionSetDof(coordSection, v, dimEmbed);
3116:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
3117:   }
3118:   PetscSectionSetUp(coordSection);
3119:   PetscSectionGetStorageSize(coordSection, &coordSize);
3120:   VecCreate(PETSC_COMM_SELF, &coordinates);
3121:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3122:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3123:   VecSetBlockSize(coordinates, dimEmbed);
3124:   VecSetType(coordinates,VECSTANDARD);
3125:   VecGetArray(coordinates, &coords);
3126:   for (v = 0; v < numPoints[0]; ++v) {
3127:     PetscInt off;

3129:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
3130:     for (d = 0; d < dimEmbed; ++d) {
3131:       coords[off+d] = vertexCoords[v*dimEmbed+d];
3132:     }
3133:   }
3134:   VecRestoreArray(coordinates, &coords);
3135:   DMSetCoordinatesLocal(dm, coordinates);
3136:   VecDestroy(&coordinates);
3137:   return(0);
3138: }

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

3143: + comm        - The MPI communicator
3144: . filename    - Name of the .dat file
3145: - interpolate - Create faces and edges in the mesh

3147:   Output Parameter:
3148: . dm  - The DM object representing the mesh

3150:   Note: The format is the simplest possible:
3151: $ Ne
3152: $ v0 v1 ... vk
3153: $ Nv
3154: $ x y z marker

3156:   Level: beginner

3158: .seealso: DMPlexCreateFromFile(), DMPlexCreateMedFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
3159: @*/
3160: PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3161: {
3162:   DMLabel         marker;
3163:   PetscViewer     viewer;
3164:   Vec             coordinates;
3165:   PetscSection    coordSection;
3166:   PetscScalar    *coords;
3167:   char            line[PETSC_MAX_PATH_LEN];
3168:   PetscInt        dim = 3, cdim = 3, coordSize, v, c, d;
3169:   PetscMPIInt     rank;
3170:   int             snum, Nv, Nc;
3171:   PetscErrorCode  ierr;

3174:   MPI_Comm_rank(comm, &rank);
3175:   PetscViewerCreate(comm, &viewer);
3176:   PetscViewerSetType(viewer, PETSCVIEWERASCII);
3177:   PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3178:   PetscViewerFileSetName(viewer, filename);
3179:   if (!rank) {
3180:     PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);
3181:     snum = sscanf(line, "%d %d", &Nc, &Nv);
3182:     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3183:   } else {
3184:     Nc = Nv = 0;
3185:   }
3186:   DMCreate(comm, dm);
3187:   DMSetType(*dm, DMPLEX);
3188:   DMPlexSetChart(*dm, 0, Nc+Nv);
3189:   DMSetDimension(*dm, dim);
3190:   DMSetCoordinateDim(*dm, cdim);
3191:   /* Read topology */
3192:   if (!rank) {
3193:     PetscInt cone[8], corners = 8;
3194:     int      vbuf[8], v;

3196:     for (c = 0; c < Nc; ++c) {DMPlexSetConeSize(*dm, c, corners);}
3197:     DMSetUp(*dm);
3198:     for (c = 0; c < Nc; ++c) {
3199:       PetscViewerRead(viewer, line, corners, NULL, PETSC_STRING);
3200:       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]);
3201:       if (snum != corners) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3202:       for (v = 0; v < corners; ++v) cone[v] = vbuf[v] + Nc;
3203:       /* Hexahedra are inverted */
3204:       {
3205:         PetscInt tmp = cone[1];
3206:         cone[1] = cone[3];
3207:         cone[3] = tmp;
3208:       }
3209:       DMPlexSetCone(*dm, c, cone);
3210:     }
3211:   }
3212:   DMPlexSymmetrize(*dm);
3213:   DMPlexStratify(*dm);
3214:   /* Read coordinates */
3215:   DMGetCoordinateSection(*dm, &coordSection);
3216:   PetscSectionSetNumFields(coordSection, 1);
3217:   PetscSectionSetFieldComponents(coordSection, 0, cdim);
3218:   PetscSectionSetChart(coordSection, Nc, Nc + Nv);
3219:   for (v = Nc; v < Nc+Nv; ++v) {
3220:     PetscSectionSetDof(coordSection, v, cdim);
3221:     PetscSectionSetFieldDof(coordSection, v, 0, cdim);
3222:   }
3223:   PetscSectionSetUp(coordSection);
3224:   PetscSectionGetStorageSize(coordSection, &coordSize);
3225:   VecCreate(PETSC_COMM_SELF, &coordinates);
3226:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
3227:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
3228:   VecSetBlockSize(coordinates, cdim);
3229:   VecSetType(coordinates, VECSTANDARD);
3230:   VecGetArray(coordinates, &coords);
3231:   if (!rank) {
3232:     double x[3];
3233:     int    val;

3235:     DMCreateLabel(*dm, "marker");
3236:     DMGetLabel(*dm, "marker", &marker);
3237:     for (v = 0; v < Nv; ++v) {
3238:       PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING);
3239:       snum = sscanf(line, "%lg %lg %lg %d", &x[0], &x[1], &x[2], &val);
3240:       if (snum != 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
3241:       for (d = 0; d < cdim; ++d) coords[v*cdim+d] = x[d];
3242:       if (val) {DMLabelSetValue(marker, v+Nc, val);}
3243:     }
3244:   }
3245:   VecRestoreArray(coordinates, &coords);
3246:   DMSetCoordinatesLocal(*dm, coordinates);
3247:   VecDestroy(&coordinates);
3248:   PetscViewerDestroy(&viewer);
3249:   if (interpolate) {
3250:     DM      idm;
3251:     DMLabel bdlabel;

3253:     DMPlexInterpolate(*dm, &idm);
3254:     DMDestroy(dm);
3255:     *dm  = idm;

3257:     DMGetLabel(*dm, "marker", &bdlabel);
3258:     DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel);
3259:     DMPlexLabelComplete(*dm, bdlabel);
3260:   }
3261:   return(0);
3262: }

3264: /*@C
3265:   DMPlexCreateFromFile - This takes a filename and produces a DM

3267:   Input Parameters:
3268: + comm - The communicator
3269: . filename - A file name
3270: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

3272:   Output Parameter:
3273: . dm - The DM

3275:   Options Database Keys:
3276: . -dm_plex_create_from_hdf5_xdmf - use the PETSC_VIEWER_HDF5_XDMF format for reading HDF5

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

3281:   Level: beginner

3283: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
3284: @*/
3285: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
3286: {
3287:   const char    *extGmsh    = ".msh";
3288:   const char    *extGmsh2   = ".msh2";
3289:   const char    *extGmsh4   = ".msh4";
3290:   const char    *extCGNS    = ".cgns";
3291:   const char    *extExodus  = ".exo";
3292:   const char    *extGenesis = ".gen";
3293:   const char    *extFluent  = ".cas";
3294:   const char    *extHDF5    = ".h5";
3295:   const char    *extMed     = ".med";
3296:   const char    *extPLY     = ".ply";
3297:   const char    *extCV      = ".dat";
3298:   size_t         len;
3299:   PetscBool      isGmsh, isGmsh2, isGmsh4, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isCV;
3300:   PetscMPIInt    rank;

3306:   DMInitializePackage();
3307:   PetscLogEventBegin(DMPLEX_CreateFromFile,0,0,0,0);
3308:   MPI_Comm_rank(comm, &rank);
3309:   PetscStrlen(filename, &len);
3310:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
3311:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,    4, &isGmsh);
3312:   PetscStrncmp(&filename[PetscMax(0,len-5)], extGmsh2,   5, &isGmsh2);
3313:   PetscStrncmp(&filename[PetscMax(0,len-5)], extGmsh4,   5, &isGmsh4);
3314:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,    5, &isCGNS);
3315:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus,  4, &isExodus);
3316:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);
3317:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent,  4, &isFluent);
3318:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,    3, &isHDF5);
3319:   PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,     4, &isMed);
3320:   PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY,     4, &isPLY);
3321:   PetscStrncmp(&filename[PetscMax(0,len-4)], extCV,      4, &isCV);
3322:   if (isGmsh || isGmsh2 || isGmsh4) {
3323:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
3324:   } else if (isCGNS) {
3325:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
3326:   } else if (isExodus || isGenesis) {
3327:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
3328:   } else if (isFluent) {
3329:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
3330:   } else if (isHDF5) {
3331:     PetscBool      load_hdf5_xdmf = PETSC_FALSE;
3332:     PetscViewer viewer;

3334:     /* PETSC_VIEWER_HDF5_XDMF is used if the filename ends with .xdmf.h5, or if -dm_plex_create_from_hdf5_xdmf option is present */
3335:     PetscStrncmp(&filename[PetscMax(0,len-8)], ".xdmf",  5, &load_hdf5_xdmf);
3336:     PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &load_hdf5_xdmf, NULL);
3337:     PetscViewerCreate(comm, &viewer);
3338:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
3339:     PetscViewerSetOptionsPrefix(viewer, "dm_plex_create_");
3340:     PetscViewerSetFromOptions(viewer);
3341:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
3342:     PetscViewerFileSetName(viewer, filename);
3343:     DMCreate(comm, dm);
3344:     DMSetType(*dm, DMPLEX);
3345:     if (load_hdf5_xdmf) {PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF);}
3346:     DMLoad(*dm, viewer);
3347:     if (load_hdf5_xdmf) {PetscViewerPopFormat(viewer);}
3348:     PetscViewerDestroy(&viewer);

3350:     if (interpolate) {
3351:       DM idm;

3353:       DMPlexInterpolate(*dm, &idm);
3354:       DMDestroy(dm);
3355:       *dm  = idm;
3356:     }
3357:   } else if (isMed) {
3358:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
3359:   } else if (isPLY) {
3360:     DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);
3361:   } else if (isCV) {
3362:     DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm);
3363:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
3364:   PetscLogEventEnd(DMPLEX_CreateFromFile,0,0,0,0);
3365:   return(0);
3366: }

3368: /*@
3369:   DMPlexCreateReferenceCellByType - Create a DMPLEX with the appropriate FEM reference cell

3371:   Collective

3373:   Input Parameters:
3374: + comm - The communicator
3375: - ct   - The cell type of the reference cell

3377:   Output Parameter:
3378: . refdm - The reference cell

3380:   Level: intermediate

3382: .seealso: DMPlexCreateReferenceCell(), DMPlexCreateBoxMesh()
3383: @*/
3384: PetscErrorCode DMPlexCreateReferenceCellByType(MPI_Comm comm, DMPolytopeType ct, DM *refdm)
3385: {
3386:   DM             rdm;
3387:   Vec            coords;

3391:   DMCreate(comm, &rdm);
3392:   DMSetType(rdm, DMPLEX);
3393:   switch (ct) {
3394:     case DM_POLYTOPE_POINT:
3395:     {
3396:       PetscInt    numPoints[1]        = {1};
3397:       PetscInt    coneSize[1]         = {0};
3398:       PetscInt    cones[1]            = {0};
3399:       PetscInt    coneOrientations[1] = {0};
3400:       PetscScalar vertexCoords[1]     = {0.0};

3402:       DMSetDimension(rdm, 0);
3403:       DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3404:     }
3405:     break;
3406:     case DM_POLYTOPE_SEGMENT:
3407:     {
3408:       PetscInt    numPoints[2]        = {2, 1};
3409:       PetscInt    coneSize[3]         = {2, 0, 0};
3410:       PetscInt    cones[2]            = {1, 2};
3411:       PetscInt    coneOrientations[2] = {0, 0};
3412:       PetscScalar vertexCoords[2]     = {-1.0,  1.0};

3414:       DMSetDimension(rdm, 1);
3415:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3416:     }
3417:     break;
3418:     case DM_POLYTOPE_TRIANGLE:
3419:     {
3420:       PetscInt    numPoints[2]        = {3, 1};
3421:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
3422:       PetscInt    cones[3]            = {1, 2, 3};
3423:       PetscInt    coneOrientations[3] = {0, 0, 0};
3424:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

3426:       DMSetDimension(rdm, 2);
3427:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3428:     }
3429:     break;
3430:     case DM_POLYTOPE_QUADRILATERAL:
3431:     {
3432:       PetscInt    numPoints[2]        = {4, 1};
3433:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3434:       PetscInt    cones[4]            = {1, 2, 3, 4};
3435:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3436:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

3438:       DMSetDimension(rdm, 2);
3439:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3440:     }
3441:     break;
3442:     case DM_POLYTOPE_SEG_PRISM_TENSOR:
3443:     {
3444:       PetscInt    numPoints[2]        = {4, 1};
3445:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3446:       PetscInt    cones[4]            = {1, 2, 3, 4};
3447:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3448:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0,  1.0, 1.0};

3450:       DMSetDimension(rdm, 2);
3451:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3452:     }
3453:     break;
3454:     case DM_POLYTOPE_TETRAHEDRON:
3455:     {
3456:       PetscInt    numPoints[2]        = {4, 1};
3457:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
3458:       PetscInt    cones[4]            = {1, 3, 2, 4};
3459:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
3460:       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};

3462:       DMSetDimension(rdm, 3);
3463:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3464:     }
3465:     break;
3466:     case DM_POLYTOPE_HEXAHEDRON:
3467:     {
3468:       PetscInt    numPoints[2]        = {8, 1};
3469:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3470:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
3471:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3472:       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,
3473:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3475:       DMSetDimension(rdm, 3);
3476:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3477:     }
3478:     break;
3479:     case DM_POLYTOPE_TRI_PRISM:
3480:     {
3481:       PetscInt    numPoints[2]        = {6, 1};
3482:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3483:       PetscInt    cones[6]            = {1, 3, 2, 4, 5, 6};
3484:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3485:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3486:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3488:       DMSetDimension(rdm, 3);
3489:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3490:     }
3491:     break;
3492:     case DM_POLYTOPE_TRI_PRISM_TENSOR:
3493:     {
3494:       PetscInt    numPoints[2]        = {6, 1};
3495:       PetscInt    coneSize[7]         = {6, 0, 0, 0, 0, 0, 0};
3496:       PetscInt    cones[6]            = {1, 2, 3, 4, 5, 6};
3497:       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3498:       PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,  1.0, -1.0, -1.0,  -1.0, 1.0, -1.0,
3499:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  -1.0, 1.0,  1.0};

3501:       DMSetDimension(rdm, 3);
3502:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3503:     }
3504:     break;
3505:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
3506:     {
3507:       PetscInt    numPoints[2]        = {8, 1};
3508:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3509:       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
3510:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3511:       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,
3512:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

3514:       DMSetDimension(rdm, 3);
3515:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
3516:     }
3517:     break;
3518:     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for cell type %s", DMPolytopeTypes[ct]);
3519:   }
3520:   {
3521:     PetscInt Nv, v;

3523:     /* Must create the celltype label here so that we do not automatically try to compute the types */
3524:     DMCreateLabel(rdm, "celltype");
3525:     DMPlexSetCellType(rdm, 0, ct);
3526:     DMPlexGetChart(rdm, NULL, &Nv);
3527:     for (v = 1; v < Nv; ++v) {DMPlexSetCellType(rdm, v, DM_POLYTOPE_POINT);}
3528:   }
3529:   DMPlexInterpolate(rdm, refdm);
3530:   if (rdm->coordinateDM) {
3531:     DM           ncdm;
3532:     PetscSection cs;
3533:     PetscInt     pEnd = -1;

3535:     DMGetLocalSection(rdm->coordinateDM, &cs);
3536:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
3537:     if (pEnd >= 0) {
3538:       DMClone(*refdm, &ncdm);
3539:       DMCopyDisc(rdm->coordinateDM, ncdm);
3540:       DMSetLocalSection(ncdm, cs);
3541:       DMSetCoordinateDM(*refdm, ncdm);
3542:       DMDestroy(&ncdm);
3543:     }
3544:   }
3545:   DMGetCoordinatesLocal(rdm, &coords);
3546:   if (coords) {
3547:     DMSetCoordinatesLocal(*refdm, coords);
3548:   } else {
3549:     DMGetCoordinates(rdm, &coords);
3550:     if (coords) {DMSetCoordinates(*refdm, coords);}
3551:   }
3552:   DMDestroy(&rdm);
3553:   return(0);
3554: }

3556: /*@
3557:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

3559:   Collective

3561:   Input Parameters:
3562: + comm    - The communicator
3563: . dim     - The spatial dimension
3564: - simplex - Flag for simplex, otherwise use a tensor-product cell

3566:   Output Parameter:
3567: . refdm - The reference cell

3569:   Level: intermediate

3571: .seealso: DMPlexCreateReferenceCellByType(), DMPlexCreateBoxMesh()
3572: @*/
3573: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
3574: {

3578:   switch (dim) {
3579:   case 0: DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_POINT, refdm);break;
3580:   case 1: DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_SEGMENT, refdm);break;
3581:   case 2:
3582:     if (simplex) {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_TRIANGLE, refdm);}
3583:     else         {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_QUADRILATERAL, refdm);}
3584:     break;
3585:   case 3:
3586:     if (simplex) {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_TETRAHEDRON, refdm);}
3587:     else         {DMPlexCreateReferenceCellByType(comm, DM_POLYTOPE_HEXAHEDRON, refdm);}
3588:     break;
3589:   default:
3590:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
3591:   }
3592:   return(0);
3593: }