Actual source code: plexcreate.c

petsc-master 2017-02-25
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: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
437: {
438:   PetscInt       markerTop      = 1, faceMarkerTop      = 1;
439:   PetscInt       markerBottom   = 1, faceMarkerBottom   = 1;
440:   PetscInt       markerFront    = 1, faceMarkerFront    = 1;
441:   PetscInt       markerBack     = 1, faceMarkerBack     = 1;
442:   PetscInt       markerRight    = 1, faceMarkerRight    = 1;
443:   PetscInt       markerLeft     = 1, faceMarkerLeft     = 1;
444:   PetscInt       dim;
445:   PetscBool      markerSeparate = PETSC_FALSE;
446:   PetscMPIInt    rank;

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

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

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

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

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

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

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

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

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

821: /*@
822:   DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice.

824:   Collective on MPI_Comm

826:   Input Parameters:
827: + comm  - The communicator for the DM object
828: . lower - The lower left corner coordinates
829: . upper - The upper right corner coordinates
830: . edges - The number of cells in each direction
831: . bdX   - The boundary type for the X direction
832: - bdY   - The boundary type for the Y direction

834:   Output Parameter:
835: . dm  - The DM object

837:   Note: Here is the numbering returned for 2 cells in each direction:
838: $ 22--8-23--9--24
839: $  |     |     |
840: $ 13  2 14  3  15
841: $  |     |     |
842: $ 19--6-20--7--21
843: $  |     |     |
844: $ 10  0 11  1 12
845: $  |     |     |
846: $ 16--4-17--5--18

848:   Level: beginner

850: .keywords: DM, create
851: .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate()
852: @*/
853: PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY)
854: {
855:   PetscReal      lower3[3], upper3[3];
856:   PetscInt       edges3[3];

860:   lower3[0] = lower[0]; lower3[1] = lower[1]; lower3[2] = 0.;
861:   upper3[0] = upper[0]; upper3[1] = upper[1]; upper3[2] = 0.;
862:   edges3[0] = edges[0]; edges3[1] = edges[1]; edges3[2] = 0;
863:   DMPlexCreateCubeMesh_Internal(dm, lower3, upper3, edges3, bdX, bdY, DM_BOUNDARY_NONE);
864:   return(0);
865: }

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

870:   Collective on MPI_Comm

872:   Input Parameters:
873: + comm - The communicator for the DM object
874: . dim - The spatial dimension
875: . numFaces - Number of faces per dimension
876: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

878:   Output Parameter:
879: . dm  - The DM object

881:   Level: beginner

883: .keywords: DM, create
884: .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate()
885: @*/
886: PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm)
887: {
888:   DM             boundary;

893:   DMCreate(comm, &boundary);
895:   DMSetType(boundary, DMPLEX);
896:   DMSetDimension(boundary, dim-1);
897:   DMSetCoordinateDim(boundary, dim);
898:   switch (dim) {
899:   case 2:
900:   {
901:     PetscReal lower[2] = {0.0, 0.0};
902:     PetscReal upper[2] = {1.0, 1.0};
903:     PetscInt  edges[2];

905:     edges[0] = numFaces; edges[1] = numFaces;
906:     DMPlexCreateSquareBoundary(boundary, lower, upper, edges);
907:     break;
908:   }
909:   case 3:
910:   {
911:     PetscReal lower[3] = {0.0, 0.0, 0.0};
912:     PetscReal upper[3] = {1.0, 1.0, 1.0};
913:     PetscInt  faces[3];

915:     faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces;
916:     DMPlexCreateCubeBoundary(boundary, lower, upper, faces);
917:     break;
918:   }
919:   default:
920:     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
921:   }
922:   DMPlexGenerate(boundary, NULL, interpolate, dm);
923:   DMDestroy(&boundary);
924:   return(0);
925: }

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

930:   Collective on MPI_Comm

932:   Input Parameters:
933: + comm  - The communicator for the DM object
934: . dim   - The spatial dimension
935: . periodicX - The boundary type for the X direction
936: . periodicY - The boundary type for the Y direction
937: . periodicZ - The boundary type for the Z direction
938: - cells - The number of cells in each direction

