Actual source code: plexcreate.c

petsc-master 2017-05-23
Report Typos and Errors
  1: #define PETSCDM_DLL
  2:  #include <petsc/private/dmpleximpl.h>
  3:  #include <petscdmda.h>
  4:  #include <petscsf.h>

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

  9:   Collective on MPI_Comm

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

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

 22:   Level: beginner

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

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

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

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

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

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

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

121:     DMPlexInterpolate(*newdm, &idm);
122:     PetscObjectGetName((PetscObject) *newdm, &name);
123:     PetscObjectSetName((PetscObject)    idm,  name);
124:     DMPlexCopyCoordinates(*newdm, idm);
125:     DMCopyLabels(*newdm, idm);
126:     DMDestroy(newdm);
127:     *newdm = idm;
128:   }
129:   {
130:     DM refinedMesh     = NULL;
131:     DM distributedMesh = NULL;

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

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

154:   Collective on MPI_Comm

156:   Input Parameters:
157: + comm  - The communicator for the DM object
158: . lower - The lower left corner coordinates
159: . upper - The upper right corner coordinates
160: - edges - The number of cells in each direction

162:   Output Parameter:
163: . dm  - The DM object

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

176:   Level: beginner

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

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

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

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

244:         cone[0] = vertex; cone[1] = vertex+1;
245:         DMPlexSetCone(dm, edge, cone);
246:         if (vy == edges[1]) {
247:           DMSetLabelValue(dm, "marker", edge,    markerTop);
248:           DMSetLabelValue(dm, "marker", cone[0], markerTop);
249:           if (ex == edges[0]-1) {
250:             DMSetLabelValue(dm, "marker", cone[1], markerTop);
251:           }
252:         } else if (vy == 0) {
253:           DMSetLabelValue(dm, "marker", edge,    markerBottom);
254:           DMSetLabelValue(dm, "marker", cone[0], markerBottom);
255:           if (ex == edges[0]-1) {
256:             DMSetLabelValue(dm, "marker", cone[1], markerBottom);
257:           }
258:         }
259:       }
260:     }
261:   }
262:   DMPlexSymmetrize(dm);
263:   DMPlexStratify(dm);
264:   /* Build coordinates */
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 the is the boundary of a cubic lattice.

296:   Collective on MPI_Comm

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

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

307:   Level: beginner

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

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

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

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

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

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

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

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

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

436: /*@
437:   DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices.

439:   Collective on MPI_Comm

441:   Input Parameters:
442: + comm - The communicator for the DM object
443: . dim - The spatial dimension
444: . numFaces - Number of faces per dimension
445: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

447:   Output Parameter:
448: . dm  - The DM object

450:   Level: beginner

452: .keywords: DM, create
453: .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
454: @*/
455: PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm)
456: {
457:   DM             boundary;

462:   DMCreate(comm, &boundary);
464:   DMSetType(boundary, DMPLEX);
465:   DMSetDimension(boundary, dim-1);
466:   DMSetCoordinateDim(boundary, dim);
467:   switch (dim) {
468:   case 2:
469:   {
470:     PetscReal lower[2] = {0.0, 0.0};
471:     PetscReal upper[2] = {1.0, 1.0};
472:     PetscInt  edges[2];

474:     edges[0] = numFaces; edges[1] = numFaces;
475:     DMPlexCreateSquareBoundary(boundary, lower, upper, edges);
476:     break;
477:   }
478:   case 3:
479:   {
480:     PetscReal lower[3] = {0.0, 0.0, 0.0};
481:     PetscReal upper[3] = {1.0, 1.0, 1.0};
482:     PetscInt  faces[3];

484:     faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces;
485:     DMPlexCreateCubeBoundary(boundary, lower, upper, faces);
486:     break;
487:   }
488:   default:
489:     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
490:   }
491:   DMPlexGenerate(boundary, NULL, interpolate, dm);
492:   DMDestroy(&boundary);
493:   return(0);
494: }

