Actual source code: plexcreate.c

petsc-master 2016-08-23
Report Typos and Errors
  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
  3: #include <petscdmda.h>
  4: #include <petscsf.h>

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

 11:   Collective on MPI_Comm

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

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

 24:   Level: beginner

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

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

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

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

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

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

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

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

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

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

158:   Collective on MPI_Comm

160:   Input Parameters:
161: + comm  - The communicator for the DM object
162: . lower - The lower left corner coordinates
163: . upper - The upper right corner coordinates
164: - edges - The number of cells in each direction

166:   Output Parameter:
167: . dm  - The DM object

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

180:   Level: beginner

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

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

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

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

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

299: /*@
300:   DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice.

302:   Collective on MPI_Comm

304:   Input Parameters:
305: + comm  - The communicator for the DM object
306: . lower - The lower left front corner coordinates
307: . upper - The upper right back corner coordinates
308: - edges - The number of cells in each direction

310:   Output Parameter:
311: . dm  - The DM object

313:   Level: beginner

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

332:   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");
333:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
334:   vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1;
335:   numVertices = vertices[0]*vertices[1]*vertices[2];
336:   if (!rank) {
337:     PetscInt f;

339:     DMPlexSetChart(dm, 0, numFaces+numVertices);
340:     for (f = 0; f < numFaces; ++f) {
341:       DMPlexSetConeSize(dm, f, 4);
342:     }
343:     DMSetUp(dm); /* Allocate space for cones */
344:     for (v = 0; v < numFaces+numVertices; ++v) {
345:       DMSetLabelValue(dm, "marker", v, 1);
346:     }

348:     /* Side 0 (Top) */
349:     for (vy = 0; vy < faces[1]; vy++) {
350:       for (vx = 0; vx < faces[0]; vx++) {
351:         voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx;
352:         cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0];
353:         DMPlexSetCone(dm, iface, cone);
354:         iface++;
355:       }
356:     }

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

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

378:     /* Side 3 (Back) */
379:     for (vz = 0; vz < faces[2]; vz++) {
380:       for (vx = 0; vx < faces[0]; vx++) {
381:         voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx;
382:         cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1;
383:         cone[2] = voffset+1; cone[3] = voffset;
384:         DMPlexSetCone(dm, iface, cone);
385:         iface++;
386:       }
387:     }

389:     /* Side 4 (Left) */
390:     for (vz = 0; vz < faces[2]; vz++) {
391:       for (vy = 0; vy < faces[1]; vy++) {
392:         voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0];
393:         cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1];
394:         cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0];
395:         DMPlexSetCone(dm, iface, cone);
396:         iface++;
397:       }
398:     }

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

444: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
445: {
446:   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
447:   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
448:   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
449:   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
450:   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
451:   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
452:   PetscInt       dim;
453:   PetscBool      markerSeparate = PETSC_FALSE;
454:   PetscMPIInt    rank;

458:   DMGetDimension(dm,&dim);
459:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
460:   DMCreateLabel(dm,"marker");
461:   DMCreateLabel(dm,"Face Sets");
462:   switch (dim) {
463:   case 2:
464:     faceMarkerTop    = 3;
465:     faceMarkerBottom = 1;
466:     faceMarkerRight  = 2;
467:     faceMarkerLeft   = 4;
468:     break;
469:   case 3:
470:     faceMarkerBottom = 1;
471:     faceMarkerTop    = 2;
472:     faceMarkerFront  = 3;
473:     faceMarkerBack   = 4;
474:     faceMarkerRight  = 5;
475:     faceMarkerLeft   = 6;
476:     break;
477:   default:
478:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim);
479:     break;
480:   }
481:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);
482:   if (markerSeparate) {
483:     markerBottom = faceMarkerBottom;
484:     markerTop    = faceMarkerTop;
485:     markerFront  = faceMarkerFront;
486:     markerBack   = faceMarkerBack;
487:     markerRight  = faceMarkerRight;
488:     markerLeft   = faceMarkerLeft;
489:   }
490:   {
491:     const PetscInt numXEdges    = !rank ? edges[0]   : 0;
492:     const PetscInt numYEdges    = !rank ? edges[1]   : 0;
493:     const PetscInt numZEdges    = !rank ? edges[2]   : 0;
494:     const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0;
495:     const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0;
496:     const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0;
497:     const PetscInt numCells     = numXEdges*numYEdges*numZEdges;
498:     const PetscInt numXFaces    = numYEdges*numZEdges;
499:     const PetscInt numYFaces    = numXEdges*numZEdges;
500:     const PetscInt numZFaces    = numXEdges*numYEdges;
501:     const PetscInt numTotXFaces = numXVertices*numXFaces;
502:     const PetscInt numTotYFaces = numYVertices*numYFaces;
503:     const PetscInt numTotZFaces = numZVertices*numZFaces;
504:     const PetscInt numFaces     = numTotXFaces + numTotYFaces + numTotZFaces;
505:     const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices;
506:     const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices;
507:     const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices;
508:     const PetscInt numVertices  = numXVertices*numYVertices*numZVertices;
509:     const PetscInt numEdges     = numTotXEdges + numTotYEdges + numTotZEdges;
510:     const PetscInt firstVertex  = (dim == 2) ? numFaces : numCells;
511:     const PetscInt firstXFace   = (dim == 2) ? 0 : numCells + numVertices;
512:     const PetscInt firstYFace   = firstXFace + numTotXFaces;
513:     const PetscInt firstZFace   = firstYFace + numTotYFaces;
514:     const PetscInt firstXEdge   = numCells + numFaces + numVertices;
515:     const PetscInt firstYEdge   = firstXEdge + numTotXEdges;
516:     const PetscInt firstZEdge   = firstYEdge + numTotYEdges;
517:     Vec            coordinates;
518:     PetscSection   coordSection;
519:     PetscScalar   *coords;
520:     PetscInt       coordSize;
521:     PetscInt       v, vx, vy, vz;
522:     PetscInt       c, f, fx, fy, fz, e, ex, ey, ez;

524:     DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);
525:     for (c = 0; c < numCells; c++) {
526:       DMPlexSetConeSize(dm, c, 6);
527:     }
528:     for (f = firstXFace; f < firstXFace+numFaces; ++f) {
529:       DMPlexSetConeSize(dm, f, 4);
530:     }
531:     for (e = firstXEdge; e < firstXEdge+numEdges; ++e) {
532:       DMPlexSetConeSize(dm, e, 2);
533:     }
534:     DMSetUp(dm); /* Allocate space for cones */
535:     /* Build cells */
536:     for (fz = 0; fz < numZEdges; ++fz) {
537:       for (fy = 0; fy < numYEdges; ++fy) {
538:         for (fx = 0; fx < numXEdges; ++fx) {
539:           PetscInt cell    = (fz*numYEdges + fy)*numXEdges + fx;
540:           PetscInt faceB   = firstZFace + (fy*numXEdges+fx)*numZVertices +   fz;
541:           PetscInt faceT   = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices);
542:           PetscInt faceF   = firstYFace + (fz*numXEdges+fx)*numYVertices +   fy;
543:           PetscInt faceK   = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices);
544:           PetscInt faceL   = firstXFace + (fz*numYEdges+fy)*numXVertices +   fx;
545:           PetscInt faceR   = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices);
546:                             /* B,  T,  F,  K,  R,  L */
547:           PetscInt ornt[8] = {-4,  0,  0, -1,  0, -4}; /* ??? */
548:           PetscInt cone[8];