940:   Output Parameter:
941: . dm  - The DM object

943:   Level: beginner

945: .keywords: DM, create
946: .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate()
947: @*/
948: PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm)
949: {
950:   PetscInt       i;

955:   DMCreate(comm, dm);
957:   DMSetType(*dm, DMPLEX);
958:   DMSetDimension(*dm, dim);
959:   switch (dim) {
960:   case 2:
961:   {
962:     PetscReal lower[2] = {0.0, 0.0};
963:     PetscReal upper[2] = {1.0, 1.0};

965:     DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);
966:     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
967:         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST) {
968:       PetscReal      L[2];
969:       PetscReal      maxCell[2];
970:       DMBoundaryType bdType[2];

972:       bdType[0] = periodicX;
973:       bdType[1] = periodicY;
974:       for (i = 0; i < dim; i++) {
975:         L[i]       = upper[i] - lower[i];
976:         maxCell[i] = 1.1 * (L[i] / cells[i]);
977:       }

979:       DMSetPeriodicity(*dm,maxCell,L,bdType);
980:     }
981:     break;
982:   }
983:   case 3:
984:   {
985:     PetscReal lower[3] = {0.0, 0.0, 0.0};
986:     PetscReal upper[3] = {1.0, 1.0, 1.0};

988:     DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);
989:     if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST ||
990:         periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST ||
991:         periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
992:       PetscReal      L[3];
993:       PetscReal      maxCell[3];
994:       DMBoundaryType bdType[3];

996:       bdType[0] = periodicX;
997:       bdType[1] = periodicY;
998:       bdType[2] = periodicZ;
999:       for (i = 0; i < dim; i++) {
1000:         L[i]       = upper[i] - lower[i];
1001:         maxCell[i] = 1.1 * (L[i] / cells[i]);
1002:       }

1004:       DMSetPeriodicity(*dm,maxCell,L,bdType);
1005:     }
1006:     break;
1007:   }
1008:   default:
1009:     SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim);
1010:   }
1011:   return(0);
1012: }

1014: /* External function declarations here */
1015: extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
1016: extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
1017: extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm);
1018: extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
1019: extern PetscErrorCode DMCreateMatrix_Plex(DM dm,  Mat *J);
1020: extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
1021: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
1022: extern PetscErrorCode DMSetUp_Plex(DM dm);
1023: extern PetscErrorCode DMDestroy_Plex(DM dm);
1024: extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
1025: extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
1026: extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm);

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

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

1058:     /* destroy the labels */
1059:     while (next) {
1060:       DMLabelLink tmp = next->next;

1062:       DMLabelDestroy(&next->label);
1063:       PetscFree(next);
1064:       next = tmp;
1065:     }
1066:     PetscFree(dm->labels);
1067:   }
1068:   dm->labels = dmNew->labels;
1069:   dm->depthLabel = dmNew->depthLabel;
1070:   DMGetCoarseDM(dmNew,&coarseDM);
1071:   DMSetCoarseDM(dm,coarseDM);
1072:   return(0);
1073: }

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

1091:   DMGetPointSF(dmA, &sfA);
1092:   DMGetPointSF(dmB, &sfB);
1093:   PetscObjectReference((PetscObject) sfA);
1094:   DMSetPointSF(dmA, sfB);
1095:   DMSetPointSF(dmB, sfA);
1096:   PetscObjectDereference((PetscObject) sfA);

1098:   DMGetCoordinateDM(dmA, &coordDMA);
1099:   DMGetCoordinateDM(dmB, &coordDMB);
1100:   PetscObjectReference((PetscObject) coordDMA);
1101:   DMSetCoordinateDM(dmA, coordDMB);
1102:   DMSetCoordinateDM(dmB, coordDMA);
1103:   PetscObjectDereference((PetscObject) coordDMA);