496: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
497: {
498:   DMLabel        cutLabel = NULL;
499:   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
500:   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
501:   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
502:   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
503:   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
504:   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
505:   PetscInt       dim;
506:   PetscBool      markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
507:   PetscMPIInt    rank;

511:   DMGetDimension(dm,&dim);
512:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
513:   DMCreateLabel(dm,"marker");
514:   DMCreateLabel(dm,"Face Sets");
515:   if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ||
516:       bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ||
517:       bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {
518:     PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);
519:     if (cutMarker) {DMCreateLabel(dm, "periodic_cut"); DMGetLabel(dm, "periodic_cut", &cutLabel);}
520:   }
521:   switch (dim) {
522:   case 2:
523:     faceMarkerTop    = 3;
524:     faceMarkerBottom = 1;
525:     faceMarkerRight  = 2;
526:     faceMarkerLeft   = 4;
527:     break;
528:   case 3:
529:     faceMarkerBottom = 1;
530:     faceMarkerTop    = 2;
531:     faceMarkerFront  = 3;
532:     faceMarkerBack   = 4;
533:     faceMarkerRight  = 5;
534:     faceMarkerLeft   = 6;
535:     break;
536:   default:
537:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim);
538:     break;
539:   }
540:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
541:   if (markerSeparate) {
542:     markerBottom = faceMarkerBottom;
543:     markerTop    = faceMarkerTop;
544:     markerFront  = faceMarkerFront;
545:     markerBack   = faceMarkerBack;
546:     markerRight  = faceMarkerRight;
547:     markerLeft   = faceMarkerLeft;
548:   }
549:   {
550:     const PetscInt numXEdges    = !rank ? edges[0] : 0;
551:     const PetscInt numYEdges    = !rank ? edges[1] : 0;
552:     const PetscInt numZEdges    = !rank ? edges[2] : 0;
553:     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
554:     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
555:     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
556:     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
557:     const PetscInt numXFaces    = numYEdges*numZEdges;
558:     const PetscInt numYFaces    = numXEdges*numZEdges;
559:     const PetscInt numZFaces    = numXEdges*numYEdges;
560:     const PetscInt numTotXFaces = numXVertices*numXFaces;
561:     const PetscInt numTotYFaces = numYVertices*numYFaces;
562:     const PetscInt numTotZFaces = numZVertices*numZFaces;
563:     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
564:     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
565:     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
566:     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
567:     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
568:     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
569:     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
570:     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
571:     const PetscInt firstYFace   = firstXFace + numTotXFaces;
572:     const PetscInt firstZFace   = firstYFace + numTotYFaces;
573:     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
574:     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
575:     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
576:     Vec            coordinates;
577:     PetscSection   coordSection;
578:     PetscScalar   *coords;
579:     PetscInt       coordSize;
580:     PetscInt       v, vx, vy, vz;
581:     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;

583:     DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);
584:     for (c = 0; c < numCells; c++) {
585:       DMPlexSetConeSize(dm, c, 6);
586:     }
587:     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
588:       DMPlexSetConeSize(dm, f, 4);
589:     }
590:     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
591:       DMPlexSetConeSize(dm, e, 2);
592:     }
593:     DMSetUp(dm); /* Allocate space for cones */
594:     /* Build cells */
595:     for (fz = 0; fz < numZEdges; ++fz) {
596:       for (fy = 0; fy < numYEdges; ++fy) {
597:         for (fx = 0; fx < numXEdges; ++fx) {
598:           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
599:           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
600:           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
601:           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
602:           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
603:           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
604:           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
605:                             /* B,  T,  F,  K,  R,  L */
606:           PetscInt ornt[6] = {-4,  0,  0, -1,  0, -4}; /* ??? */
607:           PetscInt cone[6];

609:           /* no boundary twisting in 3D */
610:           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
611:           DMPlexSetCone(dm, cell, cone);
612:           DMPlexSetConeOrientation(dm, cell, ornt);
613:           if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
614:           if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
615:           if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, cell, 2);}
616:         }
617:       }
618:     }
619:     /* Build x faces */
620:     for (fz = 0; fz < numZEdges; ++fz) {
621:       for (fy = 0; fy < numYEdges; ++fy) {
622:         for (fx = 0; fx < numXVertices; ++fx) {
623:           PetscInt face    = firstXFace + (fz*numYEdges+fy)*numXVertices + fx;
624:           PetscInt edgeL   = firstZEdge + (  fy*                 numXVertices+fx)*numZEdges + fz;
625:           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
626:           PetscInt edgeB   = firstYEdge + (  fz*                 numXVertices+fx)*numYEdges + fy;
627:           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
628:           PetscInt ornt[4] = {0, 0, -2, -2};
629:           PetscInt cone[4];

631:           if (dim == 3) {
632:             /* markers */
633:             if (bdX != DM_BOUNDARY_PERIODIC) {
634:               if (fx == numXVertices-1) {
635:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);
636:                 DMSetLabelValue(dm, "marker", face, markerRight);
637:               }
638:               else if (fx == 0) {
639:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);
640:                 DMSetLabelValue(dm, "marker", face, markerLeft);
641:               }
642:             }
643:           }
644:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
645:           DMPlexSetCone(dm, face, cone);
646:           DMPlexSetConeOrientation(dm, face, ornt);
647:         }
648:       }
649:     }
650:     /* Build y faces */
651:     for (fz = 0; fz < numZEdges; ++fz) {
652:       for (fx = 0; fx < numXEdges; ++fx) {
653:         for (fy = 0; fy < numYVertices; ++fy) {
654:           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
655:           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx                 )*numZEdges + fz;
656:           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
657:           PetscInt edgeB   = firstXEdge + (  fz                 *numYVertices+fy)*numXEdges + fx;
658:           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
659:           PetscInt ornt[4] = {0, 0, -2, -2};
660:           PetscInt cone[4];

662:           if (dim == 3) {
663:             /* markers */
664:             if (bdY != DM_BOUNDARY_PERIODIC) {
665:               if (fy == numYVertices-1) {
666:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);
667:                 DMSetLabelValue(dm, "marker", face, markerBack);
668:               }
669:               else if (fy == 0) {
670:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);
671:                 DMSetLabelValue(dm, "marker", face, markerFront);
672:               }
673:             }
674:           }
675:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
676:           DMPlexSetCone(dm, face, cone);
677:           DMPlexSetConeOrientation(dm, face, ornt);
678:         }
679:       }
680:     }
681:     /* Build z faces */
682:     for (fy = 0; fy < numYEdges; ++fy) {
683:       for (fx = 0; fx < numXEdges; ++fx) {
684:         for (fz = 0; fz < numZVertices; fz++) {
685:           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
686:           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx                 )*numYEdges + fy;
687:           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
688:           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy                 )*numXEdges + fx;
689:           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
690:           PetscInt ornt[4] = {0, 0, -2, -2};
691:           PetscInt cone[4];

693:           if (dim == 2) {
694:             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
695:             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
696:             if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, face, 2);}
697:             if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {DMLabelSetValue(cutLabel, face, 2);}
698:           } else {
699:             /* markers */
700:             if (bdZ != DM_BOUNDARY_PERIODIC) {
701:               if (fz == numZVertices-1) {
702:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);
703:                 DMSetLabelValue(dm, "marker", face, markerTop);
704:               }
705:               else if (fz == 0) {
706:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);
707:                 DMSetLabelValue(dm, "marker", face, markerBottom);
708:               }
709:             }
710:           }
711:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
712:           DMPlexSetCone(dm, face, cone);
713:           DMPlexSetConeOrientation(dm, face, ornt);
714:         }
715:       }
716:     }
717:     /* Build Z edges*/
718:     for (vy = 0; vy < numYVertices; vy++) {
719:       for (vx = 0; vx < numXVertices; vx++) {
720:         for (ez = 0; ez < numZEdges; ez++) {
721:           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
722:           const PetscInt vertexB = firstVertex + (  ez                 *numYVertices+vy)*numXVertices + vx;
723:           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
724:           PetscInt       cone[2];

726:           if (dim == 3) {
727:             if (bdX != DM_BOUNDARY_PERIODIC) {
728:               if (vx == numXVertices-1) {
729:                 DMSetLabelValue(dm, "marker", edge, markerRight);
730:               }
731:               else if (vx == 0) {
732:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
733:               }
734:             }
735:             if (bdY != DM_BOUNDARY_PERIODIC) {
736:               if (vy == numYVertices-1) {
737:                 DMSetLabelValue(dm, "marker", edge, markerBack);
738:               }
739:               else if (vy == 0) {
740:                 DMSetLabelValue(dm, "marker", edge, markerFront);
741:               }
742:             }
743:           }
744:           cone[0] = vertexB; cone[1] = vertexT;
745:           DMPlexSetCone(dm, edge, cone);
746:         }
747:       }
748:     }
749:     /* Build Y edges*/
750:     for (vz = 0; vz < numZVertices; vz++) {
751:       for (vx = 0; vx < numXVertices; vx++) {
752:         for (ey = 0; ey < numYEdges; ey++) {
753:           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
754:           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
755:           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
756:           const PetscInt vertexK = firstVertex + nextv;
757:           PetscInt       cone[2];

759:           cone[0] = vertexF; cone[1] = vertexK;
760:           DMPlexSetCone(dm, edge, cone);
761:           if (dim == 2) {
762:             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
763:               if (vx == numXVertices-1) {
764:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);
765:                 DMSetLabelValue(dm, "marker", edge,    markerRight);
766:                 DMSetLabelValue(dm, "marker", cone[0], markerRight);
767:                 if (ey == numYEdges-1) {
768:                   DMSetLabelValue(dm, "marker", cone[1], markerRight);
769:                 }
770:               } else if (vx == 0) {
771:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);
772:                 DMSetLabelValue(dm, "marker", edge,    markerLeft);
773:                 DMSetLabelValue(dm, "marker", cone[0], markerLeft);
774:                 if (ey == numYEdges-1) {
775:                   DMSetLabelValue(dm, "marker", cone[1], markerLeft);
776:                 }
777:               }
778:             } else {
779:               if (vx == 0 && cutMarker) {
780:                 DMLabelSetValue(cutLabel, edge,    1);
781:                 DMLabelSetValue(cutLabel, cone[0], 1);
782:                 if (ey == numYEdges-1) {
783:                   DMLabelSetValue(cutLabel, cone[1], 1);
784:                 }
785:               }
786:             }
787:           } else {
788:             if (bdX != DM_BOUNDARY_PERIODIC) {
789:               if (vx == numXVertices-1) {
790:                 DMSetLabelValue(dm, "marker", edge, markerRight);
791:               } else if (vx == 0) {
792:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
793:               }
794:             }
795:             if (bdZ != DM_BOUNDARY_PERIODIC) {
796:               if (vz == numZVertices-1) {
797:                 DMSetLabelValue(dm, "marker", edge, markerTop);
798:               } else if (vz == 0) {
799:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
800:               }
801:             }
802:           }
803:         }
804:       }
805:     }
806:     /* Build X edges*/
807:     for (vz = 0; vz < numZVertices; vz++) {
808:       for (vy = 0; vy < numYVertices; vy++) {
809:         for (ex = 0; ex < numXEdges; ex++) {
810:           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
811:           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
812:           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
813:           const PetscInt vertexR = firstVertex + nextv;
814:           PetscInt       cone[2];

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

902: /*@
903:   DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra.

905:   Collective on MPI_Comm

907:   Input Parameters:
908: + comm  - The communicator for the DM object
909: . dim   - The spatial dimension
910: . periodicX - The boundary type for the X direction
911: . periodicY - The boundary type for the Y direction
912: . periodicZ - The boundary type for the Z direction
913: - cells - The number of cells in each direction

915:   Output Parameter:
916: . dm  - The DM object

918:   Note: Here is the numbering returned for 2 cells in each direction:
919: $ 22--8-23--9--24
920: $  |     |     |
921: $ 13  2 14  3  15
922: $  |     |     |
923: $ 19--6-20--7--21
924: $  |     |     |
925: $ 10  0 11  1 12
926: $  |     |     |
927: $ 16--4-17--5--18

929:   Level: beginner

931: .keywords: DM, create
932: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
933: @*/
934: PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
935: {
936:   PetscInt       i;

941:   DMCreate(comm, dm);
943:   DMSetType(*dm, DMPLEX);
944:   DMSetDimension(*dm, dim);
945:   switch (dim) {
946:   case 2:
947:   {
948:     PetscReal lower[3] = {0.0, 0.0, 0.0};
949:     PetscReal upper[3] = {1.0, 1.0, 0.0};
950:     PetscInt  edges[3];

952:     edges[0] = cells[0];
953:     edges[1] = cells[1];
954:     edges[2] = 0;

956:     DMPlexCreateCubeMesh_Internal(*dm, lower, upper, edges, periodicX, periodicY, DM_BOUNDARY_NONE);
957:     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
958:         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST) {
959:       PetscReal      L[2];
960:       PetscReal      maxCell[2];
961:       DMBoundaryType bdType[2];

963:       bdType[0] = periodicX;
964:       bdType[1] = periodicY;
965:       for (i = 0; i < dim; i++) {
966:         L[i]       = upper[i] - lower[i];
967:         maxCell[i] = 1.1 * (L[i] / PetscMax(1,cells[i]));
968:       }

970:       DMSetPeriodicity(*dm,maxCell,L,bdType);
971:     }
972:     break;
973:   }
974:   case 3:
975:   {
976:     PetscReal lower[3] = {0.0, 0.0, 0.0};
977:     PetscReal upper[3] = {1.0, 1.0, 1.0};

979:     DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);
980:     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
981:         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST ||
982:         periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
983:       PetscReal      L[3];
984:       PetscReal      maxCell[3];
985:       DMBoundaryType bdType[3];

987:       bdType[0] = periodicX;
988:       bdType[1] = periodicY;
989:       bdType[2] = periodicZ;
990:       for (i = 0; i < dim; i++) {
991:         L[i]       = upper[i] - lower[i];
992:         maxCell[i] = 1.1 * (L[i] / cells[i]);
993:       }

995:       DMSetPeriodicity(*dm,maxCell,L,bdType);
996:     }
997:     break;
998:   }
999:   default:
1000:     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
1001:   }
1002:   return(0);
1003: }

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