550:           /* no boundary twisting in 3D */
551:           cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL;
552:           DMPlexSetCone(dm, cell, cone);
553:           DMPlexSetConeOrientation(dm, cell, ornt);
554:         }
555:       }
556:     }
557:     /* Build x faces */
558:     for (fz = 0; fz < numZEdges; ++fz) {
559:       for (fy = 0; fy < numYEdges; ++fy) {
560:         for (fx = 0; fx < numXVertices; ++fx) {
561:           PetscInt face    = firstXFace + (fz*numYEdges+fy)*numXVertices + fx;
562:           PetscInt edgeL   = firstZEdge + (  fy*                 numXVertices+fx)*numZEdges + fz;
563:           PetscInt edgeR   = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz;
564:           PetscInt edgeB   = firstYEdge + (  fz*                 numXVertices+fx)*numYEdges + fy;
565:           PetscInt edgeT   = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy;
566:           PetscInt ornt[4] = {0, 0, -2, -2};
567:           PetscInt cone[4];

569:           if (dim == 3) {
570:             /* markers */
571:             if (bdX != DM_BOUNDARY_PERIODIC) {
572:               if (fx == numXVertices-1) {
573:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);
574:                 DMSetLabelValue(dm, "marker", face, markerRight);
575:               }
576:               else if (fx == 0) {
577:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);
578:                 DMSetLabelValue(dm, "marker", face, markerLeft);
579:               }
580:             }
581:           }
582:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
583:           DMPlexSetCone(dm, face, cone);
584:           DMPlexSetConeOrientation(dm, face, ornt);
585:         }
586:       }
587:     }
588:     /* Build y faces */
589:     for (fz = 0; fz < numZEdges; ++fz) {
590:       for (fx = 0; fx < numYEdges; ++fx) {
591:         for (fy = 0; fy < numYVertices; ++fy) {
592:           PetscInt face    = firstYFace + (fz*numXEdges+fx)*numYVertices + fy;
593:           PetscInt edgeL   = firstZEdge + (fy*numXVertices+  fx                 )*numZEdges + fz;
594:           PetscInt edgeR   = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz;
595:           PetscInt edgeB   = firstXEdge + (  fz                 *numYVertices+fy)*numXEdges + fx;
596:           PetscInt edgeT   = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx;
597:           PetscInt ornt[4] = {0, 0, -2, -2};
598:           PetscInt cone[4];

600:           if (dim == 3) {
601:             /* markers */
602:             if (bdY != DM_BOUNDARY_PERIODIC) {
603:               if (fy == numYVertices-1) {
604:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);
605:                 DMSetLabelValue(dm, "marker", face, markerBack);
606:               }
607:               else if (fy == 0) {
608:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);
609:                 DMSetLabelValue(dm, "marker", face, markerFront);
610:               }
611:             }
612:           }
613:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
614:           DMPlexSetCone(dm, face, cone);
615:           DMPlexSetConeOrientation(dm, face, ornt);
616:         }
617:       }
618:     }
619:     /* Build z faces */
620:     for (fy = 0; fy < numYEdges; ++fy) {
621:       for (fx = 0; fx < numXEdges; ++fx) {
622:         for (fz = 0; fz < numZVertices; fz++) {
623:           PetscInt face    = firstZFace + (fy*numXEdges+fx)*numZVertices + fz;
624:           PetscInt edgeL   = firstYEdge + (fz*numXVertices+  fx                 )*numYEdges + fy;
625:           PetscInt edgeR   = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy;
626:           PetscInt edgeB   = firstXEdge + (fz*numYVertices+  fy                 )*numXEdges + fx;
627:           PetscInt edgeT   = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx;
628:           PetscInt ornt[4] = {0, 0, -2, -2};
629:           PetscInt cone[4];

631:           if (dim == 2) {
632:             if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;}
633:             if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] =  0;}
634:           }
635:           else {
636:             /* markers */
637:             if (bdZ != DM_BOUNDARY_PERIODIC) {
638:               if (fz == numZVertices-1) {
639:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);
640:                 DMSetLabelValue(dm, "marker", face, markerTop);
641:               }
642:               else if (fz == 0) {
643:                 DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);
644:                 DMSetLabelValue(dm, "marker", face, markerBottom);
645:               }
646:             }
647:           }
648:           cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL;
649:           DMPlexSetCone(dm, face, cone);
650:           DMPlexSetConeOrientation(dm, face, ornt);
651:         }
652:       }
653:     }
654:     /* Build Z edges*/
655:     for (vy = 0; vy < numYVertices; vy++) {
656:       for (vx = 0; vx < numXVertices; vx++) {
657:         for (ez = 0; ez < numZEdges; ez++) {
658:           const PetscInt edge    = firstZEdge  + (vy*numXVertices+vx)*numZEdges + ez;
659:           const PetscInt vertexB = firstVertex + (  ez                 *numYVertices+vy)*numXVertices + vx;
660:           const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx;
661:           PetscInt       cone[2];

663:           if (dim == 3) {
664:             if (bdX != DM_BOUNDARY_PERIODIC) {
665:               if (vx == numXVertices-1) {
666:                 DMSetLabelValue(dm, "marker", edge, markerRight);
667:               }
668:               else if (vx == 0) {
669:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
670:               }
671:             }
672:             if (bdY != DM_BOUNDARY_PERIODIC) {
673:               if (vy == numYVertices-1) {
674:                 DMSetLabelValue(dm, "marker", edge, markerBack);
675:               }
676:               else if (vy == 0) {
677:                 DMSetLabelValue(dm, "marker", edge, markerFront);
678:               }
679:             }
680:           }
681:           cone[0] = vertexB; cone[1] = vertexT;
682:           DMPlexSetCone(dm, edge, cone);
683:         }
684:       }
685:     }
686:     /* Build Y edges*/
687:     for (vz = 0; vz < numZVertices; vz++) {
688:       for (vx = 0; vx < numXVertices; vx++) {
689:         for (ey = 0; ey < numYEdges; ey++) {
690:           const PetscInt nextv   = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx;
691:           const PetscInt edge    = firstYEdge  + (vz*numXVertices+vx)*numYEdges + ey;
692:           const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx;
693:           const PetscInt vertexK = firstVertex + nextv;
694:           PetscInt       cone[2];

696:           cone[0] = vertexF; cone[1] = vertexK;
697:           DMPlexSetCone(dm, edge, cone);
698:           if (dim == 2) {
699:             if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
700:               if (vx == numXVertices-1) {
701:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);
702:                 DMSetLabelValue(dm, "marker", edge,    markerRight);
703:                 DMSetLabelValue(dm, "marker", cone[0], markerRight);
704:                 if (ey == numYEdges-1) {
705:                   DMSetLabelValue(dm, "marker", cone[1], markerRight);
706:                 }
707:               }
708:               else if (vx == 0) {
709:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);
710:                 DMSetLabelValue(dm, "marker", edge,    markerLeft);
711:                 DMSetLabelValue(dm, "marker", cone[0], markerLeft);
712:                 if (ey == numYEdges-1) {
713:                   DMSetLabelValue(dm, "marker", cone[1], markerLeft);
714:                 }
715:               }
716:             }
717:           }
718:           else {
719:             if (bdX != DM_BOUNDARY_PERIODIC) {
720:               if (vx == numXVertices-1) {
721:                 DMSetLabelValue(dm, "marker", edge, markerRight);
722:               }
723:               else if (vx == 0) {
724:                 DMSetLabelValue(dm, "marker", edge, markerLeft);
725:               }
726:             }
727:             if (bdZ != DM_BOUNDARY_PERIODIC) {
728:               if (vz == numZVertices-1) {
729:                 DMSetLabelValue(dm, "marker", edge, markerTop);
730:               }
731:               else if (vz == 0) {
732:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
733:               }
734:             }
735:           }
736:         }
737:       }
738:     }
739:     /* Build X edges*/
740:     for (vz = 0; vz < numZVertices; vz++) {
741:       for (vy = 0; vy < numYVertices; vy++) {
742:         for (ex = 0; ex < numXEdges; ex++) {
743:           const PetscInt nextv   = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices;
744:           const PetscInt edge    = firstXEdge  + (vz*numYVertices+vy)*numXEdges + ex;
745:           const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex;
746:           const PetscInt vertexR = firstVertex + nextv;
747:           PetscInt       cone[2];

749:           cone[0] = vertexL; cone[1] = vertexR;
750:           DMPlexSetCone(dm, edge, cone);
751:           if (dim == 2) {
752:             if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
753:               if (vy == numYVertices-1) {
754:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);
755:                 DMSetLabelValue(dm, "marker", edge,    markerTop);
756:                 DMSetLabelValue(dm, "marker", cone[0], markerTop);
757:                 if (ex == numXEdges-1) {
758:                   DMSetLabelValue(dm, "marker", cone[1], markerTop);
759:                 }
760:               }
761:               else if (vy == 0) {
762:                 DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);
763:                 DMSetLabelValue(dm, "marker", edge,    markerBottom);
764:                 DMSetLabelValue(dm, "marker", cone[0], markerBottom);
765:                 if (ex == numXEdges-1) {
766:                   DMSetLabelValue(dm, "marker", cone[1], markerBottom);
767:                 }
768:               }
769:             }
770:           }
771:           else {
772:             if (bdY != DM_BOUNDARY_PERIODIC) {
773:               if (vy == numYVertices-1) {
774:                 DMSetLabelValue(dm, "marker", edge, markerBack);
775:               }
776:               else if (vy == 0) {
777:                 DMSetLabelValue(dm, "marker", edge, markerFront);
778:               }
779:             }
780:             if (bdZ != DM_BOUNDARY_PERIODIC) {
781:               if (vz == numZVertices-1) {
782:                 DMSetLabelValue(dm, "marker", edge, markerTop);
783:               }
784:               else if (vz == 0) {
785:                 DMSetLabelValue(dm, "marker", edge, markerBottom);
786:               }
787:             }
788:           }
789:         }
790:       }
791:     }
792:     DMPlexSymmetrize(dm);
793:     DMPlexStratify(dm);
794:     /* Build coordinates */
795:     DMGetCoordinateSection(dm, &coordSection);
796:     PetscSectionSetNumFields(coordSection, 1);
797:     PetscSectionSetFieldComponents(coordSection, 0, dim);
798:     PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);
799:     for (v = firstVertex; v < firstVertex+numVertices; ++v) {
800:       PetscSectionSetDof(coordSection, v, dim);
801:       PetscSectionSetFieldDof(coordSection, v, 0, dim);
802:     }
803:     PetscSectionSetUp(coordSection);
804:     PetscSectionGetStorageSize(coordSection, &coordSize);
805:     VecCreate(PETSC_COMM_SELF, &coordinates);
806:     PetscObjectSetName((PetscObject) coordinates, "coordinates");
807:     VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
808:     VecSetBlockSize(coordinates, dim);
809:     VecSetType(coordinates,VECSTANDARD);
810:     VecGetArray(coordinates, &coords);
811:     for (vz = 0; vz < numZVertices; ++vz) {
812:       for (vy = 0; vy < numYVertices; ++vy) {
813:         for (vx = 0; vx < numXVertices; ++vx) {
814:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx;
815:           coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy;
816:           if (dim == 3) {
817:             coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz;
818:           }
819:         }
820:       }
821:     }
822:     VecRestoreArray(coordinates, &coords);
823:     DMSetCoordinatesLocal(dm, coordinates);
824:     VecDestroy(&coordinates);
825:   }
826:   return(0);
827: }