1105:   DMGetCoordinatesLocal(dmA, &coordsA);
1106:   DMGetCoordinatesLocal(dmB, &coordsB);
1107:   PetscObjectReference((PetscObject) coordsA);
1108:   DMSetCoordinatesLocal(dmA, coordsB);
1109:   DMSetCoordinatesLocal(dmB, coordsA);
1110:   PetscObjectDereference((PetscObject) coordsA);

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

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

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

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

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

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

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

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

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

1227: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec)
1228: {

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

1241: static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec)
1242: {

1246:   DMCreateLocalVector_Section_Private(dm,vec);
1247:   VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);
1248:   VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);
1249:   return(0);
1250: }

1252: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
1253: {
1254:   PetscInt       depth, d;

1258:   DMPlexGetDepth(dm, &depth);
1259:   if (depth == 1) {
1260:     DMGetDimension(dm, &d);
1261:     if (dim == 0)      {DMPlexGetDepthStratum(dm, dim, pStart, pEnd);}
1262:     else if (dim == d) {DMPlexGetDepthStratum(dm, 1, pStart, pEnd);}
1263:     else               {*pStart = 0; *pEnd = 0;}
1264:   } else {
1265:     DMPlexGetDepthStratum(dm, dim, pStart, pEnd);
1266:   }
1267:   return(0);
1268: }

1270: static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
1271: {
1272:   PetscSF        sf;

1276:   DMGetPointSF(dm, &sf);
1277:   PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);
1278:   return(0);
1279: }

1281: static PetscErrorCode DMInitialize_Plex(DM dm)
1282: {

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

1327: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
1328: {
1329:   DM_Plex        *mesh = (DM_Plex *) dm->data;

1333:   mesh->refct++;
1334:   (*newdm)->data = mesh;
1335:   PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);
1336:   DMInitialize_Plex(*newdm);
1337:   return(0);
1338: }

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

1346:   Level: intermediate

1348: .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType()
1349: M*/

1351: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
1352: {
1353:   DM_Plex        *mesh;
1354:   PetscInt       unit, d;

1359:   PetscNewLog(dm,&mesh);
1360:   dm->dim  = 0;
1361:   dm->data = mesh;

1363:   mesh->refct             = 1;
1364:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);
1365:   mesh->maxConeSize       = 0;
1366:   mesh->cones             = NULL;
1367:   mesh->coneOrientations  = NULL;
1368:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);
1369:   mesh->maxSupportSize    = 0;
1370:   mesh->supports          = NULL;
1371:   mesh->refinementUniform = PETSC_TRUE;
1372:   mesh->refinementLimit   = -1.0;

1374:   mesh->facesTmp = NULL;

1376:   mesh->tetgenOpts   = NULL;
1377:   mesh->triangleOpts = NULL;
1378:   PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);
1379:   PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);
1380:   mesh->remeshBd     = PETSC_FALSE;

1382:   mesh->subpointMap = NULL;

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

1386:   mesh->regularRefinement   = PETSC_FALSE;
1387:   mesh->depthState          = -1;
1388:   mesh->globalVertexNumbers = NULL;
1389:   mesh->globalCellNumbers   = NULL;
1390:   mesh->anchorSection       = NULL;
1391:   mesh->anchorIS            = NULL;
1392:   mesh->createanchors       = NULL;
1393:   mesh->computeanchormatrix = NULL;
1394:   mesh->parentSection       = NULL;
1395:   mesh->parents             = NULL;
1396:   mesh->childIDs            = NULL;
1397:   mesh->childSection        = NULL;
1398:   mesh->children            = NULL;
1399:   mesh->referenceTree       = NULL;
1400:   mesh->getchildsymmetry    = NULL;
1401:   for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE;
1402:   mesh->vtkCellHeight       = 0;
1403:   mesh->useCone             = PETSC_FALSE;
1404:   mesh->useClosure          = PETSC_TRUE;
1405:   mesh->useAnchors          = PETSC_FALSE;

1407:   mesh->maxProjectionHeight = 0;

1409:   mesh->printSetValues = PETSC_FALSE;
1410:   mesh->printFEM       = 0;
1411:   mesh->printTol       = 1.0e-10;

1413:   DMInitialize_Plex(dm);
1414:   return(0);
1415: }

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