1008:   Logically Collective on DM

1010:   Input Parameters:
1011: + dm - the DM context
1012: - prefix - the prefix to prepend to all option names

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

1018:   Level: advanced

1020: .seealso: SNESSetFromOptions()
1021: @*/
1022: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1023: {
1024:   DM_Plex       *mesh = (DM_Plex *) dm->data;

1029:   PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);
1030:   PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);
1031:   return(0);
1032: }

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

1037:   Collective on MPI_Comm

1039:   Input Parameters:
1040: + comm      - The communicator for the DM object
1041: . numRefine - The number of regular refinements to the basic 5 cell structure
1042: - periodicZ - The boundary type for the Z direction

1044:   Output Parameter:
1045: . dm  - The DM object

1047:   Note: Here is the output numbering looking from the bottom of the cylinder:
1048: $       17-----14
1049: $        |     |
1050: $        |  2  |
1051: $        |     |
1052: $ 17-----8-----7-----14
1053: $  |     |     |     |
1054: $  |  3  |  0  |  1  |
1055: $  |     |     |     |
1056: $ 19-----5-----6-----13
1057: $        |     |
1058: $        |  4  |
1059: $        |     |
1060: $       19-----13
1061: $
1062: $ and up through the top
1063: $
1064: $       18-----16
1065: $        |     |
1066: $        |  2  |
1067: $        |     |
1068: $ 18----10----11-----16
1069: $  |     |     |     |
1070: $  |  3  |  0  |  1  |
1071: $  |     |     |     |
1072: $ 20-----9----12-----15
1073: $        |     |
1074: $        |  4  |
1075: $        |     |
1076: $       20-----15

1078:   Level: beginner

1080: .keywords: DM, create
1081: .seealso: DMPlexCreateHexBoxMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
1082: @*/
1083: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm)
1084: {
1085:   const PetscInt dim = 3;
1086:   PetscInt       numCells, numVertices, r;

1091:   if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine);
1092:   DMCreate(comm, dm);
1093:   DMSetType(*dm, DMPLEX);
1094:   DMSetDimension(*dm, dim);
1095:   /* Create topology */
1096:   {
1097:     PetscInt cone[8], c;

1099:     numCells    = 5;
1100:     numVertices = 16;
1101:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1102:       numCells   *= 3;
1103:       numVertices = 24;
1104:     }
1105:     DMPlexSetChart(*dm, 0, numCells+numVertices);
1106:     for (c = 0; c < numCells; c++) {DMPlexSetConeSize(*dm, c, 8);}
1107:     DMSetUp(*dm);
1108:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1109:       cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16;
1110:       cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34;
1111:       DMPlexSetCone(*dm, 0, cone);
1112:       cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23;
1113:       cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */
1114:       DMPlexSetCone(*dm, 1, cone);
1115:       cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17;
1116:       cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38;
1117:       DMPlexSetCone(*dm, 2, cone);
1118:       cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15;
1119:       cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38;
1120:       DMPlexSetCone(*dm, 3, cone);
1121:       cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23;
1122:       cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31;
1123:       DMPlexSetCone(*dm, 4, cone);

1125:       cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32;
1126:       cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20;
1127:       DMPlexSetCone(*dm, 5, cone);
1128:       cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36;
1129:       cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21;
1130:       DMPlexSetCone(*dm, 6, cone);
1131:       cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33;
1132:       cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28;
1133:       DMPlexSetCone(*dm, 7, cone);
1134:       cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31;
1135:       cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28;
1136:       DMPlexSetCone(*dm, 8, cone);
1137:       cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36;
1138:       cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19;
1139:       DMPlexSetCone(*dm, 9, cone);