831: /*@
832:   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.

834:   Collective on MPI_Comm

836:   Input Parameters:
837: + comm  - The communicator for the DM object
838: . lower - The lower left corner coordinates
839: . upper - The upper right corner coordinates
840: . edges - The number of cells in each direction
841: . bdX   - The boundary type for the X direction
842: - bdY   - The boundary type for the Y direction

844:   Output Parameter:
845: . dm  - The DM object

847:   Note: Here is the numbering returned for 2 cells in each direction:
848: $ 22--8-23--9--24
849: $  |     |     |
850: $ 13  2 14  3  15
851: $  |     |     |
852: $ 19--6-20--7--21
853: $  |     |     |
854: $ 10  0 11  1 12
855: $  |     |     |
856: $ 16--4-17--5--18

858:   Level: beginner

860: .keywords: DM, create
861: .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
862: @*/
863: PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
864: {
865:   PetscReal      lower3[3], upper3[3];
866:   PetscInt       edges3[3];

870:   lower3[0] = lower[0]; lower3[1] = lower[1]; lower3[2] = 0.;
871:   upper3[0] = upper[0]; upper3[1] = upper[1]; upper3[2] = 0.;
872:   edges3[0] = edges[0]; edges3[1] = edges[1]; edges3[2] = 0;
873:   DMPlexCreateCubeMesh_Internal(dm, lower3, upper3, edges3, bdX, bdY, DM_BOUNDARY_NONE);
874:   return(0);
875: }

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

882:   Collective on MPI_Comm

884:   Input Parameters:
885: + comm - The communicator for the DM object
886: . dim - The spatial dimension
887: . numFaces - Number of faces per dimension
888: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

890:   Output Parameter:
891: . dm  - The DM object

893:   Level: beginner

895: .keywords: DM, create
896: .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
897: @*/
898: PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm)
899: {
900:   DM             boundary;

905:   DMCreate(comm, &boundary);
907:   DMSetType(boundary, DMPLEX);
908:   DMSetDimension(boundary, dim-1);
909:   DMSetCoordinateDim(boundary, dim);
910:   switch (dim) {
911:   case 2:
912:   {
913:     PetscReal lower[2] = {0.0, 0.0};
914:     PetscReal upper[2] = {1.0, 1.0};
915:     PetscInt  edges[2];

917:     edges[0] = numFaces; edges[1] = numFaces;
918:     DMPlexCreateSquareBoundary(boundary, lower, upper, edges);
919:     break;
920:   }
921:   case 3:
922:   {
923:     PetscReal lower[3] = {0.0, 0.0, 0.0};
924:     PetscReal upper[3] = {1.0, 1.0, 1.0};
925:     PetscInt  faces[3];

927:     faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces;
928:     DMPlexCreateCubeBoundary(boundary, lower, upper, faces);
929:     break;
930:   }
931:   default:
932:     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
933:   }
934:   DMPlexGenerate(boundary, NULL, interpolate, dm);
935:   DMDestroy(&boundary);
936:   return(0);
937: }

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