1420:   Collective on MPI_Comm

1422:   Input Parameter:
1423: . comm - The communicator for the DMPlex object

1425:   Output Parameter:
1426: . mesh  - The DMPlex object

1428:   Level: beginner

1430: .keywords: DMPlex, create
1431: @*/
1432: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
1433: {

1438:   DMCreate(comm, mesh);
1439:   DMSetType(*mesh, DMPLEX);
1440:   return(0);
1441: }

1443: /*
1444:   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
1445: */
1446: static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert)
1447: {
1448:   PetscSF         sfPoint;
1449:   PetscLayout     vLayout;
1450:   PetscHashI      vhash;
1451:   PetscSFNode    *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex;
1452:   const PetscInt *vrange;
1453:   PetscInt        numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g;
1454:   PETSC_UNUSED PetscHashIIter ret, iter;
1455:   PetscMPIInt     rank, size;
1456:   PetscErrorCode  ierr;

1459:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1460:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);
1461:   /* Partition vertices */
1462:   PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);
1463:   PetscLayoutSetLocalSize(vLayout, numVertices);
1464:   PetscLayoutSetBlockSize(vLayout, 1);
1465:   PetscLayoutSetUp(vLayout);
1466:   PetscLayoutGetRanges(vLayout, &vrange);
1467:   /* Count vertices and map them to procs */
1468:   PetscHashICreate(vhash);
1469:   for (c = 0; c < numCells; ++c) {
1470:     for (p = 0; p < numCorners; ++p) {
1471:       PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter);
1472:     }
1473:   }
1474:   PetscHashISize(vhash, numVerticesAdj);
1475:   PetscMalloc1(numVerticesAdj, &verticesAdj);
1476:   off = 0; PetscHashIGetKeys(vhash, &off, verticesAdj);
1477:   if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj);
1478:   PetscSortInt(numVerticesAdj, verticesAdj);
1479:   PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);
1480:   for (v = 0; v < numVerticesAdj; ++v) {
1481:     const PetscInt gv = verticesAdj[v];
1482:     PetscInt       vrank;

1484:     PetscFindInt(gv, size+1, vrange, &vrank);
1485:     vrank = vrank < 0 ? -(vrank+2) : vrank;
1486:     remoteVerticesAdj[v].index = gv - vrange[vrank];
1487:     remoteVerticesAdj[v].rank  = vrank;
1488:   }
1489:   /* Create cones */
1490:   DMPlexSetChart(dm, 0, numCells+numVerticesAdj);
1491:   for (c = 0; c < numCells; ++c) {DMPlexSetConeSize(dm, c, numCorners);}
1492:   DMSetUp(dm);
1493:   DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);
1494:   for (c = 0; c < numCells; ++c) {
1495:     for (p = 0; p < numCorners; ++p) {
1496:       const PetscInt gv = cells[c*numCorners+p];
1497:       PetscInt       lv;

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

1545: /*
1546:   This takes as input the coordinates for each owned vertex
1547: */
1548: static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscSF sfVert, const double vertexCoords[])
1549: {
1550:   PetscSection   coordSection;
1551:   Vec            coordinates;
1552:   PetscScalar   *coords;
1553:   PetscInt       numVertices, numVerticesAdj, coordSize, v;

1557:   PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);
1558:   DMGetCoordinateSection(dm, &coordSection);
1559:   PetscSectionSetNumFields(coordSection, 1);
1560:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
1561:   PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);
1562:   for (v = numCells; v < numCells+numVerticesAdj; ++v) {
1563:     PetscSectionSetDof(coordSection, v, spaceDim);
1564:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
1565:   }
1566:   PetscSectionSetUp(coordSection);
1567:   PetscSectionGetStorageSize(coordSection, &coordSize);
1568:   VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);
1569:   VecSetBlockSize(coordinates, spaceDim);
1570:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1571:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1572:   VecSetType(coordinates,VECSTANDARD);
1573:   VecGetArray(coordinates, &coords);
1574:   {
1575:     MPI_Datatype coordtype;

1577:     /* Need a temp buffer for coords if we have complex/single */
1578:     MPI_Type_contiguous(spaceDim, MPI_DOUBLE, &coordtype);
1579:     MPI_Type_commit(&coordtype);
1580:     PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);
1581:     PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);
1582:     MPI_Type_free(&coordtype);
1583:   }
1584:   VecRestoreArray(coordinates, &coords);
1585:   DMSetCoordinatesLocal(dm, coordinates);
1586:   VecDestroy(&coordinates);
1587:   return(0);
1588: }

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