1141:       cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22;
1142:       cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18;
1143:       DMPlexSetCone(*dm, 10, cone);
1144:       cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25;
1145:       cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17;
1146:       DMPlexSetCone(*dm, 11, cone);
1147:       cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21;
1148:       cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27;
1149:       DMPlexSetCone(*dm, 12, cone);
1150:       cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19;
1151:       cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27;
1152:       DMPlexSetCone(*dm, 13, cone);
1153:       cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25;
1154:       cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15;
1155:       DMPlexSetCone(*dm, 14, cone);
1156:     } else {
1157:       cone[0] =  5; cone[1] =  8; cone[2] =  7; cone[3] =  6;
1158:       cone[4] =  9; cone[5] = 12; cone[6] = 11; cone[7] = 10;
1159:       DMPlexSetCone(*dm, 0, cone);
1160:       cone[0] =  6; cone[1] =  7; cone[2] = 14; cone[3] = 13;
1161:       cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11;
1162:       DMPlexSetCone(*dm, 1, cone);
1163:       cone[0] =  8; cone[1] = 17; cone[2] = 14; cone[3] =  7;
1164:       cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18;
1165:       DMPlexSetCone(*dm, 2, cone);
1166:       cone[0] = 19; cone[1] = 17; cone[2] =  8; cone[3] =  5;
1167:       cone[4] = 20; cone[5] =  9; cone[6] = 10; cone[7] = 18;
1168:       DMPlexSetCone(*dm, 3, cone);
1169:       cone[0] = 19; cone[1] =  5; cone[2] =  6; cone[3] = 13;
1170:       cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] =  9;
1171:       DMPlexSetCone(*dm, 4, cone);
1172:     }
1173:     DMPlexSymmetrize(*dm);
1174:     DMPlexStratify(*dm);
1175:   }
1176:   /* Interpolate */
1177:   {
1178:     DM idm = NULL;

1180:     DMPlexInterpolate(*dm, &idm);
1181:     DMDestroy(dm);
1182:     *dm  = idm;
1183:   }
1184:   /* Create cube geometry */
1185:   {
1186:     Vec             coordinates;
1187:     PetscSection    coordSection;
1188:     PetscScalar    *coords;
1189:     PetscInt        coordSize, v;
1190:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1191:     const PetscReal ds2 = dis/2.0;

1193:     /* Build coordinates */
1194:     DMGetCoordinateSection(*dm, &coordSection);
1195:     PetscSectionSetNumFields(coordSection, 1);
1196:     PetscSectionSetFieldComponents(coordSection, 0, dim);
1197:     PetscSectionSetChart(coordSection, numCells, numCells+numVertices);
1198:     for (v = numCells; v < numCells+numVertices; ++v) {
1199:       PetscSectionSetDof(coordSection, v, dim);
1200:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
1201:     }
1202:     PetscSectionSetUp(coordSection);
1203:     PetscSectionGetStorageSize(coordSection, &coordSize);
1204:     VecCreate(PETSC_COMM_SELF, &coordinates);
1205:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
1206:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1207:     VecSetBlockSize(coordinates, dim);
1208:     VecSetType(coordinates,VECSTANDARD);
1209:     VecGetArray(coordinates, &coords);
1210:     coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0;
1211:     coords[1*dim+0] =  ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0;
1212:     coords[2*dim+0] =  ds2; coords[2*dim+1] =  ds2; coords[2*dim+2] = 0.0;
1213:     coords[3*dim+0] = -ds2; coords[3*dim+1] =  ds2; coords[3*dim+2] = 0.0;
1214:     coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0;
1215:     coords[5*dim+0] = -ds2; coords[5*dim+1] =  ds2; coords[5*dim+2] = 1.0;
1216:     coords[6*dim+0] =  ds2; coords[6*dim+1] =  ds2; coords[6*dim+2] = 1.0;
1217:     coords[7*dim+0] =  ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0;
1218:     coords[ 8*dim+0] =  dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0;
1219:     coords[ 9*dim+0] =  dis; coords[ 9*dim+1] =  dis; coords[ 9*dim+2] = 0.0;
1220:     coords[10*dim+0] =  dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0;
1221:     coords[11*dim+0] =  dis; coords[11*dim+1] =  dis; coords[11*dim+2] = 1.0;
1222:     coords[12*dim+0] = -dis; coords[12*dim+1] =  dis; coords[12*dim+2] = 0.0;
1223:     coords[13*dim+0] = -dis; coords[13*dim+1] =  dis; coords[13*dim+2] = 1.0;
1224:     coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0;
1225:     coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0;
1226:     if (periodicZ == DM_BOUNDARY_PERIODIC) {
1227:       /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5;
1228:       /* 16 32 22 */ coords[17*dim+0] =  ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5;
1229:       /* 17 33 21 */ coords[18*dim+0] =  ds2; coords[18*dim+1] =  ds2; coords[18*dim+2] = 0.5;
1230:       /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] =  ds2; coords[19*dim+2] = 0.5;
1231:       /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5;
1232:       /* 23 36 25 */ coords[21*dim+0] =  dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5;
1233:       /* 24 37 26 */ coords[22*dim+0] =  dis; coords[22*dim+1] =  dis; coords[22*dim+2] = 0.5;
1234:       /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] =  dis; coords[23*dim+2] = 0.5;
1235:     }
1236:     VecRestoreArray(coordinates, &coords);
1237:     DMSetCoordinatesLocal(*dm, coordinates);
1238:     VecDestroy(&coordinates);
1239:   }
1240:   /* Create periodicity */
1241:   if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1242:     PetscReal      L[3];
1243:     PetscReal      maxCell[3];
1244:     DMBoundaryType bdType[3];
1245:     PetscReal      lower[3] = {0.0, 0.0, 0.0};
1246:     PetscReal      upper[3] = {1.0, 1.0, 1.5};
1247:     PetscInt       i, numZCells = 3;

1249:     bdType[0] = DM_BOUNDARY_NONE;
1250:     bdType[1] = DM_BOUNDARY_NONE;
1251:     bdType[2] = periodicZ;
1252:     for (i = 0; i < dim; i++) {
1253:       L[i]       = upper[i] - lower[i];
1254:       maxCell[i] = 1.1 * (L[i] / numZCells);
1255:     }
1256:     DMSetPeriodicity(*dm, maxCell, L, bdType);
1257:   }
1258:   /* Refine topology */
1259:   for (r = 0; r < numRefine; ++r) {
1260:     DM rdm = NULL;

1262:     DMRefine(*dm, comm, &rdm);
1263:     DMDestroy(dm);
1264:     *dm  = rdm;
1265:   }
1266:   /* Remap geometry to cylinder
1267:        Interior square: Linear interpolation is correct
1268:        The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1269:        such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance

1271:          phi     = arctan(y/x)
1272:          d_close = sqrt(1/8 + 1/4 sin^2(phi))
1273:          d_far   = sqrt(1/2 + sin^2(phi))

1275:        so we remap them using

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

1280:        If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1281:   */
1282:   {
1283:     Vec           coordinates;
1284:     PetscSection  coordSection;
1285:     PetscScalar  *coords;
1286:     PetscInt      vStart, vEnd, v;
1287:     const PetscReal dis = 1.0/PetscSqrtReal(2.0);
1288:     const PetscReal ds2 = 0.5*dis;

1290:     DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);
1291:     DMGetCoordinateSection(*dm, &coordSection);
1292:     DMGetCoordinatesLocal(*dm, &coordinates);
1293:     VecGetArray(coordinates, &coords);
1294:     for (v = vStart; v < vEnd; ++v) {
1295:       PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1296:       PetscInt  off;

1298:       PetscSectionGetOffset(coordSection, v, &off);
1299:       if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue;
1300:       x    = PetscRealPart(coords[off]);
1301:       y    = PetscRealPart(coords[off+1]);
1302:       phi  = PetscAtan2Real(y, x);
1303:       sinp = PetscSinReal(phi);
1304:       cosp = PetscCosReal(phi);
1305:       if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) {
1306:         dc = PetscAbsReal(ds2/sinp);
1307:         df = PetscAbsReal(dis/sinp);
1308:         xc = ds2*x/PetscAbsScalar(y);
1309:         yc = ds2*PetscSignReal(y);
1310:       } else {
1311:         dc = PetscAbsReal(ds2/cosp);
1312:         df = PetscAbsReal(dis/cosp);
1313:         xc = ds2*PetscSignReal(x);
1314:         yc = ds2*y/PetscAbsScalar(x);
1315:       }
1316:       coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc);
1317:       coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc);
1318:     }
1319:     VecRestoreArray(coordinates, &coords);
1320:     if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
1321:       DMLocalizeCoordinates(*dm);
1322:     }
1323:   }
1324:   return(0);
1325: }

1327: /* External function declarations here */
1328: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
1329: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
1330: extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
1331: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
1332: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
1333: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
1334: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
1335: extern PetscErrorCode DMSetUp_Plex(DM dm);
1336: extern PetscErrorCode DMDestroy_Plex(DM dm);
1337: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
1338: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
1339: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);