944:   Collective on MPI_Comm

946:   Input Parameters:
947: + comm  - The communicator for the DM object
948: . dim   - The spatial dimension
949: . periodicX - The boundary type for the X direction
950: . periodicY - The boundary type for the Y direction
951: . periodicZ - The boundary type for the Z direction
952: - cells - The number of cells in each direction

954:   Output Parameter:
955: . dm  - The DM object

957:   Level: beginner

959: .keywords: DM, create
960: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
961: @*/
962: PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
963: {

968:   DMCreate(comm, dm);
970:   DMSetType(*dm, DMPLEX);
971:   DMSetDimension(*dm, dim);
972:   switch (dim) {
973:   case 2:
974:   {
975:     PetscReal lower[2] = {0.0, 0.0};
976:     PetscReal upper[2] = {1.0, 1.0};

978:     DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);
979:     break;
980:   }
981:   case 3:
982:   {
983:     PetscReal lower[3] = {0.0, 0.0, 0.0};
984:     PetscReal upper[3] = {1.0, 1.0, 1.0};

986:     DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);
987:     break;
988:   }
989:   default:
990:     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
991:   }
992:   return(0);
993: }

995: /* External function declarations here */
996: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
997: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
998: extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
999: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
1000: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
1001: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
1002: extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined);
1003: extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened);
1004: extern PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmRefined);
1005: extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]);
1006: extern PetscErrorCode DMCoarsenHierarchy_Plex(DM dm, PetscInt nlevels, DM dmCoarsened[]);
1007: extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
1008: extern PetscErrorCode DMSetUp_Plex(DM dm);
1009: extern PetscErrorCode DMDestroy_Plex(DM dm);
1010: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
1011: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
1012: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);
1013: extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF);
1014: extern PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
1015: extern PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
1016: extern PetscErrorCode DMProjectFieldLocal_Plex(DM,Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscScalar[]),InsertMode,Vec);
1017: extern PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
1018: extern PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
1019: extern PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);

1023: /* Replace dm with the contents of dmNew
1024:    - Share the DM_Plex structure
1025:    - Share the coordinates
1026:    - Share the SF
1027: */
1028: static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew)
1029: {
1030:   PetscSF          sf;
1031:   DM               coordDM, coarseDM;
1032:   Vec              coords;
1033:   const PetscReal *maxCell, *L;
1034:   const DMBoundaryType *bd;
1035:   PetscErrorCode   ierr;

1038:   DMGetPointSF(dmNew, &sf);
1039:   DMSetPointSF(dm, sf);
1040:   DMGetCoordinateDM(dmNew, &coordDM);
1041:   DMGetCoordinatesLocal(dmNew, &coords);
1042:   DMSetCoordinateDM(dm, coordDM);
1043:   DMSetCoordinatesLocal(dm, coords);
1044:   DMGetPeriodicity(dm, &maxCell, &L, &bd);
1045:   if (L) {DMSetPeriodicity(dmNew, maxCell, L, bd);}
1046:   DMDestroy_Plex(dm);
1047:   dm->data = dmNew->data;
1048:   ((DM_Plex *) dmNew->data)->refct++;
1049:   dmNew->labels->refct++;
1050:   if (!--(dm->labels->refct)) {
1051:     DMLabelLink next = dm->labels->next;

1053:     /* destroy the labels */
1054:     while (next) {
1055:       DMLabelLink tmp = next->next;

1057:       DMLabelDestroy(&next->label);
1058:       PetscFree(next);
1059:       next = tmp;
1060:     }
1061:     PetscFree(dm->labels);
1062:   }
1063:   dm->labels = dmNew->labels;
1064:   dm->depthLabel = dmNew->depthLabel;
1065:   DMGetCoarseDM(dmNew,&coarseDM);
1066:   DMSetCoarseDM(dm,coarseDM);
1067:   return(0);
1068: }

1072: /* Swap dm with the contents of dmNew
1073:    - Swap the DM_Plex structure
1074:    - Swap the coordinates
1075:    - Swap the point PetscSF
1076: */
1077: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
1078: {
1079:   DM              coordDMA, coordDMB;
1080:   Vec             coordsA,  coordsB;
1081:   PetscSF         sfA,      sfB;
1082:   void            *tmp;
1083:   DMLabelLinkList listTmp;
1084:   DMLabel         depthTmp;
1085:   PetscErrorCode  ierr;

1088:   DMGetPointSF(dmA, &sfA);
1089:   DMGetPointSF(dmB, &sfB);
1090:   PetscObjectReference((PetscObject) sfA);
1091:   DMSetPointSF(dmA, sfB);
1092:   DMSetPointSF(dmB, sfA);
1093:   PetscObjectDereference((PetscObject) sfA);

1095:   DMGetCoordinateDM(dmA, &coordDMA);
1096:   DMGetCoordinateDM(dmB, &coordDMB);
1097:   PetscObjectReference((PetscObject) coordDMA);
1098:   DMSetCoordinateDM(dmA, coordDMB);
1099:   DMSetCoordinateDM(dmB, coordDMA);
1100:   PetscObjectDereference((PetscObject) coordDMA);

1102:   DMGetCoordinatesLocal(dmA, &coordsA);
1103:   DMGetCoordinatesLocal(dmB, &coordsB);
1104:   PetscObjectReference((PetscObject) coordsA);
1105:   DMSetCoordinatesLocal(dmA, coordsB);
1106:   DMSetCoordinatesLocal(dmB, coordsA);
1107:   PetscObjectDereference((PetscObject) coordsA);

1109:   tmp       = dmA->data;
1110:   dmA->data = dmB->data;
1111:   dmB->data = tmp;
1112:   listTmp   = dmA->labels;
1113:   dmA->labels = dmB->labels;
1114:   dmB->labels = listTmp;
1115:   depthTmp  = dmA->depthLabel;
1116:   dmA->depthLabel = dmB->depthLabel;
1117:   dmB->depthLabel = depthTmp;
1118:   return(0);
1119: }

1123: PetscErrorCode  DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
1124: {
1125:   DM_Plex       *mesh = (DM_Plex*) dm->data;

1129:   /* Handle viewing */
1130:   PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);
1131:   PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);
1132:   PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);
1133:   /* Point Location */
1134:   PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);
1135:   /* Generation and remeshing */
1136:   PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);
1137:   /* Projection behavior */
1138:   PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);
1139:   PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);
1140:   return(0);
1141: }