1593:   Input Parameters:
1594: + comm - The communicator
1595: . dim - The topological dimension of the mesh
1596: . numCells - The number of cells owned by this process
1597: . numVertices - The number of vertices owned by this process
1598: . numCorners - The number of vertices for each cell
1599: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1600: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
1601: . spaceDim - The spatial dimension used for coordinates
1602: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

1604:   Output Parameter:
1605: + dm - The DM
1606: - vertexSF - Optional, SF describing complete vertex ownership

1608:   Note: Two triangles sharing a face
1609: $
1610: $        2
1611: $      / | \
1612: $     /  |  \
1613: $    /   |   \
1614: $   0  0 | 1  3
1615: $    \   |   /
1616: $     \  |  /
1617: $      \ | /
1618: $        1
1619: would have input
1620: $  numCells = 2, numVertices = 4
1621: $  cells = [0 1 2  1 3 2]
1622: $
1623: which would result in the DMPlex
1624: $
1625: $        4
1626: $      / | \
1627: $     /  |  \
1628: $    /   |   \
1629: $   2  0 | 1  5
1630: $    \   |   /
1631: $     \  |  /
1632: $      \ | /
1633: $        3

1635:   Level: beginner

1637: .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate()
1638: @*/
1639: PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], PetscSF *vertexSF, DM *dm)
1640: {
1641:   PetscSF        sfVert;

1645:   DMCreate(comm, dm);
1646:   DMSetType(*dm, DMPLEX);
1649:   DMSetDimension(*dm, dim);
1650:   DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);
1651:   if (interpolate) {
1652:     DM idm = NULL;

1654:     DMPlexInterpolate(*dm, &idm);
1655:     DMDestroy(dm);
1656:     *dm  = idm;
1657:   }
1658:   DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, sfVert, vertexCoords);
1659:   if (vertexSF) *vertexSF = sfVert;
1660:   else PetscSFDestroy(&sfVert);
1661:   return(0);
1662: }

1664: /*
1665:   This takes as input the common mesh generator output, a list of the vertices for each cell
1666: */
1667: static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[])
1668: {
1669:   PetscInt      *cone, c, p;

1673:   DMPlexSetChart(dm, 0, numCells+numVertices);
1674:   for (c = 0; c < numCells; ++c) {
1675:     DMPlexSetConeSize(dm, c, numCorners);
1676:   }
1677:   DMSetUp(dm);
1678:   DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);
1679:   for (c = 0; c < numCells; ++c) {
1680:     for (p = 0; p < numCorners; ++p) {
1681:       cone[p] = cells[c*numCorners+p]+numCells;
1682:     }
1683:     DMPlexSetCone(dm, c, cone);
1684:   }
1685:   DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);
1686:   DMPlexSymmetrize(dm);
1687:   DMPlexStratify(dm);
1688:   return(0);
1689: }

1691: /*
1692:   This takes as input the coordinates for each vertex
1693: */
1694: static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[])
1695: {
1696:   PetscSection   coordSection;
1697:   Vec            coordinates;
1698:   DM             cdm;
1699:   PetscScalar   *coords;
1700:   PetscInt       v, d;

1704:   DMGetCoordinateSection(dm, &coordSection);
1705:   PetscSectionSetNumFields(coordSection, 1);
1706:   PetscSectionSetFieldComponents(coordSection, 0, spaceDim);
1707:   PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
1708:   for (v = numCells; v < numCells+numVertices; ++v) {
1709:     PetscSectionSetDof(coordSection, v, spaceDim);
1710:     PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);
1711:   }
1712:   PetscSectionSetUp(coordSection);