1341: /* Replace dm with the contents of dmNew
1342:    - Share the DM_Plex structure
1343:    - Share the coordinates
1344:    - Share the SF
1345: */
1346: static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
1347: {
1348:   PetscSF          sf;
1349:   DM               coordDM, coarseDM;
1350:   Vec              coords;
1351:   const PetscReal *maxCell, *L;
1352:   const DMBoundaryType *bd;
1353:   PetscErrorCode   ierr;

1356:   DMGetPointSF(dmNew, &sf);
1357:   DMSetPointSF(dm, sf);
1358:   DMGetCoordinateDM(dmNew, &coordDM);
1359:   DMGetCoordinatesLocal(dmNew, &coords);
1360:   DMSetCoordinateDM(dm, coordDM);
1361:   DMSetCoordinatesLocal(dm, coords);
1362:   DMGetPeriodicity(dm, &maxCell, &L, &bd);
1363:   if (L) {DMSetPeriodicity(dmNew, maxCell, L, bd);}
1364:   DMDestroy_Plex(dm);
1365:   dm->data = dmNew->data;
1366:   ((DM_Plex *) dmNew->data)->refct++;
1367:   dmNew->labels->refct++;
1368:   if (!--(dm->labels->refct)) {
1369:     DMLabelLink next = dm->labels->next;

1371:     /* destroy the labels */
1372:     while (next) {
1373:       DMLabelLink tmp = next->next;

1375:       DMLabelDestroy(&next->label);
1376:       PetscFree(next);
1377:       next = tmp;
1378:     }
1379:     PetscFree(dm->labels);
1380:   }
1381:   dm->labels = dmNew->labels;
1382:   dm->depthLabel = dmNew->depthLabel;
1383:   DMGetCoarseDM(dmNew,&coarseDM);
1384:   DMSetCoarseDM(dm,coarseDM);
1385:   return(0);
1386: }

1388: /* Swap dm with the contents of dmNew
1389:    - Swap the DM_Plex structure
1390:    - Swap the coordinates
1391:    - Swap the point PetscSF
1392: */
1393: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
1394: {
1395:   DM              coordDMA, coordDMB;
1396:   Vec             coordsA,  coordsB;
1397:   PetscSF         sfA,      sfB;
1398:   void            *tmp;
1399:   DMLabelLinkList listTmp;
1400:   DMLabel         depthTmp;
1401:   PetscErrorCode  ierr;

1404:   DMGetPointSF(dmA, &sfA);
1405:   DMGetPointSF(dmB, &sfB);
1406:   PetscObjectReference((PetscObject) sfA);
1407:   DMSetPointSF(dmA, sfB);
1408:   DMSetPointSF(dmB, sfA);
1409:   PetscObjectDereference((PetscObject) sfA);

1411:   DMGetCoordinateDM(dmA, &coordDMA);
1412:   DMGetCoordinateDM(dmB, &coordDMB);
1413:   PetscObjectReference((PetscObject) coordDMA);
1414:   DMSetCoordinateDM(dmA, coordDMB);
1415:   DMSetCoordinateDM(dmB, coordDMA);
1416:   PetscObjectDereference((PetscObject) coordDMA);

1418:   DMGetCoordinatesLocal(dmA, &coordsA);
1419:   DMGetCoordinatesLocal(dmB, &coordsB);
1420:   PetscObjectReference((PetscObject) coordsA);
1421:   DMSetCoordinatesLocal(dmA, coordsB);
1422:   DMSetCoordinatesLocal(dmB, coordsA);
1423:   PetscObjectDereference((PetscObject) coordsA);

1425:   tmp       = dmA->data;
1426:   dmA->data = dmB->data;
1427:   dmB->data = tmp;
1428:   listTmp   = dmA->labels;
1429:   dmA->labels = dmB->labels;
1430:   dmB->labels = listTmp;
1431:   depthTmp  = dmA->depthLabel;
1432:   dmA->depthLabel = dmB->depthLabel;
1433:   dmB->depthLabel = depthTmp;
1434:   return(0);
1435: }

1437: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
1438: {
1439:   DM_Plex       *mesh = (DM_Plex*) dm->data;

1443:   /* Handle viewing */
1444:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);
1445:   PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);
1446:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);
1447:   /* Point Location */
1448:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);
1449:   /* Generation and remeshing */
1450:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);
1451:   /* Projection behavior */
1452:   PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);
1453:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);

1455:   PetscPartitionerSetFromOptions(mesh->partitioner);
1456:   return(0);
1457: }

1459: static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
1460: {
1461:   PetscInt       refine = 0, coarsen = 0, r;
1462:   PetscBool      isHierarchy;

1467:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
1468:   /* Handle DMPlex refinement */
1469:   PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);
1470:   PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);
1471:   if (refine) {DMPlexSetRefinementUniform(dm, PETSC_TRUE);}
1472:   if (refine && isHierarchy) {
1473:     DM *dms, coarseDM;

1475:     DMGetCoarseDM(dm, &coarseDM);
1476:     PetscObjectReference((PetscObject)coarseDM);
1477:     PetscMalloc1(refine,&dms);
1478:     DMRefineHierarchy(dm, refine, dms);
1479:     /* Total hack since we do not pass in a pointer */
1480:     DMPlexSwap_Static(dm, dms[refine-1]);
1481:     if (refine == 1) {
1482:       DMSetCoarseDM(dm, dms[0]);
1483:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
1484:     } else {
1485:       DMSetCoarseDM(dm, dms[refine-2]);
1486:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
1487:       DMSetCoarseDM(dms[0], dms[refine-1]);
1488:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
1489:     }
1490:     DMSetCoarseDM(dms[refine-1], coarseDM);
1491:     PetscObjectDereference((PetscObject)coarseDM);
1492:     /* Free DMs */
1493:     for (r = 0; r < refine; ++r) {
1494:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
1495:       DMDestroy(&dms[r]);
1496:     }
1497:     PetscFree(dms);
1498:   } else {
1499:     for (r = 0; r < refine; ++r) {
1500:       DM refinedMesh;

1502:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1503:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);
1504:       /* Total hack since we do not pass in a pointer */
1505:       DMPlexReplace_Static(dm, refinedMesh);
1506:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1507:       DMDestroy(&refinedMesh);
1508:     }
1509:   }
1510:   /* Handle DMPlex coarsening */
1511:   PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);
1512:   PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);
1513:   if (coarsen && isHierarchy) {
1514:     DM *dms;

1516:     PetscMalloc1(coarsen, &dms);
1517:     DMCoarsenHierarchy(dm, coarsen, dms);
1518:     /* Free DMs */
1519:     for (r = 0; r < coarsen; ++r) {
1520:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
1521:       DMDestroy(&dms[r]);
1522:     }
1523:     PetscFree(dms);
1524:   } else {
1525:     for (r = 0; r < coarsen; ++r) {
1526:       DM coarseMesh;

1528:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1529:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);
1530:       /* Total hack since we do not pass in a pointer */
1531:       DMPlexReplace_Static(dm, coarseMesh);
1532:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1533:       DMDestroy(&coarseMesh);
1534:     }
1535:   }
1536:   /* Handle */
1537:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1538:   PetscOptionsTail();
1539:   return(0);
1540: }

1542: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
1543: {

1547:   DMCreateGlobalVector_Section_Private(dm,vec);
1548:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
1549:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
1550:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
1551:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
1552:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
1553:   return(0);
1554: }

1556: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
1557: {

1561:   DMCreateLocalVector_Section_Private(dm,vec);
1562:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
1563:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
1564:   return(0);
1565: }

1567: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
1568: {
1569:   PetscInt       depth, d;

1573:   DMPlexGetDepth(dm, &depth);
1574:   if (depth == 1) {
1575:     DMGetDimension(dm, &d);
1576:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
1577:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
1578:     else               {*pStart = 0; *pEnd = 0;}
1579:   } else {
1580:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
1581:   }
1582:   return(0);
1583: }

1585: static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
1586: {
1587:   PetscSF        sf;

1591:   DMGetPointSF(dm, &sf);
1592:   PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);
1593:   return(0);
1594: }