1145: PetscErrorCode  DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm)
1146: {
1147:   PetscInt       refine = 0, coarsen = 0, r;
1148:   PetscBool      isHierarchy;

1153:   PetscOptionsHead(PetscOptionsObject,"DMPlex Options");
1154:   /* Handle DMPlex refinement */
1155:   PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);
1156:   PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);
1157:   if (refine) {DMPlexSetRefinementUniform(dm, PETSC_TRUE);}
1158:   if (refine && isHierarchy) {
1159:     DM *dms, coarseDM;

1161:     DMGetCoarseDM(dm, &coarseDM);
1162:     PetscObjectReference((PetscObject)coarseDM);
1163:     PetscMalloc1(refine,&dms);
1164:     DMRefineHierarchy(dm, refine, dms);
1165:     /* Total hack since we do not pass in a pointer */
1166:     DMPlexSwap_Static(dm, dms[refine-1]);
1167:     if (refine == 1) {
1168:       DMSetCoarseDM(dm, dms[0]);
1169:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
1170:     } else {
1171:       DMSetCoarseDM(dm, dms[refine-2]);
1172:       DMPlexSetRegularRefinement(dm, PETSC_TRUE);
1173:       DMSetCoarseDM(dms[0], dms[refine-1]);
1174:       DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);
1175:     }
1176:     DMSetCoarseDM(dms[refine-1], coarseDM);
1177:     PetscObjectDereference((PetscObject)coarseDM);
1178:     /* Free DMs */
1179:     for (r = 0; r < refine; ++r) {
1180:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
1181:       DMDestroy(&dms[r]);
1182:     }
1183:     PetscFree(dms);
1184:   } else {
1185:     for (r = 0; r < refine; ++r) {
1186:       DM refinedMesh;

1188:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1189:       DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);
1190:       /* Total hack since we do not pass in a pointer */
1191:       DMPlexReplace_Static(dm, refinedMesh);
1192:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1193:       DMDestroy(&refinedMesh);
1194:     }
1195:   }
1196:   /* Handle DMPlex coarsening */
1197:   PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);
1198:   PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);
1199:   if (coarsen && isHierarchy) {
1200:     DM *dms;

1202:     PetscMalloc1(coarsen, &dms);
1203:     DMCoarsenHierarchy(dm, coarsen, dms);
1204:     /* Free DMs */
1205:     for (r = 0; r < coarsen; ++r) {
1206:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);
1207:       DMDestroy(&dms[r]);
1208:     }
1209:     PetscFree(dms);
1210:   } else {
1211:     for (r = 0; r < coarsen; ++r) {
1212:       DM coarseMesh;

1214:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1215:       DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);
1216:       /* Total hack since we do not pass in a pointer */
1217:       DMPlexReplace_Static(dm, coarseMesh);
1218:       DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1219:       DMDestroy(&coarseMesh);
1220:     }
1221:   }
1222:   /* Handle */
1223:   DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);
1224:   PetscOptionsTail();
1225:   return(0);
1226: }

1230: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
1231: {

1235:   DMCreateGlobalVector_Section_Private(dm,vec);
1236:   /* VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM); */
1237:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);
1238:   VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);
1239:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);
1240:   VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);
1241:   return(0);
1242: }

1246: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
1247: {

1251:   DMCreateLocalVector_Section_Private(dm,vec);
1252:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
1253:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
1254:   return(0);
1255: }

1259: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
1260: {
1261:   PetscInt       depth, d;

1265:   DMPlexGetDepth(dm, &depth);
1266:   if (depth == 1) {
1267:     DMGetDimension(dm, &d);
1268:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
1269:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
1270:     else               {*pStart = 0; *pEnd = 0;}
1271:   } else {
1272:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
1273:   }
1274:   return(0);
1275: }

1279: static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
1280: {
1281:   PetscSF        sf;

1285:   DMGetPointSF(dm, &sf);
1286:   PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);
1287:   return(0);
1288: }

1292: PetscErrorCode DMInitialize_Plex(DM dm)
1293: {

1297:   dm->ops->view                            = DMView_Plex;
1298:   dm->ops->load                            = DMLoad_Plex;
1299:   dm->ops->setfromoptions                  = DMSetFromOptions_Plex;
1300:   dm->ops->clone                           = DMClone_Plex;
1301:   dm->ops->setup                           = DMSetUp_Plex;
1302:   dm->ops->createdefaultsection            = DMCreateDefaultSection_Plex;
1303:   dm->ops->createdefaultconstraints        = DMCreateDefaultConstraints_Plex;
1304:   dm->ops->createglobalvector              = DMCreateGlobalVector_Plex;
1305:   dm->ops->createlocalvector               = DMCreateLocalVector_Plex;
1306:   dm->ops->getlocaltoglobalmapping         = NULL;
1307:   dm->ops->createfieldis                   = NULL;
1308:   dm->ops->createcoordinatedm              = DMCreateCoordinateDM_Plex;
1309:   dm->ops->getcoloring                     = NULL;
1310:   dm->ops->creatematrix                    = DMCreateMatrix_Plex;
1311:   dm->ops->createinterpolation             = DMCreateInterpolation_Plex;
1312:   dm->ops->getaggregates                   = NULL;
1313:   dm->ops->getinjection                    = DMCreateInjection_Plex;
1314:   dm->ops->refine                          = DMRefine_Plex;
1315:   dm->ops->coarsen                         = DMCoarsen_Plex;
1316:   dm->ops->refinehierarchy                 = DMRefineHierarchy_Plex;
1317:   dm->ops->coarsenhierarchy                = DMCoarsenHierarchy_Plex;
1318:   dm->ops->globaltolocalbegin              = NULL;
1319:   dm->ops->globaltolocalend                = NULL;
1320:   dm->ops->localtoglobalbegin              = NULL;
1321:   dm->ops->localtoglobalend                = NULL;
1322:   dm->ops->destroy                         = DMDestroy_Plex;
1323:   dm->ops->createsubdm                     = DMCreateSubDM_Plex;
1324:   dm->ops->getdimpoints                    = DMGetDimPoints_Plex;
1325:   dm->ops->locatepoints                    = DMLocatePoints_Plex;
1326:   dm->ops->projectfunctionlocal            = DMProjectFunctionLocal_Plex;
1327:   dm->ops->projectfunctionlabellocal       = DMProjectFunctionLabelLocal_Plex;
1328:   dm->ops->projectfieldlocal               = DMProjectFieldLocal_Plex;
1329:   dm->ops->computel2diff                   = DMComputeL2Diff_Plex;
1330:   dm->ops->computel2gradientdiff           = DMComputeL2GradientDiff_Plex;
1331:   dm->ops->computel2fielddiff              = DMComputeL2FieldDiff_Plex;
1332:   dm->ops->getneighbors                    = DMGetNeighors_Plex;
1333:   PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Plex);
1334:   return(0);
1335: }

1339: PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
1340: {
1341:   DM_Plex        *mesh = (DM_Plex *) dm->data;

1345:   mesh->refct++;
1346:   (*newdm)->data = mesh;
1347:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
1348:   DMInitialize_Plex(*newdm);
1349:   return(0);
1350: }

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

1358:   Level: intermediate

1360: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
1361: M*/