1714:   DMGetCoordinateDM(dm, &cdm);
1715:   DMCreateLocalVector(cdm, &coordinates);
1716:   VecSetBlockSize(coordinates, spaceDim);
1717:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1718:   VecGetArray(coordinates, &coords);
1719:   for (v = 0; v < numVertices; ++v) {
1720:     for (d = 0; d < spaceDim; ++d) {
1721:       coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d];
1722:     }
1723:   }
1724:   VecRestoreArray(coordinates, &coords);
1725:   DMSetCoordinatesLocal(dm, coordinates);
1726:   VecDestroy(&coordinates);
1727:   return(0);
1728: }

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

1733:   Input Parameters:
1734: + comm - The communicator
1735: . dim - The topological dimension of the mesh
1736: . numCells - The number of cells
1737: . numVertices - The number of vertices
1738: . numCorners - The number of vertices for each cell
1739: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
1740: . cells - An array of numCells*numCorners numbers, the vertices for each cell
1741: . spaceDim - The spatial dimension used for coordinates
1742: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex

1744:   Output Parameter:
1745: . dm - The DM

1747:   Note: Two triangles sharing a face
1748: $
1749: $        2
1750: $      / | \
1751: $     /  |  \
1752: $    /   |   \
1753: $   0  0 | 1  3
1754: $    \   |   /
1755: $     \  |  /
1756: $      \ | /
1757: $        1
1758: would have input
1759: $  numCells = 2, numVertices = 4
1760: $  cells = [0 1 2  1 3 2]
1761: $
1762: which would result in the DMPlex
1763: $
1764: $        4
1765: $      / | \
1766: $     /  |  \
1767: $    /   |   \
1768: $   2  0 | 1  5
1769: $    \   |   /
1770: $     \  |  /
1771: $      \ | /
1772: $        3

1774:   Level: beginner

1776: .seealso: DMPlexCreateFromDAG(), DMPlexCreate()
1777: @*/
1778: 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)
1779: {

1783:   DMCreate(comm, dm);
1784:   DMSetType(*dm, DMPLEX);
1785:   DMSetDimension(*dm, dim);
1786:   DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);
1787:   if (interpolate) {
1788:     DM idm = NULL;

1790:     DMPlexInterpolate(*dm, &idm);
1791:     DMDestroy(dm);
1792:     *dm  = idm;
1793:   }
1794:   DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);
1795:   return(0);
1796: }

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

1801:   Input Parameters:
1802: + dm - The empty DM object, usually from DMCreate() and DMSetDimension()
1803: . depth - The depth of the DAG
1804: . numPoints - The number of points at each depth
1805: . coneSize - The cone size of each point
1806: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
1807: . coneOrientations - The orientation of each cone point
1808: - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex

1810:   Output Parameter:
1811: . dm - The DM

1813:   Note: Two triangles sharing a face would have input
1814: $  depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
1815: $  cones = [2 3 4  3 5 4], coneOrientations = [0 0 0  0 0 0]
1816: $ vertexCoords = [-1.0 0.0  0.0 -1.0  0.0 1.0  1.0 0.0]
1817: $
1818: which would result in the DMPlex
1819: $
1820: $        4
1821: $      / | \
1822: $     /  |  \
1823: $    /   |   \
1824: $   2  0 | 1  5
1825: $    \   |   /
1826: $     \  |  /
1827: $      \ | /
1828: $        3
1829: $
1830: $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList()

1832:   Level: advanced