1596: static PetscErrorCode DMInitialize_Plex(DM dm)
1597: {

1601:   dm->ops->view                            = DMView_Plex;
1602:   dm->ops->load                            = DMLoad_Plex;
1603:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
1604:   dm->ops->clone                           = DMClone_Plex;
1605:   dm->ops->setup                           = DMSetUp_Plex;
1606:   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
1607:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
1608:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
1609:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
1610:   dm->ops->getlocaltoglobalmapping         = NULL;
1611:   dm->ops->createfieldis                   = NULL;
1612:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
1613:   dm->ops->getcoloring                     = NULL;
1614:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
1615:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
1616:   dm->ops->getaggregates                   = NULL;
1617:   dm->ops->getinjection                    = DMCreateInjection_Plex;
1618:   dm->ops->refine                          = DMRefine_Plex;
1619:   dm->ops->coarsen                         = DMCoarsen_Plex;
1620:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
1621:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
1622:   dm->ops->globaltolocalbegin              = NULL;
1623:   dm->ops->globaltolocalend                = NULL;
1624:   dm->ops->localtoglobalbegin              = NULL;
1625:   dm->ops->localtoglobalend                = NULL;
1626:   dm->ops->destroy                         = DMDestroy_Plex;
1627:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
1628:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
1629:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
1630:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
1631:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
1632:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
1633:   dm->ops->projectfieldlabellocal          = DMProjectFieldLabelLocal_Plex;
1634:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
1635:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
1636:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
1637:   dm->ops->getneighbors                    = DMGetNeighors_Plex;
1638:   PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Plex);
1639:   return(0);
1640: }

1642: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
1643: {
1644:   DM_Plex        *mesh = (DM_Plex *) dm->data;

1648:   mesh->refct++;
1649:   (*newdm)->data = mesh;
1650:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
1651:   DMInitialize_Plex(*newdm);
1652:   return(0);
1653: }

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

1661:   Level: intermediate

1663: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
1664: M*/

1666: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1667: {
1668:   DM_Plex        *mesh;
1669:   PetscInt       unit, d;

1674:   PetscNewLog(dm,&mesh);
1675:   dm->dim  = 0;
1676:   dm->data = mesh;

1678:   mesh->refct             = 1;
1679:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
1680:   mesh->maxConeSize       = 0;
1681:   mesh->cones             = NULL;
1682:   mesh->coneOrientations  = NULL;
1683:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
1684:   mesh->maxSupportSize    = 0;
1685:   mesh->supports          = NULL;
1686:   mesh->refinementUniform = PETSC_TRUE;
1687:   mesh->refinementLimit   = -1.0;

1689:   mesh->facesTmp = NULL;

1691:   mesh->tetgenOpts   = NULL;
1692:   mesh->triangleOpts = NULL;
1693:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
1694:   mesh->remeshBd     = PETSC_FALSE;

1696:   mesh->subpointMap = NULL;

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

1700:   mesh->regularRefinement   = PETSC_FALSE;
1701:   mesh->depthState          = -1;
1702:   mesh->globalVertexNumbers = NULL;
1703:   mesh->globalCellNumbers   = NULL;
1704:   mesh->anchorSection       = NULL;
1705:   mesh->anchorIS            = NULL;
1706:   mesh->createanchors       = NULL;
1707:   mesh->computeanchormatrix = NULL;
1708:   mesh->parentSection       = NULL;
1709:   mesh->parents             = NULL;
1710:   mesh->childIDs            = NULL;
1711:   mesh->childSection        = NULL;
1712:   mesh->children            = NULL;
1713:   mesh->referenceTree       = NULL;
1714:   mesh->getchildsymmetry    = NULL;
1715:   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1716:   mesh->vtkCellHeight       = 0;
1717:   mesh->useCone             = PETSC_FALSE;
1718:   mesh->useClosure          = PETSC_TRUE;
1719:   mesh->useAnchors          = PETSC_FALSE;

1721:   mesh->maxProjectionHeight = 0;

1723:   mesh->printSetValues = PETSC_FALSE;
1724:   mesh->printFEM       = 0;
1725:   mesh->printTol       = 1.0e-10;

1727:   DMInitialize_Plex(dm);
1728:   return(0);
1729: }

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

1734:   Collective on MPI_Comm

1736:   Input Parameter:
1737: . comm - The communicator for the DMPlex object

1739:   Output Parameter:
1740: . mesh  - The DMPlex object

1742:   Level: beginner

1744: .keywords: DMPlex, create
1745: @*/
1746: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1747: {

1752:   DMCreate(comm, mesh);
1753:   DMSetType(*mesh, DMPLEX);
1754:   return(0);
1755: }

1757: /*
1758:   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
1759: */
1760: static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert)
1761: {
1762:   PetscSF         sfPoint;
1763:   PetscLayout     vLayout;
1764:   PetscHashI      vhash;
1765:   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
1766:   const PetscInt *vrange;
1767:   PetscInt        numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
1768:   PETSC_UNUSED PetscHashIIter ret, iter;
1769:   PetscMPIInt     rank, size;
1770:   PetscErrorCode  ierr;

1773:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1774:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);
1775:   /* Partition vertices */
1776:   PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);
1777:   PetscLayoutSetLocalSize(vLayout, numVertices);
1778:   PetscLayoutSetBlockSize(vLayout, 1);
1779:   PetscLayoutSetUp(vLayout);
1780:   PetscLayoutGetRanges(vLayout, &vrange);
1781:   /* Count vertices and map them to procs */
1782:   PetscHashICreate(vhash);
1783:   for (c = 0; c < numCells; ++c) {
1784:     for (p = 0; p < numCorners; ++p) {
1785:       PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter);
1786:     }
1787:   }
1788:   PetscHashISize(vhash, numVerticesAdj);
1789:   PetscMalloc1(numVerticesAdj, &verticesAdj);
1790:   off = 0; PetscHashIGetKeys(vhash, &off, verticesAdj);
1791:   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
1792:   PetscSortInt(numVerticesAdj, verticesAdj);
1793:   PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);
1794:   for (v = 0; v < numVerticesAdj; ++v) {
1795:     const PetscInt gv = verticesAdj[v];
1796:     PetscInt       vrank;

1798:     PetscFindInt(gv, size+1, vrange, &vrank);
1799:     vrank = vrank < 0 ? -(vrank+2) : vrank;
1800:     remoteVerticesAdj[v].index = gv - vrange[vrank];
1801:     remoteVerticesAdj[v].rank  = vrank;
1802:   }
1803:   /* Create cones */
1804:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
1805:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
1806:   DMSetUp(dm);
1807:   DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);
1808:   for (c = 0; c < numCells; ++c) {
1809:     for (p = 0; p < numCorners; ++p) {
1810:       const PetscInt gv = cells[c*numCorners+p];
1811:       PetscInt       lv;

1813:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
1814:       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
1815:       cone[p] = lv+numCells;
1816:     }
1817:     DMPlexSetCone(dm, c, cone);
1818:   }
1819:   DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);
1820:   /* Create SF for vertices */
1821:   PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);
1822:   PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");
1823:   PetscSFSetFromOptions(*sfVert);
1824:   PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);
1825:   PetscFree(verticesAdj);
1826:   /* Build pointSF */
1827:   PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);
1828:   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
1829:   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
1830:   PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
1831:   PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
1832:   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);
1833:   PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);
1834:   PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);
1835:   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
1836:   PetscMalloc1(numVerticesGhost, &localVertex);
1837:   PetscMalloc1(numVerticesGhost, &remoteVertex);
1838:   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
1839:     if (vertexLocal[v].rank != rank) {
1840:       localVertex[g]        = v+numCells;
1841:       remoteVertex[g].index = vertexLocal[v].index;
1842:       remoteVertex[g].rank  = vertexLocal[v].rank;
1843:       ++g;
1844:     }
1845:   }
1846:   PetscFree2(vertexLocal, vertexOwner);
1847:   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
1848:   DMGetPointSF(dm, &sfPoint);
1849:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
1850:   PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);
1851:   PetscLayoutDestroy(&vLayout);
1852:   PetscHashIDestroy(vhash);
1853:   /* Fill in the rest of the topology structure */
1854:   DMPlexSymmetrize(dm);
1855:   DMPlexStratify(dm);
1856:   return(0);
1857: }