1365: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1366: {
1367:   DM_Plex        *mesh;
1368:   PetscInt       unit, d;

1373:   PetscNewLog(dm,&mesh);
1374:   dm->dim  = 0;
1375:   dm->data = mesh;

1377:   mesh->refct             = 1;
1378:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
1379:   mesh->maxConeSize       = 0;
1380:   mesh->cones             = NULL;
1381:   mesh->coneOrientations  = NULL;
1382:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
1383:   mesh->maxSupportSize    = 0;
1384:   mesh->supports          = NULL;
1385:   mesh->refinementUniform = PETSC_TRUE;
1386:   mesh->refinementLimit   = -1.0;

1388:   mesh->facesTmp = NULL;

1390:   mesh->tetgenOpts   = NULL;
1391:   mesh->triangleOpts = NULL;
1392:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
1393:   PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);
1394:   mesh->remeshBd     = PETSC_FALSE;

1396:   mesh->subpointMap = NULL;

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

1400:   mesh->regularRefinement   = PETSC_FALSE;
1401:   mesh->depthState          = -1;
1402:   mesh->globalVertexNumbers = NULL;
1403:   mesh->globalCellNumbers   = NULL;
1404:   mesh->anchorSection       = NULL;
1405:   mesh->anchorIS            = NULL;
1406:   mesh->createanchors       = NULL;
1407:   mesh->computeanchormatrix = NULL;
1408:   mesh->parentSection       = NULL;
1409:   mesh->parents             = NULL;
1410:   mesh->childIDs            = NULL;
1411:   mesh->childSection        = NULL;
1412:   mesh->children            = NULL;
1413:   mesh->referenceTree       = NULL;
1414:   mesh->getchildsymmetry    = NULL;
1415:   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1416:   mesh->vtkCellHeight       = 0;
1417:   mesh->useCone             = PETSC_FALSE;
1418:   mesh->useClosure          = PETSC_TRUE;
1419:   mesh->useAnchors          = PETSC_FALSE;

1421:   mesh->maxProjectionHeight = 0;

1423:   mesh->printSetValues = PETSC_FALSE;
1424:   mesh->printFEM       = 0;
1425:   mesh->printTol       = 1.0e-10;

1427:   DMInitialize_Plex(dm);
1428:   return(0);
1429: }

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

1436:   Collective on MPI_Comm

1438:   Input Parameter:
1439: . comm - The communicator for the DMPlex object

1441:   Output Parameter:
1442: . mesh  - The DMPlex object

1444:   Level: beginner

1446: .keywords: DMPlex, create
1447: @*/
1448: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1449: {

1454:   DMCreate(comm, mesh);
1455:   DMSetType(*mesh, DMPLEX);
1456:   return(0);
1457: }

1461: /*
1462:   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
1463: */
1464: static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert)
1465: {
1466:   PetscSF         sfPoint;
1467:   PetscLayout     vLayout;
1468:   PetscHashI      vhash;
1469:   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
1470:   const PetscInt *vrange;
1471:   PetscInt        numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
1472:   PETSC_UNUSED PetscHashIIter ret, iter;
1473:   PetscMPIInt     rank, numProcs;
1474:   PetscErrorCode  ierr;

1477:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1478:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);
1479:   /* Partition vertices */
1480:   PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);
1481:   PetscLayoutSetLocalSize(vLayout, numVertices);
1482:   PetscLayoutSetBlockSize(vLayout, 1);
1483:   PetscLayoutSetUp(vLayout);
1484:   PetscLayoutGetRanges(vLayout, &vrange);
1485:   /* Count vertices and map them to procs */
1486:   PetscHashICreate(vhash);
1487:   for (c = 0; c < numCells; ++c) {
1488:     for (p = 0; p < numCorners; ++p) {
1489:       PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter);
1490:     }
1491:   }
1492:   PetscHashISize(vhash, numVerticesAdj);
1493:   PetscMalloc1(numVerticesAdj, &verticesAdj);
1494:   off = 0; PetscHashIGetKeys(vhash, &off, verticesAdj);
1495:   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
1496:   PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);
1497:   for (v = 0; v < numVerticesAdj; ++v) {
1498:     const PetscInt gv = verticesAdj[v];
1499:     PetscInt       vrank;

1501:     PetscFindInt(gv, numProcs+1, vrange, &vrank);
1502:     vrank = vrank < 0 ? -(vrank+2) : vrank;
1503:     remoteVerticesAdj[v].index = gv - vrange[vrank];
1504:     remoteVerticesAdj[v].rank  = vrank;
1505:   }
1506:   /* Create cones */
1507:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
1508:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
1509:   DMSetUp(dm);
1510:   DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);
1511:   for (c = 0; c < numCells; ++c) {
1512:     for (p = 0; p < numCorners; ++p) {
1513:       const PetscInt gv = cells[c*numCorners+p];
1514:       PetscInt       lv;

1516:       PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);
1517:       if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv);
1518:       cone[p] = lv+numCells;
1519:     }
1520:     DMPlexSetCone(dm, c, cone);
1521:   }
1522:   DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);
1523:   /* Create SF for vertices */
1524:   PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);
1525:   PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");
1526:   PetscSFSetFromOptions(*sfVert);
1527:   PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);
1528:   PetscFree(verticesAdj);
1529:   /* Build pointSF */
1530:   PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);
1531:   for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;}
1532:   for (v = 0; v < numVertices;    ++v) {vertexOwner[v].index = -1;         vertexOwner[v].rank = -1;}
1533:   PetscSFReduceBegin(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
1534:   PetscSFReduceEnd(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);
1535:   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);
1536:   PetscSFBcastBegin(*sfVert, MPI_2INT, vertexOwner, vertexLocal);
1537:   PetscSFBcastEnd(*sfVert, MPI_2INT, vertexOwner, vertexLocal);
1538:   for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost;
1539:   PetscMalloc1(numVerticesGhost, &localVertex);
1540:   PetscMalloc1(numVerticesGhost, &remoteVertex);
1541:   for (v = 0, g = 0; v < numVerticesAdj; ++v) {
1542:     if (vertexLocal[v].rank != rank) {
1543:       localVertex[g]        = v+numCells;
1544:       remoteVertex[g].index = vertexLocal[v].index;
1545:       remoteVertex[g].rank  = vertexLocal[v].rank;
1546:       ++g;
1547:     }
1548:   }
1549:   PetscFree2(vertexLocal, vertexOwner);
1550:   if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost);
1551:   DMGetPointSF(dm, &sfPoint);
1552:   PetscObjectSetName((PetscObject) sfPoint, "point SF");
1553:   PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);
1554:   PetscLayoutDestroy(&vLayout);
1555:   PetscHashIDestroy(vhash);
1556:   /* Fill in the rest of the topology structure */
1557:   DMPlexSymmetrize(dm);
1558:   DMPlexStratify(dm);
1559:   return(0);
1560: }