1834: .seealso: DMPlexCreateFromCellList(), DMPlexCreate()
1835: @*/
1836: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
1837: {
1838:   Vec            coordinates;
1839:   PetscSection   coordSection;
1840:   PetscScalar    *coords;
1841:   PetscInt       coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;

1845:   DMGetDimension(dm, &dim);
1846:   DMGetCoordinateDim(dm, &dimEmbed);
1847:   if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim);
1848:   for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
1849:   DMPlexSetChart(dm, pStart, pEnd);
1850:   for (p = pStart; p < pEnd; ++p) {
1851:     DMPlexSetConeSize(dm, p, coneSize[p-pStart]);
1852:     if (firstVertex < 0 && !coneSize[p - pStart]) {
1853:       firstVertex = p - pStart;
1854:     }
1855:   }
1856:   if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]);
1857:   DMSetUp(dm); /* Allocate space for cones */
1858:   for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) {
1859:     DMPlexSetCone(dm, p, &cones[off]);
1860:     DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);
1861:   }
1862:   DMPlexSymmetrize(dm);
1863:   DMPlexStratify(dm);
1864:   /* Build coordinates */
1865:   DMGetCoordinateSection(dm, &coordSection);
1866:   PetscSectionSetNumFields(coordSection, 1);
1867:   PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);
1868:   PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);
1869:   for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) {
1870:     PetscSectionSetDof(coordSection, v, dimEmbed);
1871:     PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);
1872:   }
1873:   PetscSectionSetUp(coordSection);
1874:   PetscSectionGetStorageSize(coordSection, &coordSize);
1875:   VecCreate(PETSC_COMM_SELF, &coordinates);
1876:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
1877:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
1878:   VecSetBlockSize(coordinates, dimEmbed);
1879:   VecSetType(coordinates,VECSTANDARD);
1880:   VecGetArray(coordinates, &coords);
1881:   for (v = 0; v < numPoints[0]; ++v) {
1882:     PetscInt off;

1884:     PetscSectionGetOffset(coordSection, v+firstVertex, &off);
1885:     for (d = 0; d < dimEmbed; ++d) {
1886:       coords[off+d] = vertexCoords[v*dimEmbed+d];
1887:     }
1888:   }
1889:   VecRestoreArray(coordinates, &coords);
1890:   DMSetCoordinatesLocal(dm, coordinates);
1891:   VecDestroy(&coordinates);
1892:   return(0);
1893: }

1895: /*@C
1896:   DMPlexCreateFromFile - This takes a filename and produces a DM

1898:   Input Parameters:
1899: + comm - The communicator
1900: . filename - A file name
1901: - interpolate - Flag to create intermediate mesh pieces (edges, faces)

1903:   Output Parameter:
1904: . dm - The DM

1906:   Level: beginner

1908: .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate()
1909: @*/
1910: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1911: {
1912:   const char    *extGmsh   = ".msh";
1913:   const char    *extCGNS   = ".cgns";
1914:   const char    *extExodus = ".exo";
1915:   const char    *extFluent = ".cas";
1916:   const char    *extHDF5   = ".h5";
1917:   const char    *extMed    = ".med";
1918:   size_t         len;
1919:   PetscBool      isGmsh, isCGNS, isExodus, isFluent, isHDF5, isMed;
1920:   PetscMPIInt    rank;

1926:   MPI_Comm_rank(comm, &rank);
1927:   PetscStrlen(filename, &len);
1928:   if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
1929:   PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh,   4, &isGmsh);
1930:   PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS,   5, &isCGNS);
1931:   PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);
1932:   PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);
1933:   PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5,   3, &isHDF5);
1934:   PetscStrncmp(&filename[PetscMax(0,len-4)], extMed,    4, &isMed);
1935:   if (isGmsh) {
1936:     DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);
1937:   } else if (isCGNS) {
1938:     DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);
1939:   } else if (isExodus) {
1940:     DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);
1941:   } else if (isFluent) {
1942:     DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);
1943:   } else if (isHDF5) {
1944:     PetscViewer viewer;

1946:     PetscViewerCreate(comm, &viewer);
1947:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
1948:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1949:     PetscViewerFileSetName(viewer, filename);
1950:     DMCreate(comm, dm);
1951:     DMSetType(*dm, DMPLEX);
1952:     DMLoad(*dm, viewer);
1953:     PetscViewerDestroy(&viewer);
1954:   } else if (isMed) {
1955:     DMPlexCreateMedFromFile(comm, filename, interpolate, dm);
1956:   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
1957:   return(0);
1958: }