1859: /*
1860:   This takes as input the coordinates for each owned vertex
1861: */
1862: static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[])
1863: {
1864:   PetscSection   coordSection;
1865:   Vec            coordinates;
1866:   PetscScalar   *coords;
1867:   PetscInt       numVertices, numVerticesAdj, coordSize, v;

1871:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
1872:   DMGetCoordinateSection(dm, &coordSection);
1873:   PetscSectionSetNumFields(coordSection, 1);
1874:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
1875:   PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);
1876:   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
1877:     PetscSectionSetDof(coordSection, v, spaceDim);
1878:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
1879:   }
1880:   PetscSectionSetUp(coordSection);
1881:   PetscSectionGetStorageSize(coordSection, &coordSize);
1882:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
1883:   VecSetBlockSize(coordinates, spaceDim);
1884:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1885:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1886:   VecSetType(coordinates,VECSTANDARD);
1887:   VecGetArray(coordinates, &coords);
1888:   {
1889:     MPI_Datatype coordtype;

1891:     /* Need a temp buffer for coords if we have complex/single */
1892:     MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);
1893:     MPI_Type_commit(&coordtype);
1894: #if defined(PETSC_USE_COMPLEX)
1895:     {
1896:     PetscScalar *svertexCoords;
1897:     PetscInt    i;
1898:     PetscMalloc1(numV*spaceDim,&svertexCoords);
1899:     for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i];
1900:     PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);
1901:     PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);
1902:     PetscFree(svertexCoords);
1903:     }
1904: #else
1905:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);
1906:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);
1907: #endif
1908:     MPI_Type_free(&coordtype);
1909:   }
1910:   VecRestoreArray(coordinates, &coords);
1911:   DMSetCoordinatesLocal(dm, coordinates);
1912:   VecDestroy(&coordinates);
1913:   return(0);
1914: }

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

1919:   Input Parameters:
1920: + comm - The communicator
1921: . dim - The topological dimension of the mesh
1922: . numCells - The number of cells owned by this process
1923: . numVertices - The number of vertices owned by this process
1924: . numCorners - The number of vertices for each cell
1925: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1926: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
1927: . spaceDim - The spatial dimension used for coordinates
1928: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

1930:   Output Parameter:
1931: + dm - The DM
1932: - vertexSF - Optional, SF describing complete vertex ownership

1934:   Note: Two triangles sharing a face
1935: $
1936: $        2
1937: $      / | \
1938: $     /  |  \
1939: $    /   |   \
1940: $   0  0 | 1  3
1941: $    \   |   /
1942: $     \  |  /
1943: $      \ | /
1944: $        1
1945: would have input
1946: $  numCells = 2, numVertices = 4
1947: $  cells = [0 1 2  1 3 2]
1948: $
1949: which would result in the DMPlex
1950: $
1951: $        4
1952: $      / | \
1953: $     /  |  \
1954: $    /   |   \
1955: $   2  0 | 1  5
1956: $    \   |   /
1957: $     \  |  /
1958: $      \ | /
1959: $        3

1961:   Level: beginner

1963: .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
1964: @*/
1965: 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)
1966: {
1967:   PetscSF        sfVert;

1971:   DMCreate(comm, dm);
1972:   DMSetType(*dm, DMPLEX);
1975:   DMSetDimension(*dm, dim);
1976:   DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);
1977:   if (interpolate) {
1978:     DM idm = NULL;

1980:     DMPlexInterpolate(*dm, &idm);
1981:     DMDestroy(dm);
1982:     *dm  = idm;
1983:   }
1984:   DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);
1985:   if (vertexSF) *vertexSF = sfVert;
1986:   else {PetscSFDestroy(&sfVert);}
1987:   return(0);
1988: }

1990: /*
1991:   This takes as input the common mesh generator output, a list of the vertices for each cell
1992: */
1993: static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
1994: {
1995:   PetscInt      *cone, c, p;

1999:   DMPlexSetChart(dm, 0, numCells+numVertices);
2000:   for (c = 0; c < numCells; ++c) {
2001:     DMPlexSetConeSize(dm, c, numCorners);
2002:   }
2003:   DMSetUp(dm);
2004:   DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);
2005:   for (c = 0; c < numCells; ++c) {
2006:     for (p = 0; p < numCorners; ++p) {
2007:       cone[p] = cells[c*numCorners+p]+numCells;
2008:     }
2009:     DMPlexSetCone(dm, c, cone);
2010:   }
2011:   DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);
2012:   DMPlexSymmetrize(dm);
2013:   DMPlexStratify(dm);
2014:   return(0);
2015: }

2017: /*
2018:   This takes as input the coordinates for each vertex
2019: */
2020: static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
2021: {
2022:   PetscSection   coordSection;
2023:   Vec            coordinates;
2024:   DM             cdm;
2025:   PetscScalar   *coords;
2026:   PetscInt       v, d;

2030:   DMGetCoordinateSection(dm, &coordSection);
2031:   PetscSectionSetNumFields(coordSection, 1);
2032:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
2033:   PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
2034:   for (v = numCells; v < numCells+numVertices; ++v) {
2035:     PetscSectionSetDof(coordSection, v, spaceDim);
2036:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
2037:   }
2038:   PetscSectionSetUp(coordSection);

2040:   DMGetCoordinateDM(dm, &cdm);
2041:   DMCreateLocalVector(cdm, &coordinates);
2042:   VecSetBlockSize(coordinates, spaceDim);
2043:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2044:   VecGetArray(coordinates, &coords);
2045:   for (v = 0; v < numVertices; ++v) {
2046:     for (d = 0; d < spaceDim; ++d) {
2047:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
2048:     }
2049:   }
2050:   VecRestoreArray(coordinates, &coords);
2051:   DMSetCoordinatesLocal(dm, coordinates);
2052:   VecDestroy(&coordinates);
2053:   return(0);
2054: }

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

2059:   Input Parameters:
2060: + comm - The communicator
2061: . dim - The topological dimension of the mesh
2062: . numCells - The number of cells
2063: . numVertices - The number of vertices
2064: . numCorners - The number of vertices for each cell
2065: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
2066: . cells - An array of numCells*numCorners numbers, the vertices for each cell
2067: . spaceDim - The spatial dimension used for coordinates
2068: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

2070:   Output Parameter:
2071: . dm - The DM

2073:   Note: Two triangles sharing a face
2074: $
2075: $        2
2076: $      / | \
2077: $     /  |  \
2078: $    /   |   \
2079: $   0  0 | 1  3
2080: $    \   |   /
2081: $     \  |  /
2082: $      \ | /
2083: $        1
2084: would have input
2085: $  numCells = 2, numVertices = 4
2086: $  cells = [0 1 2  1 3 2]
2087: $
2088: which would result in the DMPlex
2089: $
2090: $        4
2091: $      / | \
2092: $     /  |  \
2093: $    /   |   \
2094: $   2  0 | 1  5
2095: $    \   |   /
2096: $     \  |  /
2097: $      \ | /
2098: $        3

2100:   Level: beginner

2102: .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
2103: @*/
2104: 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)
2105: {

2109:   DMCreate(comm, dm);
2110:   DMSetType(*dm, DMPLEX);
2111:   DMSetDimension(*dm, dim);
2112:   DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);
2113:   if (interpolate) {
2114:     DM idm = NULL;

2116:     DMPlexInterpolate(*dm, &idm);
2117:     DMDestroy(dm);
2118:     *dm  = idm;
2119:   }
2120:   DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);
2121:   return(0);
2122: }

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

2127:   Input Parameters:
2128: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
2129: . depth - The depth of the DAG
2130: . numPoints - The number of points at each depth
2131: . coneSize - The cone size of each point
2132: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
2133: . coneOrientations - The orientation of each cone point
2134: - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex

2136:   Output Parameter:
2137: . dm - The DM