1564: /*
1565:   This takes as input the coordinates for each owned vertex
1566: */
1567: static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscSF sfVert, const double vertexCoords[])
1568: {
1569:   PetscSection   coordSection;
1570:   Vec            coordinates;
1571:   PetscScalar   *coords;
1572:   PetscInt       numVertices, numVerticesAdj, coordSize, v;

1576:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
1577:   DMGetCoordinateSection(dm, &coordSection);
1578:   PetscSectionSetNumFields(coordSection, 1);
1579:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
1580:   PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);
1581:   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
1582:     PetscSectionSetDof(coordSection, v, spaceDim);
1583:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
1584:   }
1585:   PetscSectionSetUp(coordSection);
1586:   PetscSectionGetStorageSize(coordSection, &coordSize);
1587:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
1588:   VecSetBlockSize(coordinates, spaceDim);
1589:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1590:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1591:   VecSetType(coordinates,VECSTANDARD);
1592:   VecGetArray(coordinates, &coords);
1593:   {
1594:     MPI_Datatype coordtype;

1596:     /* Need a temp buffer for coords if we have complex/single */
1597:     MPI_Type_contiguous(spaceDim, MPI_DOUBLE, &coordtype);
1598:     MPI_Type_commit(&coordtype);
1599:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);
1600:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);
1601:     MPI_Type_free(&coordtype);
1602:   }
1603:   VecRestoreArray(coordinates, &coords);
1604:   DMSetCoordinatesLocal(dm, coordinates);
1605:   VecDestroy(&coordinates);
1606:   return(0);
1607: }

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

1614:   Input Parameters:
1615: + comm - The communicator
1616: . dim - The topological dimension of the mesh
1617: . numCells - The number of cells owned by this process
1618: . numVertices - The number of vertices owned by this process
1619: . numCorners - The number of vertices for each cell
1620: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1621: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
1622: . spaceDim - The spatial dimension used for coordinates
1623: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

1625:   Output Parameter:
1626: . dm - The DM

1628:   Note: Two triangles sharing a face
1629: $
1630: $        2
1631: $      / | \
1632: $     /  |  \
1633: $    /   |   \
1634: $   0  0 | 1  3
1635: $    \   |   /
1636: $     \  |  /
1637: $      \ | /
1638: $        1
1639: would have input
1640: $  numCells = 2, numVertices = 4
1641: $  cells = [0 1 2  1 3 2]
1642: $
1643: which would result in the DMPlex
1644: $
1645: $        4
1646: $      / | \
1647: $     /  |  \
1648: $    /   |   \
1649: $   2  0 | 1  5
1650: $    \   |   /
1651: $     \  |  /
1652: $      \ | /
1653: $        3

1655:   Level: beginner

1657: .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
1658: @*/
1659: PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm)
1660: {
1661:   PetscSF        sfVert;

1665:   DMCreate(comm, dm);
1666:   DMSetType(*dm, DMPLEX);
1669:   DMSetDimension(*dm, dim);
1670:   DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);
1671:   if (interpolate) {
1672:     DM idm = NULL;

1674:     DMPlexInterpolate(*dm, &idm);
1675:     DMDestroy(dm);
1676:     *dm  = idm;
1677:   }
1678:   DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, sfVert, vertexCoords);
1679:   PetscSFDestroy(&sfVert);
1680:   return(0);
1681: }

1685: /*
1686:   This takes as input the common mesh generator output, a list of the vertices for each cell
1687: */
1688: static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
1689: {
1690:   PetscInt      *cone, c, p;

1694:   DMPlexSetChart(dm, 0, numCells+numVertices);
1695:   for (c = 0; c < numCells; ++c) {
1696:     DMPlexSetConeSize(dm, c, numCorners);
1697:   }
1698:   DMSetUp(dm);
1699:   DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);
1700:   for (c = 0; c < numCells; ++c) {
1701:     for (p = 0; p < numCorners; ++p) {
1702:       cone[p] = cells[c*numCorners+p]+numCells;
1703:     }
1704:     DMPlexSetCone(dm, c, cone);
1705:   }
1706:   DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);
1707:   DMPlexSymmetrize(dm);
1708:   DMPlexStratify(dm);
1709:   return(0);
1710: }

1714: /*
1715:   This takes as input the coordinates for each vertex
1716: */
1717: static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
1718: {
1719:   PetscSection   coordSection;
1720:   Vec            coordinates;
1721:   DM             cdm;
1722:   PetscScalar   *coords;
1723:   PetscInt       v, d;

1727:   DMGetCoordinateSection(dm, &coordSection);
1728:   PetscSectionSetNumFields(coordSection, 1);
1729:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
1730:   PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
1731:   for (v = numCells; v < numCells+numVertices; ++v) {
1732:     PetscSectionSetDof(coordSection, v, spaceDim);
1733:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
1734:   }
1735:   PetscSectionSetUp(coordSection);

1737:   DMGetCoordinateDM(dm, &cdm);
1738:   DMCreateLocalVector(cdm, &coordinates);
1739:   VecSetBlockSize(coordinates, spaceDim);
1740:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1741:   VecGetArray(coordinates, &coords);
1742:   for (v = 0; v < numVertices; ++v) {
1743:     for (d = 0; d < spaceDim; ++d) {
1744:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
1745:     }
1746:   }
1747:   VecRestoreArray(coordinates, &coords);
1748:   DMSetCoordinatesLocal(dm, coordinates);
1749:   VecDestroy(&coordinates);
1750:   return(0);
1751: }

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

1758:   Input Parameters:
1759: + comm - The communicator
1760: . dim - The topological dimension of the mesh
1761: . numCells - The number of cells
1762: . numVertices - The number of vertices
1763: . numCorners - The number of vertices for each cell
1764: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1765: . cells - An array of numCells*numCorners numbers, the vertices for each cell
1766: . spaceDim - The spatial dimension used for coordinates
1767: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

1769:   Output Parameter:
1770: . dm - The DM

1772:   Note: Two triangles sharing a face
1773: $
1774: $        2
1775: $      / | \
1776: $     /  |  \
1777: $    /   |   \
1778: $   0  0 | 1  3
1779: $    \   |   /
1780: $     \  |  /
1781: $      \ | /
1782: $        1
1783: would have input
1784: $  numCells = 2, numVertices = 4
1785: $  cells = [0 1 2  1 3 2]
1786: $
1787: which would result in the DMPlex
1788: $
1789: $        4
1790: $      / | \
1791: $     /  |  \
1792: $    /   |   \
1793: $   2  0 | 1  5
1794: $    \   |   /
1795: $     \  |  /
1796: $      \ | /
1797: $        3

1799:   Level: beginner

1801: .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
1802: @*/
1803: 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)
1804: {

1808:   DMCreate(comm, dm);
1809:   DMSetType(*dm, DMPLEX);
1810:   DMSetDimension(*dm, dim);
1811:   DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);
1812:   if (interpolate) {
1813:     DM idm = NULL;

1815:     DMPlexInterpolate(*dm, &idm);
1816:     DMDestroy(dm);
1817:     *dm  = idm;
1818:   }
1819:   DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);
1820:   return(0);
1821: }

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

1828:   Input Parameters:
1829: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
1830: . depth - The depth of the DAG
1831: . numPoints - The number of points at each depth
1832: . coneSize - The cone size of each point
1833: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1834: . coneOrientations - The orientation of each cone point
1835: - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex

1837:   Output Parameter:
1838: . dm - The DM

1840:   Note: Two triangles sharing a face would have input
1841: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1842: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1843: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1844: $
1845: which would result in the DMPlex
1846: $
1847: $        4
1848: $      / | \
1849: $     /  |  \
1850: $    /   |   \
1851: $   2  0 | 1  5
1852: $    \   |   /
1853: $     \  |  /
1854: $      \ | /
1855: $        3
1856: $
1857: $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()

1859:   Level: advanced