1960: /*@
1961:   DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell

1963:   Collective on comm

1965:   Input Parameters:
1966: + comm    - The communicator
1967: . dim     - The spatial dimension
1968: - simplex - Flag for simplex, otherwise use a tensor-product cell

1970:   Output Parameter:
1971: . refdm - The reference cell

1973:   Level: intermediate

1975: .keywords: reference cell
1976: .seealso:
1977: @*/
1978: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm)
1979: {
1980:   DM             rdm;
1981:   Vec            coords;

1985:   DMCreate(comm, &rdm);
1986:   DMSetType(rdm, DMPLEX);
1987:   DMSetDimension(rdm, dim);
1988:   switch (dim) {
1989:   case 0:
1990:   {
1991:     PetscInt    numPoints[1]        = {1};
1992:     PetscInt    coneSize[1]         = {0};
1993:     PetscInt    cones[1]            = {0};
1994:     PetscInt    coneOrientations[1] = {0};
1995:     PetscScalar vertexCoords[1]     = {0.0};

1997:     DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);
1998:   }
1999:   break;
2000:   case 1:
2001:   {
2002:     PetscInt    numPoints[2]        = {2, 1};
2003:     PetscInt    coneSize[3]         = {2, 0, 0};
2004:     PetscInt    cones[2]            = {1, 2};
2005:     PetscInt    coneOrientations[2] = {0, 0};
2006:     PetscScalar vertexCoords[2]     = {-1.0,  1.0};

2008:     DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2009:   }
2010:   break;
2011:   case 2:
2012:     if (simplex) {
2013:       PetscInt    numPoints[2]        = {3, 1};
2014:       PetscInt    coneSize[4]         = {3, 0, 0, 0};
2015:       PetscInt    cones[3]            = {1, 2, 3};
2016:       PetscInt    coneOrientations[3] = {0, 0, 0};
2017:       PetscScalar vertexCoords[6]     = {-1.0, -1.0,  1.0, -1.0,  -1.0, 1.0};

2019:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2020:     } else {
2021:       PetscInt    numPoints[2]        = {4, 1};
2022:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2023:       PetscInt    cones[4]            = {1, 2, 3, 4};
2024:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2025:       PetscScalar vertexCoords[8]     = {-1.0, -1.0,  1.0, -1.0,  1.0, 1.0,  -1.0, 1.0};

2027:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2028:     }
2029:   break;
2030:   case 3:
2031:     if (simplex) {
2032:       PetscInt    numPoints[2]        = {4, 1};
2033:       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
2034:       PetscInt    cones[4]            = {1, 3, 2, 4};
2035:       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
2036:       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};

2038:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2039:     } else {
2040:       PetscInt    numPoints[2]        = {8, 1};
2041:       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
2042:       PetscInt    cones[8]            = {1, 4, 3, 2, 5, 6, 7, 8};
2043:       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
2044:       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,
2045:                                          -1.0, -1.0,  1.0,  1.0, -1.0,  1.0,  1.0, 1.0,  1.0,  -1.0, 1.0,  1.0};

2047:       DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);
2048:     }
2049:   break;
2050:   default:
2051:     SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim);
2052:   }
2053:   *refdm = NULL;
2054:   DMPlexInterpolate(rdm, refdm);
2055:   if (rdm->coordinateDM) {
2056:     DM           ncdm;
2057:     PetscSection cs;
2058:     PetscInt     pEnd = -1;

2060:     DMGetDefaultSection(rdm->coordinateDM, &cs);
2061:     if (cs) {PetscSectionGetChart(cs, NULL, &pEnd);}
2062:     if (pEnd >= 0) {
2063:       DMClone(rdm->coordinateDM, &ncdm);
2064:       DMSetDefaultSection(ncdm, cs);
2065:       DMSetCoordinateDM(*refdm, ncdm);
2066:       DMDestroy(&ncdm);
2067:     }
2068:   }
2069:   DMGetCoordinatesLocal(rdm, &coords);
2070:   if (coords) {
2071:     DMSetCoordinatesLocal(*refdm, coords);
2072:   } else {
2073:     DMGetCoordinates(rdm, &coords);
2074:     if (coords) {DMSetCoordinates(*refdm, coords);}
2075:   }
2076:   DMDestroy(&rdm);
2077:   return(0);
2078: }