2139:   Note: Two triangles sharing a face would have input
2140: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
2141: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
2142: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
2143: $
2144: which would result in the DMPlex
2145: $
2146: $        4
2147: $      / | \
2148: $     /  |  \
2149: $    /   |   \
2150: $   2  0 | 1  5
2151: $    \   |   /
2152: $     \  |  /
2153: $      \ | /
2154: $        3
2155: $
2156: $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()

2158:   Level: advanced

2160: .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
2161: @*/
2162: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
2163: {
2164:   Vec            coordinates;
2165:   PetscSection   coordSection;
2166:   PetscScalar    *coords;
2167:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

2171:   DMGetDimension(dm, &dim);
2172:   DMGetCoordinateDim(dm, &dimEmbed);
2173:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
2174:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
2175:   DMPlexSetChart(dm, pStart, pEnd);
2176:   for (p = pStart; p < pEnd; ++p) {
2177:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
2178:     if (firstVertex < 0 && !coneSize[p - pStart]) {
2179:       firstVertex = p - pStart;
2180:     }
2181:   }
2182:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
2183:   DMSetUp(dm); /* Allocate space for cones */
2184:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
2185:     DMPlexSetCone(dm, p, &cones[off]);
2186:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
2187:   }
2188:   DMPlexSymmetrize(dm);
2189:   DMPlexStratify(dm);
2190:   /* Build coordinates */
2191:   DMGetCoordinateSection(dm, &coordSection);
2192:   PetscSectionSetNumFields(coordSection, 1);
2193:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
2194:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
2195:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
2196:     PetscSectionSetDof(coordSection, v, dimEmbed);
2197:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
2198:   }
2199:   PetscSectionSetUp(coordSection);
2200:   PetscSectionGetStorageSize(coordSection, &coordSize);
2201:   VecCreate(PETSC_COMM_SELF, &coordinates);
2202:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
2203:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
2204:   VecSetBlockSize(coordinates, dimEmbed);
2205:   VecSetType(coordinates,VECSTANDARD);
2206:   VecGetArray(coordinates, &coords);
2207:   for (v = 0; v < numPoints[0]; ++v) {
2208:     PetscInt off;

2210:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
2211:     for (d = 0; d < dimEmbed; ++d) {
2212:       coords[off+d] = vertexCoords[v*dimEmbed+d];
2213:     }
2214:   }
2215:   VecRestoreArray(coordinates, &coords);
2216:   DMSetCoordinatesLocal(dm, coordinates);
2217:   VecDestroy(&coordinates);
2218:   return(0);
2219: }

2221: /*@C
2222:   DMPlexCreateFromFile - This takes a filename and produces a DM

2224:   Input Parameters:
2225: + comm - The communicator
2226: . filename - A file name
2227: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

2229:   Output Parameter:
2230: . dm - The DM

2232:   Level: beginner

2234: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
2235: @*/
2236: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
2237: {
2238:   const char    *extGmsh   = ".msh";
2239:   const char    *extCGNS   = ".cgns";
2240:   const char    *extExodus = ".exo";
2241:   const char    *extFluent = ".cas";
2242:   const char    *extHDF5   = ".h5";
2243:   const char    *extMed    = ".med";
2244:   const char    *extPLY    = ".ply";
2245:   size_t         len;
2246:   PetscBool      isGmsh, isCGNS, isExodus, isFluent, isHDF5, isMed, isPLY;
2247:   PetscMPIInt    rank;

2253:   MPI_Comm_rank(comm, &rank);
2254:   PetscStrlen(filename, &len);
2255:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
2256:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,   4, &isGmsh);
2257:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,   5, &isCGNS);
2258:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);
2259:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);
2260:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,   3, &isHDF5);
2261:   PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,    4, &isMed);
2262:   PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY,    4, &isPLY);
2263:   if (isGmsh) {
2264:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
2265:   } else if (isCGNS) {
2266:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
2267:   } else if (isExodus) {
2268:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
2269:   } else if (isFluent) {
2270:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
2271:   } else if (isHDF5) {
2272:     PetscViewer viewer;

2274:     PetscViewerCreate(comm, &viewer);
2275:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
2276:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
2277:     PetscViewerFileSetName(viewer, filename);
2278:     DMCreate(comm, dm);
2279:     DMSetType(*dm, DMPLEX);
2280:     DMLoad(*dm, viewer);
2281:     PetscViewerDestroy(&viewer);
2282:   } else if (isMed) {
2283:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
2284:   } else if (isPLY) {
2285:     DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);
2286:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
2287:   return(0);
2288: }

2290: /*@
2291:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

2293:   Collective on comm

2295:   Input Parameters:
2296: + comm    - The communicator
2297: . dim     - The spatial dimension
2298: - simplex - Flag for simplex, otherwise use a tensor-product cell

2300:   Output Parameter:
2301: . refdm - The reference cell

2303:   Level: intermediate

2305: .keywords: reference cell
2306: .seealso:
2307: @*/
2308: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
2309: {
2310:   DM             rdm;
2311:   Vec            coords;

2315:   DMCreate(comm, &rdm);
2316:   DMSetType(rdm, DMPLEX);
2317:   DMSetDimension(rdm, dim);
2318:   switch (dim) {
2319:   case 0:
2320:   {
2321:     PetscInt    numPoints[1]        = {1};
2322:     PetscInt    coneSize[1]         = {0};
2323:     PetscInt    cones[1]            = {0};
2324:     PetscInt    coneOrientations[1] = {0};
2325:     PetscScalar vertexCoords[1]     = {0.0};

2327:     DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2328:   }
2329:   break;
2330:   case 1:
2331:   {
2332:     PetscInt    numPoints[2]        = {2, 1};
2333:     PetscInt    coneSize[3]         = {2, 0, 0};
2334:     PetscInt    cones[2]            = {1, 2};
2335:     PetscInt    coneOrientations[2] = {0, 0};
2336:     PetscScalar vertexCoords[2]     = {-1.0,  1.0};

2338:     DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2339:   }
2340:   break;
2341:   case 2:
2342:     if (simplex) {
2343:       PetscInt    numPoints[2]        = {3, 1};
2344:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
2345:       PetscInt    cones[3]            = {1, 2, 3};
2346:       PetscInt    coneOrientations[3] = {0, 0, 0};
2347:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

2349:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2350:     } else {
2351:       PetscInt    numPoints[2]        = {4, 1};
2352:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2353:       PetscInt    cones[4]            = {1, 2, 3, 4};
2354:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2355:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

2357:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2358:     }
2359:   break;
2360:   case 3:
2361:     if (simplex) {
2362:       PetscInt    numPoints[2]        = {4, 1};
2363:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2364:       PetscInt    cones[4]            = {1, 3, 2, 4};
2365:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2366:       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};

2368:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2369:     } else {
2370:       PetscInt    numPoints[2]        = {8, 1};
2371:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
2372:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
2373:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
2374:       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,
2375:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

2377:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2378:     }
2379:   break;
2380:   default:
2381:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
2382:   }
2383:   *refdm = NULL;
2384:   DMPlexInterpolate(rdm, refdm);
2385:   if (rdm->coordinateDM) {
2386:     DM           ncdm;
2387:     PetscSection cs;
2388:     PetscInt     pEnd = -1;

2390:     DMGetDefaultSection(rdm->coordinateDM, &cs);
2391:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
2392:     if (pEnd >= 0) {
2393:       DMClone(rdm->coordinateDM, &ncdm);
2394:       DMSetDefaultSection(ncdm, cs);
2395:       DMSetCoordinateDM(*refdm, ncdm);
2396:       DMDestroy(&ncdm);
2397:     }
2398:   }
2399:   DMGetCoordinatesLocal(rdm, &coords);
2400:   if (coords) {
2401:     DMSetCoordinatesLocal(*refdm, coords);
2402:   } else {
2403:     DMGetCoordinates(rdm, &coords);
2404:     if (coords) {DMSetCoordinates(*refdm, coords);}
2405:   }
2406:   DMDestroy(&rdm);
2407:   return(0);
2408: }