1861: .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1862: @*/
1863: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
1864: {
1865:   Vec            coordinates;
1866:   PetscSection   coordSection;
1867:   PetscScalar    *coords;
1868:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

1872:   DMGetDimension(dm, &dim);
1873:   DMGetCoordinateDim(dm, &dimEmbed);
1874:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
1875:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
1876:   DMPlexSetChart(dm, pStart, pEnd);
1877:   for (p = pStart; p < pEnd; ++p) {
1878:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
1879:     if (firstVertex < 0 && !coneSize[p - pStart]) {
1880:       firstVertex = p - pStart;
1881:     }
1882:   }
1883:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
1884:   DMSetUp(dm); /* Allocate space for cones */
1885:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
1886:     DMPlexSetCone(dm, p, &cones[off]);
1887:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
1888:   }
1889:   DMPlexSymmetrize(dm);
1890:   DMPlexStratify(dm);
1891:   /* Build coordinates */
1892:   DMGetCoordinateSection(dm, &coordSection);
1893:   PetscSectionSetNumFields(coordSection, 1);
1894:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
1895:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
1896:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
1897:     PetscSectionSetDof(coordSection, v, dimEmbed);
1898:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
1899:   }
1900:   PetscSectionSetUp(coordSection);
1901:   PetscSectionGetStorageSize(coordSection, &coordSize);
1902:   VecCreate(PETSC_COMM_SELF, &coordinates);
1903:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1904:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1905:   VecSetBlockSize(coordinates, dimEmbed);
1906:   VecSetType(coordinates,VECSTANDARD);
1907:   VecGetArray(coordinates, &coords);
1908:   for (v = 0; v < numPoints[0]; ++v) {
1909:     PetscInt off;

1911:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
1912:     for (d = 0; d < dimEmbed; ++d) {
1913:       coords[off+d] = vertexCoords[v*dimEmbed+d];
1914:     }
1915:   }
1916:   VecRestoreArray(coordinates, &coords);
1917:   DMSetCoordinatesLocal(dm, coordinates);
1918:   VecDestroy(&coordinates);
1919:   return(0);
1920: }

1924: /*@C
1925:   DMPlexCreateFromFile - This takes a filename and produces a DM

1927:   Input Parameters:
1928: + comm - The communicator
1929: . filename - A file name
1930: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1932:   Output Parameter:
1933: . dm - The DM

1935:   Level: beginner

1937: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
1938: @*/
1939: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1940: {
1941:   const char    *extGmsh   = ".msh";
1942:   const char    *extCGNS   = ".cgns";
1943:   const char    *extExodus = ".exo";
1944:   const char    *extFluent = ".cas";
1945:   const char    *extHDF5   = ".h5";
1946:   size_t         len;
1947:   PetscBool      isGmsh, isCGNS, isExodus, isFluent, isHDF5;
1948:   PetscMPIInt    rank;

1954:   MPI_Comm_rank(comm, &rank);
1955:   PetscStrlen(filename, &len);
1956:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
1957:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,   4, &isGmsh);
1958:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,   5, &isCGNS);
1959:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);
1960:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);
1961:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,   3, &isHDF5);
1962:   if (isGmsh) {
1963:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
1964:   } else if (isCGNS) {
1965:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
1966:   } else if (isExodus) {
1967:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
1968:   } else if (isFluent) {
1969:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
1970:   } else if (isHDF5) {
1971:     PetscViewer viewer;

1973:     PetscViewerCreate(comm, &viewer);
1974:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
1975:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1976:     PetscViewerFileSetName(viewer, filename);
1977:     DMCreate(comm, dm);
1978:     DMSetType(*dm, DMPLEX);
1979:     DMLoad(*dm, viewer);
1980:     PetscViewerDestroy(&viewer);
1981:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
1982:   return(0);
1983: }

1987: /*@
1988:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

1990:   Collective on comm

1992:   Input Parameters:
1993: + comm    - The communicator
1994: . dim     - The spatial dimension
1995: - simplex - Flag for simplex, otherwise use a tensor-product cell

1997:   Output Parameter:
1998: . refdm - The reference cell

2000:   Level: intermediate

2002: .keywords: reference cell
2003: .seealso:
2004: @*/
2005: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
2006: {
2007:   DM             rdm;
2008:   Vec            coords;

2012:   DMCreate(comm, &rdm);
2013:   DMSetType(rdm, DMPLEX);
2014:   DMSetDimension(rdm, dim);
2015:   switch (dim) {
2016:   case 0:
2017:   {
2018:     PetscInt    numPoints[1]        = {1};
2019:     PetscInt    coneSize[1]         = {0};
2020:     PetscInt    cones[1]            = {0};
2021:     PetscInt    coneOrientations[1] = {0};
2022:     PetscScalar vertexCoords[1]     = {0.0};

2024:     DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2025:   }
2026:   break;
2027:   case 1:
2028:   {
2029:     PetscInt    numPoints[2]        = {2, 1};
2030:     PetscInt    coneSize[3]         = {2, 0, 0};
2031:     PetscInt    cones[2]            = {1, 2};
2032:     PetscInt    coneOrientations[2] = {0, 0};
2033:     PetscScalar vertexCoords[2]     = {-1.0,  1.0};

2035:     DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2036:   }
2037:   break;
2038:   case 2:
2039:     if (simplex) {
2040:       PetscInt    numPoints[2]        = {3, 1};
2041:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
2042:       PetscInt    cones[3]            = {1, 2, 3};
2043:       PetscInt    coneOrientations[3] = {0, 0, 0};
2044:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

2046:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2047:     } else {
2048:       PetscInt    numPoints[2]        = {4, 1};
2049:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2050:       PetscInt    cones[4]            = {1, 2, 3, 4};
2051:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2052:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

2054:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2055:     }
2056:   break;
2057:   case 3:
2058:     if (simplex) {
2059:       PetscInt    numPoints[2]        = {4, 1};
2060:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2061:       PetscInt    cones[4]            = {1, 3, 2, 4};
2062:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2063:       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};

2065:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2066:     } else {
2067:       PetscInt    numPoints[2]        = {8, 1};
2068:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
2069:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
2070:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
2071:       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,
2072:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

2074:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2075:     }
2076:   break;
2077:   default:
2078:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
2079:   }
2080:   *refdm = NULL;
2081:   DMPlexInterpolate(rdm, refdm);
2082:   if (rdm->coordinateDM) {
2083:     DM           ncdm;
2084:     PetscSection cs;
2085:     PetscInt     pEnd = -1;

2087:     DMGetDefaultSection(rdm->coordinateDM, &cs);
2088:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
2089:     if (pEnd >= 0) {
2090:       DMClone(rdm->coordinateDM, &ncdm);
2091:       DMSetDefaultSection(ncdm, cs);
2092:       DMSetCoordinateDM(*refdm, ncdm);
2093:       DMDestroy(&ncdm);
2094:     }
2095:   }
2096:   DMGetCoordinatesLocal(rdm, &coords);
2097:   if (coords) {
2098:     DMSetCoordinatesLocal(*refdm, coords);
2099:   } else {
2100:     DMGetCoordinates(rdm, &coords);
2101:     if (coords) {DMSetCoordinates(*refdm, coords);}
2102:   }
2103:   DMDestroy(&rdm);
2104:   return(0);
2105: }