Actual source code: plexcreate.c
1: #define PETSCDM_DLL
2: #include <petsc/private/dmpleximpl.h>
3: #include <petsc/private/hashseti.h>
4: #include <petscsf.h>
5: #include <petscdmplextransform.h>
6: #include <petscdmlabelephemeral.h>
7: #include <petsc/private/kernels/blockmatmult.h>
8: #include <petsc/private/kernels/blockinvert.h>
10: PetscLogEvent DMPLEX_CreateFromFile, DMPLEX_CreateFromOptions, DMPLEX_BuildFromCellList, DMPLEX_BuildCoordinatesFromCellList;
12: /* External function declarations here */
13: static PetscErrorCode DMInitialize_Plex(DM dm);
15: /* This copies internal things in the Plex structure that we generally want when making a new, related Plex */
16: PetscErrorCode DMPlexCopy_Internal(DM dmin, PetscBool copyPeriodicity, PetscBool copyOverlap, DM dmout)
17: {
18: const PetscReal *maxCell, *Lstart, *L;
19: PetscBool dist;
20: DMPlexReorderDefaultFlag reorder;
22: PetscFunctionBegin;
23: if (copyPeriodicity) {
24: PetscCall(DMGetPeriodicity(dmin, &maxCell, &Lstart, &L));
25: PetscCall(DMSetPeriodicity(dmout, maxCell, Lstart, L));
26: }
27: PetscCall(DMPlexDistributeGetDefault(dmin, &dist));
28: PetscCall(DMPlexDistributeSetDefault(dmout, dist));
29: PetscCall(DMPlexReorderGetDefault(dmin, &reorder));
30: PetscCall(DMPlexReorderSetDefault(dmout, reorder));
31: ((DM_Plex *)dmout->data)->useHashLocation = ((DM_Plex *)dmin->data)->useHashLocation;
32: if (copyOverlap) PetscCall(DMPlexSetOverlap_Plex(dmout, dmin, 0));
33: PetscFunctionReturn(PETSC_SUCCESS);
34: }
36: /* Replace dm with the contents of ndm, and then destroy ndm
37: - Share the DM_Plex structure
38: - Share the coordinates
39: - Share the SF
40: */
41: PetscErrorCode DMPlexReplace_Internal(DM dm, DM *ndm)
42: {
43: PetscSF sf;
44: DM dmNew = *ndm, coordDM, coarseDM;
45: Vec coords;
46: const PetscReal *maxCell, *Lstart, *L;
47: PetscInt dim, cdim;
49: PetscFunctionBegin;
50: if (dm == dmNew) {
51: PetscCall(DMDestroy(ndm));
52: PetscFunctionReturn(PETSC_SUCCESS);
53: }
54: dm->setupcalled = dmNew->setupcalled;
55: PetscCall(DMGetDimension(dmNew, &dim));
56: PetscCall(DMSetDimension(dm, dim));
57: PetscCall(DMGetCoordinateDim(dmNew, &cdim));
58: PetscCall(DMSetCoordinateDim(dm, cdim));
59: PetscCall(DMGetPointSF(dmNew, &sf));
60: PetscCall(DMSetPointSF(dm, sf));
61: PetscCall(DMGetCoordinateDM(dmNew, &coordDM));
62: PetscCall(DMGetCoordinatesLocal(dmNew, &coords));
63: PetscCall(DMSetCoordinateDM(dm, coordDM));
64: PetscCall(DMSetCoordinatesLocal(dm, coords));
65: PetscCall(DMGetCellCoordinateDM(dmNew, &coordDM));
66: PetscCall(DMGetCellCoordinatesLocal(dmNew, &coords));
67: PetscCall(DMSetCellCoordinateDM(dm, coordDM));
68: PetscCall(DMSetCellCoordinatesLocal(dm, coords));
69: /* Do not want to create the coordinate field if it does not already exist, so do not call DMGetCoordinateField() */
70: PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
71: dm->coordinates[0].field = dmNew->coordinates[0].field;
72: ((DM_Plex *)dmNew->data)->coordFunc = ((DM_Plex *)dm->data)->coordFunc;
73: PetscCall(DMGetPeriodicity(dmNew, &maxCell, &Lstart, &L));
74: PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L));
75: PetscCall(DMDestroy_Plex(dm));
76: PetscCall(DMInitialize_Plex(dm));
77: dm->data = dmNew->data;
78: ((DM_Plex *)dmNew->data)->refct++;
79: PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &sf));
80: PetscCall(DMPlexSetIsoperiodicFaceSF(dm, sf)); // for the compose function effect on dm
81: PetscCall(DMDestroyLabelLinkList_Internal(dm));
82: PetscCall(DMCopyLabels(dmNew, dm, PETSC_OWN_POINTER, PETSC_TRUE, DM_COPY_LABELS_FAIL));
83: PetscCall(DMGetCoarseDM(dmNew, &coarseDM));
84: PetscCall(DMSetCoarseDM(dm, coarseDM));
85: PetscCall(DMDestroy(ndm));
86: PetscFunctionReturn(PETSC_SUCCESS);
87: }
89: /* Swap dm with the contents of dmNew
90: - Swap the DM_Plex structure
91: - Swap the coordinates
92: - Swap the point PetscSF
93: */
94: static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB)
95: {
96: DM coordDMA, coordDMB;
97: Vec coordsA, coordsB;
98: PetscSF sfA, sfB;
99: DMField fieldTmp;
100: void *tmp;
101: DMLabelLink listTmp;
102: DMLabel depthTmp;
103: PetscInt tmpI;
105: PetscFunctionBegin;
106: if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS);
107: PetscCall(DMGetPointSF(dmA, &sfA));
108: PetscCall(DMGetPointSF(dmB, &sfB));
109: PetscCall(PetscObjectReference((PetscObject)sfA));
110: PetscCall(DMSetPointSF(dmA, sfB));
111: PetscCall(DMSetPointSF(dmB, sfA));
112: PetscCall(PetscObjectDereference((PetscObject)sfA));
114: PetscCall(DMGetCoordinateDM(dmA, &coordDMA));
115: PetscCall(DMGetCoordinateDM(dmB, &coordDMB));
116: PetscCall(PetscObjectReference((PetscObject)coordDMA));
117: PetscCall(DMSetCoordinateDM(dmA, coordDMB));
118: PetscCall(DMSetCoordinateDM(dmB, coordDMA));
119: PetscCall(PetscObjectDereference((PetscObject)coordDMA));
121: PetscCall(DMGetCoordinatesLocal(dmA, &coordsA));
122: PetscCall(DMGetCoordinatesLocal(dmB, &coordsB));
123: PetscCall(PetscObjectReference((PetscObject)coordsA));
124: PetscCall(DMSetCoordinatesLocal(dmA, coordsB));
125: PetscCall(DMSetCoordinatesLocal(dmB, coordsA));
126: PetscCall(PetscObjectDereference((PetscObject)coordsA));
128: PetscCall(DMGetCellCoordinateDM(dmA, &coordDMA));
129: PetscCall(DMGetCellCoordinateDM(dmB, &coordDMB));
130: PetscCall(PetscObjectReference((PetscObject)coordDMA));
131: PetscCall(DMSetCellCoordinateDM(dmA, coordDMB));
132: PetscCall(DMSetCellCoordinateDM(dmB, coordDMA));
133: PetscCall(PetscObjectDereference((PetscObject)coordDMA));
135: PetscCall(DMGetCellCoordinatesLocal(dmA, &coordsA));
136: PetscCall(DMGetCellCoordinatesLocal(dmB, &coordsB));
137: PetscCall(PetscObjectReference((PetscObject)coordsA));
138: PetscCall(DMSetCellCoordinatesLocal(dmA, coordsB));
139: PetscCall(DMSetCellCoordinatesLocal(dmB, coordsA));
140: PetscCall(PetscObjectDereference((PetscObject)coordsA));
142: fieldTmp = dmA->coordinates[0].field;
143: dmA->coordinates[0].field = dmB->coordinates[0].field;
144: dmB->coordinates[0].field = fieldTmp;
145: fieldTmp = dmA->coordinates[1].field;
146: dmA->coordinates[1].field = dmB->coordinates[1].field;
147: dmB->coordinates[1].field = fieldTmp;
148: tmp = dmA->data;
149: dmA->data = dmB->data;
150: dmB->data = tmp;
151: listTmp = dmA->labels;
152: dmA->labels = dmB->labels;
153: dmB->labels = listTmp;
154: depthTmp = dmA->depthLabel;
155: dmA->depthLabel = dmB->depthLabel;
156: dmB->depthLabel = depthTmp;
157: depthTmp = dmA->celltypeLabel;
158: dmA->celltypeLabel = dmB->celltypeLabel;
159: dmB->celltypeLabel = depthTmp;
160: tmpI = dmA->levelup;
161: dmA->levelup = dmB->levelup;
162: dmB->levelup = tmpI;
163: PetscFunctionReturn(PETSC_SUCCESS);
164: }
166: PetscErrorCode DMPlexInterpolateInPlace_Internal(DM dm)
167: {
168: DM idm;
170: PetscFunctionBegin;
171: PetscCall(DMPlexInterpolate(dm, &idm));
172: PetscCall(DMPlexCopyCoordinates(dm, idm));
173: PetscCall(DMPlexReplace_Internal(dm, &idm));
174: PetscFunctionReturn(PETSC_SUCCESS);
175: }
177: /*@C
178: DMPlexCreateCoordinateSpace - Creates a finite element space for the coordinates
180: Collective
182: Input Parameters:
183: + dm - The `DMPLEX`
184: . degree - The degree of the finite element or `PETSC_DECIDE`
185: - coordFunc - An optional function to map new points from refinement to the surface
187: Level: advanced
189: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscPointFunc`, `PetscFECreateLagrange()`, `DMGetCoordinateDM()`
190: @*/
191: PetscErrorCode DMPlexCreateCoordinateSpace(DM dm, PetscInt degree, PetscPointFunc coordFunc)
192: {
193: DM_Plex *mesh = (DM_Plex *)dm->data;
194: DM cdm;
195: PetscDS cds;
196: PetscFE fe;
197: PetscClassId id;
199: PetscFunctionBegin;
200: PetscCall(DMGetCoordinateDM(dm, &cdm));
201: PetscCall(DMGetDS(cdm, &cds));
202: PetscCall(PetscDSGetDiscretization(cds, 0, (PetscObject *)&fe));
203: PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
204: if (id != PETSCFE_CLASSID) {
205: PetscInt dim, dE, qorder;
207: PetscCall(DMGetDimension(dm, &dim));
208: PetscCall(DMGetCoordinateDim(dm, &dE));
209: qorder = degree;
210: PetscObjectOptionsBegin((PetscObject)cdm);
211: PetscCall(PetscOptionsBoundedInt("-default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "DMPlexCreateCoordinateSpace", qorder, &qorder, NULL, 0));
212: PetscOptionsEnd();
213: if (degree == PETSC_DECIDE) fe = NULL;
214: else {
215: DMPolytopeType ct;
216: PetscInt cStart, cEnd, ctTmp;
218: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
219: // Want to match cell types
220: if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct));
221: else ct = DM_POLYTOPE_UNKNOWN;
222: ctTmp = (PetscInt)ct;
223: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &ctTmp, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
224: ct = (DMPolytopeType)ctTmp;
225: // Work around current bug
226: if (ct == DM_POLYTOPE_SEG_PRISM_TENSOR || ct == DM_POLYTOPE_TRI_PRISM_TENSOR || ct == DM_POLYTOPE_QUAD_PRISM_TENSOR) fe = NULL;
227: else PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, degree, qorder, &fe));
228: }
229: if (fe) PetscCall(DMProjectCoordinates(dm, fe));
230: PetscCall(PetscFEDestroy(&fe));
231: }
232: mesh->coordFunc = coordFunc;
233: PetscFunctionReturn(PETSC_SUCCESS);
234: }
236: /*@
237: DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement.
239: Collective
241: Input Parameters:
242: + comm - The communicator for the `DM` object
243: . dim - The spatial dimension
244: . simplex - Flag for simplicial cells, otherwise they are tensor product cells
245: . interpolate - Flag to create intermediate mesh pieces (edges, faces)
246: - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
248: Output Parameter:
249: . newdm - The `DM` object
251: Level: beginner
253: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMSetType()`, `DMCreate()`
254: @*/
255: PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscReal refinementLimit, DM *newdm)
256: {
257: DM dm;
258: PetscMPIInt rank;
260: PetscFunctionBegin;
261: PetscCall(DMCreate(comm, &dm));
262: PetscCall(DMSetType(dm, DMPLEX));
263: PetscCall(DMSetDimension(dm, dim));
264: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
265: PetscCallMPI(MPI_Comm_rank(comm, &rank));
266: switch (dim) {
267: case 2:
268: if (simplex) PetscCall(PetscObjectSetName((PetscObject)dm, "triangular"));
269: else PetscCall(PetscObjectSetName((PetscObject)dm, "quadrilateral"));
270: break;
271: case 3:
272: if (simplex) PetscCall(PetscObjectSetName((PetscObject)dm, "tetrahedral"));
273: else PetscCall(PetscObjectSetName((PetscObject)dm, "hexahedral"));
274: break;
275: default:
276: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
277: }
278: if (rank) {
279: PetscInt numPoints[2] = {0, 0};
280: PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL));
281: } else {
282: switch (dim) {
283: case 2:
284: if (simplex) {
285: PetscInt numPoints[2] = {4, 2};
286: PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0};
287: PetscInt cones[6] = {2, 3, 4, 5, 4, 3};
288: PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0};
289: PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
291: PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
292: } else {
293: PetscInt numPoints[2] = {6, 2};
294: PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0};
295: PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4};
296: PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
297: 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};
299: PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
300: }
301: break;
302: case 3:
303: if (simplex) {
304: PetscInt numPoints[2] = {5, 2};
305: PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0};
306: PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6};
307: PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
308: 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};
310: PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
311: } else {
312: PetscInt numPoints[2] = {12, 2};
313: PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
314: PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8};
315: PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
316: 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, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5};
318: PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
319: }
320: break;
321: default:
322: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
323: }
324: }
325: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
326: *newdm = dm;
327: if (refinementLimit > 0.0) {
328: DM rdm;
329: const char *name;
331: PetscCall(DMPlexSetRefinementUniform(*newdm, PETSC_FALSE));
332: PetscCall(DMPlexSetRefinementLimit(*newdm, refinementLimit));
333: PetscCall(DMRefine(*newdm, comm, &rdm));
334: PetscCall(PetscObjectGetName((PetscObject)*newdm, &name));
335: PetscCall(PetscObjectSetName((PetscObject)rdm, name));
336: PetscCall(DMDestroy(newdm));
337: *newdm = rdm;
338: }
339: if (interpolate) {
340: DM idm;
342: PetscCall(DMPlexInterpolate(*newdm, &idm));
343: PetscCall(DMDestroy(newdm));
344: *newdm = idm;
345: }
346: PetscFunctionReturn(PETSC_SUCCESS);
347: }
349: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_1D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
350: {
351: const PetscInt numVertices = 2;
352: PetscInt markerRight = 1;
353: PetscInt markerLeft = 1;
354: PetscBool markerSeparate = PETSC_FALSE;
355: Vec coordinates;
356: PetscSection coordSection;
357: PetscScalar *coords;
358: PetscInt coordSize;
359: PetscMPIInt rank;
360: PetscInt cdim = 1, v;
362: PetscFunctionBegin;
363: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL));
364: if (markerSeparate) {
365: markerRight = 2;
366: markerLeft = 1;
367: }
368: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
369: if (rank == 0) {
370: PetscCall(DMPlexSetChart(dm, 0, numVertices));
371: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
372: PetscCall(DMSetLabelValue(dm, "marker", 0, markerLeft));
373: PetscCall(DMSetLabelValue(dm, "marker", 1, markerRight));
374: }
375: PetscCall(DMPlexSymmetrize(dm));
376: PetscCall(DMPlexStratify(dm));
377: /* Build coordinates */
378: PetscCall(DMSetCoordinateDim(dm, cdim));
379: PetscCall(DMGetCoordinateSection(dm, &coordSection));
380: PetscCall(PetscSectionSetNumFields(coordSection, 1));
381: PetscCall(PetscSectionSetChart(coordSection, 0, numVertices));
382: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim));
383: for (v = 0; v < numVertices; ++v) {
384: PetscCall(PetscSectionSetDof(coordSection, v, cdim));
385: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim));
386: }
387: PetscCall(PetscSectionSetUp(coordSection));
388: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
389: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
390: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
391: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
392: PetscCall(VecSetBlockSize(coordinates, cdim));
393: PetscCall(VecSetType(coordinates, VECSTANDARD));
394: PetscCall(VecGetArray(coordinates, &coords));
395: coords[0] = lower[0];
396: coords[1] = upper[0];
397: PetscCall(VecRestoreArray(coordinates, &coords));
398: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
399: PetscCall(VecDestroy(&coordinates));
400: PetscFunctionReturn(PETSC_SUCCESS);
401: }
403: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_2D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[])
404: {
405: const PetscInt numVertices = (edges[0] + 1) * (edges[1] + 1);
406: const PetscInt numEdges = edges[0] * (edges[1] + 1) + (edges[0] + 1) * edges[1];
407: PetscInt markerTop = 1;
408: PetscInt markerBottom = 1;
409: PetscInt markerRight = 1;
410: PetscInt markerLeft = 1;
411: PetscBool markerSeparate = PETSC_FALSE;
412: Vec coordinates;
413: PetscSection coordSection;
414: PetscScalar *coords;
415: PetscInt coordSize;
416: PetscMPIInt rank;
417: PetscInt v, vx, vy;
419: PetscFunctionBegin;
420: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL));
421: if (markerSeparate) {
422: markerTop = 3;
423: markerBottom = 1;
424: markerRight = 2;
425: markerLeft = 4;
426: }
427: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
428: if (rank == 0) {
429: PetscInt e, ex, ey;
431: PetscCall(DMPlexSetChart(dm, 0, numEdges + numVertices));
432: for (e = 0; e < numEdges; ++e) PetscCall(DMPlexSetConeSize(dm, e, 2));
433: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
434: for (vx = 0; vx <= edges[0]; vx++) {
435: for (ey = 0; ey < edges[1]; ey++) {
436: PetscInt edge = vx * edges[1] + ey + edges[0] * (edges[1] + 1);
437: PetscInt vertex = ey * (edges[0] + 1) + vx + numEdges;
438: PetscInt cone[2];
440: cone[0] = vertex;
441: cone[1] = vertex + edges[0] + 1;
442: PetscCall(DMPlexSetCone(dm, edge, cone));
443: if (vx == edges[0]) {
444: PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight));
445: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight));
446: if (ey == edges[1] - 1) {
447: PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight));
448: PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerRight));
449: }
450: } else if (vx == 0) {
451: PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft));
452: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft));
453: if (ey == edges[1] - 1) {
454: PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft));
455: PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft));
456: }
457: }
458: }
459: }
460: for (vy = 0; vy <= edges[1]; vy++) {
461: for (ex = 0; ex < edges[0]; ex++) {
462: PetscInt edge = vy * edges[0] + ex;
463: PetscInt vertex = vy * (edges[0] + 1) + ex + numEdges;
464: PetscInt cone[2];
466: cone[0] = vertex;
467: cone[1] = vertex + 1;
468: PetscCall(DMPlexSetCone(dm, edge, cone));
469: if (vy == edges[1]) {
470: PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop));
471: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop));
472: if (ex == edges[0] - 1) {
473: PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop));
474: PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerTop));
475: }
476: } else if (vy == 0) {
477: PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom));
478: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom));
479: if (ex == edges[0] - 1) {
480: PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom));
481: PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom));
482: }
483: }
484: }
485: }
486: }
487: PetscCall(DMPlexSymmetrize(dm));
488: PetscCall(DMPlexStratify(dm));
489: /* Build coordinates */
490: PetscCall(DMSetCoordinateDim(dm, 2));
491: PetscCall(DMGetCoordinateSection(dm, &coordSection));
492: PetscCall(PetscSectionSetNumFields(coordSection, 1));
493: PetscCall(PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices));
494: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, 2));
495: for (v = numEdges; v < numEdges + numVertices; ++v) {
496: PetscCall(PetscSectionSetDof(coordSection, v, 2));
497: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, 2));
498: }
499: PetscCall(PetscSectionSetUp(coordSection));
500: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
501: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
502: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
503: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
504: PetscCall(VecSetBlockSize(coordinates, 2));
505: PetscCall(VecSetType(coordinates, VECSTANDARD));
506: PetscCall(VecGetArray(coordinates, &coords));
507: for (vy = 0; vy <= edges[1]; ++vy) {
508: for (vx = 0; vx <= edges[0]; ++vx) {
509: coords[(vy * (edges[0] + 1) + vx) * 2 + 0] = lower[0] + ((upper[0] - lower[0]) / edges[0]) * vx;
510: coords[(vy * (edges[0] + 1) + vx) * 2 + 1] = lower[1] + ((upper[1] - lower[1]) / edges[1]) * vy;
511: }
512: }
513: PetscCall(VecRestoreArray(coordinates, &coords));
514: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
515: PetscCall(VecDestroy(&coordinates));
516: PetscFunctionReturn(PETSC_SUCCESS);
517: }
519: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_3D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[])
520: {
521: PetscInt vertices[3], numVertices;
522: PetscInt numFaces = 2 * faces[0] * faces[1] + 2 * faces[1] * faces[2] + 2 * faces[0] * faces[2];
523: PetscInt markerTop = 1;
524: PetscInt markerBottom = 1;
525: PetscInt markerFront = 1;
526: PetscInt markerBack = 1;
527: PetscInt markerRight = 1;
528: PetscInt markerLeft = 1;
529: PetscBool markerSeparate = PETSC_FALSE;
530: Vec coordinates;
531: PetscSection coordSection;
532: PetscScalar *coords;
533: PetscInt coordSize;
534: PetscMPIInt rank;
535: PetscInt v, vx, vy, vz;
536: PetscInt voffset, iface = 0, cone[4];
538: PetscFunctionBegin;
539: PetscCheck(faces[0] >= 1 && faces[1] >= 1 && faces[2] >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side");
540: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
541: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL));
542: if (markerSeparate) {
543: markerBottom = 1;
544: markerTop = 2;
545: markerFront = 3;
546: markerBack = 4;
547: markerRight = 5;
548: markerLeft = 6;
549: }
550: vertices[0] = faces[0] + 1;
551: vertices[1] = faces[1] + 1;
552: vertices[2] = faces[2] + 1;
553: numVertices = vertices[0] * vertices[1] * vertices[2];
554: if (rank == 0) {
555: PetscInt f;
557: PetscCall(DMPlexSetChart(dm, 0, numFaces + numVertices));
558: for (f = 0; f < numFaces; ++f) PetscCall(DMPlexSetConeSize(dm, f, 4));
559: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
561: /* Side 0 (Top) */
562: for (vy = 0; vy < faces[1]; vy++) {
563: for (vx = 0; vx < faces[0]; vx++) {
564: voffset = numFaces + vertices[0] * vertices[1] * (vertices[2] - 1) + vy * vertices[0] + vx;
565: cone[0] = voffset;
566: cone[1] = voffset + 1;
567: cone[2] = voffset + vertices[0] + 1;
568: cone[3] = voffset + vertices[0];
569: PetscCall(DMPlexSetCone(dm, iface, cone));
570: PetscCall(DMSetLabelValue(dm, "marker", iface, markerTop));
571: PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerTop));
572: PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerTop));
573: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerTop));
574: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 1, markerTop));
575: iface++;
576: }
577: }
579: /* Side 1 (Bottom) */
580: for (vy = 0; vy < faces[1]; vy++) {
581: for (vx = 0; vx < faces[0]; vx++) {
582: voffset = numFaces + vy * (faces[0] + 1) + vx;
583: cone[0] = voffset + 1;
584: cone[1] = voffset;
585: cone[2] = voffset + vertices[0];
586: cone[3] = voffset + vertices[0] + 1;
587: PetscCall(DMPlexSetCone(dm, iface, cone));
588: PetscCall(DMSetLabelValue(dm, "marker", iface, markerBottom));
589: PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerBottom));
590: PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerBottom));
591: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerBottom));
592: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 1, markerBottom));
593: iface++;
594: }
595: }
597: /* Side 2 (Front) */
598: for (vz = 0; vz < faces[2]; vz++) {
599: for (vx = 0; vx < faces[0]; vx++) {
600: voffset = numFaces + vz * vertices[0] * vertices[1] + vx;
601: cone[0] = voffset;
602: cone[1] = voffset + 1;
603: cone[2] = voffset + vertices[0] * vertices[1] + 1;
604: cone[3] = voffset + vertices[0] * vertices[1];
605: PetscCall(DMPlexSetCone(dm, iface, cone));
606: PetscCall(DMSetLabelValue(dm, "marker", iface, markerFront));
607: PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerFront));
608: PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerFront));
609: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 0, markerFront));
610: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 1, markerFront));
611: iface++;
612: }
613: }
615: /* Side 3 (Back) */
616: for (vz = 0; vz < faces[2]; vz++) {
617: for (vx = 0; vx < faces[0]; vx++) {
618: voffset = numFaces + vz * vertices[0] * vertices[1] + vertices[0] * (vertices[1] - 1) + vx;
619: cone[0] = voffset + vertices[0] * vertices[1];
620: cone[1] = voffset + vertices[0] * vertices[1] + 1;
621: cone[2] = voffset + 1;
622: cone[3] = voffset;
623: PetscCall(DMPlexSetCone(dm, iface, cone));
624: PetscCall(DMSetLabelValue(dm, "marker", iface, markerBack));
625: PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerBack));
626: PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerBack));
627: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 0, markerBack));
628: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 1, markerBack));
629: iface++;
630: }
631: }
633: /* Side 4 (Left) */
634: for (vz = 0; vz < faces[2]; vz++) {
635: for (vy = 0; vy < faces[1]; vy++) {
636: voffset = numFaces + vz * vertices[0] * vertices[1] + vy * vertices[0];
637: cone[0] = voffset;
638: cone[1] = voffset + vertices[0] * vertices[1];
639: cone[2] = voffset + vertices[0] * vertices[1] + vertices[0];
640: cone[3] = voffset + vertices[0];
641: PetscCall(DMPlexSetCone(dm, iface, cone));
642: PetscCall(DMSetLabelValue(dm, "marker", iface, markerLeft));
643: PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerLeft));
644: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerLeft));
645: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[1] + 0, markerLeft));
646: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + vertices[0], markerLeft));
647: iface++;
648: }
649: }
651: /* Side 5 (Right) */
652: for (vz = 0; vz < faces[2]; vz++) {
653: for (vy = 0; vy < faces[1]; vy++) {
654: voffset = numFaces + vz * vertices[0] * vertices[1] + vy * vertices[0] + faces[0];
655: cone[0] = voffset + vertices[0] * vertices[1];
656: cone[1] = voffset;
657: cone[2] = voffset + vertices[0];
658: cone[3] = voffset + vertices[0] * vertices[1] + vertices[0];
659: PetscCall(DMPlexSetCone(dm, iface, cone));
660: PetscCall(DMSetLabelValue(dm, "marker", iface, markerRight));
661: PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerRight));
662: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerRight));
663: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 0, markerRight));
664: PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + vertices[0], markerRight));
665: iface++;
666: }
667: }
668: }
669: PetscCall(DMPlexSymmetrize(dm));
670: PetscCall(DMPlexStratify(dm));
671: /* Build coordinates */
672: PetscCall(DMSetCoordinateDim(dm, 3));
673: PetscCall(DMGetCoordinateSection(dm, &coordSection));
674: PetscCall(PetscSectionSetNumFields(coordSection, 1));
675: PetscCall(PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices));
676: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, 3));
677: for (v = numFaces; v < numFaces + numVertices; ++v) {
678: PetscCall(PetscSectionSetDof(coordSection, v, 3));
679: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, 3));
680: }
681: PetscCall(PetscSectionSetUp(coordSection));
682: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
683: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
684: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
685: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
686: PetscCall(VecSetBlockSize(coordinates, 3));
687: PetscCall(VecSetType(coordinates, VECSTANDARD));
688: PetscCall(VecGetArray(coordinates, &coords));
689: for (vz = 0; vz <= faces[2]; ++vz) {
690: for (vy = 0; vy <= faces[1]; ++vy) {
691: for (vx = 0; vx <= faces[0]; ++vx) {
692: coords[((vz * (faces[1] + 1) + vy) * (faces[0] + 1) + vx) * 3 + 0] = lower[0] + ((upper[0] - lower[0]) / faces[0]) * vx;
693: coords[((vz * (faces[1] + 1) + vy) * (faces[0] + 1) + vx) * 3 + 1] = lower[1] + ((upper[1] - lower[1]) / faces[1]) * vy;
694: coords[((vz * (faces[1] + 1) + vy) * (faces[0] + 1) + vx) * 3 + 2] = lower[2] + ((upper[2] - lower[2]) / faces[2]) * vz;
695: }
696: }
697: }
698: PetscCall(VecRestoreArray(coordinates, &coords));
699: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
700: PetscCall(VecDestroy(&coordinates));
701: PetscFunctionReturn(PETSC_SUCCESS);
702: }
704: static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], PetscBool interpolate)
705: {
706: PetscFunctionBegin;
708: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
709: PetscCall(DMSetDimension(dm, dim - 1));
710: PetscCall(DMSetCoordinateDim(dm, dim));
711: switch (dim) {
712: case 1:
713: PetscCall(DMPlexCreateBoxSurfaceMesh_Tensor_1D_Internal(dm, lower, upper, faces));
714: break;
715: case 2:
716: PetscCall(DMPlexCreateBoxSurfaceMesh_Tensor_2D_Internal(dm, lower, upper, faces));
717: break;
718: case 3:
719: PetscCall(DMPlexCreateBoxSurfaceMesh_Tensor_3D_Internal(dm, lower, upper, faces));
720: break;
721: default:
722: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Dimension not supported: %" PetscInt_FMT, dim);
723: }
724: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
725: if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(dm));
726: PetscFunctionReturn(PETSC_SUCCESS);
727: }
729: /*@C
730: DMPlexCreateBoxSurfaceMesh - Creates a mesh on the surface of the tensor product of unit intervals (box) using tensor cells (hexahedra).
732: Collective
734: Input Parameters:
735: + comm - The communicator for the `DM` object
736: . dim - The spatial dimension of the box, so the resulting mesh is has dimension `dim`-1
737: . faces - Number of faces per dimension, or `NULL` for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
738: . lower - The lower left corner, or `NULL` for (0, 0, 0)
739: . upper - The upper right corner, or `NULL` for (1, 1, 1)
740: - interpolate - Flag to create intermediate mesh pieces (edges, faces)
742: Output Parameter:
743: . dm - The `DM` object
745: Level: beginner
747: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMSetFromOptions()`, `DMPlexCreateBoxMesh()`, `DMPlexCreateFromFile()`, `DMSetType()`, `DMCreate()`
748: @*/
749: PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], PetscBool interpolate, DM *dm)
750: {
751: PetscInt fac[3] = {1, 1, 1};
752: PetscReal low[3] = {0, 0, 0};
753: PetscReal upp[3] = {1, 1, 1};
755: PetscFunctionBegin;
756: PetscCall(DMCreate(comm, dm));
757: PetscCall(DMSetType(*dm, DMPLEX));
758: PetscCall(DMPlexCreateBoxSurfaceMesh_Internal(*dm, dim, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, interpolate));
759: PetscFunctionReturn(PETSC_SUCCESS);
760: }
762: static PetscErrorCode DMPlexCreateLineMesh_Internal(DM dm, PetscInt segments, PetscReal lower, PetscReal upper, DMBoundaryType bd)
763: {
764: PetscInt i, fStart, fEnd, numCells = 0, numVerts = 0;
765: PetscInt numPoints[2], *coneSize, *cones, *coneOrientations;
766: PetscScalar *vertexCoords;
767: PetscReal L, maxCell;
768: PetscBool markerSeparate = PETSC_FALSE;
769: PetscInt markerLeft = 1, faceMarkerLeft = 1;
770: PetscInt markerRight = 1, faceMarkerRight = 2;
771: PetscBool wrap = (bd == DM_BOUNDARY_PERIODIC || bd == DM_BOUNDARY_TWIST) ? PETSC_TRUE : PETSC_FALSE;
772: PetscMPIInt rank;
774: PetscFunctionBegin;
775: PetscAssertPointer(dm, 1);
777: PetscCall(DMSetDimension(dm, 1));
778: PetscCall(DMCreateLabel(dm, "marker"));
779: PetscCall(DMCreateLabel(dm, "Face Sets"));
781: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
782: if (rank == 0) numCells = segments;
783: if (rank == 0) numVerts = segments + (wrap ? 0 : 1);
785: numPoints[0] = numVerts;
786: numPoints[1] = numCells;
787: PetscCall(PetscMalloc4(numCells + numVerts, &coneSize, numCells * 2, &cones, numCells + numVerts, &coneOrientations, numVerts, &vertexCoords));
788: PetscCall(PetscArrayzero(coneOrientations, numCells + numVerts));
789: for (i = 0; i < numCells; ++i) coneSize[i] = 2;
790: for (i = 0; i < numVerts; ++i) coneSize[numCells + i] = 0;
791: for (i = 0; i < numCells; ++i) {
792: cones[2 * i] = numCells + i % numVerts;
793: cones[2 * i + 1] = numCells + (i + 1) % numVerts;
794: }
795: for (i = 0; i < numVerts; ++i) vertexCoords[i] = lower + (upper - lower) * ((PetscReal)i / (PetscReal)numCells);
796: PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
797: PetscCall(PetscFree4(coneSize, cones, coneOrientations, vertexCoords));
799: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL));
800: if (markerSeparate) {
801: markerLeft = faceMarkerLeft;
802: markerRight = faceMarkerRight;
803: }
804: if (!wrap && rank == 0) {
805: PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
806: PetscCall(DMSetLabelValue(dm, "marker", fStart, markerLeft));
807: PetscCall(DMSetLabelValue(dm, "marker", fEnd - 1, markerRight));
808: PetscCall(DMSetLabelValue(dm, "Face Sets", fStart, faceMarkerLeft));
809: PetscCall(DMSetLabelValue(dm, "Face Sets", fEnd - 1, faceMarkerRight));
810: }
811: if (wrap) {
812: L = upper - lower;
813: maxCell = (PetscReal)1.1 * (L / (PetscReal)PetscMax(1, segments));
814: PetscCall(DMSetPeriodicity(dm, &maxCell, &lower, &L));
815: }
816: PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
817: PetscFunctionReturn(PETSC_SUCCESS);
818: }
820: static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
821: {
822: DM boundary, vol;
823: DMLabel bdlabel;
825: PetscFunctionBegin;
826: PetscAssertPointer(dm, 1);
827: for (PetscInt i = 0; i < dim; ++i) PetscCheck(periodicity[i] == DM_BOUNDARY_NONE, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Periodicity is not supported for simplex meshes");
828: PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &boundary));
829: PetscCall(DMSetType(boundary, DMPLEX));
830: PetscCall(DMPlexCreateBoxSurfaceMesh_Internal(boundary, dim, faces, lower, upper, PETSC_FALSE));
831: PetscCall(DMPlexGenerate(boundary, NULL, interpolate, &vol));
832: PetscCall(DMGetLabel(vol, "marker", &bdlabel));
833: if (bdlabel) PetscCall(DMPlexLabelComplete(vol, bdlabel));
834: PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_FALSE, vol));
835: PetscCall(DMPlexReplace_Internal(dm, &vol));
836: PetscCall(DMDestroy(&boundary));
837: PetscFunctionReturn(PETSC_SUCCESS);
838: }
840: static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ)
841: {
842: DMLabel cutLabel = NULL;
843: PetscInt markerTop = 1, faceMarkerTop = 1;
844: PetscInt markerBottom = 1, faceMarkerBottom = 1;
845: PetscInt markerFront = 1, faceMarkerFront = 1;
846: PetscInt markerBack = 1, faceMarkerBack = 1;
847: PetscInt markerRight = 1, faceMarkerRight = 1;
848: PetscInt markerLeft = 1, faceMarkerLeft = 1;
849: PetscInt dim;
850: PetscBool markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE;
851: PetscMPIInt rank;
853: PetscFunctionBegin;
854: PetscCall(DMGetDimension(dm, &dim));
855: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
856: PetscCall(DMCreateLabel(dm, "marker"));
857: PetscCall(DMCreateLabel(dm, "Face Sets"));
858: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL));
859: if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST || bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST || bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) {
860: if (cutMarker) {
861: PetscCall(DMCreateLabel(dm, "periodic_cut"));
862: PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel));
863: }
864: }
865: switch (dim) {
866: case 2:
867: faceMarkerTop = 3;
868: faceMarkerBottom = 1;
869: faceMarkerRight = 2;
870: faceMarkerLeft = 4;
871: break;
872: case 3:
873: faceMarkerBottom = 1;
874: faceMarkerTop = 2;
875: faceMarkerFront = 3;
876: faceMarkerBack = 4;
877: faceMarkerRight = 5;
878: faceMarkerLeft = 6;
879: break;
880: default:
881: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim);
882: }
883: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL));
884: if (markerSeparate) {
885: markerBottom = faceMarkerBottom;
886: markerTop = faceMarkerTop;
887: markerFront = faceMarkerFront;
888: markerBack = faceMarkerBack;
889: markerRight = faceMarkerRight;
890: markerLeft = faceMarkerLeft;
891: }
892: {
893: const PetscInt numXEdges = rank == 0 ? edges[0] : 0;
894: const PetscInt numYEdges = rank == 0 ? edges[1] : 0;
895: const PetscInt numZEdges = rank == 0 ? edges[2] : 0;
896: const PetscInt numXVertices = rank == 0 ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0] + 1) : 0;
897: const PetscInt numYVertices = rank == 0 ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1] + 1) : 0;
898: const PetscInt numZVertices = rank == 0 ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2] + 1) : 0;
899: const PetscInt numCells = numXEdges * numYEdges * numZEdges;
900: const PetscInt numXFaces = numYEdges * numZEdges;
901: const PetscInt numYFaces = numXEdges * numZEdges;
902: const PetscInt numZFaces = numXEdges * numYEdges;
903: const PetscInt numTotXFaces = numXVertices * numXFaces;
904: const PetscInt numTotYFaces = numYVertices * numYFaces;
905: const PetscInt numTotZFaces = numZVertices * numZFaces;
906: const PetscInt numFaces = numTotXFaces + numTotYFaces + numTotZFaces;
907: const PetscInt numTotXEdges = numXEdges * numYVertices * numZVertices;
908: const PetscInt numTotYEdges = numYEdges * numXVertices * numZVertices;
909: const PetscInt numTotZEdges = numZEdges * numXVertices * numYVertices;
910: const PetscInt numVertices = numXVertices * numYVertices * numZVertices;
911: const PetscInt numEdges = numTotXEdges + numTotYEdges + numTotZEdges;
912: const PetscInt firstVertex = (dim == 2) ? numFaces : numCells;
913: const PetscInt firstXFace = (dim == 2) ? 0 : numCells + numVertices;
914: const PetscInt firstYFace = firstXFace + numTotXFaces;
915: const PetscInt firstZFace = firstYFace + numTotYFaces;
916: const PetscInt firstXEdge = numCells + numFaces + numVertices;
917: const PetscInt firstYEdge = firstXEdge + numTotXEdges;
918: const PetscInt firstZEdge = firstYEdge + numTotYEdges;
919: Vec coordinates;
920: PetscSection coordSection;
921: PetscScalar *coords;
922: PetscInt coordSize;
923: PetscInt v, vx, vy, vz;
924: PetscInt c, f, fx, fy, fz, e, ex, ey, ez;
926: PetscCall(DMPlexSetChart(dm, 0, numCells + numFaces + numEdges + numVertices));
927: for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 6));
928: for (f = firstXFace; f < firstXFace + numFaces; ++f) PetscCall(DMPlexSetConeSize(dm, f, 4));
929: for (e = firstXEdge; e < firstXEdge + numEdges; ++e) PetscCall(DMPlexSetConeSize(dm, e, 2));
930: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
931: /* Build cells */
932: for (fz = 0; fz < numZEdges; ++fz) {
933: for (fy = 0; fy < numYEdges; ++fy) {
934: for (fx = 0; fx < numXEdges; ++fx) {
935: PetscInt cell = (fz * numYEdges + fy) * numXEdges + fx;
936: PetscInt faceB = firstZFace + (fy * numXEdges + fx) * numZVertices + fz;
937: PetscInt faceT = firstZFace + (fy * numXEdges + fx) * numZVertices + ((fz + 1) % numZVertices);
938: PetscInt faceF = firstYFace + (fz * numXEdges + fx) * numYVertices + fy;
939: PetscInt faceK = firstYFace + (fz * numXEdges + fx) * numYVertices + ((fy + 1) % numYVertices);
940: PetscInt faceL = firstXFace + (fz * numYEdges + fy) * numXVertices + fx;
941: PetscInt faceR = firstXFace + (fz * numYEdges + fy) * numXVertices + ((fx + 1) % numXVertices);
942: /* B, T, F, K, R, L */
943: PetscInt ornt[6] = {-2, 0, 0, -3, 0, -2}; /* ??? */
944: PetscInt cone[6];
946: /* no boundary twisting in 3D */
947: cone[0] = faceB;
948: cone[1] = faceT;
949: cone[2] = faceF;
950: cone[3] = faceK;
951: cone[4] = faceR;
952: cone[5] = faceL;
953: PetscCall(DMPlexSetCone(dm, cell, cone));
954: PetscCall(DMPlexSetConeOrientation(dm, cell, ornt));
955: if (bdX != DM_BOUNDARY_NONE && fx == numXEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, cell, 2));
956: if (bdY != DM_BOUNDARY_NONE && fy == numYEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, cell, 2));
957: if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, cell, 2));
958: }
959: }
960: }
961: /* Build x faces */
962: for (fz = 0; fz < numZEdges; ++fz) {
963: for (fy = 0; fy < numYEdges; ++fy) {
964: for (fx = 0; fx < numXVertices; ++fx) {
965: PetscInt face = firstXFace + (fz * numYEdges + fy) * numXVertices + fx;
966: PetscInt edgeL = firstZEdge + (fy * numXVertices + fx) * numZEdges + fz;
967: PetscInt edgeR = firstZEdge + (((fy + 1) % numYVertices) * numXVertices + fx) * numZEdges + fz;
968: PetscInt edgeB = firstYEdge + (fz * numXVertices + fx) * numYEdges + fy;
969: PetscInt edgeT = firstYEdge + (((fz + 1) % numZVertices) * numXVertices + fx) * numYEdges + fy;
970: PetscInt ornt[4] = {0, 0, -1, -1};
971: PetscInt cone[4];
973: if (dim == 3) {
974: /* markers */
975: if (bdX != DM_BOUNDARY_PERIODIC) {
976: if (fx == numXVertices - 1) {
977: PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight));
978: PetscCall(DMSetLabelValue(dm, "marker", face, markerRight));
979: } else if (fx == 0) {
980: PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft));
981: PetscCall(DMSetLabelValue(dm, "marker", face, markerLeft));
982: }
983: }
984: }
985: cone[0] = edgeB;
986: cone[1] = edgeR;
987: cone[2] = edgeT;
988: cone[3] = edgeL;
989: PetscCall(DMPlexSetCone(dm, face, cone));
990: PetscCall(DMPlexSetConeOrientation(dm, face, ornt));
991: }
992: }
993: }
994: /* Build y faces */
995: for (fz = 0; fz < numZEdges; ++fz) {
996: for (fx = 0; fx < numXEdges; ++fx) {
997: for (fy = 0; fy < numYVertices; ++fy) {
998: PetscInt face = firstYFace + (fz * numXEdges + fx) * numYVertices + fy;
999: PetscInt edgeL = firstZEdge + (fy * numXVertices + fx) * numZEdges + fz;
1000: PetscInt edgeR = firstZEdge + (fy * numXVertices + ((fx + 1) % numXVertices)) * numZEdges + fz;
1001: PetscInt edgeB = firstXEdge + (fz * numYVertices + fy) * numXEdges + fx;
1002: PetscInt edgeT = firstXEdge + (((fz + 1) % numZVertices) * numYVertices + fy) * numXEdges + fx;
1003: PetscInt ornt[4] = {0, 0, -1, -1};
1004: PetscInt cone[4];
1006: if (dim == 3) {
1007: /* markers */
1008: if (bdY != DM_BOUNDARY_PERIODIC) {
1009: if (fy == numYVertices - 1) {
1010: PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack));
1011: PetscCall(DMSetLabelValue(dm, "marker", face, markerBack));
1012: } else if (fy == 0) {
1013: PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront));
1014: PetscCall(DMSetLabelValue(dm, "marker", face, markerFront));
1015: }
1016: }
1017: }
1018: cone[0] = edgeB;
1019: cone[1] = edgeR;
1020: cone[2] = edgeT;
1021: cone[3] = edgeL;
1022: PetscCall(DMPlexSetCone(dm, face, cone));
1023: PetscCall(DMPlexSetConeOrientation(dm, face, ornt));
1024: }
1025: }
1026: }
1027: /* Build z faces */
1028: for (fy = 0; fy < numYEdges; ++fy) {
1029: for (fx = 0; fx < numXEdges; ++fx) {
1030: for (fz = 0; fz < numZVertices; fz++) {
1031: PetscInt face = firstZFace + (fy * numXEdges + fx) * numZVertices + fz;
1032: PetscInt edgeL = firstYEdge + (fz * numXVertices + fx) * numYEdges + fy;
1033: PetscInt edgeR = firstYEdge + (fz * numXVertices + ((fx + 1) % numXVertices)) * numYEdges + fy;
1034: PetscInt edgeB = firstXEdge + (fz * numYVertices + fy) * numXEdges + fx;
1035: PetscInt edgeT = firstXEdge + (fz * numYVertices + ((fy + 1) % numYVertices)) * numXEdges + fx;
1036: PetscInt ornt[4] = {0, 0, -1, -1};
1037: PetscInt cone[4];
1039: if (dim == 2) {
1040: if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges - 1) {
1041: edgeR += numYEdges - 1 - 2 * fy;
1042: ornt[1] = -1;
1043: }
1044: if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges - 1) {
1045: edgeT += numXEdges - 1 - 2 * fx;
1046: ornt[2] = 0;
1047: }
1048: if (bdX != DM_BOUNDARY_NONE && fx == numXEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, face, 2));
1049: if (bdY != DM_BOUNDARY_NONE && fy == numYEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, face, 2));
1050: } else {
1051: /* markers */
1052: if (bdZ != DM_BOUNDARY_PERIODIC) {
1053: if (fz == numZVertices - 1) {
1054: PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop));
1055: PetscCall(DMSetLabelValue(dm, "marker", face, markerTop));
1056: } else if (fz == 0) {
1057: PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom));
1058: PetscCall(DMSetLabelValue(dm, "marker", face, markerBottom));
1059: }
1060: }
1061: }
1062: cone[0] = edgeB;
1063: cone[1] = edgeR;
1064: cone[2] = edgeT;
1065: cone[3] = edgeL;
1066: PetscCall(DMPlexSetCone(dm, face, cone));
1067: PetscCall(DMPlexSetConeOrientation(dm, face, ornt));
1068: }
1069: }
1070: }
1071: /* Build Z edges*/
1072: for (vy = 0; vy < numYVertices; vy++) {
1073: for (vx = 0; vx < numXVertices; vx++) {
1074: for (ez = 0; ez < numZEdges; ez++) {
1075: const PetscInt edge = firstZEdge + (vy * numXVertices + vx) * numZEdges + ez;
1076: const PetscInt vertexB = firstVertex + (ez * numYVertices + vy) * numXVertices + vx;
1077: const PetscInt vertexT = firstVertex + (((ez + 1) % numZVertices) * numYVertices + vy) * numXVertices + vx;
1078: PetscInt cone[2];
1080: cone[0] = vertexB;
1081: cone[1] = vertexT;
1082: PetscCall(DMPlexSetCone(dm, edge, cone));
1083: if (dim == 3) {
1084: if (bdX != DM_BOUNDARY_PERIODIC) {
1085: if (vx == numXVertices - 1) {
1086: PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight));
1087: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight));
1088: if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight));
1089: } else if (vx == 0) {
1090: PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft));
1091: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft));
1092: if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft));
1093: }
1094: }
1095: if (bdY != DM_BOUNDARY_PERIODIC) {
1096: if (vy == numYVertices - 1) {
1097: PetscCall(DMSetLabelValue(dm, "marker", edge, markerBack));
1098: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBack));
1099: if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBack));
1100: } else if (vy == 0) {
1101: PetscCall(DMSetLabelValue(dm, "marker", edge, markerFront));
1102: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerFront));
1103: if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerFront));
1104: }
1105: }
1106: }
1107: }
1108: }
1109: }
1110: /* Build Y edges*/
1111: for (vz = 0; vz < numZVertices; vz++) {
1112: for (vx = 0; vx < numXVertices; vx++) {
1113: for (ey = 0; ey < numYEdges; ey++) {
1114: const PetscInt nextv = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges - 1) ? (numXVertices - vx - 1) : (vz * numYVertices + ((ey + 1) % numYVertices)) * numXVertices + vx;
1115: const PetscInt edge = firstYEdge + (vz * numXVertices + vx) * numYEdges + ey;
1116: const PetscInt vertexF = firstVertex + (vz * numYVertices + ey) * numXVertices + vx;
1117: const PetscInt vertexK = firstVertex + nextv;
1118: PetscInt cone[2];
1120: cone[0] = vertexF;
1121: cone[1] = vertexK;
1122: PetscCall(DMPlexSetCone(dm, edge, cone));
1123: if (dim == 2) {
1124: if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) {
1125: if (vx == numXVertices - 1) {
1126: PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight));
1127: PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight));
1128: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight));
1129: if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight));
1130: } else if (vx == 0) {
1131: PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft));
1132: PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft));
1133: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft));
1134: if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft));
1135: }
1136: } else {
1137: if (vx == 0 && cutLabel) {
1138: PetscCall(DMLabelSetValue(cutLabel, edge, 1));
1139: PetscCall(DMLabelSetValue(cutLabel, cone[0], 1));
1140: if (ey == numYEdges - 1) PetscCall(DMLabelSetValue(cutLabel, cone[1], 1));
1141: }
1142: }
1143: } else {
1144: if (bdX != DM_BOUNDARY_PERIODIC) {
1145: if (vx == numXVertices - 1) {
1146: PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight));
1147: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight));
1148: if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight));
1149: } else if (vx == 0) {
1150: PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft));
1151: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft));
1152: if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft));
1153: }
1154: }
1155: if (bdZ != DM_BOUNDARY_PERIODIC) {
1156: if (vz == numZVertices - 1) {
1157: PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop));
1158: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop));
1159: if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop));
1160: } else if (vz == 0) {
1161: PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom));
1162: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom));
1163: if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom));
1164: }
1165: }
1166: }
1167: }
1168: }
1169: }
1170: /* Build X edges*/
1171: for (vz = 0; vz < numZVertices; vz++) {
1172: for (vy = 0; vy < numYVertices; vy++) {
1173: for (ex = 0; ex < numXEdges; ex++) {
1174: const PetscInt nextv = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges - 1) ? (numYVertices - vy - 1) * numXVertices : (vz * numYVertices + vy) * numXVertices + (ex + 1) % numXVertices;
1175: const PetscInt edge = firstXEdge + (vz * numYVertices + vy) * numXEdges + ex;
1176: const PetscInt vertexL = firstVertex + (vz * numYVertices + vy) * numXVertices + ex;
1177: const PetscInt vertexR = firstVertex + nextv;
1178: PetscInt cone[2];
1180: cone[0] = vertexL;
1181: cone[1] = vertexR;
1182: PetscCall(DMPlexSetCone(dm, edge, cone));
1183: if (dim == 2) {
1184: if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) {
1185: if (vy == numYVertices - 1) {
1186: PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop));
1187: PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop));
1188: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop));
1189: if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop));
1190: } else if (vy == 0) {
1191: PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom));
1192: PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom));
1193: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom));
1194: if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom));
1195: }
1196: } else {
1197: if (vy == 0 && cutLabel) {
1198: PetscCall(DMLabelSetValue(cutLabel, edge, 1));
1199: PetscCall(DMLabelSetValue(cutLabel, cone[0], 1));
1200: if (ex == numXEdges - 1) PetscCall(DMLabelSetValue(cutLabel, cone[1], 1));
1201: }
1202: }
1203: } else {
1204: if (bdY != DM_BOUNDARY_PERIODIC) {
1205: if (vy == numYVertices - 1) {
1206: PetscCall(DMSetLabelValue(dm, "marker", edge, markerBack));
1207: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBack));
1208: if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBack));
1209: } else if (vy == 0) {
1210: PetscCall(DMSetLabelValue(dm, "marker", edge, markerFront));
1211: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerFront));
1212: if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerFront));
1213: }
1214: }
1215: if (bdZ != DM_BOUNDARY_PERIODIC) {
1216: if (vz == numZVertices - 1) {
1217: PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop));
1218: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop));
1219: if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop));
1220: } else if (vz == 0) {
1221: PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom));
1222: PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom));
1223: if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom));
1224: }
1225: }
1226: }
1227: }
1228: }
1229: }
1230: PetscCall(DMPlexSymmetrize(dm));
1231: PetscCall(DMPlexStratify(dm));
1232: /* Build coordinates */
1233: PetscCall(DMGetCoordinateSection(dm, &coordSection));
1234: PetscCall(PetscSectionSetNumFields(coordSection, 1));
1235: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim));
1236: PetscCall(PetscSectionSetChart(coordSection, firstVertex, firstVertex + numVertices));
1237: for (v = firstVertex; v < firstVertex + numVertices; ++v) {
1238: PetscCall(PetscSectionSetDof(coordSection, v, dim));
1239: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim));
1240: }
1241: PetscCall(PetscSectionSetUp(coordSection));
1242: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
1243: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
1244: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
1245: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
1246: PetscCall(VecSetBlockSize(coordinates, dim));
1247: PetscCall(VecSetType(coordinates, VECSTANDARD));
1248: PetscCall(VecGetArray(coordinates, &coords));
1249: for (vz = 0; vz < numZVertices; ++vz) {
1250: for (vy = 0; vy < numYVertices; ++vy) {
1251: for (vx = 0; vx < numXVertices; ++vx) {
1252: coords[((vz * numYVertices + vy) * numXVertices + vx) * dim + 0] = lower[0] + ((upper[0] - lower[0]) / numXEdges) * vx;
1253: coords[((vz * numYVertices + vy) * numXVertices + vx) * dim + 1] = lower[1] + ((upper[1] - lower[1]) / numYEdges) * vy;
1254: if (dim == 3) coords[((vz * numYVertices + vy) * numXVertices + vx) * dim + 2] = lower[2] + ((upper[2] - lower[2]) / numZEdges) * vz;
1255: }
1256: }
1257: }
1258: PetscCall(VecRestoreArray(coordinates, &coords));
1259: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
1260: PetscCall(VecDestroy(&coordinates));
1261: }
1262: PetscFunctionReturn(PETSC_SUCCESS);
1263: }
1265: static PetscErrorCode DMPlexCreateBoxMesh_Tensor_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[])
1266: {
1267: DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1268: PetscInt fac[3] = {0, 0, 0}, d;
1270: PetscFunctionBegin;
1271: PetscAssertPointer(dm, 1);
1273: PetscCall(DMSetDimension(dm, dim));
1274: for (d = 0; d < dim; ++d) {
1275: fac[d] = faces[d];
1276: bdt[d] = periodicity[d];
1277: }
1278: PetscCall(DMPlexCreateCubeMesh_Internal(dm, lower, upper, fac, bdt[0], bdt[1], bdt[2]));
1279: if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST || periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST || (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) {
1280: PetscReal L[3] = {-1., -1., 0.};
1281: PetscReal maxCell[3] = {-1., -1., 0.};
1283: for (d = 0; d < dim; ++d) {
1284: if (periodicity[d] != DM_BOUNDARY_NONE) {
1285: L[d] = upper[d] - lower[d];
1286: maxCell[d] = 1.1 * (L[d] / PetscMax(1, faces[d]));
1287: }
1288: }
1289: PetscCall(DMSetPeriodicity(dm, maxCell, lower, L));
1290: }
1291: PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
1292: PetscFunctionReturn(PETSC_SUCCESS);
1293: }
1295: static PetscErrorCode DMPlexCreateBoxMesh_Internal(DM dm, DMPlexShape shape, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
1296: {
1297: PetscFunctionBegin;
1298: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
1299: if (shape == DM_SHAPE_ZBOX) PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Internal(dm, dim, faces, lower, upper, periodicity, interpolate));
1300: else if (dim == 1) PetscCall(DMPlexCreateLineMesh_Internal(dm, faces[0], lower[0], upper[0], periodicity[0]));
1301: else if (simplex) PetscCall(DMPlexCreateBoxMesh_Simplex_Internal(dm, dim, faces, lower, upper, periodicity, interpolate));
1302: else PetscCall(DMPlexCreateBoxMesh_Tensor_Internal(dm, dim, faces, lower, upper, periodicity));
1303: if (!interpolate && dim > 1 && !simplex) {
1304: DM udm;
1306: PetscCall(DMPlexUninterpolate(dm, &udm));
1307: PetscCall(DMPlexCopyCoordinates(dm, udm));
1308: PetscCall(DMPlexReplace_Internal(dm, &udm));
1309: }
1310: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
1311: PetscFunctionReturn(PETSC_SUCCESS);
1312: }
1314: /*@C
1315: DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices or tensor cells (hexahedra).
1317: Collective
1319: Input Parameters:
1320: + comm - The communicator for the `DM` object
1321: . dim - The spatial dimension
1322: . simplex - `PETSC_TRUE` for simplices, `PETSC_FALSE` for tensor cells
1323: . faces - Number of faces per dimension, or `NULL` for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
1324: . lower - The lower left corner, or `NULL` for (0, 0, 0)
1325: . upper - The upper right corner, or `NULL` for (1, 1, 1)
1326: . periodicity - The boundary type for the X,Y,Z direction, or `NULL` for `DM_BOUNDARY_NONE`
1327: - interpolate - Flag to create intermediate mesh pieces (edges, faces)
1329: Output Parameter:
1330: . dm - The `DM` object
1332: Level: beginner
1334: Note:
1335: To customize this mesh using options, use
1336: .vb
1337: DMCreate(comm, &dm);
1338: DMSetType(dm, DMPLEX);
1339: DMSetFromOptions(dm);
1340: .ve
1341: and use the options in `DMSetFromOptions()`.
1343: Here is the numbering returned for 2 faces in each direction for tensor cells\:
1344: .vb
1345: 10---17---11---18----12
1346: | | |
1347: | | |
1348: 20 2 22 3 24
1349: | | |
1350: | | |
1351: 7---15----8---16----9
1352: | | |
1353: | | |
1354: 19 0 21 1 23
1355: | | |
1356: | | |
1357: 4---13----5---14----6
1358: .ve
1359: and for simplicial cells
1360: .vb
1361: 14----8---15----9----16
1362: |\ 5 |\ 7 |
1363: | \ | \ |
1364: 13 2 14 3 15
1365: | 4 \ | 6 \ |
1366: | \ | \ |
1367: 11----6---12----7----13
1368: |\ |\ |
1369: | \ 1 | \ 3 |
1370: 10 0 11 1 12
1371: | 0 \ | 2 \ |
1372: | \ | \ |
1373: 8----4----9----5----10
1374: .ve
1376: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMSetFromOptions()`, `DMPlexCreateFromFile()`, `DMPlexCreateHexCylinderMesh()`, `DMSetType()`, `DMCreate()`
1377: @*/
1378: PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)
1379: {
1380: PetscInt fac[3] = {1, 1, 1};
1381: PetscReal low[3] = {0, 0, 0};
1382: PetscReal upp[3] = {1, 1, 1};
1383: DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1385: PetscFunctionBegin;
1386: PetscCall(DMCreate(comm, dm));
1387: PetscCall(DMSetType(*dm, DMPLEX));
1388: PetscCall(DMPlexCreateBoxMesh_Internal(*dm, DM_SHAPE_BOX, dim, simplex, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, periodicity ? periodicity : bdt, interpolate));
1389: if (periodicity) PetscCall(DMLocalizeCoordinates(*dm));
1390: PetscFunctionReturn(PETSC_SUCCESS);
1391: }
1393: static PetscErrorCode DMPlexCreateWedgeBoxMesh_Internal(DM dm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[])
1394: {
1395: DM bdm, vol;
1396: PetscInt i;
1398: PetscFunctionBegin;
1399: // TODO Now we can support periodicity
1400: for (i = 0; i < 3; ++i) PetscCheck(periodicity[i] == DM_BOUNDARY_NONE, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Periodicity not yet supported");
1401: PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &bdm));
1402: PetscCall(DMSetType(bdm, DMPLEX));
1403: PetscCall(DMSetDimension(bdm, 2));
1404: PetscCall(PetscLogEventBegin(DMPLEX_Generate, bdm, 0, 0, 0));
1405: PetscCall(DMPlexCreateBoxMesh_Simplex_Internal(bdm, 2, faces, lower, upper, periodicity, PETSC_TRUE));
1406: PetscCall(DMPlexExtrude(bdm, faces[2], upper[2] - lower[2], PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, NULL, NULL, &vol));
1407: PetscCall(PetscLogEventEnd(DMPLEX_Generate, bdm, 0, 0, 0));
1408: PetscCall(DMDestroy(&bdm));
1409: PetscCall(DMPlexReplace_Internal(dm, &vol));
1410: if (lower[2] != 0.0) {
1411: Vec v;
1412: PetscScalar *x;
1413: PetscInt cDim, n;
1415: PetscCall(DMGetCoordinatesLocal(dm, &v));
1416: PetscCall(VecGetBlockSize(v, &cDim));
1417: PetscCall(VecGetLocalSize(v, &n));
1418: PetscCall(VecGetArray(v, &x));
1419: x += cDim;
1420: for (i = 0; i < n; i += cDim) x[i] += lower[2];
1421: PetscCall(VecRestoreArray(v, &x));
1422: PetscCall(DMSetCoordinatesLocal(dm, v));
1423: }
1424: PetscFunctionReturn(PETSC_SUCCESS);
1425: }
1427: /*@
1428: DMPlexCreateWedgeBoxMesh - Creates a 3-D mesh tesselating the (x,y) plane and extruding in the third direction using wedge cells.
1430: Collective
1432: Input Parameters:
1433: + comm - The communicator for the `DM` object
1434: . faces - Number of faces per dimension, or `NULL` for (1, 1, 1)
1435: . lower - The lower left corner, or `NULL` for (0, 0, 0)
1436: . upper - The upper right corner, or `NULL` for (1, 1, 1)
1437: . periodicity - The boundary type for the X,Y,Z direction, or `NULL` for `DM_BOUNDARY_NONE`
1438: . orderHeight - If `PETSC_TRUE`, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
1439: - interpolate - Flag to create intermediate mesh pieces (edges, faces)
1441: Output Parameter:
1442: . dm - The `DM` object
1444: Level: beginner
1446: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateHexCylinderMesh()`, `DMPlexCreateWedgeCylinderMesh()`, `DMExtrude()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()`
1447: @*/
1448: PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool orderHeight, PetscBool interpolate, DM *dm)
1449: {
1450: PetscInt fac[3] = {1, 1, 1};
1451: PetscReal low[3] = {0, 0, 0};
1452: PetscReal upp[3] = {1, 1, 1};
1453: DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1455: PetscFunctionBegin;
1456: PetscCall(DMCreate(comm, dm));
1457: PetscCall(DMSetType(*dm, DMPLEX));
1458: PetscCall(DMPlexCreateWedgeBoxMesh_Internal(*dm, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, periodicity ? periodicity : bdt));
1459: if (!interpolate) {
1460: DM udm;
1462: PetscCall(DMPlexUninterpolate(*dm, &udm));
1463: PetscCall(DMPlexReplace_Internal(*dm, &udm));
1464: }
1465: if (periodicity) PetscCall(DMLocalizeCoordinates(*dm));
1466: PetscFunctionReturn(PETSC_SUCCESS);
1467: }
1469: /*
1470: DMPlexTensorPointLexicographic_Private - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max' for that dimension.
1472: Input Parameters:
1473: + len - The length of the tuple
1474: . max - The maximum for each dimension, so values are in [0, max)
1475: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
1477: Output Parameter:
1478: . tup - A tuple of `len` integers whose entries are at most `max`
1480: Level: developer
1482: Note:
1483: Ordering is lexicographic with lowest index as least significant in ordering.
1484: e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.
1486: .seealso: PetscDualSpaceTensorPointLexicographic_Internal(), PetscDualSpaceLatticePointLexicographic_Internal()
1487: */
1488: static PetscErrorCode DMPlexTensorPointLexicographic_Private(PetscInt len, const PetscInt max[], PetscInt tup[])
1489: {
1490: PetscInt i;
1492: PetscFunctionBegin;
1493: for (i = 0; i < len; ++i) {
1494: if (tup[i] < max[i] - 1) {
1495: break;
1496: } else {
1497: tup[i] = 0;
1498: }
1499: }
1500: if (i == len) tup[i - 1] = max[i - 1];
1501: else ++tup[i];
1502: PetscFunctionReturn(PETSC_SUCCESS);
1503: }
1505: static PetscInt TupleToIndex_Private(PetscInt len, const PetscInt max[], const PetscInt tup[])
1506: {
1507: PetscInt i, idx = tup[len - 1];
1509: for (i = len - 2; i >= 0; --i) {
1510: idx *= max[i];
1511: idx += tup[i];
1512: }
1513: return idx;
1514: }
1516: static PetscErrorCode DestroyExtent_Private(void *extent)
1517: {
1518: return PetscFree(extent);
1519: }
1521: static PetscErrorCode DMPlexCreateHypercubicMesh_Internal(DM dm, PetscInt dim, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], const DMBoundaryType bd[])
1522: {
1523: Vec coordinates;
1524: PetscSection coordSection;
1525: DMLabel cutLabel = NULL;
1526: PetscBool cutMarker = PETSC_FALSE;
1527: PetscBool periodic = PETSC_FALSE;
1528: PetscInt numCells = 1, c;
1529: PetscInt numVertices = 1, v;
1530: PetscScalar *coords;
1531: PetscInt *vertices, *vert, *vtmp, *supp, cone[2];
1532: PetscInt d, e, cell = 0, coordSize;
1533: PetscMPIInt rank;
1535: PetscFunctionBegin;
1536: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1537: PetscCall(DMSetDimension(dm, dim));
1538: PetscCall(PetscCalloc4(dim, &vertices, dim, &vert, dim, &vtmp, 2 * dim, &supp));
1539: PetscCall(DMCreateLabel(dm, "marker"));
1540: PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL));
1541: for (d = 0; d < dim; ++d) periodic = (periodic || bd[d] == DM_BOUNDARY_PERIODIC) ? PETSC_TRUE : PETSC_FALSE;
1542: if (periodic && cutMarker) {
1543: PetscCall(DMCreateLabel(dm, "periodic_cut"));
1544: PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel));
1545: }
1546: for (d = 0; d < dim; ++d) PetscCheck(bd[d] == DM_BOUNDARY_PERIODIC, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Hypercubic mesh must be periodic now");
1547: for (d = 0; d < dim; ++d) {
1548: vertices[d] = edges[d];
1549: numVertices *= vertices[d];
1550: }
1551: numCells = numVertices * dim;
1552: PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices));
1553: for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(dm, c, 2));
1554: for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetSupportSize(dm, v, 2 * dim));
1555: /* TODO Loop over boundary and reset support sizes */
1556: PetscCall(DMSetUp(dm)); /* Allocate space for cones and supports */
1557: /* Build cell cones and vertex supports */
1558: PetscCall(DMCreateLabel(dm, "celltype"));
1559: while (vert[dim - 1] < vertices[dim - 1]) {
1560: const PetscInt vertex = TupleToIndex_Private(dim, vertices, vert) + numCells;
1561: PetscInt s = 0;
1563: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT ":", vertex));
1564: for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, vert[d]));
1565: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1566: PetscCall(DMPlexSetCellType(dm, vertex, DM_POLYTOPE_POINT));
1567: for (d = 0; d < dim; ++d) {
1568: for (e = 0; e < dim; ++e) vtmp[e] = vert[e];
1569: vtmp[d] = (vert[d] + 1) % vertices[d];
1570: cone[0] = vertex;
1571: cone[1] = TupleToIndex_Private(dim, vertices, vtmp) + numCells;
1572: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Vertex %" PetscInt_FMT ":", cone[1]));
1573: for (e = 0; e < dim; ++e) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, vtmp[e]));
1574: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1575: PetscCall(DMPlexSetCone(dm, cell, cone));
1576: PetscCall(DMPlexSetCellType(dm, cell, DM_POLYTOPE_SEGMENT));
1577: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Edge %" PetscInt_FMT " (%" PetscInt_FMT " %" PetscInt_FMT ")\n", cell, cone[0], cone[1]));
1578: ++cell;
1579: }
1580: for (d = 0; d < dim; ++d) {
1581: for (e = 0; e < dim; ++e) vtmp[e] = vert[e];
1582: vtmp[d] = (vert[d] + vertices[d] - 1) % vertices[d];
1583: supp[s++] = TupleToIndex_Private(dim, vertices, vtmp) * dim + d;
1584: supp[s++] = (vertex - numCells) * dim + d;
1585: PetscCall(DMPlexSetSupport(dm, vertex, supp));
1586: }
1587: PetscCall(DMPlexTensorPointLexicographic_Private(dim, vertices, vert));
1588: }
1589: PetscCall(DMPlexStratify(dm));
1590: /* Build coordinates */
1591: PetscCall(DMGetCoordinateSection(dm, &coordSection));
1592: PetscCall(PetscSectionSetNumFields(coordSection, 1));
1593: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim));
1594: PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
1595: for (v = numCells; v < numCells + numVertices; ++v) {
1596: PetscCall(PetscSectionSetDof(coordSection, v, dim));
1597: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim));
1598: }
1599: PetscCall(PetscSectionSetUp(coordSection));
1600: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
1601: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
1602: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
1603: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
1604: PetscCall(VecSetBlockSize(coordinates, dim));
1605: PetscCall(VecSetType(coordinates, VECSTANDARD));
1606: PetscCall(VecGetArray(coordinates, &coords));
1607: for (d = 0; d < dim; ++d) vert[d] = 0;
1608: while (vert[dim - 1] < vertices[dim - 1]) {
1609: const PetscInt vertex = TupleToIndex_Private(dim, vertices, vert);
1611: for (d = 0; d < dim; ++d) coords[vertex * dim + d] = lower[d] + ((upper[d] - lower[d]) / vertices[d]) * vert[d];
1612: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT ":", vertex));
1613: for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, vert[d]));
1614: for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(coords[vertex * dim + d])));
1615: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1616: PetscCall(DMPlexTensorPointLexicographic_Private(dim, vertices, vert));
1617: }
1618: PetscCall(VecRestoreArray(coordinates, &coords));
1619: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
1620: PetscCall(VecDestroy(&coordinates));
1621: PetscCall(PetscFree4(vertices, vert, vtmp, supp));
1622: //PetscCall(DMSetPeriodicity(dm, NULL, lower, upper));
1623: // Attach the extent
1624: {
1625: PetscContainer c;
1626: PetscInt *extent;
1628: PetscCall(PetscMalloc1(dim, &extent));
1629: for (PetscInt d = 0; d < dim; ++d) extent[d] = edges[d];
1630: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &c));
1631: PetscCall(PetscContainerSetUserDestroy(c, DestroyExtent_Private));
1632: PetscCall(PetscContainerSetPointer(c, extent));
1633: PetscCall(PetscObjectCompose((PetscObject)dm, "_extent", (PetscObject)c));
1634: PetscCall(PetscContainerDestroy(&c));
1635: }
1636: PetscFunctionReturn(PETSC_SUCCESS);
1637: }
1639: /*@C
1640: DMPlexCreateHypercubicMesh - Creates a periodic mesh on the tensor product of unit intervals using only vertices and edges.
1642: Collective
1644: Input Parameters:
1645: + comm - The communicator for the DM object
1646: . dim - The spatial dimension
1647: . edges - Number of edges per dimension, or `NULL` for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D
1648: . lower - The lower left corner, or `NULL` for (0, 0, 0)
1649: - upper - The upper right corner, or `NULL` for (1, 1, 1)
1651: Output Parameter:
1652: . dm - The DM object
1654: Level: beginner
1656: Note:
1657: If you want to customize this mesh using options, you just need to
1658: .vb
1659: DMCreate(comm, &dm);
1660: DMSetType(dm, DMPLEX);
1661: DMSetFromOptions(dm);
1662: .ve
1663: and use the options on the `DMSetFromOptions()` page.
1665: The vertices are numbered is lexicographic order, and the dim edges exiting a vertex in the positive orthant are number consecutively,
1666: .vb
1667: 18--0-19--2-20--4-18
1668: | | | |
1669: 13 15 17 13
1670: | | | |
1671: 24-12-25-14-26-16-24
1672: | | | |
1673: 7 9 11 7
1674: | | | |
1675: 21--6-22--8-23-10-21
1676: | | | |
1677: 1 3 5 1
1678: | | | |
1679: 18--0-19--2-20--4-18
1680: .ve
1682: .seealso: `DMSetFromOptions()`, `DMPlexCreateFromFile()`, `DMPlexCreateHexCylinderMesh()`, `DMSetType()`, `DMCreate()`
1683: @*/
1684: PetscErrorCode DMPlexCreateHypercubicMesh(MPI_Comm comm, PetscInt dim, const PetscInt edges[], const PetscReal lower[], const PetscReal upper[], DM *dm)
1685: {
1686: PetscInt *edg;
1687: PetscReal *low, *upp;
1688: DMBoundaryType *bdt;
1689: PetscInt d;
1691: PetscFunctionBegin;
1692: PetscCall(DMCreate(comm, dm));
1693: PetscCall(DMSetType(*dm, DMPLEX));
1694: PetscCall(PetscMalloc4(dim, &edg, dim, &low, dim, &upp, dim, &bdt));
1695: for (d = 0; d < dim; ++d) {
1696: edg[d] = edges ? edges[d] : 1;
1697: low[d] = lower ? lower[d] : 0.;
1698: upp[d] = upper ? upper[d] : 1.;
1699: bdt[d] = DM_BOUNDARY_PERIODIC;
1700: }
1701: PetscCall(DMPlexCreateHypercubicMesh_Internal(*dm, dim, low, upp, edg, bdt));
1702: PetscCall(PetscFree4(edg, low, upp, bdt));
1703: PetscFunctionReturn(PETSC_SUCCESS);
1704: }
1706: /*@C
1707: DMPlexSetOptionsPrefix - Sets the prefix used for searching for all `DM` options in the database.
1709: Logically Collective
1711: Input Parameters:
1712: + dm - the `DM` context
1713: - prefix - the prefix to prepend to all option names
1715: Level: advanced
1717: Note:
1718: A hyphen (-) must NOT be given at the beginning of the prefix name.
1719: The first character of all runtime options is AUTOMATICALLY the hyphen.
1721: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `SNESSetFromOptions()`
1722: @*/
1723: PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[])
1724: {
1725: DM_Plex *mesh = (DM_Plex *)dm->data;
1727: PetscFunctionBegin;
1729: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, prefix));
1730: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mesh->partitioner, prefix));
1731: PetscFunctionReturn(PETSC_SUCCESS);
1732: }
1734: /* Remap geometry to cylinder
1735: TODO: This only works for a single refinement, then it is broken
1737: Interior square: Linear interpolation is correct
1738: The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing
1739: such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance
1741: phi = arctan(y/x)
1742: d_close = sqrt(1/8 + 1/4 sin^2(phi))
1743: d_far = sqrt(1/2 + sin^2(phi))
1745: so we remap them using
1747: x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close)
1748: y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close)
1750: If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y.
1751: */
1752: static void snapToCylinder(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1753: {
1754: const PetscReal dis = 1.0 / PetscSqrtReal(2.0);
1755: const PetscReal ds2 = 0.5 * dis;
1757: if ((PetscAbsScalar(u[0]) <= ds2) && (PetscAbsScalar(u[1]) <= ds2)) {
1758: f0[0] = u[0];
1759: f0[1] = u[1];
1760: } else {
1761: PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc;
1763: x = PetscRealPart(u[0]);
1764: y = PetscRealPart(u[1]);
1765: phi = PetscAtan2Real(y, x);
1766: sinp = PetscSinReal(phi);
1767: cosp = PetscCosReal(phi);
1768: if ((PetscAbsReal(phi) > PETSC_PI / 4.0) && (PetscAbsReal(phi) < 3.0 * PETSC_PI / 4.0)) {
1769: dc = PetscAbsReal(ds2 / sinp);
1770: df = PetscAbsReal(dis / sinp);
1771: xc = ds2 * x / PetscAbsReal(y);
1772: yc = ds2 * PetscSignReal(y);
1773: } else {
1774: dc = PetscAbsReal(ds2 / cosp);
1775: df = PetscAbsReal(dis / cosp);
1776: xc = ds2 * PetscSignReal(x);
1777: yc = ds2 * y / PetscAbsReal(x);
1778: }
1779: f0[0] = xc + (u[0] - xc) * (1.0 - dc) / (df - dc);
1780: f0[1] = yc + (u[1] - yc) * (1.0 - dc) / (df - dc);
1781: }
1782: f0[2] = u[2];
1783: }
1785: static PetscErrorCode DMPlexCreateHexCylinderMesh_Internal(DM dm, DMBoundaryType periodicZ)
1786: {
1787: const PetscInt dim = 3;
1788: PetscInt numCells, numVertices;
1789: PetscMPIInt rank;
1791: PetscFunctionBegin;
1792: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
1793: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1794: PetscCall(DMSetDimension(dm, dim));
1795: /* Create topology */
1796: {
1797: PetscInt cone[8], c;
1799: numCells = rank == 0 ? 5 : 0;
1800: numVertices = rank == 0 ? 16 : 0;
1801: if (periodicZ == DM_BOUNDARY_PERIODIC) {
1802: numCells *= 3;
1803: numVertices = rank == 0 ? 24 : 0;
1804: }
1805: PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices));
1806: for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 8));
1807: PetscCall(DMSetUp(dm));
1808: if (rank == 0) {
1809: if (periodicZ == DM_BOUNDARY_PERIODIC) {
1810: cone[0] = 15;
1811: cone[1] = 18;
1812: cone[2] = 17;
1813: cone[3] = 16;
1814: cone[4] = 31;
1815: cone[5] = 32;
1816: cone[6] = 33;
1817: cone[7] = 34;
1818: PetscCall(DMPlexSetCone(dm, 0, cone));
1819: cone[0] = 16;
1820: cone[1] = 17;
1821: cone[2] = 24;
1822: cone[3] = 23;
1823: cone[4] = 32;
1824: cone[5] = 36;
1825: cone[6] = 37;
1826: cone[7] = 33; /* 22 25 26 21 */
1827: PetscCall(DMPlexSetCone(dm, 1, cone));
1828: cone[0] = 18;
1829: cone[1] = 27;
1830: cone[2] = 24;
1831: cone[3] = 17;
1832: cone[4] = 34;
1833: cone[5] = 33;
1834: cone[6] = 37;
1835: cone[7] = 38;
1836: PetscCall(DMPlexSetCone(dm, 2, cone));
1837: cone[0] = 29;
1838: cone[1] = 27;
1839: cone[2] = 18;
1840: cone[3] = 15;
1841: cone[4] = 35;
1842: cone[5] = 31;
1843: cone[6] = 34;
1844: cone[7] = 38;
1845: PetscCall(DMPlexSetCone(dm, 3, cone));
1846: cone[0] = 29;
1847: cone[1] = 15;
1848: cone[2] = 16;
1849: cone[3] = 23;
1850: cone[4] = 35;
1851: cone[5] = 36;
1852: cone[6] = 32;
1853: cone[7] = 31;
1854: PetscCall(DMPlexSetCone(dm, 4, cone));
1856: cone[0] = 31;
1857: cone[1] = 34;
1858: cone[2] = 33;
1859: cone[3] = 32;
1860: cone[4] = 19;
1861: cone[5] = 22;
1862: cone[6] = 21;
1863: cone[7] = 20;
1864: PetscCall(DMPlexSetCone(dm, 5, cone));
1865: cone[0] = 32;
1866: cone[1] = 33;
1867: cone[2] = 37;
1868: cone[3] = 36;
1869: cone[4] = 22;
1870: cone[5] = 25;
1871: cone[6] = 26;
1872: cone[7] = 21;
1873: PetscCall(DMPlexSetCone(dm, 6, cone));
1874: cone[0] = 34;
1875: cone[1] = 38;
1876: cone[2] = 37;
1877: cone[3] = 33;
1878: cone[4] = 20;
1879: cone[5] = 21;
1880: cone[6] = 26;
1881: cone[7] = 28;
1882: PetscCall(DMPlexSetCone(dm, 7, cone));
1883: cone[0] = 35;
1884: cone[1] = 38;
1885: cone[2] = 34;
1886: cone[3] = 31;
1887: cone[4] = 30;
1888: cone[5] = 19;
1889: cone[6] = 20;
1890: cone[7] = 28;
1891: PetscCall(DMPlexSetCone(dm, 8, cone));
1892: cone[0] = 35;
1893: cone[1] = 31;
1894: cone[2] = 32;
1895: cone[3] = 36;
1896: cone[4] = 30;
1897: cone[5] = 25;
1898: cone[6] = 22;
1899: cone[7] = 19;
1900: PetscCall(DMPlexSetCone(dm, 9, cone));
1902: cone[0] = 19;
1903: cone[1] = 20;
1904: cone[2] = 21;
1905: cone[3] = 22;
1906: cone[4] = 15;
1907: cone[5] = 16;
1908: cone[6] = 17;
1909: cone[7] = 18;
1910: PetscCall(DMPlexSetCone(dm, 10, cone));
1911: cone[0] = 22;
1912: cone[1] = 21;
1913: cone[2] = 26;
1914: cone[3] = 25;
1915: cone[4] = 16;
1916: cone[5] = 23;
1917: cone[6] = 24;
1918: cone[7] = 17;
1919: PetscCall(DMPlexSetCone(dm, 11, cone));
1920: cone[0] = 20;
1921: cone[1] = 28;
1922: cone[2] = 26;
1923: cone[3] = 21;
1924: cone[4] = 18;
1925: cone[5] = 17;
1926: cone[6] = 24;
1927: cone[7] = 27;
1928: PetscCall(DMPlexSetCone(dm, 12, cone));
1929: cone[0] = 30;
1930: cone[1] = 28;
1931: cone[2] = 20;
1932: cone[3] = 19;
1933: cone[4] = 29;
1934: cone[5] = 15;
1935: cone[6] = 18;
1936: cone[7] = 27;
1937: PetscCall(DMPlexSetCone(dm, 13, cone));
1938: cone[0] = 30;
1939: cone[1] = 19;
1940: cone[2] = 22;
1941: cone[3] = 25;
1942: cone[4] = 29;
1943: cone[5] = 23;
1944: cone[6] = 16;
1945: cone[7] = 15;
1946: PetscCall(DMPlexSetCone(dm, 14, cone));
1947: } else {
1948: cone[0] = 5;
1949: cone[1] = 8;
1950: cone[2] = 7;
1951: cone[3] = 6;
1952: cone[4] = 9;
1953: cone[5] = 12;
1954: cone[6] = 11;
1955: cone[7] = 10;
1956: PetscCall(DMPlexSetCone(dm, 0, cone));
1957: cone[0] = 6;
1958: cone[1] = 7;
1959: cone[2] = 14;
1960: cone[3] = 13;
1961: cone[4] = 12;
1962: cone[5] = 15;
1963: cone[6] = 16;
1964: cone[7] = 11;
1965: PetscCall(DMPlexSetCone(dm, 1, cone));
1966: cone[0] = 8;
1967: cone[1] = 17;
1968: cone[2] = 14;
1969: cone[3] = 7;
1970: cone[4] = 10;
1971: cone[5] = 11;
1972: cone[6] = 16;
1973: cone[7] = 18;
1974: PetscCall(DMPlexSetCone(dm, 2, cone));
1975: cone[0] = 19;
1976: cone[1] = 17;
1977: cone[2] = 8;
1978: cone[3] = 5;
1979: cone[4] = 20;
1980: cone[5] = 9;
1981: cone[6] = 10;
1982: cone[7] = 18;
1983: PetscCall(DMPlexSetCone(dm, 3, cone));
1984: cone[0] = 19;
1985: cone[1] = 5;
1986: cone[2] = 6;
1987: cone[3] = 13;
1988: cone[4] = 20;
1989: cone[5] = 15;
1990: cone[6] = 12;
1991: cone[7] = 9;
1992: PetscCall(DMPlexSetCone(dm, 4, cone));
1993: }
1994: }
1995: PetscCall(DMPlexSymmetrize(dm));
1996: PetscCall(DMPlexStratify(dm));
1997: }
1998: /* Create cube geometry */
1999: {
2000: Vec coordinates;
2001: PetscSection coordSection;
2002: PetscScalar *coords;
2003: PetscInt coordSize, v;
2004: const PetscReal dis = 1.0 / PetscSqrtReal(2.0);
2005: const PetscReal ds2 = dis / 2.0;
2007: /* Build coordinates */
2008: PetscCall(DMGetCoordinateSection(dm, &coordSection));
2009: PetscCall(PetscSectionSetNumFields(coordSection, 1));
2010: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim));
2011: PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
2012: for (v = numCells; v < numCells + numVertices; ++v) {
2013: PetscCall(PetscSectionSetDof(coordSection, v, dim));
2014: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim));
2015: }
2016: PetscCall(PetscSectionSetUp(coordSection));
2017: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
2018: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
2019: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
2020: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
2021: PetscCall(VecSetBlockSize(coordinates, dim));
2022: PetscCall(VecSetType(coordinates, VECSTANDARD));
2023: PetscCall(VecGetArray(coordinates, &coords));
2024: if (rank == 0) {
2025: coords[0 * dim + 0] = -ds2;
2026: coords[0 * dim + 1] = -ds2;
2027: coords[0 * dim + 2] = 0.0;
2028: coords[1 * dim + 0] = ds2;
2029: coords[1 * dim + 1] = -ds2;
2030: coords[1 * dim + 2] = 0.0;
2031: coords[2 * dim + 0] = ds2;
2032: coords[2 * dim + 1] = ds2;
2033: coords[2 * dim + 2] = 0.0;
2034: coords[3 * dim + 0] = -ds2;
2035: coords[3 * dim + 1] = ds2;
2036: coords[3 * dim + 2] = 0.0;
2037: coords[4 * dim + 0] = -ds2;
2038: coords[4 * dim + 1] = -ds2;
2039: coords[4 * dim + 2] = 1.0;
2040: coords[5 * dim + 0] = -ds2;
2041: coords[5 * dim + 1] = ds2;
2042: coords[5 * dim + 2] = 1.0;
2043: coords[6 * dim + 0] = ds2;
2044: coords[6 * dim + 1] = ds2;
2045: coords[6 * dim + 2] = 1.0;
2046: coords[7 * dim + 0] = ds2;
2047: coords[7 * dim + 1] = -ds2;
2048: coords[7 * dim + 2] = 1.0;
2049: coords[8 * dim + 0] = dis;
2050: coords[8 * dim + 1] = -dis;
2051: coords[8 * dim + 2] = 0.0;
2052: coords[9 * dim + 0] = dis;
2053: coords[9 * dim + 1] = dis;
2054: coords[9 * dim + 2] = 0.0;
2055: coords[10 * dim + 0] = dis;
2056: coords[10 * dim + 1] = -dis;
2057: coords[10 * dim + 2] = 1.0;
2058: coords[11 * dim + 0] = dis;
2059: coords[11 * dim + 1] = dis;
2060: coords[11 * dim + 2] = 1.0;
2061: coords[12 * dim + 0] = -dis;
2062: coords[12 * dim + 1] = dis;
2063: coords[12 * dim + 2] = 0.0;
2064: coords[13 * dim + 0] = -dis;
2065: coords[13 * dim + 1] = dis;
2066: coords[13 * dim + 2] = 1.0;
2067: coords[14 * dim + 0] = -dis;
2068: coords[14 * dim + 1] = -dis;
2069: coords[14 * dim + 2] = 0.0;
2070: coords[15 * dim + 0] = -dis;
2071: coords[15 * dim + 1] = -dis;
2072: coords[15 * dim + 2] = 1.0;
2073: if (periodicZ == DM_BOUNDARY_PERIODIC) {
2074: /* 15 31 19 */ coords[16 * dim + 0] = -ds2;
2075: coords[16 * dim + 1] = -ds2;
2076: coords[16 * dim + 2] = 0.5;
2077: /* 16 32 22 */ coords[17 * dim + 0] = ds2;
2078: coords[17 * dim + 1] = -ds2;
2079: coords[17 * dim + 2] = 0.5;
2080: /* 17 33 21 */ coords[18 * dim + 0] = ds2;
2081: coords[18 * dim + 1] = ds2;
2082: coords[18 * dim + 2] = 0.5;
2083: /* 18 34 20 */ coords[19 * dim + 0] = -ds2;
2084: coords[19 * dim + 1] = ds2;
2085: coords[19 * dim + 2] = 0.5;
2086: /* 29 35 30 */ coords[20 * dim + 0] = -dis;
2087: coords[20 * dim + 1] = -dis;
2088: coords[20 * dim + 2] = 0.5;
2089: /* 23 36 25 */ coords[21 * dim + 0] = dis;
2090: coords[21 * dim + 1] = -dis;
2091: coords[21 * dim + 2] = 0.5;
2092: /* 24 37 26 */ coords[22 * dim + 0] = dis;
2093: coords[22 * dim + 1] = dis;
2094: coords[22 * dim + 2] = 0.5;
2095: /* 27 38 28 */ coords[23 * dim + 0] = -dis;
2096: coords[23 * dim + 1] = dis;
2097: coords[23 * dim + 2] = 0.5;
2098: }
2099: }
2100: PetscCall(VecRestoreArray(coordinates, &coords));
2101: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
2102: PetscCall(VecDestroy(&coordinates));
2103: }
2104: /* Create periodicity */
2105: if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) {
2106: PetscReal L[3] = {-1., -1., 0.};
2107: PetscReal maxCell[3] = {-1., -1., 0.};
2108: PetscReal lower[3] = {0.0, 0.0, 0.0};
2109: PetscReal upper[3] = {1.0, 1.0, 1.5};
2110: PetscInt numZCells = 3;
2112: L[2] = upper[2] - lower[2];
2113: maxCell[2] = 1.1 * (L[2] / numZCells);
2114: PetscCall(DMSetPeriodicity(dm, maxCell, lower, L));
2115: }
2116: {
2117: DM cdm;
2118: PetscDS cds;
2119: PetscScalar c[2] = {1.0, 1.0};
2121: PetscCall(DMPlexCreateCoordinateSpace(dm, 1, snapToCylinder));
2122: PetscCall(DMGetCoordinateDM(dm, &cdm));
2123: PetscCall(DMGetDS(cdm, &cds));
2124: PetscCall(PetscDSSetConstants(cds, 2, c));
2125: }
2126: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
2128: /* Wait for coordinate creation before doing in-place modification */
2129: PetscCall(DMPlexInterpolateInPlace_Internal(dm));
2130: PetscFunctionReturn(PETSC_SUCCESS);
2131: }
2133: /*@
2134: DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra.
2136: Collective
2138: Input Parameters:
2139: + comm - The communicator for the `DM` object
2140: - periodicZ - The boundary type for the Z direction
2142: Output Parameter:
2143: . dm - The `DM` object
2145: Level: beginner
2147: Note:
2148: Here is the output numbering looking from the bottom of the cylinder\:
2149: .vb
2150: 17-----14
2151: | |
2152: | 2 |
2153: | |
2154: 17-----8-----7-----14
2155: | | | |
2156: | 3 | 0 | 1 |
2157: | | | |
2158: 19-----5-----6-----13
2159: | |
2160: | 4 |
2161: | |
2162: 19-----13
2164: and up through the top
2166: 18-----16
2167: | |
2168: | 2 |
2169: | |
2170: 18----10----11-----16
2171: | | | |
2172: | 3 | 0 | 1 |
2173: | | | |
2174: 20-----9----12-----15
2175: | |
2176: | 4 |
2177: | |
2178: 20-----15
2179: .ve
2181: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()`
2182: @*/
2183: PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, DMBoundaryType periodicZ, DM *dm)
2184: {
2185: PetscFunctionBegin;
2186: PetscAssertPointer(dm, 3);
2187: PetscCall(DMCreate(comm, dm));
2188: PetscCall(DMSetType(*dm, DMPLEX));
2189: PetscCall(DMPlexCreateHexCylinderMesh_Internal(*dm, periodicZ));
2190: PetscFunctionReturn(PETSC_SUCCESS);
2191: }
2193: static PetscErrorCode DMPlexCreateWedgeCylinderMesh_Internal(DM dm, PetscInt n, PetscBool interpolate)
2194: {
2195: const PetscInt dim = 3;
2196: PetscInt numCells, numVertices, v;
2197: PetscMPIInt rank;
2199: PetscFunctionBegin;
2200: PetscCheck(n >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %" PetscInt_FMT " cannot be negative", n);
2201: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
2202: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
2203: PetscCall(DMSetDimension(dm, dim));
2204: /* Must create the celltype label here so that we do not automatically try to compute the types */
2205: PetscCall(DMCreateLabel(dm, "celltype"));
2206: /* Create topology */
2207: {
2208: PetscInt cone[6], c;
2210: numCells = rank == 0 ? n : 0;
2211: numVertices = rank == 0 ? 2 * (n + 1) : 0;
2212: PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices));
2213: for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 6));
2214: PetscCall(DMSetUp(dm));
2215: for (c = 0; c < numCells; c++) {
2216: cone[0] = c + n * 1;
2217: cone[1] = (c + 1) % n + n * 1;
2218: cone[2] = 0 + 3 * n;
2219: cone[3] = c + n * 2;
2220: cone[4] = (c + 1) % n + n * 2;
2221: cone[5] = 1 + 3 * n;
2222: PetscCall(DMPlexSetCone(dm, c, cone));
2223: PetscCall(DMPlexSetCellType(dm, c, DM_POLYTOPE_TRI_PRISM_TENSOR));
2224: }
2225: PetscCall(DMPlexSymmetrize(dm));
2226: PetscCall(DMPlexStratify(dm));
2227: }
2228: for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(dm, v, DM_POLYTOPE_POINT));
2229: /* Create cylinder geometry */
2230: {
2231: Vec coordinates;
2232: PetscSection coordSection;
2233: PetscScalar *coords;
2234: PetscInt coordSize, c;
2236: /* Build coordinates */
2237: PetscCall(DMGetCoordinateSection(dm, &coordSection));
2238: PetscCall(PetscSectionSetNumFields(coordSection, 1));
2239: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim));
2240: PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
2241: for (v = numCells; v < numCells + numVertices; ++v) {
2242: PetscCall(PetscSectionSetDof(coordSection, v, dim));
2243: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim));
2244: }
2245: PetscCall(PetscSectionSetUp(coordSection));
2246: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
2247: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
2248: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
2249: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
2250: PetscCall(VecSetBlockSize(coordinates, dim));
2251: PetscCall(VecSetType(coordinates, VECSTANDARD));
2252: PetscCall(VecGetArray(coordinates, &coords));
2253: for (c = 0; c < numCells; c++) {
2254: coords[(c + 0 * n) * dim + 0] = PetscCosReal(2.0 * c * PETSC_PI / n);
2255: coords[(c + 0 * n) * dim + 1] = PetscSinReal(2.0 * c * PETSC_PI / n);
2256: coords[(c + 0 * n) * dim + 2] = 1.0;
2257: coords[(c + 1 * n) * dim + 0] = PetscCosReal(2.0 * c * PETSC_PI / n);
2258: coords[(c + 1 * n) * dim + 1] = PetscSinReal(2.0 * c * PETSC_PI / n);
2259: coords[(c + 1 * n) * dim + 2] = 0.0;
2260: }
2261: if (rank == 0) {
2262: coords[(2 * n + 0) * dim + 0] = 0.0;
2263: coords[(2 * n + 0) * dim + 1] = 0.0;
2264: coords[(2 * n + 0) * dim + 2] = 1.0;
2265: coords[(2 * n + 1) * dim + 0] = 0.0;
2266: coords[(2 * n + 1) * dim + 1] = 0.0;
2267: coords[(2 * n + 1) * dim + 2] = 0.0;
2268: }
2269: PetscCall(VecRestoreArray(coordinates, &coords));
2270: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
2271: PetscCall(VecDestroy(&coordinates));
2272: }
2273: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
2274: /* Interpolate */
2275: if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(dm));
2276: PetscFunctionReturn(PETSC_SUCCESS);
2277: }
2279: /*@
2280: DMPlexCreateWedgeCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using wedges.
2282: Collective
2284: Input Parameters:
2285: + comm - The communicator for the `DM` object
2286: . n - The number of wedges around the origin
2287: - interpolate - Create edges and faces
2289: Output Parameter:
2290: . dm - The `DM` object
2292: Level: beginner
2294: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateHexCylinderMesh()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()`
2295: @*/
2296: PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm)
2297: {
2298: PetscFunctionBegin;
2299: PetscAssertPointer(dm, 4);
2300: PetscCall(DMCreate(comm, dm));
2301: PetscCall(DMSetType(*dm, DMPLEX));
2302: PetscCall(DMPlexCreateWedgeCylinderMesh_Internal(*dm, n, interpolate));
2303: PetscFunctionReturn(PETSC_SUCCESS);
2304: }
2306: static inline PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
2307: {
2308: PetscReal prod = 0.0;
2309: PetscInt i;
2310: for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]);
2311: return PetscSqrtReal(prod);
2312: }
2313: static inline PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[])
2314: {
2315: PetscReal prod = 0.0;
2316: PetscInt i;
2317: for (i = 0; i < dim; ++i) prod += x[i] * y[i];
2318: return prod;
2319: }
2321: /* The first constant is the sphere radius */
2322: static void snapToSphere(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
2323: {
2324: PetscReal r = PetscRealPart(constants[0]);
2325: PetscReal norm2 = 0.0, fac;
2326: PetscInt n = uOff[1] - uOff[0], d;
2328: for (d = 0; d < n; ++d) norm2 += PetscSqr(PetscRealPart(u[d]));
2329: fac = r / PetscSqrtReal(norm2);
2330: for (d = 0; d < n; ++d) f0[d] = u[d] * fac;
2331: }
2333: static PetscErrorCode DMPlexCreateSphereMesh_Internal(DM dm, PetscInt dim, PetscBool simplex, PetscReal R)
2334: {
2335: const PetscInt embedDim = dim + 1;
2336: PetscSection coordSection;
2337: Vec coordinates;
2338: PetscScalar *coords;
2339: PetscReal *coordsIn;
2340: PetscInt numCells, numEdges, numVerts = 0, firstVertex = 0, v, firstEdge, coordSize, d, c, e;
2341: PetscMPIInt rank;
2343: PetscFunctionBegin;
2345: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
2346: PetscCall(DMSetDimension(dm, dim));
2347: PetscCall(DMSetCoordinateDim(dm, dim + 1));
2348: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
2349: switch (dim) {
2350: case 2:
2351: if (simplex) {
2352: const PetscReal radius = PetscSqrtReal(1 + PETSC_PHI * PETSC_PHI) / (1.0 + PETSC_PHI);
2353: const PetscReal edgeLen = 2.0 / (1.0 + PETSC_PHI) * (R / radius);
2354: const PetscInt degree = 5;
2355: PetscReal vertex[3] = {0.0, 1.0 / (1.0 + PETSC_PHI), PETSC_PHI / (1.0 + PETSC_PHI)};
2356: PetscInt s[3] = {1, 1, 1};
2357: PetscInt cone[3];
2358: PetscInt *graph, p, i, j, k;
2360: vertex[0] *= R / radius;
2361: vertex[1] *= R / radius;
2362: vertex[2] *= R / radius;
2363: numCells = rank == 0 ? 20 : 0;
2364: numVerts = rank == 0 ? 12 : 0;
2365: firstVertex = numCells;
2366: /* Use icosahedron, which for a R-sphere has coordinates which are all cyclic permutations of
2368: (0, \pm 1/\phi+1, \pm \phi/\phi+1)
2370: where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
2371: length is then given by 2/(1+\phi) = 2 * 0.38197 = 0.76393.
2372: */
2373: /* Construct vertices */
2374: PetscCall(PetscCalloc1(numVerts * embedDim, &coordsIn));
2375: if (rank == 0) {
2376: for (p = 0, i = 0; p < embedDim; ++p) {
2377: for (s[1] = -1; s[1] < 2; s[1] += 2) {
2378: for (s[2] = -1; s[2] < 2; s[2] += 2) {
2379: for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[(d + p) % embedDim] * vertex[(d + p) % embedDim];
2380: ++i;
2381: }
2382: }
2383: }
2384: }
2385: /* Construct graph */
2386: PetscCall(PetscCalloc1(numVerts * numVerts, &graph));
2387: for (i = 0; i < numVerts; ++i) {
2388: for (j = 0, k = 0; j < numVerts; ++j) {
2389: if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i * embedDim], &coordsIn[j * embedDim]) - edgeLen) < PETSC_SMALL) {
2390: graph[i * numVerts + j] = 1;
2391: ++k;
2392: }
2393: }
2394: PetscCheck(k == degree, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid icosahedron, vertex %" PetscInt_FMT " degree %" PetscInt_FMT " != %" PetscInt_FMT, i, k, degree);
2395: }
2396: /* Build Topology */
2397: PetscCall(DMPlexSetChart(dm, 0, numCells + numVerts));
2398: for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, embedDim));
2399: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
2400: /* Cells */
2401: for (i = 0, c = 0; i < numVerts; ++i) {
2402: for (j = 0; j < i; ++j) {
2403: for (k = 0; k < j; ++k) {
2404: if (graph[i * numVerts + j] && graph[j * numVerts + k] && graph[k * numVerts + i]) {
2405: cone[0] = firstVertex + i;
2406: cone[1] = firstVertex + j;
2407: cone[2] = firstVertex + k;
2408: /* Check orientation */
2409: {
2410: const PetscInt epsilon[3][3][3] = {
2411: {{0, 0, 0}, {0, 0, 1}, {0, -1, 0}},
2412: {{0, 0, -1}, {0, 0, 0}, {1, 0, 0} },
2413: {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0} }
2414: };
2415: PetscReal normal[3];
2416: PetscInt e, f;
2418: for (d = 0; d < embedDim; ++d) {
2419: normal[d] = 0.0;
2420: for (e = 0; e < embedDim; ++e) {
2421: for (f = 0; f < embedDim; ++f) normal[d] += epsilon[d][e][f] * (coordsIn[j * embedDim + e] - coordsIn[i * embedDim + e]) * (coordsIn[k * embedDim + f] - coordsIn[i * embedDim + f]);
2422: }
2423: }
2424: if (DotReal(embedDim, normal, &coordsIn[i * embedDim]) < 0) {
2425: PetscInt tmp = cone[1];
2426: cone[1] = cone[2];
2427: cone[2] = tmp;
2428: }
2429: }
2430: PetscCall(DMPlexSetCone(dm, c++, cone));
2431: }
2432: }
2433: }
2434: }
2435: PetscCall(DMPlexSymmetrize(dm));
2436: PetscCall(DMPlexStratify(dm));
2437: PetscCall(PetscFree(graph));
2438: } else {
2439: /*
2440: 12-21--13
2441: | |
2442: 25 4 24
2443: | |
2444: 12-25--9-16--8-24--13
2445: | | | |
2446: 23 5 17 0 15 3 22
2447: | | | |
2448: 10-20--6-14--7-19--11
2449: | |
2450: 20 1 19
2451: | |
2452: 10-18--11
2453: | |
2454: 23 2 22
2455: | |
2456: 12-21--13
2457: */
2458: PetscInt cone[4], ornt[4];
2460: numCells = rank == 0 ? 6 : 0;
2461: numEdges = rank == 0 ? 12 : 0;
2462: numVerts = rank == 0 ? 8 : 0;
2463: firstVertex = numCells;
2464: firstEdge = numCells + numVerts;
2465: /* Build Topology */
2466: PetscCall(DMPlexSetChart(dm, 0, numCells + numEdges + numVerts));
2467: for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 4));
2468: for (e = firstEdge; e < firstEdge + numEdges; ++e) PetscCall(DMPlexSetConeSize(dm, e, 2));
2469: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
2470: if (rank == 0) {
2471: /* Cell 0 */
2472: cone[0] = 14;
2473: cone[1] = 15;
2474: cone[2] = 16;
2475: cone[3] = 17;
2476: PetscCall(DMPlexSetCone(dm, 0, cone));
2477: ornt[0] = 0;
2478: ornt[1] = 0;
2479: ornt[2] = 0;
2480: ornt[3] = 0;
2481: PetscCall(DMPlexSetConeOrientation(dm, 0, ornt));
2482: /* Cell 1 */
2483: cone[0] = 18;
2484: cone[1] = 19;
2485: cone[2] = 14;
2486: cone[3] = 20;
2487: PetscCall(DMPlexSetCone(dm, 1, cone));
2488: ornt[0] = 0;
2489: ornt[1] = 0;
2490: ornt[2] = -1;
2491: ornt[3] = 0;
2492: PetscCall(DMPlexSetConeOrientation(dm, 1, ornt));
2493: /* Cell 2 */
2494: cone[0] = 21;
2495: cone[1] = 22;
2496: cone[2] = 18;
2497: cone[3] = 23;
2498: PetscCall(DMPlexSetCone(dm, 2, cone));
2499: ornt[0] = 0;
2500: ornt[1] = 0;
2501: ornt[2] = -1;
2502: ornt[3] = 0;
2503: PetscCall(DMPlexSetConeOrientation(dm, 2, ornt));
2504: /* Cell 3 */
2505: cone[0] = 19;
2506: cone[1] = 22;
2507: cone[2] = 24;
2508: cone[3] = 15;
2509: PetscCall(DMPlexSetCone(dm, 3, cone));
2510: ornt[0] = -1;
2511: ornt[1] = -1;
2512: ornt[2] = 0;
2513: ornt[3] = -1;
2514: PetscCall(DMPlexSetConeOrientation(dm, 3, ornt));
2515: /* Cell 4 */
2516: cone[0] = 16;
2517: cone[1] = 24;
2518: cone[2] = 21;
2519: cone[3] = 25;
2520: PetscCall(DMPlexSetCone(dm, 4, cone));
2521: ornt[0] = -1;
2522: ornt[1] = -1;
2523: ornt[2] = -1;
2524: ornt[3] = 0;
2525: PetscCall(DMPlexSetConeOrientation(dm, 4, ornt));
2526: /* Cell 5 */
2527: cone[0] = 20;
2528: cone[1] = 17;
2529: cone[2] = 25;
2530: cone[3] = 23;
2531: PetscCall(DMPlexSetCone(dm, 5, cone));
2532: ornt[0] = -1;
2533: ornt[1] = -1;
2534: ornt[2] = -1;
2535: ornt[3] = -1;
2536: PetscCall(DMPlexSetConeOrientation(dm, 5, ornt));
2537: /* Edges */
2538: cone[0] = 6;
2539: cone[1] = 7;
2540: PetscCall(DMPlexSetCone(dm, 14, cone));
2541: cone[0] = 7;
2542: cone[1] = 8;
2543: PetscCall(DMPlexSetCone(dm, 15, cone));
2544: cone[0] = 8;
2545: cone[1] = 9;
2546: PetscCall(DMPlexSetCone(dm, 16, cone));
2547: cone[0] = 9;
2548: cone[1] = 6;
2549: PetscCall(DMPlexSetCone(dm, 17, cone));
2550: cone[0] = 10;
2551: cone[1] = 11;
2552: PetscCall(DMPlexSetCone(dm, 18, cone));
2553: cone[0] = 11;
2554: cone[1] = 7;
2555: PetscCall(DMPlexSetCone(dm, 19, cone));
2556: cone[0] = 6;
2557: cone[1] = 10;
2558: PetscCall(DMPlexSetCone(dm, 20, cone));
2559: cone[0] = 12;
2560: cone[1] = 13;
2561: PetscCall(DMPlexSetCone(dm, 21, cone));
2562: cone[0] = 13;
2563: cone[1] = 11;
2564: PetscCall(DMPlexSetCone(dm, 22, cone));
2565: cone[0] = 10;
2566: cone[1] = 12;
2567: PetscCall(DMPlexSetCone(dm, 23, cone));
2568: cone[0] = 13;
2569: cone[1] = 8;
2570: PetscCall(DMPlexSetCone(dm, 24, cone));
2571: cone[0] = 12;
2572: cone[1] = 9;
2573: PetscCall(DMPlexSetCone(dm, 25, cone));
2574: }
2575: PetscCall(DMPlexSymmetrize(dm));
2576: PetscCall(DMPlexStratify(dm));
2577: /* Build coordinates */
2578: PetscCall(PetscCalloc1(numVerts * embedDim, &coordsIn));
2579: if (rank == 0) {
2580: coordsIn[0 * embedDim + 0] = -R;
2581: coordsIn[0 * embedDim + 1] = R;
2582: coordsIn[0 * embedDim + 2] = -R;
2583: coordsIn[1 * embedDim + 0] = R;
2584: coordsIn[1 * embedDim + 1] = R;
2585: coordsIn[1 * embedDim + 2] = -R;
2586: coordsIn[2 * embedDim + 0] = R;
2587: coordsIn[2 * embedDim + 1] = -R;
2588: coordsIn[2 * embedDim + 2] = -R;
2589: coordsIn[3 * embedDim + 0] = -R;
2590: coordsIn[3 * embedDim + 1] = -R;
2591: coordsIn[3 * embedDim + 2] = -R;
2592: coordsIn[4 * embedDim + 0] = -R;
2593: coordsIn[4 * embedDim + 1] = R;
2594: coordsIn[4 * embedDim + 2] = R;
2595: coordsIn[5 * embedDim + 0] = R;
2596: coordsIn[5 * embedDim + 1] = R;
2597: coordsIn[5 * embedDim + 2] = R;
2598: coordsIn[6 * embedDim + 0] = -R;
2599: coordsIn[6 * embedDim + 1] = -R;
2600: coordsIn[6 * embedDim + 2] = R;
2601: coordsIn[7 * embedDim + 0] = R;
2602: coordsIn[7 * embedDim + 1] = -R;
2603: coordsIn[7 * embedDim + 2] = R;
2604: }
2605: }
2606: break;
2607: case 3:
2608: if (simplex) {
2609: const PetscReal edgeLen = 1.0 / PETSC_PHI;
2610: PetscReal vertexA[4] = {0.5, 0.5, 0.5, 0.5};
2611: PetscReal vertexB[4] = {1.0, 0.0, 0.0, 0.0};
2612: PetscReal vertexC[4] = {0.5, 0.5 * PETSC_PHI, 0.5 / PETSC_PHI, 0.0};
2613: const PetscInt degree = 12;
2614: PetscInt s[4] = {1, 1, 1};
2615: PetscInt evenPerm[12][4] = {
2616: {0, 1, 2, 3},
2617: {0, 2, 3, 1},
2618: {0, 3, 1, 2},
2619: {1, 0, 3, 2},
2620: {1, 2, 0, 3},
2621: {1, 3, 2, 0},
2622: {2, 0, 1, 3},
2623: {2, 1, 3, 0},
2624: {2, 3, 0, 1},
2625: {3, 0, 2, 1},
2626: {3, 1, 0, 2},
2627: {3, 2, 1, 0}
2628: };
2629: PetscInt cone[4];
2630: PetscInt *graph, p, i, j, k, l;
2632: vertexA[0] *= R;
2633: vertexA[1] *= R;
2634: vertexA[2] *= R;
2635: vertexA[3] *= R;
2636: vertexB[0] *= R;
2637: vertexB[1] *= R;
2638: vertexB[2] *= R;
2639: vertexB[3] *= R;
2640: vertexC[0] *= R;
2641: vertexC[1] *= R;
2642: vertexC[2] *= R;
2643: vertexC[3] *= R;
2644: numCells = rank == 0 ? 600 : 0;
2645: numVerts = rank == 0 ? 120 : 0;
2646: firstVertex = numCells;
2647: /* Use the 600-cell, which for a unit sphere has coordinates which are
2649: 1/2 (\pm 1, \pm 1, \pm 1, \pm 1) 16
2650: (\pm 1, 0, 0, 0) all cyclic permutations 8
2651: 1/2 (\pm 1, \pm phi, \pm 1/phi, 0) all even permutations 96
2653: where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge
2654: length is then given by 1/\phi = 0.61803.
2656: http://buzzard.pugetsound.edu/sage-practice/ch03s03.html
2657: http://mathworld.wolfram.com/600-Cell.html
2658: */
2659: /* Construct vertices */
2660: PetscCall(PetscCalloc1(numVerts * embedDim, &coordsIn));
2661: i = 0;
2662: if (rank == 0) {
2663: for (s[0] = -1; s[0] < 2; s[0] += 2) {
2664: for (s[1] = -1; s[1] < 2; s[1] += 2) {
2665: for (s[2] = -1; s[2] < 2; s[2] += 2) {
2666: for (s[3] = -1; s[3] < 2; s[3] += 2) {
2667: for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[d] * vertexA[d];
2668: ++i;
2669: }
2670: }
2671: }
2672: }
2673: for (p = 0; p < embedDim; ++p) {
2674: s[1] = s[2] = s[3] = 1;
2675: for (s[0] = -1; s[0] < 2; s[0] += 2) {
2676: for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[(d + p) % embedDim] * vertexB[(d + p) % embedDim];
2677: ++i;
2678: }
2679: }
2680: for (p = 0; p < 12; ++p) {
2681: s[3] = 1;
2682: for (s[0] = -1; s[0] < 2; s[0] += 2) {
2683: for (s[1] = -1; s[1] < 2; s[1] += 2) {
2684: for (s[2] = -1; s[2] < 2; s[2] += 2) {
2685: for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[evenPerm[p][d]] * vertexC[evenPerm[p][d]];
2686: ++i;
2687: }
2688: }
2689: }
2690: }
2691: }
2692: PetscCheck(i == numVerts, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid 600-cell, vertices %" PetscInt_FMT " != %" PetscInt_FMT, i, numVerts);
2693: /* Construct graph */
2694: PetscCall(PetscCalloc1(numVerts * numVerts, &graph));
2695: for (i = 0; i < numVerts; ++i) {
2696: for (j = 0, k = 0; j < numVerts; ++j) {
2697: if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i * embedDim], &coordsIn[j * embedDim]) - edgeLen) < PETSC_SMALL) {
2698: graph[i * numVerts + j] = 1;
2699: ++k;
2700: }
2701: }
2702: PetscCheck(k == degree, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid 600-cell, vertex %" PetscInt_FMT " degree %" PetscInt_FMT " != %" PetscInt_FMT, i, k, degree);
2703: }
2704: /* Build Topology */
2705: PetscCall(DMPlexSetChart(dm, 0, numCells + numVerts));
2706: for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, embedDim));
2707: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
2708: /* Cells */
2709: if (rank == 0) {
2710: for (i = 0, c = 0; i < numVerts; ++i) {
2711: for (j = 0; j < i; ++j) {
2712: for (k = 0; k < j; ++k) {
2713: for (l = 0; l < k; ++l) {
2714: if (graph[i * numVerts + j] && graph[j * numVerts + k] && graph[k * numVerts + i] && graph[l * numVerts + i] && graph[l * numVerts + j] && graph[l * numVerts + k]) {
2715: cone[0] = firstVertex + i;
2716: cone[1] = firstVertex + j;
2717: cone[2] = firstVertex + k;
2718: cone[3] = firstVertex + l;
2719: /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */
2720: {
2721: const PetscInt epsilon[4][4][4][4] = {
2722: {{{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 1}, {0, 0, -1, 0}}, {{0, 0, 0, 0}, {0, 0, 0, -1}, {0, 0, 0, 0}, {0, 1, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 1, 0}, {0, -1, 0, 0}, {0, 0, 0, 0}}},
2724: {{{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, -1}, {0, 0, 1, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 0, 1}, {0, 0, 0, 0}, {0, 0, 0, 0}, {-1, 0, 0, 0}}, {{0, 0, -1, 0}, {0, 0, 0, 0}, {1, 0, 0, 0}, {0, 0, 0, 0}}},
2726: {{{0, 0, 0, 0}, {0, 0, 0, 1}, {0, 0, 0, 0}, {0, -1, 0, 0}}, {{0, 0, 0, -1}, {0, 0, 0, 0}, {0, 0, 0, 0}, {1, 0, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 1, 0, 0}, {-1, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}},
2728: {{{0, 0, 0, 0}, {0, 0, -1, 0}, {0, 1, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 1, 0}, {0, 0, 0, 0}, {-1, 0, 0, 0}, {0, 0, 0, 0}}, {{0, -1, 0, 0}, {1, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}} }
2729: };
2730: PetscReal normal[4];
2731: PetscInt e, f, g;
2733: for (d = 0; d < embedDim; ++d) {
2734: normal[d] = 0.0;
2735: for (e = 0; e < embedDim; ++e) {
2736: for (f = 0; f < embedDim; ++f) {
2737: for (g = 0; g < embedDim; ++g) {
2738: normal[d] += epsilon[d][e][f][g] * (coordsIn[j * embedDim + e] - coordsIn[i * embedDim + e]) * (coordsIn[k * embedDim + f] - coordsIn[i * embedDim + f]) * (coordsIn[l * embedDim + f] - coordsIn[i * embedDim + f]);
2739: }
2740: }
2741: }
2742: }
2743: if (DotReal(embedDim, normal, &coordsIn[i * embedDim]) < 0) {
2744: PetscInt tmp = cone[1];
2745: cone[1] = cone[2];
2746: cone[2] = tmp;
2747: }
2748: }
2749: PetscCall(DMPlexSetCone(dm, c++, cone));
2750: }
2751: }
2752: }
2753: }
2754: }
2755: }
2756: PetscCall(DMPlexSymmetrize(dm));
2757: PetscCall(DMPlexStratify(dm));
2758: PetscCall(PetscFree(graph));
2759: }
2760: break;
2761: default:
2762: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension for sphere: %" PetscInt_FMT, dim);
2763: }
2764: /* Create coordinates */
2765: PetscCall(DMGetCoordinateSection(dm, &coordSection));
2766: PetscCall(PetscSectionSetNumFields(coordSection, 1));
2767: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, embedDim));
2768: PetscCall(PetscSectionSetChart(coordSection, firstVertex, firstVertex + numVerts));
2769: for (v = firstVertex; v < firstVertex + numVerts; ++v) {
2770: PetscCall(PetscSectionSetDof(coordSection, v, embedDim));
2771: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, embedDim));
2772: }
2773: PetscCall(PetscSectionSetUp(coordSection));
2774: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
2775: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
2776: PetscCall(VecSetBlockSize(coordinates, embedDim));
2777: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
2778: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
2779: PetscCall(VecSetType(coordinates, VECSTANDARD));
2780: PetscCall(VecGetArray(coordinates, &coords));
2781: for (v = 0; v < numVerts; ++v)
2782: for (d = 0; d < embedDim; ++d) coords[v * embedDim + d] = coordsIn[v * embedDim + d];
2783: PetscCall(VecRestoreArray(coordinates, &coords));
2784: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
2785: PetscCall(VecDestroy(&coordinates));
2786: PetscCall(PetscFree(coordsIn));
2787: {
2788: DM cdm;
2789: PetscDS cds;
2790: PetscScalar c = R;
2792: PetscCall(DMPlexCreateCoordinateSpace(dm, 1, snapToSphere));
2793: PetscCall(DMGetCoordinateDM(dm, &cdm));
2794: PetscCall(DMGetDS(cdm, &cds));
2795: PetscCall(PetscDSSetConstants(cds, 1, &c));
2796: }
2797: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
2798: /* Wait for coordinate creation before doing in-place modification */
2799: if (simplex) PetscCall(DMPlexInterpolateInPlace_Internal(dm));
2800: PetscFunctionReturn(PETSC_SUCCESS);
2801: }
2803: typedef void (*TPSEvaluateFunc)(const PetscReal[], PetscReal *, PetscReal[], PetscReal (*)[3]);
2805: /*
2806: The Schwarz P implicit surface is
2808: f(x) = cos(x0) + cos(x1) + cos(x2) = 0
2809: */
2810: static void TPSEvaluate_SchwarzP(const PetscReal y[3], PetscReal *f, PetscReal grad[], PetscReal (*hess)[3])
2811: {
2812: PetscReal c[3] = {PetscCosReal(y[0] * PETSC_PI), PetscCosReal(y[1] * PETSC_PI), PetscCosReal(y[2] * PETSC_PI)};
2813: PetscReal g[3] = {-PetscSinReal(y[0] * PETSC_PI), -PetscSinReal(y[1] * PETSC_PI), -PetscSinReal(y[2] * PETSC_PI)};
2814: f[0] = c[0] + c[1] + c[2];
2815: for (PetscInt i = 0; i < 3; i++) {
2816: grad[i] = PETSC_PI * g[i];
2817: for (PetscInt j = 0; j < 3; j++) hess[i][j] = (i == j) ? -PetscSqr(PETSC_PI) * c[i] : 0.;
2818: }
2819: }
2821: // u[] is a tentative normal on input. Replace with the implicit function gradient in the same direction
2822: static PetscErrorCode TPSExtrudeNormalFunc_SchwarzP(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx)
2823: {
2824: for (PetscInt i = 0; i < 3; i++) u[i] = -PETSC_PI * PetscSinReal(x[i] * PETSC_PI);
2825: return PETSC_SUCCESS;
2826: }
2828: /*
2829: The Gyroid implicit surface is
2831: f(x,y,z) = sin(pi * x) * cos (pi * (y + 1/2)) + sin(pi * (y + 1/2)) * cos(pi * (z + 1/4)) + sin(pi * (z + 1/4)) * cos(pi * x)
2833: */
2834: static void TPSEvaluate_Gyroid(const PetscReal y[3], PetscReal *f, PetscReal grad[], PetscReal (*hess)[3])
2835: {
2836: PetscReal s[3] = {PetscSinReal(PETSC_PI * y[0]), PetscSinReal(PETSC_PI * (y[1] + .5)), PetscSinReal(PETSC_PI * (y[2] + .25))};
2837: PetscReal c[3] = {PetscCosReal(PETSC_PI * y[0]), PetscCosReal(PETSC_PI * (y[1] + .5)), PetscCosReal(PETSC_PI * (y[2] + .25))};
2838: f[0] = s[0] * c[1] + s[1] * c[2] + s[2] * c[0];
2839: grad[0] = PETSC_PI * (c[0] * c[1] - s[2] * s[0]);
2840: grad[1] = PETSC_PI * (c[1] * c[2] - s[0] * s[1]);
2841: grad[2] = PETSC_PI * (c[2] * c[0] - s[1] * s[2]);
2842: hess[0][0] = -PetscSqr(PETSC_PI) * (s[0] * c[1] + s[2] * c[0]);
2843: hess[0][1] = -PetscSqr(PETSC_PI) * (c[0] * s[1]);
2844: hess[0][2] = -PetscSqr(PETSC_PI) * (c[2] * s[0]);
2845: hess[1][0] = -PetscSqr(PETSC_PI) * (s[1] * c[2] + s[0] * c[1]);
2846: hess[1][1] = -PetscSqr(PETSC_PI) * (c[1] * s[2]);
2847: hess[2][2] = -PetscSqr(PETSC_PI) * (c[0] * s[1]);
2848: hess[2][0] = -PetscSqr(PETSC_PI) * (s[2] * c[0] + s[1] * c[2]);
2849: hess[2][1] = -PetscSqr(PETSC_PI) * (c[2] * s[0]);
2850: hess[2][2] = -PetscSqr(PETSC_PI) * (c[1] * s[2]);
2851: }
2853: // u[] is a tentative normal on input. Replace with the implicit function gradient in the same direction
2854: static PetscErrorCode TPSExtrudeNormalFunc_Gyroid(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx)
2855: {
2856: PetscReal s[3] = {PetscSinReal(PETSC_PI * x[0]), PetscSinReal(PETSC_PI * (x[1] + .5)), PetscSinReal(PETSC_PI * (x[2] + .25))};
2857: PetscReal c[3] = {PetscCosReal(PETSC_PI * x[0]), PetscCosReal(PETSC_PI * (x[1] + .5)), PetscCosReal(PETSC_PI * (x[2] + .25))};
2858: u[0] = PETSC_PI * (c[0] * c[1] - s[2] * s[0]);
2859: u[1] = PETSC_PI * (c[1] * c[2] - s[0] * s[1]);
2860: u[2] = PETSC_PI * (c[2] * c[0] - s[1] * s[2]);
2861: return PETSC_SUCCESS;
2862: }
2864: /*
2865: We wish to solve
2867: min_y || y - x ||^2 subject to f(y) = 0
2869: Let g(y) = grad(f). The minimization problem is equivalent to asking to satisfy
2870: f(y) = 0 and (y-x) is parallel to g(y). We do this by using Householder QR to obtain a basis for the
2871: tangent space and ask for both components in the tangent space to be zero.
2873: Take g to be a column vector and compute the "full QR" factorization Q R = g,
2874: where Q = I - 2 n n^T is a symmetric orthogonal matrix.
2875: The first column of Q is parallel to g so the remaining two columns span the null space.
2876: Let Qn = Q[:,1:] be those remaining columns. Then Qn Qn^T is an orthogonal projector into the tangent space.
2877: Since Q is symmetric, this is equivalent to multiplying by Q and taking the last two entries.
2878: In total, we have a system of 3 equations in 3 unknowns:
2880: f(y) = 0 1 equation
2881: Qn^T (y - x) = 0 2 equations
2883: Here, we compute the residual and Jacobian of this system.
2884: */
2885: static void TPSNearestPointResJac(TPSEvaluateFunc feval, const PetscScalar x[], const PetscScalar y[], PetscScalar res[], PetscScalar J[])
2886: {
2887: PetscReal yreal[3] = {PetscRealPart(y[0]), PetscRealPart(y[1]), PetscRealPart(y[2])};
2888: PetscReal d[3] = {PetscRealPart(y[0] - x[0]), PetscRealPart(y[1] - x[1]), PetscRealPart(y[2] - x[2])};
2889: PetscReal f, grad[3], n[3], norm, norm_y[3], nd, nd_y[3], sign;
2890: PetscReal n_y[3][3] = {
2891: {0, 0, 0},
2892: {0, 0, 0},
2893: {0, 0, 0}
2894: };
2896: feval(yreal, &f, grad, n_y);
2898: for (PetscInt i = 0; i < 3; i++) n[i] = grad[i];
2899: norm = PetscSqrtReal(PetscSqr(n[0]) + PetscSqr(n[1]) + PetscSqr(n[2]));
2900: for (PetscInt i = 0; i < 3; i++) norm_y[i] = 1. / norm * n[i] * n_y[i][i];
2902: // Define the Householder reflector
2903: sign = n[0] >= 0 ? 1. : -1.;
2904: n[0] += norm * sign;
2905: for (PetscInt i = 0; i < 3; i++) n_y[0][i] += norm_y[i] * sign;
2907: norm = PetscSqrtReal(PetscSqr(n[0]) + PetscSqr(n[1]) + PetscSqr(n[2]));
2908: norm_y[0] = 1. / norm * (n[0] * n_y[0][0]);
2909: norm_y[1] = 1. / norm * (n[0] * n_y[0][1] + n[1] * n_y[1][1]);
2910: norm_y[2] = 1. / norm * (n[0] * n_y[0][2] + n[2] * n_y[2][2]);
2912: for (PetscInt i = 0; i < 3; i++) {
2913: n[i] /= norm;
2914: for (PetscInt j = 0; j < 3; j++) {
2915: // note that n[i] is n_old[i]/norm when executing the code below
2916: n_y[i][j] = n_y[i][j] / norm - n[i] / norm * norm_y[j];
2917: }
2918: }
2920: nd = n[0] * d[0] + n[1] * d[1] + n[2] * d[2];
2921: for (PetscInt i = 0; i < 3; i++) nd_y[i] = n[i] + n_y[0][i] * d[0] + n_y[1][i] * d[1] + n_y[2][i] * d[2];
2923: res[0] = f;
2924: res[1] = d[1] - 2 * n[1] * nd;
2925: res[2] = d[2] - 2 * n[2] * nd;
2926: // J[j][i] is J_{ij} (column major)
2927: for (PetscInt j = 0; j < 3; j++) {
2928: J[0 + j * 3] = grad[j];
2929: J[1 + j * 3] = (j == 1) * 1. - 2 * (n_y[1][j] * nd + n[1] * nd_y[j]);
2930: J[2 + j * 3] = (j == 2) * 1. - 2 * (n_y[2][j] * nd + n[2] * nd_y[j]);
2931: }
2932: }
2934: /*
2935: Project x to the nearest point on the implicit surface using Newton's method.
2936: */
2937: static PetscErrorCode TPSNearestPoint(TPSEvaluateFunc feval, PetscScalar x[])
2938: {
2939: PetscScalar y[3] = {x[0], x[1], x[2]}; // Initial guess
2941: PetscFunctionBegin;
2942: for (PetscInt iter = 0; iter < 10; iter++) {
2943: PetscScalar res[3], J[9];
2944: PetscReal resnorm;
2945: TPSNearestPointResJac(feval, x, y, res, J);
2946: resnorm = PetscSqrtReal(PetscSqr(PetscRealPart(res[0])) + PetscSqr(PetscRealPart(res[1])) + PetscSqr(PetscRealPart(res[2])));
2947: if (0) { // Turn on this monitor if you need to confirm quadratic convergence
2948: PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT "] res [%g %g %g]\n", iter, (double)PetscRealPart(res[0]), (double)PetscRealPart(res[1]), (double)PetscRealPart(res[2])));
2949: }
2950: if (resnorm < PETSC_SMALL) break;
2952: // Take the Newton step
2953: PetscCall(PetscKernel_A_gets_inverse_A_3(J, 0., PETSC_FALSE, NULL));
2954: PetscKernel_v_gets_v_minus_A_times_w_3(y, J, res);
2955: }
2956: for (PetscInt i = 0; i < 3; i++) x[i] = y[i];
2957: PetscFunctionReturn(PETSC_SUCCESS);
2958: }
2960: const char *const DMPlexTPSTypes[] = {"SCHWARZ_P", "GYROID", "DMPlexTPSType", "DMPLEX_TPS_", NULL};
2962: static PetscErrorCode DMPlexCreateTPSMesh_Internal(DM dm, DMPlexTPSType tpstype, const PetscInt extent[], const DMBoundaryType periodic[], PetscBool tps_distribute, PetscInt refinements, PetscInt layers, PetscReal thickness)
2963: {
2964: PetscMPIInt rank;
2965: PetscInt topoDim = 2, spaceDim = 3, numFaces = 0, numVertices = 0, numEdges = 0;
2966: PetscInt(*edges)[2] = NULL, *edgeSets = NULL;
2967: PetscInt *cells_flat = NULL;
2968: PetscReal *vtxCoords = NULL;
2969: TPSEvaluateFunc evalFunc = NULL;
2970: PetscSimplePointFunc normalFunc = NULL;
2971: DMLabel label;
2973: PetscFunctionBegin;
2974: PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0));
2975: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
2976: PetscCheck((layers != 0) ^ (thickness == 0.), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Layers %" PetscInt_FMT " must be nonzero iff thickness %g is nonzero", layers, (double)thickness);
2977: switch (tpstype) {
2978: case DMPLEX_TPS_SCHWARZ_P:
2979: PetscCheck(!periodic || (periodic[0] == DM_BOUNDARY_NONE && periodic[1] == DM_BOUNDARY_NONE && periodic[2] == DM_BOUNDARY_NONE), PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Schwarz P does not support periodic meshes");
2980: if (rank == 0) {
2981: PetscInt(*cells)[6][4][4] = NULL; // [junction, junction-face, cell, conn]
2982: PetscInt Njunctions = 0, Ncuts = 0, Npipes[3], vcount;
2983: PetscReal L = 1;
2985: Npipes[0] = (extent[0] + 1) * extent[1] * extent[2];
2986: Npipes[1] = extent[0] * (extent[1] + 1) * extent[2];
2987: Npipes[2] = extent[0] * extent[1] * (extent[2] + 1);
2988: Njunctions = extent[0] * extent[1] * extent[2];
2989: Ncuts = 2 * (extent[0] * extent[1] + extent[1] * extent[2] + extent[2] * extent[0]);
2990: numVertices = 4 * (Npipes[0] + Npipes[1] + Npipes[2]) + 8 * Njunctions;
2991: PetscCall(PetscMalloc1(3 * numVertices, &vtxCoords));
2992: PetscCall(PetscMalloc1(Njunctions, &cells));
2993: PetscCall(PetscMalloc1(Ncuts * 4, &edges));
2994: PetscCall(PetscMalloc1(Ncuts * 4, &edgeSets));
2995: // x-normal pipes
2996: vcount = 0;
2997: for (PetscInt i = 0; i < extent[0] + 1; i++) {
2998: for (PetscInt j = 0; j < extent[1]; j++) {
2999: for (PetscInt k = 0; k < extent[2]; k++) {
3000: for (PetscInt l = 0; l < 4; l++) {
3001: vtxCoords[vcount++] = (2 * i - 1) * L;
3002: vtxCoords[vcount++] = 2 * j * L + PetscCosReal((2 * l + 1) * PETSC_PI / 4) * L / 2;
3003: vtxCoords[vcount++] = 2 * k * L + PetscSinReal((2 * l + 1) * PETSC_PI / 4) * L / 2;
3004: }
3005: }
3006: }
3007: }
3008: // y-normal pipes
3009: for (PetscInt i = 0; i < extent[0]; i++) {
3010: for (PetscInt j = 0; j < extent[1] + 1; j++) {
3011: for (PetscInt k = 0; k < extent[2]; k++) {
3012: for (PetscInt l = 0; l < 4; l++) {
3013: vtxCoords[vcount++] = 2 * i * L + PetscSinReal((2 * l + 1) * PETSC_PI / 4) * L / 2;
3014: vtxCoords[vcount++] = (2 * j - 1) * L;
3015: vtxCoords[vcount++] = 2 * k * L + PetscCosReal((2 * l + 1) * PETSC_PI / 4) * L / 2;
3016: }
3017: }
3018: }
3019: }
3020: // z-normal pipes
3021: for (PetscInt i = 0; i < extent[0]; i++) {
3022: for (PetscInt j = 0; j < extent[1]; j++) {
3023: for (PetscInt k = 0; k < extent[2] + 1; k++) {
3024: for (PetscInt l = 0; l < 4; l++) {
3025: vtxCoords[vcount++] = 2 * i * L + PetscCosReal((2 * l + 1) * PETSC_PI / 4) * L / 2;
3026: vtxCoords[vcount++] = 2 * j * L + PetscSinReal((2 * l + 1) * PETSC_PI / 4) * L / 2;
3027: vtxCoords[vcount++] = (2 * k - 1) * L;
3028: }
3029: }
3030: }
3031: }
3032: // junctions
3033: for (PetscInt i = 0; i < extent[0]; i++) {
3034: for (PetscInt j = 0; j < extent[1]; j++) {
3035: for (PetscInt k = 0; k < extent[2]; k++) {
3036: const PetscInt J = (i * extent[1] + j) * extent[2] + k, Jvoff = (Npipes[0] + Npipes[1] + Npipes[2]) * 4 + J * 8;
3037: PetscCheck(vcount / 3 == Jvoff, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected vertex count");
3038: for (PetscInt ii = 0; ii < 2; ii++) {
3039: for (PetscInt jj = 0; jj < 2; jj++) {
3040: for (PetscInt kk = 0; kk < 2; kk++) {
3041: double Ls = (1 - sqrt(2) / 4) * L;
3042: vtxCoords[vcount++] = 2 * i * L + (2 * ii - 1) * Ls;
3043: vtxCoords[vcount++] = 2 * j * L + (2 * jj - 1) * Ls;
3044: vtxCoords[vcount++] = 2 * k * L + (2 * kk - 1) * Ls;
3045: }
3046: }
3047: }
3048: const PetscInt jfaces[3][2][4] = {
3049: {{3, 1, 0, 2}, {7, 5, 4, 6}}, // x-aligned
3050: {{5, 4, 0, 1}, {7, 6, 2, 3}}, // y-aligned
3051: {{6, 2, 0, 4}, {7, 3, 1, 5}} // z-aligned
3052: };
3053: const PetscInt pipe_lo[3] = {// vertex numbers of pipes
3054: ((i * extent[1] + j) * extent[2] + k) * 4, ((i * (extent[1] + 1) + j) * extent[2] + k + Npipes[0]) * 4, ((i * extent[1] + j) * (extent[2] + 1) + k + Npipes[0] + Npipes[1]) * 4};
3055: const PetscInt pipe_hi[3] = {// vertex numbers of pipes
3056: (((i + 1) * extent[1] + j) * extent[2] + k) * 4, ((i * (extent[1] + 1) + j + 1) * extent[2] + k + Npipes[0]) * 4, ((i * extent[1] + j) * (extent[2] + 1) + k + 1 + Npipes[0] + Npipes[1]) * 4};
3057: for (PetscInt dir = 0; dir < 3; dir++) { // x,y,z
3058: const PetscInt ijk[3] = {i, j, k};
3059: for (PetscInt l = 0; l < 4; l++) { // rotations
3060: cells[J][dir * 2 + 0][l][0] = pipe_lo[dir] + l;
3061: cells[J][dir * 2 + 0][l][1] = Jvoff + jfaces[dir][0][l];
3062: cells[J][dir * 2 + 0][l][2] = Jvoff + jfaces[dir][0][(l - 1 + 4) % 4];
3063: cells[J][dir * 2 + 0][l][3] = pipe_lo[dir] + (l - 1 + 4) % 4;
3064: cells[J][dir * 2 + 1][l][0] = Jvoff + jfaces[dir][1][l];
3065: cells[J][dir * 2 + 1][l][1] = pipe_hi[dir] + l;
3066: cells[J][dir * 2 + 1][l][2] = pipe_hi[dir] + (l - 1 + 4) % 4;
3067: cells[J][dir * 2 + 1][l][3] = Jvoff + jfaces[dir][1][(l - 1 + 4) % 4];
3068: if (ijk[dir] == 0) {
3069: edges[numEdges][0] = pipe_lo[dir] + l;
3070: edges[numEdges][1] = pipe_lo[dir] + (l + 1) % 4;
3071: edgeSets[numEdges] = dir * 2 + 1;
3072: numEdges++;
3073: }
3074: if (ijk[dir] + 1 == extent[dir]) {
3075: edges[numEdges][0] = pipe_hi[dir] + l;
3076: edges[numEdges][1] = pipe_hi[dir] + (l + 1) % 4;
3077: edgeSets[numEdges] = dir * 2 + 2;
3078: numEdges++;
3079: }
3080: }
3081: }
3082: }
3083: }
3084: }
3085: PetscCheck(numEdges == Ncuts * 4, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Edge count %" PetscInt_FMT " incompatible with number of cuts %" PetscInt_FMT, numEdges, Ncuts);
3086: numFaces = 24 * Njunctions;
3087: cells_flat = cells[0][0][0];
3088: }
3089: evalFunc = TPSEvaluate_SchwarzP;
3090: normalFunc = TPSExtrudeNormalFunc_SchwarzP;
3091: break;
3092: case DMPLEX_TPS_GYROID:
3093: if (rank == 0) {
3094: // This is a coarse mesh approximation of the gyroid shifted to being the zero of the level set
3095: //
3096: // sin(pi*x)*cos(pi*(y+1/2)) + sin(pi*(y+1/2))*cos(pi*(z+1/4)) + sin(pi*(z+1/4))*cos(x)
3097: //
3098: // on the cell [0,2]^3.
3099: //
3100: // Think about dividing that cell into four columns, and focus on the column [0,1]x[0,1]x[0,2].
3101: // If you looked at the gyroid in that column at different slices of z you would see that it kind of spins
3102: // like a boomerang:
3103: //
3104: // z = 0 z = 1/4 z = 1/2 z = 3/4 //
3105: // ----- ------- ------- ------- //
3106: // //
3107: // + + + + + + + \ + //
3108: // \ / \ //
3109: // \ `-_ _-' / } //
3110: // *-_ `-' _-' / //
3111: // + `-+ + + +-' + + / + //
3112: // //
3113: // //
3114: // z = 1 z = 5/4 z = 3/2 z = 7/4 //
3115: // ----- ------- ------- ------- //
3116: // //
3117: // +-_ + + + + _-+ + / + //
3118: // `-_ _-_ _-` / //
3119: // \ _-' `-_ / { //
3120: // \ / \ //
3121: // + + + + + + + \ + //
3122: //
3123: //
3124: // This course mesh approximates each of these slices by two line segments,
3125: // and then connects the segments in consecutive layers with quadrilateral faces.
3126: // All of the end points of the segments are multiples of 1/4 except for the
3127: // point * in the picture for z = 0 above and the similar points in other layers.
3128: // That point is at (gamma, gamma, 0), where gamma is calculated below.
3129: //
3130: // The column [1,2]x[1,2]x[0,2] looks the same as this column;
3131: // The columns [1,2]x[0,1]x[0,2] and [0,1]x[1,2]x[0,2] are mirror images.
3132: //
3133: // As for how this method turned into the names given to the vertices:
3134: // that was not systematic, it was just the way it worked out in my handwritten notes.
3136: PetscInt facesPerBlock = 64;
3137: PetscInt vertsPerBlock = 56;
3138: PetscInt extentPlus[3];
3139: PetscInt numBlocks, numBlocksPlus;
3140: const PetscInt A = 0, B = 1, C = 2, D = 3, E = 4, F = 5, G = 6, H = 7, II = 8, J = 9, K = 10, L = 11, M = 12, N = 13, O = 14, P = 15, Q = 16, R = 17, S = 18, T = 19, U = 20, V = 21, W = 22, X = 23, Y = 24, Z = 25, Ap = 26, Bp = 27, Cp = 28, Dp = 29, Ep = 30, Fp = 31, Gp = 32, Hp = 33, Ip = 34, Jp = 35, Kp = 36, Lp = 37, Mp = 38, Np = 39, Op = 40, Pp = 41, Qp = 42, Rp = 43, Sp = 44, Tp = 45, Up = 46, Vp = 47, Wp = 48, Xp = 49, Yp = 50, Zp = 51, Aq = 52, Bq = 53, Cq = 54, Dq = 55;
3141: const PetscInt pattern[64][4] = {
3142: /* face to vertex within the coarse discretization of a single gyroid block */
3143: /* layer 0 */
3144: {A, C, K, G },
3145: {C, B, II, K },
3146: {D, A, H, L },
3147: {B + 56 * 1, D, L, J },
3148: {E, B + 56 * 1, J, N },
3149: {A + 56 * 2, E, N, H + 56 * 2 },
3150: {F, A + 56 * 2, G + 56 * 2, M },
3151: {B, F, M, II },
3152: /* layer 1 */
3153: {G, K, Q, O },
3154: {K, II, P, Q },
3155: {L, H, O + 56 * 1, R },
3156: {J, L, R, P },
3157: {N, J, P, S },
3158: {H + 56 * 2, N, S, O + 56 * 3 },
3159: {M, G + 56 * 2, O + 56 * 2, T },
3160: {II, M, T, P },
3161: /* layer 2 */
3162: {O, Q, Y, U },
3163: {Q, P, W, Y },
3164: {R, O + 56 * 1, U + 56 * 1, Ap },
3165: {P, R, Ap, W },
3166: {S, P, X, Bp },
3167: {O + 56 * 3, S, Bp, V + 56 * 1 },
3168: {T, O + 56 * 2, V, Z },
3169: {P, T, Z, X },
3170: /* layer 3 */
3171: {U, Y, Ep, Dp },
3172: {Y, W, Cp, Ep },
3173: {Ap, U + 56 * 1, Dp + 56 * 1, Gp },
3174: {W, Ap, Gp, Cp },
3175: {Bp, X, Cp + 56 * 2, Fp },
3176: {V + 56 * 1, Bp, Fp, Dp + 56 * 1},
3177: {Z, V, Dp, Hp },
3178: {X, Z, Hp, Cp + 56 * 2},
3179: /* layer 4 */
3180: {Dp, Ep, Mp, Kp },
3181: {Ep, Cp, Ip, Mp },
3182: {Gp, Dp + 56 * 1, Lp, Np },
3183: {Cp, Gp, Np, Jp },
3184: {Fp, Cp + 56 * 2, Jp + 56 * 2, Pp },
3185: {Dp + 56 * 1, Fp, Pp, Lp },
3186: {Hp, Dp, Kp, Op },
3187: {Cp + 56 * 2, Hp, Op, Ip + 56 * 2},
3188: /* layer 5 */
3189: {Kp, Mp, Sp, Rp },
3190: {Mp, Ip, Qp, Sp },
3191: {Np, Lp, Rp, Tp },
3192: {Jp, Np, Tp, Qp + 56 * 1},
3193: {Pp, Jp + 56 * 2, Qp + 56 * 3, Up },
3194: {Lp, Pp, Up, Rp },
3195: {Op, Kp, Rp, Vp },
3196: {Ip + 56 * 2, Op, Vp, Qp + 56 * 2},
3197: /* layer 6 */
3198: {Rp, Sp, Aq, Yp },
3199: {Sp, Qp, Wp, Aq },
3200: {Tp, Rp, Yp, Cq },
3201: {Qp + 56 * 1, Tp, Cq, Wp + 56 * 1},
3202: {Up, Qp + 56 * 3, Xp + 56 * 1, Dq },
3203: {Rp, Up, Dq, Zp },
3204: {Vp, Rp, Zp, Bq },
3205: {Qp + 56 * 2, Vp, Bq, Xp },
3206: /* layer 7 (the top is the periodic image of the bottom of layer 0) */
3207: {Yp, Aq, C + 56 * 4, A + 56 * 4 },
3208: {Aq, Wp, B + 56 * 4, C + 56 * 4 },
3209: {Cq, Yp, A + 56 * 4, D + 56 * 4 },
3210: {Wp + 56 * 1, Cq, D + 56 * 4, B + 56 * 5 },
3211: {Dq, Xp + 56 * 1, B + 56 * 5, E + 56 * 4 },
3212: {Zp, Dq, E + 56 * 4, A + 56 * 6 },
3213: {Bq, Zp, A + 56 * 6, F + 56 * 4 },
3214: {Xp, Bq, F + 56 * 4, B + 56 * 4 }
3215: };
3216: const PetscReal gamma = PetscAcosReal((PetscSqrtReal(3.) - 1.) / PetscSqrtReal(2.)) / PETSC_PI;
3217: const PetscReal patternCoords[56][3] = {
3218: {1., 0., 0. }, /* A */
3219: {0., 1., 0. }, /* B */
3220: {gamma, gamma, 0. }, /* C */
3221: {1 + gamma, 1 - gamma, 0. }, /* D */
3222: {2 - gamma, 2 - gamma, 0. }, /* E */
3223: {1 - gamma, 1 + gamma, 0. }, /* F */
3225: {.5, 0, .25 }, /* G */
3226: {1.5, 0., .25 }, /* H */
3227: {.5, 1., .25 }, /* II */
3228: {1.5, 1., .25 }, /* J */
3229: {.25, .5, .25 }, /* K */
3230: {1.25, .5, .25 }, /* L */
3231: {.75, 1.5, .25 }, /* M */
3232: {1.75, 1.5, .25 }, /* N */
3234: {0., 0., .5 }, /* O */
3235: {1., 1., .5 }, /* P */
3236: {gamma, 1 - gamma, .5 }, /* Q */
3237: {1 + gamma, gamma, .5 }, /* R */
3238: {2 - gamma, 1 + gamma, .5 }, /* S */
3239: {1 - gamma, 2 - gamma, .5 }, /* T */
3241: {0., .5, .75 }, /* U */
3242: {0., 1.5, .75 }, /* V */
3243: {1., .5, .75 }, /* W */
3244: {1., 1.5, .75 }, /* X */
3245: {.5, .75, .75 }, /* Y */
3246: {.5, 1.75, .75 }, /* Z */
3247: {1.5, .25, .75 }, /* Ap */
3248: {1.5, 1.25, .75 }, /* Bp */
3250: {1., 0., 1. }, /* Cp */
3251: {0., 1., 1. }, /* Dp */
3252: {1 - gamma, 1 - gamma, 1. }, /* Ep */
3253: {1 + gamma, 1 + gamma, 1. }, /* Fp */
3254: {2 - gamma, gamma, 1. }, /* Gp */
3255: {gamma, 2 - gamma, 1. }, /* Hp */
3257: {.5, 0., 1.25}, /* Ip */
3258: {1.5, 0., 1.25}, /* Jp */
3259: {.5, 1., 1.25}, /* Kp */
3260: {1.5, 1., 1.25}, /* Lp */
3261: {.75, .5, 1.25}, /* Mp */
3262: {1.75, .5, 1.25}, /* Np */
3263: {.25, 1.5, 1.25}, /* Op */
3264: {1.25, 1.5, 1.25}, /* Pp */
3266: {0., 0., 1.5 }, /* Qp */
3267: {1., 1., 1.5 }, /* Rp */
3268: {1 - gamma, gamma, 1.5 }, /* Sp */
3269: {2 - gamma, 1 - gamma, 1.5 }, /* Tp */
3270: {1 + gamma, 2 - gamma, 1.5 }, /* Up */
3271: {gamma, 1 + gamma, 1.5 }, /* Vp */
3273: {0., .5, 1.75}, /* Wp */
3274: {0., 1.5, 1.75}, /* Xp */
3275: {1., .5, 1.75}, /* Yp */
3276: {1., 1.5, 1.75}, /* Zp */
3277: {.5, .25, 1.75}, /* Aq */
3278: {.5, 1.25, 1.75}, /* Bq */
3279: {1.5, .75, 1.75}, /* Cq */
3280: {1.5, 1.75, 1.75}, /* Dq */
3281: };
3282: PetscInt(*cells)[64][4] = NULL;
3283: PetscBool *seen;
3284: PetscInt *vertToTrueVert;
3285: PetscInt count;
3287: for (PetscInt i = 0; i < 3; i++) extentPlus[i] = extent[i] + 1;
3288: numBlocks = 1;
3289: for (PetscInt i = 0; i < 3; i++) numBlocks *= extent[i];
3290: numBlocksPlus = 1;
3291: for (PetscInt i = 0; i < 3; i++) numBlocksPlus *= extentPlus[i];
3292: numFaces = numBlocks * facesPerBlock;
3293: PetscCall(PetscMalloc1(numBlocks, &cells));
3294: PetscCall(PetscCalloc1(numBlocksPlus * vertsPerBlock, &seen));
3295: for (PetscInt k = 0; k < extent[2]; k++) {
3296: for (PetscInt j = 0; j < extent[1]; j++) {
3297: for (PetscInt i = 0; i < extent[0]; i++) {
3298: for (PetscInt f = 0; f < facesPerBlock; f++) {
3299: for (PetscInt v = 0; v < 4; v++) {
3300: PetscInt vertRaw = pattern[f][v];
3301: PetscInt blockidx = vertRaw / 56;
3302: PetscInt patternvert = vertRaw % 56;
3303: PetscInt xplus = (blockidx & 1);
3304: PetscInt yplus = (blockidx & 2) >> 1;
3305: PetscInt zplus = (blockidx & 4) >> 2;
3306: PetscInt zcoord = (periodic && periodic[2] == DM_BOUNDARY_PERIODIC) ? ((k + zplus) % extent[2]) : (k + zplus);
3307: PetscInt ycoord = (periodic && periodic[1] == DM_BOUNDARY_PERIODIC) ? ((j + yplus) % extent[1]) : (j + yplus);
3308: PetscInt xcoord = (periodic && periodic[0] == DM_BOUNDARY_PERIODIC) ? ((i + xplus) % extent[0]) : (i + xplus);
3309: PetscInt vert = ((zcoord * extentPlus[1] + ycoord) * extentPlus[0] + xcoord) * 56 + patternvert;
3311: cells[(k * extent[1] + j) * extent[0] + i][f][v] = vert;
3312: seen[vert] = PETSC_TRUE;
3313: }
3314: }
3315: }
3316: }
3317: }
3318: for (PetscInt i = 0; i < numBlocksPlus * vertsPerBlock; i++)
3319: if (seen[i]) numVertices++;
3320: count = 0;
3321: PetscCall(PetscMalloc1(numBlocksPlus * vertsPerBlock, &vertToTrueVert));
3322: PetscCall(PetscMalloc1(numVertices * 3, &vtxCoords));
3323: for (PetscInt i = 0; i < numBlocksPlus * vertsPerBlock; i++) vertToTrueVert[i] = -1;
3324: for (PetscInt k = 0; k < extentPlus[2]; k++) {
3325: for (PetscInt j = 0; j < extentPlus[1]; j++) {
3326: for (PetscInt i = 0; i < extentPlus[0]; i++) {
3327: for (PetscInt v = 0; v < vertsPerBlock; v++) {
3328: PetscInt vIdx = ((k * extentPlus[1] + j) * extentPlus[0] + i) * vertsPerBlock + v;
3330: if (seen[vIdx]) {
3331: PetscInt thisVert;
3333: vertToTrueVert[vIdx] = thisVert = count++;
3335: for (PetscInt d = 0; d < 3; d++) vtxCoords[3 * thisVert + d] = patternCoords[v][d];
3336: vtxCoords[3 * thisVert + 0] += i * 2;
3337: vtxCoords[3 * thisVert + 1] += j * 2;
3338: vtxCoords[3 * thisVert + 2] += k * 2;
3339: }
3340: }
3341: }
3342: }
3343: }
3344: for (PetscInt i = 0; i < numBlocks; i++) {
3345: for (PetscInt f = 0; f < facesPerBlock; f++) {
3346: for (PetscInt v = 0; v < 4; v++) cells[i][f][v] = vertToTrueVert[cells[i][f][v]];
3347: }
3348: }
3349: PetscCall(PetscFree(vertToTrueVert));
3350: PetscCall(PetscFree(seen));
3351: cells_flat = cells[0][0];
3352: numEdges = 0;
3353: for (PetscInt i = 0; i < numFaces; i++) {
3354: for (PetscInt e = 0; e < 4; e++) {
3355: PetscInt ev[] = {cells_flat[i * 4 + e], cells_flat[i * 4 + ((e + 1) % 4)]};
3356: const PetscReal *evCoords[] = {&vtxCoords[3 * ev[0]], &vtxCoords[3 * ev[1]]};
3358: for (PetscInt d = 0; d < 3; d++) {
3359: if (!periodic || periodic[0] != DM_BOUNDARY_PERIODIC) {
3360: if (evCoords[0][d] == 0. && evCoords[1][d] == 0.) numEdges++;
3361: if (evCoords[0][d] == 2. * extent[d] && evCoords[1][d] == 2. * extent[d]) numEdges++;
3362: }
3363: }
3364: }
3365: }
3366: PetscCall(PetscMalloc1(numEdges, &edges));
3367: PetscCall(PetscMalloc1(numEdges, &edgeSets));
3368: for (PetscInt edge = 0, i = 0; i < numFaces; i++) {
3369: for (PetscInt e = 0; e < 4; e++) {
3370: PetscInt ev[] = {cells_flat[i * 4 + e], cells_flat[i * 4 + ((e + 1) % 4)]};
3371: const PetscReal *evCoords[] = {&vtxCoords[3 * ev[0]], &vtxCoords[3 * ev[1]]};
3373: for (PetscInt d = 0; d < 3; d++) {
3374: if (!periodic || periodic[d] != DM_BOUNDARY_PERIODIC) {
3375: if (evCoords[0][d] == 0. && evCoords[1][d] == 0.) {
3376: edges[edge][0] = ev[0];
3377: edges[edge][1] = ev[1];
3378: edgeSets[edge++] = 2 * d;
3379: }
3380: if (evCoords[0][d] == 2. * extent[d] && evCoords[1][d] == 2. * extent[d]) {
3381: edges[edge][0] = ev[0];
3382: edges[edge][1] = ev[1];
3383: edgeSets[edge++] = 2 * d + 1;
3384: }
3385: }
3386: }
3387: }
3388: }
3389: }
3390: evalFunc = TPSEvaluate_Gyroid;
3391: normalFunc = TPSExtrudeNormalFunc_Gyroid;
3392: break;
3393: }
3395: PetscCall(DMSetDimension(dm, topoDim));
3396: if (rank == 0) PetscCall(DMPlexBuildFromCellList(dm, numFaces, numVertices, 4, cells_flat));
3397: else PetscCall(DMPlexBuildFromCellList(dm, 0, 0, 0, NULL));
3398: PetscCall(PetscFree(cells_flat));
3399: {
3400: DM idm;
3401: PetscCall(DMPlexInterpolate(dm, &idm));
3402: PetscCall(DMPlexReplace_Internal(dm, &idm));
3403: }
3404: if (rank == 0) PetscCall(DMPlexBuildCoordinatesFromCellList(dm, spaceDim, vtxCoords));
3405: else PetscCall(DMPlexBuildCoordinatesFromCellList(dm, spaceDim, NULL));
3406: PetscCall(PetscFree(vtxCoords));
3408: PetscCall(DMCreateLabel(dm, "Face Sets"));
3409: PetscCall(DMGetLabel(dm, "Face Sets", &label));
3410: for (PetscInt e = 0; e < numEdges; e++) {
3411: PetscInt njoin;
3412: const PetscInt *join, verts[] = {numFaces + edges[e][0], numFaces + edges[e][1]};
3413: PetscCall(DMPlexGetJoin(dm, 2, verts, &njoin, &join));
3414: PetscCheck(njoin == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected unique join of vertices %" PetscInt_FMT " and %" PetscInt_FMT, edges[e][0], edges[e][1]);
3415: PetscCall(DMLabelSetValue(label, join[0], edgeSets[e]));
3416: PetscCall(DMPlexRestoreJoin(dm, 2, verts, &njoin, &join));
3417: }
3418: PetscCall(PetscFree(edges));
3419: PetscCall(PetscFree(edgeSets));
3420: if (tps_distribute) {
3421: DM pdm = NULL;
3422: PetscPartitioner part;
3424: PetscCall(DMPlexGetPartitioner(dm, &part));
3425: PetscCall(PetscPartitionerSetFromOptions(part));
3426: PetscCall(DMPlexDistribute(dm, 0, NULL, &pdm));
3427: if (pdm) PetscCall(DMPlexReplace_Internal(dm, &pdm));
3428: // Do not auto-distribute again
3429: PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
3430: }
3432: PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
3433: for (PetscInt refine = 0; refine < refinements; refine++) {
3434: PetscInt m;
3435: DM dmf;
3436: Vec X;
3437: PetscScalar *x;
3438: PetscCall(DMRefine(dm, MPI_COMM_NULL, &dmf));
3439: PetscCall(DMPlexReplace_Internal(dm, &dmf));
3441: PetscCall(DMGetCoordinatesLocal(dm, &X));
3442: PetscCall(VecGetLocalSize(X, &m));
3443: PetscCall(VecGetArray(X, &x));
3444: for (PetscInt i = 0; i < m; i += 3) PetscCall(TPSNearestPoint(evalFunc, &x[i]));
3445: PetscCall(VecRestoreArray(X, &x));
3446: }
3448: // Face Sets has already been propagated to new vertices during refinement; this propagates to the initial vertices.
3449: PetscCall(DMGetLabel(dm, "Face Sets", &label));
3450: PetscCall(DMPlexLabelComplete(dm, label));
3452: PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0));
3454: if (thickness > 0) {
3455: DM edm, cdm, ecdm;
3456: DMPlexTransform tr;
3457: const char *prefix;
3458: PetscOptions options;
3459: // Code from DMPlexExtrude
3460: PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
3461: PetscCall(DMPlexTransformSetDM(tr, dm));
3462: PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE));
3463: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
3464: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
3465: PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
3466: PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
3467: PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers));
3468: PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness));
3469: PetscCall(DMPlexTransformExtrudeSetTensor(tr, PETSC_FALSE));
3470: PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, PETSC_TRUE));
3471: PetscCall(DMPlexTransformExtrudeSetNormalFunction(tr, normalFunc));
3472: PetscCall(DMPlexTransformSetFromOptions(tr));
3473: PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
3474: PetscCall(DMPlexTransformSetUp(tr));
3475: PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_tps_transform_view"));
3476: PetscCall(DMPlexTransformApply(tr, dm, &edm));
3477: PetscCall(DMCopyDisc(dm, edm));
3478: PetscCall(DMGetCoordinateDM(dm, &cdm));
3479: PetscCall(DMGetCoordinateDM(edm, &ecdm));
3480: PetscCall(DMCopyDisc(cdm, ecdm));
3481: PetscCall(DMPlexTransformCreateDiscLabels(tr, edm));
3482: PetscCall(DMPlexTransformDestroy(&tr));
3483: if (edm) {
3484: ((DM_Plex *)edm->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
3485: ((DM_Plex *)edm->data)->printL2 = ((DM_Plex *)dm->data)->printL2;
3486: ((DM_Plex *)edm->data)->printLocate = ((DM_Plex *)dm->data)->printLocate;
3487: }
3488: PetscCall(DMPlexReplace_Internal(dm, &edm));
3489: }
3490: PetscFunctionReturn(PETSC_SUCCESS);
3491: }
3493: /*@
3494: DMPlexCreateTPSMesh - Create a distributed, interpolated mesh of a triply-periodic surface
3496: Collective
3498: Input Parameters:
3499: + comm - The communicator for the `DM` object
3500: . tpstype - Type of triply-periodic surface
3501: . extent - Array of length 3 containing number of periods in each direction
3502: . periodic - array of length 3 with periodicity, or `NULL` for non-periodic
3503: . tps_distribute - Distribute 2D manifold mesh prior to refinement and extrusion (more scalable)
3504: . refinements - Number of factor-of-2 refinements of 2D manifold mesh
3505: . layers - Number of cell layers extruded in normal direction
3506: - thickness - Thickness in normal direction
3508: Output Parameter:
3509: . dm - The `DM` object
3511: Level: beginner
3513: Notes:
3514: This meshes the surface of the Schwarz P or Gyroid surfaces. Schwarz P is is the simplest member of the triply-periodic minimal surfaces.
3515: <https://en.wikipedia.org/wiki/Schwarz_minimal_surface#Schwarz_P_(%22Primitive%22)> and can be cut with "clean" boundaries.
3516: The Gyroid <https://en.wikipedia.org/wiki/Gyroid> is another triply-periodic minimal surface with applications in additive manufacturing; it is much more difficult to "cut" since there are no planes of symmetry.
3517: Our implementation creates a very coarse mesh of the surface and refines (by 4-way splitting) as many times as requested.
3518: On each refinement, all vertices are projected to their nearest point on the surface.
3519: This projection could readily be extended to related surfaces.
3521: See {cite}`maskery2018insights`
3523: The face (edge) sets for the Schwarz P surface are numbered $1(-x), 2(+x), 3(-y), 4(+y), 5(-z), 6(+z)$.
3524: When the mesh is refined, "Face Sets" contain the new vertices (created during refinement).
3525: Use `DMPlexLabelComplete()` to propagate to coarse-level vertices.
3527: Developer Notes:
3528: The Gyroid mesh does not currently mark boundary sets.
3530: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateSphereMesh()`, `DMSetType()`, `DMCreate()`
3531: @*/
3532: PetscErrorCode DMPlexCreateTPSMesh(MPI_Comm comm, DMPlexTPSType tpstype, const PetscInt extent[], const DMBoundaryType periodic[], PetscBool tps_distribute, PetscInt refinements, PetscInt layers, PetscReal thickness, DM *dm)
3533: {
3534: PetscFunctionBegin;
3535: PetscCall(DMCreate(comm, dm));
3536: PetscCall(DMSetType(*dm, DMPLEX));
3537: PetscCall(DMPlexCreateTPSMesh_Internal(*dm, tpstype, extent, periodic, tps_distribute, refinements, layers, thickness));
3538: PetscFunctionReturn(PETSC_SUCCESS);
3539: }
3541: /*@
3542: DMPlexCreateSphereMesh - Creates a mesh on the d-dimensional sphere, S^d.
3544: Collective
3546: Input Parameters:
3547: + comm - The communicator for the `DM` object
3548: . dim - The dimension
3549: . simplex - Use simplices, or tensor product cells
3550: - R - The radius
3552: Output Parameter:
3553: . dm - The `DM` object
3555: Level: beginner
3557: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateBallMesh()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()`
3558: @*/
3559: PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscReal R, DM *dm)
3560: {
3561: PetscFunctionBegin;
3562: PetscAssertPointer(dm, 5);
3563: PetscCall(DMCreate(comm, dm));
3564: PetscCall(DMSetType(*dm, DMPLEX));
3565: PetscCall(DMPlexCreateSphereMesh_Internal(*dm, dim, simplex, R));
3566: PetscFunctionReturn(PETSC_SUCCESS);
3567: }
3569: static PetscErrorCode DMPlexCreateBallMesh_Internal(DM dm, PetscInt dim, PetscReal R)
3570: {
3571: DM sdm, vol;
3572: DMLabel bdlabel;
3574: PetscFunctionBegin;
3575: PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &sdm));
3576: PetscCall(DMSetType(sdm, DMPLEX));
3577: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sdm, "bd_"));
3578: PetscCall(DMPlexCreateSphereMesh_Internal(sdm, dim - 1, PETSC_TRUE, R));
3579: PetscCall(DMSetFromOptions(sdm));
3580: PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view"));
3581: PetscCall(DMPlexGenerate(sdm, NULL, PETSC_TRUE, &vol));
3582: PetscCall(DMDestroy(&sdm));
3583: PetscCall(DMPlexReplace_Internal(dm, &vol));
3584: PetscCall(DMCreateLabel(dm, "marker"));
3585: PetscCall(DMGetLabel(dm, "marker", &bdlabel));
3586: PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, bdlabel));
3587: PetscCall(DMPlexLabelComplete(dm, bdlabel));
3588: PetscFunctionReturn(PETSC_SUCCESS);
3589: }
3591: /*@
3592: DMPlexCreateBallMesh - Creates a simplex mesh on the d-dimensional ball, B^d.
3594: Collective
3596: Input Parameters:
3597: + comm - The communicator for the `DM` object
3598: . dim - The dimension
3599: - R - The radius
3601: Output Parameter:
3602: . dm - The `DM` object
3604: Options Database Key:
3605: . bd_dm_refine - This will refine the surface mesh preserving the sphere geometry
3607: Level: beginner
3609: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateSphereMesh()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()`
3610: @*/
3611: PetscErrorCode DMPlexCreateBallMesh(MPI_Comm comm, PetscInt dim, PetscReal R, DM *dm)
3612: {
3613: PetscFunctionBegin;
3614: PetscCall(DMCreate(comm, dm));
3615: PetscCall(DMSetType(*dm, DMPLEX));
3616: PetscCall(DMPlexCreateBallMesh_Internal(*dm, dim, R));
3617: PetscFunctionReturn(PETSC_SUCCESS);
3618: }
3620: static PetscErrorCode DMPlexCreateReferenceCell_Internal(DM rdm, DMPolytopeType ct)
3621: {
3622: PetscFunctionBegin;
3623: switch (ct) {
3624: case DM_POLYTOPE_POINT: {
3625: PetscInt numPoints[1] = {1};
3626: PetscInt coneSize[1] = {0};
3627: PetscInt cones[1] = {0};
3628: PetscInt coneOrientations[1] = {0};
3629: PetscScalar vertexCoords[1] = {0.0};
3631: PetscCall(DMSetDimension(rdm, 0));
3632: PetscCall(DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3633: } break;
3634: case DM_POLYTOPE_SEGMENT: {
3635: PetscInt numPoints[2] = {2, 1};
3636: PetscInt coneSize[3] = {2, 0, 0};
3637: PetscInt cones[2] = {1, 2};
3638: PetscInt coneOrientations[2] = {0, 0};
3639: PetscScalar vertexCoords[2] = {-1.0, 1.0};
3641: PetscCall(DMSetDimension(rdm, 1));
3642: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3643: } break;
3644: case DM_POLYTOPE_POINT_PRISM_TENSOR: {
3645: PetscInt numPoints[2] = {2, 1};
3646: PetscInt coneSize[3] = {2, 0, 0};
3647: PetscInt cones[2] = {1, 2};
3648: PetscInt coneOrientations[2] = {0, 0};
3649: PetscScalar vertexCoords[2] = {-1.0, 1.0};
3651: PetscCall(DMSetDimension(rdm, 1));
3652: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3653: } break;
3654: case DM_POLYTOPE_TRIANGLE: {
3655: PetscInt numPoints[2] = {3, 1};
3656: PetscInt coneSize[4] = {3, 0, 0, 0};
3657: PetscInt cones[3] = {1, 2, 3};
3658: PetscInt coneOrientations[3] = {0, 0, 0};
3659: PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0};
3661: PetscCall(DMSetDimension(rdm, 2));
3662: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3663: } break;
3664: case DM_POLYTOPE_QUADRILATERAL: {
3665: PetscInt numPoints[2] = {4, 1};
3666: PetscInt coneSize[5] = {4, 0, 0, 0, 0};
3667: PetscInt cones[4] = {1, 2, 3, 4};
3668: PetscInt coneOrientations[4] = {0, 0, 0, 0};
3669: PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0};
3671: PetscCall(DMSetDimension(rdm, 2));
3672: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3673: } break;
3674: case DM_POLYTOPE_SEG_PRISM_TENSOR: {
3675: PetscInt numPoints[2] = {4, 1};
3676: PetscInt coneSize[5] = {4, 0, 0, 0, 0};
3677: PetscInt cones[4] = {1, 2, 3, 4};
3678: PetscInt coneOrientations[4] = {0, 0, 0, 0};
3679: PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0, 1.0, 1.0};
3681: PetscCall(DMSetDimension(rdm, 2));
3682: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3683: } break;
3684: case DM_POLYTOPE_TETRAHEDRON: {
3685: PetscInt numPoints[2] = {4, 1};
3686: PetscInt coneSize[5] = {4, 0, 0, 0, 0};
3687: PetscInt cones[4] = {1, 2, 3, 4};
3688: PetscInt coneOrientations[4] = {0, 0, 0, 0};
3689: 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};
3691: PetscCall(DMSetDimension(rdm, 3));
3692: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3693: } break;
3694: case DM_POLYTOPE_HEXAHEDRON: {
3695: PetscInt numPoints[2] = {8, 1};
3696: PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3697: PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8};
3698: PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3699: 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, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0};
3701: PetscCall(DMSetDimension(rdm, 3));
3702: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3703: } break;
3704: case DM_POLYTOPE_TRI_PRISM: {
3705: PetscInt numPoints[2] = {6, 1};
3706: PetscInt coneSize[7] = {6, 0, 0, 0, 0, 0, 0};
3707: PetscInt cones[6] = {1, 2, 3, 4, 5, 6};
3708: PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3709: PetscScalar vertexCoords[18] = {-1.0, -1.0, -1.0, -1.0, 1.0, -1.0, 1.0, -1.0, -1.0, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, 1.0, 1.0};
3711: PetscCall(DMSetDimension(rdm, 3));
3712: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3713: } break;
3714: case DM_POLYTOPE_TRI_PRISM_TENSOR: {
3715: PetscInt numPoints[2] = {6, 1};
3716: PetscInt coneSize[7] = {6, 0, 0, 0, 0, 0, 0};
3717: PetscInt cones[6] = {1, 2, 3, 4, 5, 6};
3718: PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0};
3719: PetscScalar vertexCoords[18] = {-1.0, -1.0, -1.0, 1.0, -1.0, -1.0, -1.0, 1.0, -1.0, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, 1.0, 1.0};
3721: PetscCall(DMSetDimension(rdm, 3));
3722: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3723: } break;
3724: case DM_POLYTOPE_QUAD_PRISM_TENSOR: {
3725: PetscInt numPoints[2] = {8, 1};
3726: PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0};
3727: PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8};
3728: PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3729: 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, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0};
3731: PetscCall(DMSetDimension(rdm, 3));
3732: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3733: } break;
3734: case DM_POLYTOPE_PYRAMID: {
3735: PetscInt numPoints[2] = {5, 1};
3736: PetscInt coneSize[6] = {5, 0, 0, 0, 0, 0};
3737: PetscInt cones[5] = {1, 2, 3, 4, 5};
3738: PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
3739: 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, 0.0, 0.0, 1.0};
3741: PetscCall(DMSetDimension(rdm, 3));
3742: PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3743: } break;
3744: default:
3745: SETERRQ(PetscObjectComm((PetscObject)rdm), PETSC_ERR_ARG_WRONG, "Cannot create reference cell for cell type %s", DMPolytopeTypes[ct]);
3746: }
3747: {
3748: PetscInt Nv, v;
3750: /* Must create the celltype label here so that we do not automatically try to compute the types */
3751: PetscCall(DMCreateLabel(rdm, "celltype"));
3752: PetscCall(DMPlexSetCellType(rdm, 0, ct));
3753: PetscCall(DMPlexGetChart(rdm, NULL, &Nv));
3754: for (v = 1; v < Nv; ++v) PetscCall(DMPlexSetCellType(rdm, v, DM_POLYTOPE_POINT));
3755: }
3756: PetscCall(DMPlexInterpolateInPlace_Internal(rdm));
3757: PetscCall(PetscObjectSetName((PetscObject)rdm, DMPolytopeTypes[ct]));
3758: PetscFunctionReturn(PETSC_SUCCESS);
3759: }
3761: /*@
3762: DMPlexCreateReferenceCell - Create a `DMPLEX` with the appropriate FEM reference cell
3764: Collective
3766: Input Parameters:
3767: + comm - The communicator
3768: - ct - The cell type of the reference cell
3770: Output Parameter:
3771: . refdm - The reference cell
3773: Level: intermediate
3775: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateBoxMesh()`
3776: @*/
3777: PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, DMPolytopeType ct, DM *refdm)
3778: {
3779: PetscFunctionBegin;
3780: PetscCall(DMCreate(comm, refdm));
3781: PetscCall(DMSetType(*refdm, DMPLEX));
3782: PetscCall(DMPlexCreateReferenceCell_Internal(*refdm, ct));
3783: PetscFunctionReturn(PETSC_SUCCESS);
3784: }
3786: static PetscErrorCode DMPlexCreateBoundaryLabel_Private(DM dm, const char name[])
3787: {
3788: DM plex;
3789: DMLabel label;
3790: PetscBool hasLabel;
3792: PetscFunctionBegin;
3793: PetscCall(DMHasLabel(dm, name, &hasLabel));
3794: if (hasLabel) PetscFunctionReturn(PETSC_SUCCESS);
3795: PetscCall(DMCreateLabel(dm, name));
3796: PetscCall(DMGetLabel(dm, name, &label));
3797: PetscCall(DMConvert(dm, DMPLEX, &plex));
3798: PetscCall(DMPlexMarkBoundaryFaces(plex, 1, label));
3799: PetscCall(DMPlexLabelComplete(plex, label));
3800: PetscCall(DMDestroy(&plex));
3801: PetscFunctionReturn(PETSC_SUCCESS);
3802: }
3804: /*
3805: We use the last coordinate as the radius, the inner radius is lower[dim-1] and the outer radius is upper[dim-1]. Then we map the first coordinate around the circle.
3807: (x, y) -> (r, theta) = (x[1], (x[0] - lower[0]) * 2\pi/(upper[0] - lower[0]))
3808: */
3809: static void boxToAnnulus(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
3810: {
3811: const PetscReal low = PetscRealPart(constants[0]);
3812: const PetscReal upp = PetscRealPart(constants[1]);
3813: const PetscReal r = PetscRealPart(u[1]);
3814: const PetscReal th = 2. * PETSC_PI * (PetscRealPart(u[0]) - low) / (upp - low);
3816: f0[0] = r * PetscCosReal(th);
3817: f0[1] = r * PetscSinReal(th);
3818: }
3820: const char *const DMPlexShapes[] = {"box", "box_surface", "ball", "sphere", "cylinder", "schwarz_p", "gyroid", "doublet", "annulus", "hypercubic", "zbox", "unknown", "DMPlexShape", "DM_SHAPE_", NULL};
3822: static PetscErrorCode DMPlexCreateFromOptions_Internal(PetscOptionItems *PetscOptionsObject, PetscBool *useCoordSpace, DM dm)
3823: {
3824: DMPlexShape shape = DM_SHAPE_BOX;
3825: DMPolytopeType cell = DM_POLYTOPE_TRIANGLE;
3826: PetscInt dim = 2;
3827: PetscBool simplex = PETSC_TRUE, interpolate = PETSC_TRUE, adjCone = PETSC_FALSE, adjClosure = PETSC_TRUE, refDomain = PETSC_FALSE;
3828: PetscBool flg, flg2, fflg, bdfflg, nameflg;
3829: MPI_Comm comm;
3830: char filename[PETSC_MAX_PATH_LEN] = "<unspecified>";
3831: char bdFilename[PETSC_MAX_PATH_LEN] = "<unspecified>";
3832: char plexname[PETSC_MAX_PATH_LEN] = "";
3834: PetscFunctionBegin;
3835: PetscCall(PetscLogEventBegin(DMPLEX_CreateFromOptions, dm, 0, 0, 0));
3836: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3837: /* TODO Turn this into a registration interface */
3838: PetscCall(PetscOptionsString("-dm_plex_filename", "File containing a mesh", "DMPlexCreateFromFile", filename, filename, sizeof(filename), &fflg));
3839: PetscCall(PetscOptionsString("-dm_plex_boundary_filename", "File containing a mesh boundary", "DMPlexCreateFromFile", bdFilename, bdFilename, sizeof(bdFilename), &bdfflg));
3840: PetscCall(PetscOptionsString("-dm_plex_name", "Name of the mesh in the file", "DMPlexCreateFromFile", plexname, plexname, sizeof(plexname), &nameflg));
3841: PetscCall(PetscOptionsEnum("-dm_plex_cell", "Cell shape", "", DMPolytopeTypes, (PetscEnum)cell, (PetscEnum *)&cell, NULL));
3842: PetscCall(PetscOptionsBool("-dm_plex_reference_cell_domain", "Use a reference cell domain", "", refDomain, &refDomain, NULL));
3843: PetscCall(PetscOptionsEnum("-dm_plex_shape", "Shape for built-in mesh", "", DMPlexShapes, (PetscEnum)shape, (PetscEnum *)&shape, &flg));
3844: PetscCall(PetscOptionsBoundedInt("-dm_plex_dim", "Topological dimension of the mesh", "DMGetDimension", dim, &dim, &flg, 0));
3845: PetscCheck(dim >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " should be in [0, infinity)", dim);
3846: PetscCall(PetscOptionsBool("-dm_plex_simplex", "Mesh cell shape", "", simplex, &simplex, &flg));
3847: PetscCall(PetscOptionsBool("-dm_plex_interpolate", "Flag to create edges and faces automatically", "", interpolate, &interpolate, &flg));
3848: PetscCall(PetscOptionsBool("-dm_plex_adj_cone", "Set adjacency direction", "DMSetBasicAdjacency", adjCone, &adjCone, &flg));
3849: PetscCall(PetscOptionsBool("-dm_plex_adj_closure", "Set adjacency size", "DMSetBasicAdjacency", adjClosure, &adjClosure, &flg2));
3850: if (flg || flg2) PetscCall(DMSetBasicAdjacency(dm, adjCone, adjClosure));
3852: switch (cell) {
3853: case DM_POLYTOPE_POINT:
3854: case DM_POLYTOPE_SEGMENT:
3855: case DM_POLYTOPE_POINT_PRISM_TENSOR:
3856: case DM_POLYTOPE_TRIANGLE:
3857: case DM_POLYTOPE_QUADRILATERAL:
3858: case DM_POLYTOPE_TETRAHEDRON:
3859: case DM_POLYTOPE_HEXAHEDRON:
3860: *useCoordSpace = PETSC_TRUE;
3861: break;
3862: default:
3863: *useCoordSpace = PETSC_FALSE;
3864: break;
3865: }
3867: if (fflg) {
3868: DM dmnew;
3870: PetscCall(DMPlexCreateFromFile(PetscObjectComm((PetscObject)dm), filename, plexname, interpolate, &dmnew));
3871: PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_FALSE, dmnew));
3872: PetscCall(DMPlexReplace_Internal(dm, &dmnew));
3873: } else if (refDomain) {
3874: PetscCall(DMPlexCreateReferenceCell_Internal(dm, cell));
3875: } else if (bdfflg) {
3876: DM bdm, dmnew;
3878: PetscCall(DMPlexCreateFromFile(PetscObjectComm((PetscObject)dm), bdFilename, plexname, interpolate, &bdm));
3879: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)bdm, "bd_"));
3880: PetscCall(DMSetFromOptions(bdm));
3881: PetscCall(DMPlexGenerate(bdm, NULL, interpolate, &dmnew));
3882: PetscCall(DMDestroy(&bdm));
3883: PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_FALSE, dmnew));
3884: PetscCall(DMPlexReplace_Internal(dm, &dmnew));
3885: } else {
3886: PetscCall(PetscObjectSetName((PetscObject)dm, DMPlexShapes[shape]));
3887: switch (shape) {
3888: case DM_SHAPE_BOX:
3889: case DM_SHAPE_ZBOX:
3890: case DM_SHAPE_ANNULUS: {
3891: PetscInt faces[3] = {0, 0, 0};
3892: PetscReal lower[3] = {0, 0, 0};
3893: PetscReal upper[3] = {1, 1, 1};
3894: DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
3895: PetscBool isAnnular = shape == DM_SHAPE_ANNULUS ? PETSC_TRUE : PETSC_FALSE;
3896: PetscInt i, n;
3898: n = dim;
3899: for (i = 0; i < dim; ++i) faces[i] = (dim == 1 ? 1 : 4 - dim);
3900: PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", faces, &n, &flg));
3901: n = 3;
3902: PetscCall(PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg));
3903: PetscCheck(!flg || !(n != dim), comm, PETSC_ERR_ARG_SIZ, "Lower box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim);
3904: n = 3;
3905: PetscCall(PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg));
3906: PetscCheck(!flg || !(n != dim), comm, PETSC_ERR_ARG_SIZ, "Upper box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim);
3907: n = 3;
3908: PetscCall(PetscOptionsEnumArray("-dm_plex_box_bd", "Boundary type for each dimension", "", DMBoundaryTypes, (PetscEnum *)bdt, &n, &flg));
3909: PetscCheck(!flg || !(n != dim), comm, PETSC_ERR_ARG_SIZ, "Box boundary types had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim);
3911: PetscCheck(!isAnnular || dim == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Only two dimensional annuli have been implemented");
3912: if (isAnnular)
3913: for (i = 0; i < dim - 1; ++i) bdt[i] = DM_BOUNDARY_PERIODIC;
3915: switch (cell) {
3916: case DM_POLYTOPE_TRI_PRISM_TENSOR:
3917: PetscCall(DMPlexCreateWedgeBoxMesh_Internal(dm, faces, lower, upper, bdt));
3918: if (!interpolate) {
3919: DM udm;
3921: PetscCall(DMPlexUninterpolate(dm, &udm));
3922: PetscCall(DMPlexReplace_Internal(dm, &udm));
3923: }
3924: break;
3925: default:
3926: PetscCall(DMPlexCreateBoxMesh_Internal(dm, shape, dim, simplex, faces, lower, upper, bdt, interpolate));
3927: break;
3928: }
3929: if (isAnnular) {
3930: DM cdm;
3931: PetscDS cds;
3932: PetscScalar bounds[2] = {lower[0], upper[0]};
3934: // Fix coordinates for annular region
3935: PetscCall(DMSetPeriodicity(dm, NULL, NULL, NULL));
3936: PetscCall(DMSetCellCoordinatesLocal(dm, NULL));
3937: PetscCall(DMSetCellCoordinates(dm, NULL));
3938: PetscCall(DMPlexCreateCoordinateSpace(dm, 1, NULL));
3939: PetscCall(DMGetCoordinateDM(dm, &cdm));
3940: PetscCall(DMGetDS(cdm, &cds));
3941: PetscCall(PetscDSSetConstants(cds, 2, bounds));
3942: PetscCall(DMPlexRemapGeometry(dm, 0.0, boxToAnnulus));
3943: }
3944: } break;
3945: case DM_SHAPE_BOX_SURFACE: {
3946: PetscInt faces[3] = {0, 0, 0};
3947: PetscReal lower[3] = {0, 0, 0};
3948: PetscReal upper[3] = {1, 1, 1};
3949: PetscInt i, n;
3951: n = dim + 1;
3952: for (i = 0; i < dim + 1; ++i) faces[i] = (dim + 1 == 1 ? 1 : 4 - (dim + 1));
3953: PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", faces, &n, &flg));
3954: n = 3;
3955: PetscCall(PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg));
3956: PetscCheck(!flg || !(n != dim + 1), comm, PETSC_ERR_ARG_SIZ, "Lower box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim + 1);
3957: n = 3;
3958: PetscCall(PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg));
3959: PetscCheck(!flg || !(n != dim + 1), comm, PETSC_ERR_ARG_SIZ, "Upper box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim + 1);
3960: PetscCall(DMPlexCreateBoxSurfaceMesh_Internal(dm, dim + 1, faces, lower, upper, interpolate));
3961: } break;
3962: case DM_SHAPE_SPHERE: {
3963: PetscReal R = 1.0;
3965: PetscCall(PetscOptionsReal("-dm_plex_sphere_radius", "Radius of the sphere", "", R, &R, &flg));
3966: PetscCall(DMPlexCreateSphereMesh_Internal(dm, dim, simplex, R));
3967: } break;
3968: case DM_SHAPE_BALL: {
3969: PetscReal R = 1.0;
3971: PetscCall(PetscOptionsReal("-dm_plex_ball_radius", "Radius of the ball", "", R, &R, &flg));
3972: PetscCall(DMPlexCreateBallMesh_Internal(dm, dim, R));
3973: } break;
3974: case DM_SHAPE_CYLINDER: {
3975: DMBoundaryType bdt = DM_BOUNDARY_NONE;
3976: PetscInt Nw = 6;
3978: PetscCall(PetscOptionsEnum("-dm_plex_cylinder_bd", "Boundary type in the z direction", "", DMBoundaryTypes, (PetscEnum)bdt, (PetscEnum *)&bdt, NULL));
3979: PetscCall(PetscOptionsInt("-dm_plex_cylinder_num_wedges", "Number of wedges around the cylinder", "", Nw, &Nw, NULL));
3980: switch (cell) {
3981: case DM_POLYTOPE_TRI_PRISM_TENSOR:
3982: PetscCall(DMPlexCreateWedgeCylinderMesh_Internal(dm, Nw, interpolate));
3983: break;
3984: default:
3985: PetscCall(DMPlexCreateHexCylinderMesh_Internal(dm, bdt));
3986: break;
3987: }
3988: } break;
3989: case DM_SHAPE_SCHWARZ_P: // fallthrough
3990: case DM_SHAPE_GYROID: {
3991: PetscInt extent[3] = {1, 1, 1}, refine = 0, layers = 0, three;
3992: PetscReal thickness = 0.;
3993: DMBoundaryType periodic[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
3994: DMPlexTPSType tps_type = shape == DM_SHAPE_SCHWARZ_P ? DMPLEX_TPS_SCHWARZ_P : DMPLEX_TPS_GYROID;
3995: PetscBool tps_distribute;
3996: PetscCall(PetscOptionsIntArray("-dm_plex_tps_extent", "Number of replicas for each of three dimensions", NULL, extent, (three = 3, &three), NULL));
3997: PetscCall(PetscOptionsInt("-dm_plex_tps_refine", "Number of refinements", NULL, refine, &refine, NULL));
3998: PetscCall(PetscOptionsEnumArray("-dm_plex_tps_periodic", "Periodicity in each of three dimensions", NULL, DMBoundaryTypes, (PetscEnum *)periodic, (three = 3, &three), NULL));
3999: PetscCall(PetscOptionsInt("-dm_plex_tps_layers", "Number of layers in volumetric extrusion (or zero to not extrude)", NULL, layers, &layers, NULL));
4000: PetscCall(PetscOptionsReal("-dm_plex_tps_thickness", "Thickness of volumetric extrusion", NULL, thickness, &thickness, NULL));
4001: PetscCall(DMPlexDistributeGetDefault(dm, &tps_distribute));
4002: PetscCall(PetscOptionsBool("-dm_plex_tps_distribute", "Distribute the 2D mesh prior to refinement and extrusion", NULL, tps_distribute, &tps_distribute, NULL));
4003: PetscCall(DMPlexCreateTPSMesh_Internal(dm, tps_type, extent, periodic, tps_distribute, refine, layers, thickness));
4004: } break;
4005: case DM_SHAPE_DOUBLET: {
4006: DM dmnew;
4007: PetscReal rl = 0.0;
4009: PetscCall(PetscOptionsReal("-dm_plex_doublet_refinementlimit", "Refinement limit", NULL, rl, &rl, NULL));
4010: PetscCall(DMPlexCreateDoublet(PetscObjectComm((PetscObject)dm), dim, simplex, interpolate, rl, &dmnew));
4011: PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_FALSE, dmnew));
4012: PetscCall(DMPlexReplace_Internal(dm, &dmnew));
4013: } break;
4014: case DM_SHAPE_HYPERCUBIC: {
4015: PetscInt *edges;
4016: PetscReal *lower, *upper;
4017: DMBoundaryType *bdt;
4018: PetscInt n, d;
4020: *useCoordSpace = PETSC_FALSE;
4021: PetscCall(PetscMalloc4(dim, &edges, dim, &lower, dim, &upper, dim, &bdt));
4022: for (d = 0; d < dim; ++d) {
4023: edges[d] = 1;
4024: lower[d] = 0.;
4025: upper[d] = 1.;
4026: bdt[d] = DM_BOUNDARY_PERIODIC;
4027: }
4028: n = dim;
4029: PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", edges, &n, &flg));
4030: n = dim;
4031: PetscCall(PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg));
4032: PetscCheck(!flg || n == dim, comm, PETSC_ERR_ARG_SIZ, "Lower box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim);
4033: n = dim;
4034: PetscCall(PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg));
4035: PetscCheck(!flg || n == dim, comm, PETSC_ERR_ARG_SIZ, "Upper box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim);
4036: n = dim;
4037: PetscCall(PetscOptionsEnumArray("-dm_plex_box_bd", "Boundary type for each dimension", "", DMBoundaryTypes, (PetscEnum *)bdt, &n, &flg));
4038: PetscCheck(!flg || n == dim, comm, PETSC_ERR_ARG_SIZ, "Box boundary types had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim);
4039: PetscCall(DMPlexCreateHypercubicMesh_Internal(dm, dim, lower, upper, edges, bdt));
4040: PetscCall(PetscFree4(edges, lower, upper, bdt));
4041: } break;
4042: default:
4043: SETERRQ(comm, PETSC_ERR_SUP, "Domain shape %s is unsupported", DMPlexShapes[shape]);
4044: }
4045: }
4046: PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
4047: if (!((PetscObject)dm)->name && nameflg) PetscCall(PetscObjectSetName((PetscObject)dm, plexname));
4048: PetscCall(PetscLogEventEnd(DMPLEX_CreateFromOptions, dm, 0, 0, 0));
4049: PetscFunctionReturn(PETSC_SUCCESS);
4050: }
4052: PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM dm, PetscOptionItems *PetscOptionsObject)
4053: {
4054: DM_Plex *mesh = (DM_Plex *)dm->data;
4055: PetscBool flg, flg2;
4056: char bdLabel[PETSC_MAX_PATH_LEN];
4058: PetscFunctionBegin;
4059: /* Handle viewing */
4060: PetscCall(PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL));
4061: PetscCall(PetscOptionsBoundedInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL, 0));
4062: PetscCall(PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL));
4063: PetscCall(PetscOptionsBoundedInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL, 0));
4064: PetscCall(PetscOptionsBoundedInt("-dm_plex_print_locate", "Debug output level all point location computations", "DMLocatePoints", 0, &mesh->printLocate, NULL, 0));
4065: PetscCall(DMMonitorSetFromOptions(dm, "-dm_plex_monitor_throughput", "Monitor the simulation throughput", "DMPlexMonitorThroughput", DMPlexMonitorThroughput, NULL, &flg));
4066: if (flg) PetscCall(PetscLogDefaultBegin());
4067: /* Labeling */
4068: PetscCall(PetscOptionsString("-dm_plex_boundary_label", "Label to mark the mesh boundary", "", bdLabel, bdLabel, sizeof(bdLabel), &flg));
4069: if (flg) PetscCall(DMPlexCreateBoundaryLabel_Private(dm, bdLabel));
4070: /* Point Location */
4071: PetscCall(PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL));
4072: /* Partitioning and distribution */
4073: PetscCall(PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL));
4074: /* Generation and remeshing */
4075: PetscCall(PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL));
4076: /* Projection behavior */
4077: PetscCall(PetscOptionsBoundedInt("-dm_plex_max_projection_height", "Maximum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL, 0));
4078: PetscCall(PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL));
4079: /* Checking structure */
4080: {
4081: PetscBool all = PETSC_FALSE;
4083: PetscCall(PetscOptionsBool("-dm_plex_check_all", "Perform all basic checks", "DMPlexCheck", PETSC_FALSE, &all, NULL));
4084: if (all) {
4085: PetscCall(DMPlexCheck(dm));
4086: } else {
4087: PetscCall(PetscOptionsBool("-dm_plex_check_symmetry", "Check that the adjacency information in the mesh is symmetric", "DMPlexCheckSymmetry", PETSC_FALSE, &flg, &flg2));
4088: if (flg && flg2) PetscCall(DMPlexCheckSymmetry(dm));
4089: PetscCall(PetscOptionsBool("-dm_plex_check_skeleton", "Check that each cell has the correct number of vertices (only for homogeneous simplex or tensor meshes)", "DMPlexCheckSkeleton", PETSC_FALSE, &flg, &flg2));
4090: if (flg && flg2) PetscCall(DMPlexCheckSkeleton(dm, 0));
4091: PetscCall(PetscOptionsBool("-dm_plex_check_faces", "Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type", "DMPlexCheckFaces", PETSC_FALSE, &flg, &flg2));
4092: if (flg && flg2) PetscCall(DMPlexCheckFaces(dm, 0));
4093: PetscCall(PetscOptionsBool("-dm_plex_check_geometry", "Check that cells have positive volume", "DMPlexCheckGeometry", PETSC_FALSE, &flg, &flg2));
4094: if (flg && flg2) PetscCall(DMPlexCheckGeometry(dm));
4095: PetscCall(PetscOptionsBool("-dm_plex_check_pointsf", "Check some necessary conditions for PointSF", "DMPlexCheckPointSF", PETSC_FALSE, &flg, &flg2));
4096: if (flg && flg2) PetscCall(DMPlexCheckPointSF(dm, NULL, PETSC_FALSE));
4097: PetscCall(PetscOptionsBool("-dm_plex_check_interface_cones", "Check points on inter-partition interfaces have conforming order of cone points", "DMPlexCheckInterfaceCones", PETSC_FALSE, &flg, &flg2));
4098: if (flg && flg2) PetscCall(DMPlexCheckInterfaceCones(dm));
4099: }
4100: PetscCall(PetscOptionsBool("-dm_plex_check_cell_shape", "Check cell shape", "DMPlexCheckCellShape", PETSC_FALSE, &flg, &flg2));
4101: if (flg && flg2) PetscCall(DMPlexCheckCellShape(dm, PETSC_TRUE, PETSC_DETERMINE));
4102: }
4103: {
4104: PetscReal scale = 1.0;
4106: PetscCall(PetscOptionsReal("-dm_plex_scale", "Scale factor for mesh coordinates", "DMPlexScale", scale, &scale, &flg));
4107: if (flg) {
4108: Vec coordinates, coordinatesLocal;
4110: PetscCall(DMGetCoordinates(dm, &coordinates));
4111: PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal));
4112: PetscCall(VecScale(coordinates, scale));
4113: PetscCall(VecScale(coordinatesLocal, scale));
4114: }
4115: }
4116: PetscCall(PetscPartitionerSetFromOptions(mesh->partitioner));
4117: PetscFunctionReturn(PETSC_SUCCESS);
4118: }
4120: PetscErrorCode DMSetFromOptions_Overlap_Plex(DM dm, PetscOptionItems *PetscOptionsObject, PetscInt *overlap)
4121: {
4122: PetscInt numOvLabels = 16, numOvExLabels = 16;
4123: char *ovLabelNames[16], *ovExLabelNames[16];
4124: PetscInt numOvValues = 16, numOvExValues = 16, l;
4125: PetscBool flg;
4127: PetscFunctionBegin;
4128: PetscCall(PetscOptionsBoundedInt("-dm_distribute_overlap", "The size of the overlap halo", "DMPlexDistribute", *overlap, overlap, NULL, 0));
4129: PetscCall(PetscOptionsStringArray("-dm_distribute_overlap_labels", "List of overlap label names", "DMPlexDistribute", ovLabelNames, &numOvLabels, &flg));
4130: if (!flg) numOvLabels = 0;
4131: if (numOvLabels) {
4132: ((DM_Plex *)dm->data)->numOvLabels = numOvLabels;
4133: for (l = 0; l < numOvLabels; ++l) {
4134: PetscCall(DMGetLabel(dm, ovLabelNames[l], &((DM_Plex *)dm->data)->ovLabels[l]));
4135: PetscCheck(((DM_Plex *)dm->data)->ovLabels[l], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid label name %s", ovLabelNames[l]);
4136: PetscCall(PetscFree(ovLabelNames[l]));
4137: }
4138: PetscCall(PetscOptionsIntArray("-dm_distribute_overlap_values", "List of overlap label values", "DMPlexDistribute", ((DM_Plex *)dm->data)->ovValues, &numOvValues, &flg));
4139: if (!flg) numOvValues = 0;
4140: PetscCheck(numOvLabels == numOvValues, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The number of labels %" PetscInt_FMT " must match the number of values %" PetscInt_FMT, numOvLabels, numOvValues);
4142: PetscCall(PetscOptionsStringArray("-dm_distribute_overlap_exclude_labels", "List of overlap exclude label names", "DMPlexDistribute", ovExLabelNames, &numOvExLabels, &flg));
4143: if (!flg) numOvExLabels = 0;
4144: ((DM_Plex *)dm->data)->numOvExLabels = numOvExLabels;
4145: for (l = 0; l < numOvExLabels; ++l) {
4146: PetscCall(DMGetLabel(dm, ovExLabelNames[l], &((DM_Plex *)dm->data)->ovExLabels[l]));
4147: PetscCheck(((DM_Plex *)dm->data)->ovExLabels[l], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid label name %s", ovExLabelNames[l]);
4148: PetscCall(PetscFree(ovExLabelNames[l]));
4149: }
4150: PetscCall(PetscOptionsIntArray("-dm_distribute_overlap_exclude_values", "List of overlap exclude label values", "DMPlexDistribute", ((DM_Plex *)dm->data)->ovExValues, &numOvExValues, &flg));
4151: if (!flg) numOvExValues = 0;
4152: PetscCheck(numOvExLabels == numOvExValues, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The number of exclude labels %" PetscInt_FMT " must match the number of values %" PetscInt_FMT, numOvExLabels, numOvExValues);
4153: }
4154: PetscFunctionReturn(PETSC_SUCCESS);
4155: }
4157: static PetscErrorCode DMSetFromOptions_Plex(DM dm, PetscOptionItems *PetscOptionsObject)
4158: {
4159: PetscFunctionList ordlist;
4160: char oname[256];
4161: PetscReal volume = -1.0;
4162: PetscInt prerefine = 0, refine = 0, r, coarsen = 0, overlap = 0, extLayers = 0, dim;
4163: PetscBool uniformOrig, created = PETSC_FALSE, uniform = PETSC_TRUE, distribute, interpolate = PETSC_TRUE, coordSpace = PETSC_TRUE, remap = PETSC_TRUE, ghostCells = PETSC_FALSE, isHierarchy, ignoreModel = PETSC_FALSE, flg;
4164: DMPlexReorderDefaultFlag reorder;
4166: PetscFunctionBegin;
4167: PetscOptionsHeadBegin(PetscOptionsObject, "DMPlex Options");
4168: if (dm->cloneOpts) goto non_refine;
4169: /* Handle automatic creation */
4170: PetscCall(DMGetDimension(dm, &dim));
4171: if (dim < 0) {
4172: PetscCall(DMPlexCreateFromOptions_Internal(PetscOptionsObject, &coordSpace, dm));
4173: created = PETSC_TRUE;
4174: }
4175: PetscCall(DMGetDimension(dm, &dim));
4176: /* Handle interpolation before distribution */
4177: PetscCall(PetscOptionsBool("-dm_plex_interpolate_pre", "Flag to interpolate mesh before distribution", "", interpolate, &interpolate, &flg));
4178: if (flg) {
4179: DMPlexInterpolatedFlag interpolated;
4181: PetscCall(DMPlexIsInterpolated(dm, &interpolated));
4182: if (interpolated == DMPLEX_INTERPOLATED_FULL && !interpolate) {
4183: DM udm;
4185: PetscCall(DMPlexUninterpolate(dm, &udm));
4186: PetscCall(DMPlexReplace_Internal(dm, &udm));
4187: } else if (interpolated != DMPLEX_INTERPOLATED_FULL && interpolate) {
4188: DM idm;
4190: PetscCall(DMPlexInterpolate(dm, &idm));
4191: PetscCall(DMPlexReplace_Internal(dm, &idm));
4192: }
4193: }
4194: /* Handle DMPlex refinement before distribution */
4195: PetscCall(PetscOptionsBool("-dm_refine_ignore_model", "Flag to ignore the geometry model when refining", "DMCreate", ignoreModel, &ignoreModel, &flg));
4196: if (flg) ((DM_Plex *)dm->data)->ignoreModel = ignoreModel;
4197: PetscCall(DMPlexGetRefinementUniform(dm, &uniformOrig));
4198: PetscCall(PetscOptionsBoundedInt("-dm_refine_pre", "The number of refinements before distribution", "DMCreate", prerefine, &prerefine, NULL, 0));
4199: PetscCall(PetscOptionsBool("-dm_refine_remap_pre", "Flag to control coordinate remapping", "DMCreate", remap, &remap, NULL));
4200: PetscCall(PetscOptionsBool("-dm_refine_uniform_pre", "Flag for uniform refinement before distribution", "DMCreate", uniform, &uniform, &flg));
4201: if (flg) PetscCall(DMPlexSetRefinementUniform(dm, uniform));
4202: PetscCall(PetscOptionsReal("-dm_refine_volume_limit_pre", "The maximum cell volume after refinement before distribution", "DMCreate", volume, &volume, &flg));
4203: if (flg) {
4204: PetscCall(DMPlexSetRefinementUniform(dm, PETSC_FALSE));
4205: PetscCall(DMPlexSetRefinementLimit(dm, volume));
4206: prerefine = PetscMax(prerefine, 1);
4207: }
4208: if (prerefine) PetscCall(DMLocalizeCoordinates(dm));
4209: for (r = 0; r < prerefine; ++r) {
4210: DM rdm;
4211: PetscPointFunc coordFunc = ((DM_Plex *)dm->data)->coordFunc;
4213: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4214: PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &rdm));
4215: PetscCall(DMPlexReplace_Internal(dm, &rdm));
4216: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4217: if (coordFunc && remap) {
4218: PetscCall(DMPlexRemapGeometry(dm, 0.0, coordFunc));
4219: ((DM_Plex *)dm->data)->coordFunc = coordFunc;
4220: }
4221: }
4222: PetscCall(DMPlexSetRefinementUniform(dm, uniformOrig));
4223: /* Handle DMPlex extrusion before distribution */
4224: PetscCall(PetscOptionsBoundedInt("-dm_extrude", "The number of layers to extrude", "", extLayers, &extLayers, NULL, 0));
4225: if (extLayers) {
4226: DM edm;
4228: PetscCall(DMExtrude(dm, extLayers, &edm));
4229: PetscCall(DMPlexReplace_Internal(dm, &edm));
4230: ((DM_Plex *)dm->data)->coordFunc = NULL;
4231: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4232: extLayers = 0;
4233: PetscCall(DMGetDimension(dm, &dim));
4234: }
4235: /* Handle DMPlex reordering before distribution */
4236: PetscCall(DMPlexReorderGetDefault(dm, &reorder));
4237: PetscCall(MatGetOrderingList(&ordlist));
4238: PetscCall(PetscStrncpy(oname, MATORDERINGNATURAL, sizeof(oname)));
4239: PetscCall(PetscOptionsFList("-dm_plex_reorder", "Set mesh reordering type", "DMPlexGetOrdering", ordlist, MATORDERINGNATURAL, oname, sizeof(oname), &flg));
4240: if (reorder == DMPLEX_REORDER_DEFAULT_TRUE || flg) {
4241: DM pdm;
4242: IS perm;
4244: PetscCall(DMPlexGetOrdering(dm, oname, NULL, &perm));
4245: PetscCall(DMPlexPermute(dm, perm, &pdm));
4246: PetscCall(ISDestroy(&perm));
4247: PetscCall(DMPlexReplace_Internal(dm, &pdm));
4248: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4249: }
4250: /* Handle DMPlex distribution */
4251: PetscCall(DMPlexDistributeGetDefault(dm, &distribute));
4252: PetscCall(PetscOptionsBool("-dm_distribute", "Flag to redistribute a mesh among processes", "DMPlexDistribute", distribute, &distribute, NULL));
4253: PetscCall(DMSetFromOptions_Overlap_Plex(dm, PetscOptionsObject, &overlap));
4254: if (distribute) {
4255: DM pdm = NULL;
4256: PetscPartitioner part;
4258: PetscCall(DMPlexGetPartitioner(dm, &part));
4259: PetscCall(PetscPartitionerSetFromOptions(part));
4260: PetscCall(DMPlexDistribute(dm, overlap, NULL, &pdm));
4261: if (pdm) PetscCall(DMPlexReplace_Internal(dm, &pdm));
4262: }
4263: /* Must check CEED options before creating function space for coordinates */
4264: {
4265: PetscBool useCeed = PETSC_FALSE, flg;
4267: PetscCall(PetscOptionsBool("-dm_plex_use_ceed", "Use LibCEED as the FEM backend", "DMPlexSetUseCeed", useCeed, &useCeed, &flg));
4268: if (flg) PetscCall(DMPlexSetUseCeed(dm, useCeed));
4269: }
4270: /* Create coordinate space */
4271: if (created) {
4272: DM_Plex *mesh = (DM_Plex *)dm->data;
4273: PetscInt degree = 1;
4274: PetscInt height = 0;
4275: DM cdm;
4276: PetscBool flg;
4278: PetscCall(PetscOptionsBool("-dm_coord_space", "Use an FEM space for coordinates", "", coordSpace, &coordSpace, &flg));
4279: PetscCall(PetscOptionsInt("-dm_coord_petscspace_degree", "FEM degree for coordinate space", "", degree, °ree, NULL));
4280: if (coordSpace) PetscCall(DMPlexCreateCoordinateSpace(dm, degree, mesh->coordFunc));
4281: PetscCall(DMGetCoordinateDM(dm, &cdm));
4282: if (flg && !coordSpace) {
4283: PetscDS cds;
4284: PetscObject obj;
4285: PetscClassId id;
4287: PetscCall(DMGetDS(cdm, &cds));
4288: PetscCall(PetscDSGetDiscretization(cds, 0, &obj));
4289: PetscCall(PetscObjectGetClassId(obj, &id));
4290: if (id == PETSCFE_CLASSID) {
4291: PetscContainer dummy;
4293: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &dummy));
4294: PetscCall(PetscObjectSetName((PetscObject)dummy, "coordinates"));
4295: PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)dummy));
4296: PetscCall(PetscContainerDestroy(&dummy));
4297: PetscCall(DMClearDS(cdm));
4298: }
4299: mesh->coordFunc = NULL;
4300: }
4301: PetscCall(PetscOptionsBool("-dm_sparse_localize", "Localize only necessary cells", "", dm->sparseLocalize, &dm->sparseLocalize, &flg));
4302: PetscCall(PetscOptionsInt("-dm_localize_height", "Localize edges and faces in addition to cells", "", height, &height, &flg));
4303: if (flg) PetscCall(DMPlexSetMaxProjectionHeight(cdm, height));
4304: PetscCall(DMLocalizeCoordinates(dm));
4305: }
4306: /* Handle DMPlex refinement */
4307: remap = PETSC_TRUE;
4308: PetscCall(PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL, 0));
4309: PetscCall(PetscOptionsBool("-dm_refine_remap", "Flag to control coordinate remapping", "DMCreate", remap, &remap, NULL));
4310: PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy, 0));
4311: if (refine) PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
4312: if (refine && isHierarchy) {
4313: DM *dms, coarseDM;
4315: PetscCall(DMGetCoarseDM(dm, &coarseDM));
4316: PetscCall(PetscObjectReference((PetscObject)coarseDM));
4317: PetscCall(PetscMalloc1(refine, &dms));
4318: PetscCall(DMRefineHierarchy(dm, refine, dms));
4319: /* Total hack since we do not pass in a pointer */
4320: PetscCall(DMPlexSwap_Static(dm, dms[refine - 1]));
4321: if (refine == 1) {
4322: PetscCall(DMSetCoarseDM(dm, dms[0]));
4323: PetscCall(DMPlexSetRegularRefinement(dm, PETSC_TRUE));
4324: } else {
4325: PetscCall(DMSetCoarseDM(dm, dms[refine - 2]));
4326: PetscCall(DMPlexSetRegularRefinement(dm, PETSC_TRUE));
4327: PetscCall(DMSetCoarseDM(dms[0], dms[refine - 1]));
4328: PetscCall(DMPlexSetRegularRefinement(dms[0], PETSC_TRUE));
4329: }
4330: PetscCall(DMSetCoarseDM(dms[refine - 1], coarseDM));
4331: PetscCall(PetscObjectDereference((PetscObject)coarseDM));
4332: /* Free DMs */
4333: for (r = 0; r < refine; ++r) {
4334: PetscCall(DMSetFromOptions_NonRefinement_Plex(dms[r], PetscOptionsObject));
4335: PetscCall(DMDestroy(&dms[r]));
4336: }
4337: PetscCall(PetscFree(dms));
4338: } else {
4339: for (r = 0; r < refine; ++r) {
4340: DM rdm;
4341: PetscPointFunc coordFunc = ((DM_Plex *)dm->data)->coordFunc;
4343: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4344: PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &rdm));
4345: /* Total hack since we do not pass in a pointer */
4346: PetscCall(DMPlexReplace_Internal(dm, &rdm));
4347: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4348: if (coordFunc && remap) {
4349: PetscCall(DMPlexRemapGeometry(dm, 0.0, coordFunc));
4350: ((DM_Plex *)dm->data)->coordFunc = coordFunc;
4351: }
4352: }
4353: }
4354: /* Handle DMPlex coarsening */
4355: PetscCall(PetscOptionsBoundedInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL, 0));
4356: PetscCall(PetscOptionsBoundedInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy, 0));
4357: if (coarsen && isHierarchy) {
4358: DM *dms;
4360: PetscCall(PetscMalloc1(coarsen, &dms));
4361: PetscCall(DMCoarsenHierarchy(dm, coarsen, dms));
4362: /* Free DMs */
4363: for (r = 0; r < coarsen; ++r) {
4364: PetscCall(DMSetFromOptions_NonRefinement_Plex(dms[r], PetscOptionsObject));
4365: PetscCall(DMDestroy(&dms[r]));
4366: }
4367: PetscCall(PetscFree(dms));
4368: } else {
4369: for (r = 0; r < coarsen; ++r) {
4370: DM cdm;
4371: PetscPointFunc coordFunc = ((DM_Plex *)dm->data)->coordFunc;
4373: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4374: PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &cdm));
4375: /* Total hack since we do not pass in a pointer */
4376: PetscCall(DMPlexReplace_Internal(dm, &cdm));
4377: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4378: if (coordFunc) {
4379: PetscCall(DMPlexRemapGeometry(dm, 0.0, coordFunc));
4380: ((DM_Plex *)dm->data)->coordFunc = coordFunc;
4381: }
4382: }
4383: }
4384: /* Handle ghost cells */
4385: PetscCall(PetscOptionsBool("-dm_plex_create_fv_ghost_cells", "Flag to create finite volume ghost cells on the boundary", "DMCreate", ghostCells, &ghostCells, NULL));
4386: if (ghostCells) {
4387: DM gdm;
4388: char lname[PETSC_MAX_PATH_LEN];
4390: lname[0] = '\0';
4391: PetscCall(PetscOptionsString("-dm_plex_fv_ghost_cells_label", "Label name for ghost cells boundary", "DMCreate", lname, lname, sizeof(lname), &flg));
4392: PetscCall(DMPlexConstructGhostCells(dm, flg ? lname : NULL, NULL, &gdm));
4393: PetscCall(DMPlexReplace_Internal(dm, &gdm));
4394: }
4395: /* Handle 1D order */
4396: if (reorder != DMPLEX_REORDER_DEFAULT_FALSE && dim == 1) {
4397: DM cdm, rdm;
4398: PetscDS cds;
4399: PetscObject obj;
4400: PetscClassId id = PETSC_OBJECT_CLASSID;
4401: IS perm;
4402: PetscInt Nf;
4403: PetscBool distributed;
4405: PetscCall(DMPlexIsDistributed(dm, &distributed));
4406: PetscCall(DMGetCoordinateDM(dm, &cdm));
4407: PetscCall(DMGetDS(cdm, &cds));
4408: PetscCall(PetscDSGetNumFields(cds, &Nf));
4409: if (Nf) {
4410: PetscCall(PetscDSGetDiscretization(cds, 0, &obj));
4411: PetscCall(PetscObjectGetClassId(obj, &id));
4412: }
4413: if (!distributed && id != PETSCFE_CLASSID) {
4414: PetscCall(DMPlexGetOrdering1D(dm, &perm));
4415: PetscCall(DMPlexPermute(dm, perm, &rdm));
4416: PetscCall(DMPlexReplace_Internal(dm, &rdm));
4417: PetscCall(ISDestroy(&perm));
4418: }
4419: }
4420: /* Handle */
4421: non_refine:
4422: PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject));
4423: PetscOptionsHeadEnd();
4424: PetscFunctionReturn(PETSC_SUCCESS);
4425: }
4427: static PetscErrorCode DMCreateGlobalVector_Plex(DM dm, Vec *vec)
4428: {
4429: PetscFunctionBegin;
4430: PetscCall(DMCreateGlobalVector_Section_Private(dm, vec));
4431: /* PetscCall(VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM)); */
4432: PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex));
4433: PetscCall(VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void))VecView_Plex_Native));
4434: PetscCall(VecSetOperation(*vec, VECOP_LOAD, (void (*)(void))VecLoad_Plex));
4435: PetscCall(VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void))VecLoad_Plex_Native));
4436: PetscFunctionReturn(PETSC_SUCCESS);
4437: }
4439: static PetscErrorCode DMCreateLocalVector_Plex(DM dm, Vec *vec)
4440: {
4441: PetscFunctionBegin;
4442: PetscCall(DMCreateLocalVector_Section_Private(dm, vec));
4443: PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex_Local));
4444: PetscCall(VecSetOperation(*vec, VECOP_LOAD, (void (*)(void))VecLoad_Plex_Local));
4445: PetscFunctionReturn(PETSC_SUCCESS);
4446: }
4448: static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
4449: {
4450: PetscInt depth, d;
4452: PetscFunctionBegin;
4453: PetscCall(DMPlexGetDepth(dm, &depth));
4454: if (depth == 1) {
4455: PetscCall(DMGetDimension(dm, &d));
4456: if (dim == 0) PetscCall(DMPlexGetDepthStratum(dm, dim, pStart, pEnd));
4457: else if (dim == d) PetscCall(DMPlexGetDepthStratum(dm, 1, pStart, pEnd));
4458: else {
4459: *pStart = 0;
4460: *pEnd = 0;
4461: }
4462: } else {
4463: PetscCall(DMPlexGetDepthStratum(dm, dim, pStart, pEnd));
4464: }
4465: PetscFunctionReturn(PETSC_SUCCESS);
4466: }
4468: static PetscErrorCode DMGetNeighbors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
4469: {
4470: PetscSF sf;
4471: PetscInt niranks, njranks, n;
4472: const PetscMPIInt *iranks, *jranks;
4473: DM_Plex *data = (DM_Plex *)dm->data;
4475: PetscFunctionBegin;
4476: PetscCall(DMGetPointSF(dm, &sf));
4477: if (!data->neighbors) {
4478: PetscCall(PetscSFSetUp(sf));
4479: PetscCall(PetscSFGetRootRanks(sf, &njranks, &jranks, NULL, NULL, NULL));
4480: PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, NULL, NULL));
4481: PetscCall(PetscMalloc1(njranks + niranks + 1, &data->neighbors));
4482: PetscCall(PetscArraycpy(data->neighbors + 1, jranks, njranks));
4483: PetscCall(PetscArraycpy(data->neighbors + njranks + 1, iranks, niranks));
4484: n = njranks + niranks;
4485: PetscCall(PetscSortRemoveDupsMPIInt(&n, data->neighbors + 1));
4486: /* The following cast should never fail: can't have more neighbors than PETSC_MPI_INT_MAX */
4487: PetscCall(PetscMPIIntCast(n, data->neighbors));
4488: }
4489: if (nranks) *nranks = data->neighbors[0];
4490: if (ranks) {
4491: if (data->neighbors[0]) *ranks = data->neighbors + 1;
4492: else *ranks = NULL;
4493: }
4494: PetscFunctionReturn(PETSC_SUCCESS);
4495: }
4497: PETSC_INTERN PetscErrorCode DMInterpolateSolution_Plex(DM, DM, Mat, Vec, Vec);
4499: static PetscErrorCode DMInitialize_Plex(DM dm)
4500: {
4501: PetscFunctionBegin;
4502: dm->ops->view = DMView_Plex;
4503: dm->ops->load = DMLoad_Plex;
4504: dm->ops->setfromoptions = DMSetFromOptions_Plex;
4505: dm->ops->clone = DMClone_Plex;
4506: dm->ops->setup = DMSetUp_Plex;
4507: dm->ops->createlocalsection = DMCreateLocalSection_Plex;
4508: dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex;
4509: dm->ops->createglobalvector = DMCreateGlobalVector_Plex;
4510: dm->ops->createlocalvector = DMCreateLocalVector_Plex;
4511: dm->ops->getlocaltoglobalmapping = NULL;
4512: dm->ops->createfieldis = NULL;
4513: dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex;
4514: dm->ops->createcoordinatefield = DMCreateCoordinateField_Plex;
4515: dm->ops->getcoloring = NULL;
4516: dm->ops->creatematrix = DMCreateMatrix_Plex;
4517: dm->ops->createinterpolation = DMCreateInterpolation_Plex;
4518: dm->ops->createmassmatrix = DMCreateMassMatrix_Plex;
4519: dm->ops->createmassmatrixlumped = DMCreateMassMatrixLumped_Plex;
4520: dm->ops->createinjection = DMCreateInjection_Plex;
4521: dm->ops->refine = DMRefine_Plex;
4522: dm->ops->coarsen = DMCoarsen_Plex;
4523: dm->ops->refinehierarchy = DMRefineHierarchy_Plex;
4524: dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex;
4525: dm->ops->extrude = DMExtrude_Plex;
4526: dm->ops->globaltolocalbegin = NULL;
4527: dm->ops->globaltolocalend = NULL;
4528: dm->ops->localtoglobalbegin = NULL;
4529: dm->ops->localtoglobalend = NULL;
4530: dm->ops->destroy = DMDestroy_Plex;
4531: dm->ops->createsubdm = DMCreateSubDM_Plex;
4532: dm->ops->createsuperdm = DMCreateSuperDM_Plex;
4533: dm->ops->getdimpoints = DMGetDimPoints_Plex;
4534: dm->ops->locatepoints = DMLocatePoints_Plex;
4535: dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex;
4536: dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex;
4537: dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex;
4538: dm->ops->projectfieldlabellocal = DMProjectFieldLabelLocal_Plex;
4539: dm->ops->projectbdfieldlabellocal = DMProjectBdFieldLabelLocal_Plex;
4540: dm->ops->computel2diff = DMComputeL2Diff_Plex;
4541: dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex;
4542: dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex;
4543: dm->ops->getneighbors = DMGetNeighbors_Plex;
4544: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_Plex));
4545: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertTimeDerviativeBoundaryValues_C", DMPlexInsertTimeDerivativeBoundaryValues_Plex));
4546: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Plex));
4547: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", DMCreateNeumannOverlap_Plex));
4548: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", DMPlexGetOverlap_Plex));
4549: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeGetDefault_C", DMPlexDistributeGetDefault_Plex));
4550: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeSetDefault_C", DMPlexDistributeSetDefault_Plex));
4551: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderGetDefault_C", DMPlexReorderGetDefault_Plex));
4552: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderSetDefault_C", DMPlexReorderSetDefault_Plex));
4553: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMInterpolateSolution_C", DMInterpolateSolution_Plex));
4554: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", DMPlexGetOverlap_Plex));
4555: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetOverlap_C", DMPlexSetOverlap_Plex));
4556: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetUseCeed_C", DMPlexGetUseCeed_Plex));
4557: PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetUseCeed_C", DMPlexSetUseCeed_Plex));
4558: PetscFunctionReturn(PETSC_SUCCESS);
4559: }
4561: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm)
4562: {
4563: DM_Plex *mesh = (DM_Plex *)dm->data;
4564: PetscSF face_sf;
4566: PetscFunctionBegin;
4567: mesh->refct++;
4568: (*newdm)->data = mesh;
4569: PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &face_sf));
4570: PetscCall(DMPlexSetIsoperiodicFaceSF(*newdm, face_sf));
4571: PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMPLEX));
4572: PetscCall(DMInitialize_Plex(*newdm));
4573: PetscFunctionReturn(PETSC_SUCCESS);
4574: }
4576: /*MC
4577: DMPLEX = "plex" - A `DM` object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram.
4578: In the local representation, `Vec`s contain all unknowns in the interior and shared boundary. This is
4579: specified by a PetscSection object. Ownership in the global representation is determined by
4580: ownership of the underlying `DMPLEX` points. This is specified by another `PetscSection` object.
4582: Options Database Keys:
4583: + -dm_refine_pre - Refine mesh before distribution
4584: + -dm_refine_uniform_pre - Choose uniform or generator-based refinement
4585: + -dm_refine_volume_limit_pre - Cell volume limit after pre-refinement using generator
4586: . -dm_distribute - Distribute mesh across processes
4587: . -dm_distribute_overlap - Number of cells to overlap for distribution
4588: . -dm_refine - Refine mesh after distribution
4589: . -dm_plex_hash_location - Use grid hashing for point location
4590: . -dm_plex_hash_box_faces <n,m,p> - The number of divisions in each direction of the grid hash
4591: . -dm_plex_partition_balance - Attempt to evenly divide points on partition boundary between processes
4592: . -dm_plex_remesh_bd - Allow changes to the boundary on remeshing
4593: . -dm_plex_max_projection_height - Maximum mesh point height used to project locally
4594: . -dm_plex_regular_refinement - Use special nested projection algorithm for regular refinement
4595: . -dm_plex_check_all - Perform all checks below
4596: . -dm_plex_check_symmetry - Check that the adjacency information in the mesh is symmetric
4597: . -dm_plex_check_skeleton <celltype> - Check that each cell has the correct number of vertices
4598: . -dm_plex_check_faces <celltype> - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type
4599: . -dm_plex_check_geometry - Check that cells have positive volume
4600: . -dm_view :mesh.tex:ascii_latex - View the mesh in LaTeX/TikZ
4601: . -dm_plex_view_scale <num> - Scale the TikZ
4602: - -dm_plex_print_fem <num> - View FEM assembly information, such as element vectors and matrices
4604: Level: intermediate
4606: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMType`, `DMPlexCreate()`, `DMCreate()`, `DMSetType()`, `PetscSection`
4607: M*/
4609: PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm)
4610: {
4611: DM_Plex *mesh;
4612: PetscInt unit;
4614: PetscFunctionBegin;
4615: PetscCall(PetscCitationsRegister(PlexCitation, &Plexcite));
4617: PetscCall(PetscNew(&mesh));
4618: dm->data = mesh;
4620: mesh->refct = 1;
4621: PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection));
4622: PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection));
4623: mesh->refinementUniform = PETSC_TRUE;
4624: mesh->refinementLimit = -1.0;
4625: mesh->distDefault = PETSC_TRUE;
4626: mesh->reorderDefault = DMPLEX_REORDER_DEFAULT_NOTSET;
4627: mesh->distributionName = NULL;
4628: mesh->interpolated = DMPLEX_INTERPOLATED_INVALID;
4629: mesh->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID;
4631: PetscCall(PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner));
4632: mesh->remeshBd = PETSC_FALSE;
4634: for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0;
4636: mesh->depthState = -1;
4637: mesh->celltypeState = -1;
4638: mesh->printTol = 1.0e-10;
4640: PetscCall(DMInitialize_Plex(dm));
4641: PetscFunctionReturn(PETSC_SUCCESS);
4642: }
4644: /*@
4645: DMPlexCreate - Creates a `DMPLEX` object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram.
4647: Collective
4649: Input Parameter:
4650: . comm - The communicator for the `DMPLEX` object
4652: Output Parameter:
4653: . mesh - The `DMPLEX` object
4655: Level: beginner
4657: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMType`, `DMCreate()`, `DMSetType()`
4658: @*/
4659: PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh)
4660: {
4661: PetscFunctionBegin;
4662: PetscAssertPointer(mesh, 2);
4663: PetscCall(DMCreate(comm, mesh));
4664: PetscCall(DMSetType(*mesh, DMPLEX));
4665: PetscFunctionReturn(PETSC_SUCCESS);
4666: }
4668: /*@C
4669: DMPlexBuildFromCellListParallel - Build distributed `DMPLEX` topology from a list of vertices for each cell (common mesh generator output)
4671: Collective; No Fortran Support
4673: Input Parameters:
4674: + dm - The `DM`
4675: . numCells - The number of cells owned by this process
4676: . numVertices - The number of vertices to be owned by this process, or `PETSC_DECIDE`
4677: . NVertices - The global number of vertices, or `PETSC_DETERMINE`
4678: . numCorners - The number of vertices for each cell
4679: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
4681: Output Parameters:
4682: + vertexSF - (Optional) `PetscSF` describing complete vertex ownership
4683: - verticesAdjSaved - (Optional) vertex adjacency array
4685: Level: advanced
4687: Notes:
4688: Two triangles sharing a face
4689: .vb
4691: 2
4692: / | \
4693: / | \
4694: / | \
4695: 0 0 | 1 3
4696: \ | /
4697: \ | /
4698: \ | /
4699: 1
4700: .ve
4701: would have input
4702: .vb
4703: numCells = 2, numVertices = 4
4704: cells = [0 1 2 1 3 2]
4705: .ve
4706: which would result in the `DMPLEX`
4707: .vb
4709: 4
4710: / | \
4711: / | \
4712: / | \
4713: 2 0 | 1 5
4714: \ | /
4715: \ | /
4716: \ | /
4717: 3
4718: .ve
4720: Vertices are implicitly numbered consecutively 0,...,NVertices.
4721: Each rank owns a chunk of numVertices consecutive vertices.
4722: If numVertices is `PETSC_DECIDE`, PETSc will distribute them as evenly as possible using PetscLayout.
4723: If NVertices is `PETSC_DETERMINE` and numVertices is PETSC_DECIDE, NVertices is computed by PETSc as the maximum vertex index in cells + 1.
4724: If only NVertices is `PETSC_DETERMINE`, it is computed as the sum of numVertices over all ranks.
4726: The cell distribution is arbitrary non-overlapping, independent of the vertex distribution.
4728: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildFromCellList()`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildCoordinatesFromCellListParallel()`,
4729: `PetscSF`
4730: @*/
4731: PetscErrorCode DMPlexBuildFromCellListParallel(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, const PetscInt cells[], PetscSF *vertexSF, PetscInt **verticesAdjSaved)
4732: {
4733: PetscSF sfPoint;
4734: PetscLayout layout;
4735: PetscInt numVerticesAdj, *verticesAdj, *cones, c, p;
4737: PetscFunctionBegin;
4739: PetscCall(PetscLogEventBegin(DMPLEX_BuildFromCellList, dm, 0, 0, 0));
4740: /* Get/check global number of vertices */
4741: {
4742: PetscInt NVerticesInCells, i;
4743: const PetscInt len = numCells * numCorners;
4745: /* NVerticesInCells = max(cells) + 1 */
4746: NVerticesInCells = PETSC_MIN_INT;
4747: for (i = 0; i < len; i++)
4748: if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
4749: ++NVerticesInCells;
4750: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &NVerticesInCells, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
4752: if (numVertices == PETSC_DECIDE && NVertices == PETSC_DECIDE) NVertices = NVerticesInCells;
4753: else
4754: PetscCheck(NVertices == PETSC_DECIDE || NVertices >= NVerticesInCells, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Specified global number of vertices %" PetscInt_FMT " must be greater than or equal to the number of vertices in cells %" PetscInt_FMT, NVertices, NVerticesInCells);
4755: }
4756: /* Count locally unique vertices */
4757: {
4758: PetscHSetI vhash;
4759: PetscInt off = 0;
4761: PetscCall(PetscHSetICreate(&vhash));
4762: for (c = 0; c < numCells; ++c) {
4763: for (p = 0; p < numCorners; ++p) PetscCall(PetscHSetIAdd(vhash, cells[c * numCorners + p]));
4764: }
4765: PetscCall(PetscHSetIGetSize(vhash, &numVerticesAdj));
4766: if (!verticesAdjSaved) PetscCall(PetscMalloc1(numVerticesAdj, &verticesAdj));
4767: else verticesAdj = *verticesAdjSaved;
4768: PetscCall(PetscHSetIGetElems(vhash, &off, verticesAdj));
4769: PetscCall(PetscHSetIDestroy(&vhash));
4770: PetscCheck(off == numVerticesAdj, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %" PetscInt_FMT " should be %" PetscInt_FMT, off, numVerticesAdj);
4771: }
4772: PetscCall(PetscSortInt(numVerticesAdj, verticesAdj));
4773: /* Create cones */
4774: PetscCall(DMPlexSetChart(dm, 0, numCells + numVerticesAdj));
4775: for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(dm, c, numCorners));
4776: PetscCall(DMSetUp(dm));
4777: PetscCall(DMPlexGetCones(dm, &cones));
4778: for (c = 0; c < numCells; ++c) {
4779: for (p = 0; p < numCorners; ++p) {
4780: const PetscInt gv = cells[c * numCorners + p];
4781: PetscInt lv;
4783: /* Positions within verticesAdj form 0-based local vertex numbering;
4784: we need to shift it by numCells to get correct DAG points (cells go first) */
4785: PetscCall(PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv));
4786: PetscCheck(lv >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %" PetscInt_FMT " in local connectivity", gv);
4787: cones[c * numCorners + p] = lv + numCells;
4788: }
4789: }
4790: /* Build point sf */
4791: PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &layout));
4792: PetscCall(PetscLayoutSetSize(layout, NVertices));
4793: PetscCall(PetscLayoutSetLocalSize(layout, numVertices));
4794: PetscCall(PetscLayoutSetBlockSize(layout, 1));
4795: PetscCall(PetscSFCreateByMatchingIndices(layout, numVerticesAdj, verticesAdj, NULL, numCells, numVerticesAdj, verticesAdj, NULL, numCells, vertexSF, &sfPoint));
4796: PetscCall(PetscLayoutDestroy(&layout));
4797: if (!verticesAdjSaved) PetscCall(PetscFree(verticesAdj));
4798: PetscCall(PetscObjectSetName((PetscObject)sfPoint, "point SF"));
4799: if (dm->sf) {
4800: const char *prefix;
4802: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm->sf, &prefix));
4803: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sfPoint, prefix));
4804: }
4805: PetscCall(DMSetPointSF(dm, sfPoint));
4806: PetscCall(PetscSFDestroy(&sfPoint));
4807: if (vertexSF) PetscCall(PetscObjectSetName((PetscObject)(*vertexSF), "Vertex Ownership SF"));
4808: /* Fill in the rest of the topology structure */
4809: PetscCall(DMPlexSymmetrize(dm));
4810: PetscCall(DMPlexStratify(dm));
4811: PetscCall(PetscLogEventEnd(DMPLEX_BuildFromCellList, dm, 0, 0, 0));
4812: PetscFunctionReturn(PETSC_SUCCESS);
4813: }
4815: /*@C
4816: DMPlexBuildCoordinatesFromCellListParallel - Build `DM` coordinates from a list of coordinates for each owned vertex (common mesh generator output)
4818: Collective; No Fortran Support
4820: Input Parameters:
4821: + dm - The `DM`
4822: . spaceDim - The spatial dimension used for coordinates
4823: . sfVert - `PetscSF` describing complete vertex ownership
4824: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
4826: Level: advanced
4828: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildFromCellListParallel()`
4829: @*/
4830: PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM dm, PetscInt spaceDim, PetscSF sfVert, const PetscReal vertexCoords[])
4831: {
4832: PetscSection coordSection;
4833: Vec coordinates;
4834: PetscScalar *coords;
4835: PetscInt numVertices, numVerticesAdj, coordSize, v, vStart, vEnd;
4837: PetscFunctionBegin;
4838: PetscCall(PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0));
4839: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4840: PetscCheck(vStart >= 0 && vEnd >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first.");
4841: PetscCall(DMSetCoordinateDim(dm, spaceDim));
4842: PetscCall(PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL));
4843: PetscCheck(vEnd - vStart == numVerticesAdj, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Supplied sfVert has wrong number of leaves = %" PetscInt_FMT " != %" PetscInt_FMT " = vEnd - vStart", numVerticesAdj, vEnd - vStart);
4844: PetscCall(DMGetCoordinateSection(dm, &coordSection));
4845: PetscCall(PetscSectionSetNumFields(coordSection, 1));
4846: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim));
4847: PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd));
4848: for (v = vStart; v < vEnd; ++v) {
4849: PetscCall(PetscSectionSetDof(coordSection, v, spaceDim));
4850: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim));
4851: }
4852: PetscCall(PetscSectionSetUp(coordSection));
4853: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
4854: PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &coordinates));
4855: PetscCall(VecSetBlockSize(coordinates, spaceDim));
4856: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
4857: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
4858: PetscCall(VecSetType(coordinates, VECSTANDARD));
4859: PetscCall(VecGetArray(coordinates, &coords));
4860: {
4861: MPI_Datatype coordtype;
4863: /* Need a temp buffer for coords if we have complex/single */
4864: PetscCallMPI(MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype));
4865: PetscCallMPI(MPI_Type_commit(&coordtype));
4866: #if defined(PETSC_USE_COMPLEX)
4867: {
4868: PetscScalar *svertexCoords;
4869: PetscInt i;
4870: PetscCall(PetscMalloc1(numVertices * spaceDim, &svertexCoords));
4871: for (i = 0; i < numVertices * spaceDim; i++) svertexCoords[i] = vertexCoords[i];
4872: PetscCall(PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords, MPI_REPLACE));
4873: PetscCall(PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords, MPI_REPLACE));
4874: PetscCall(PetscFree(svertexCoords));
4875: }
4876: #else
4877: PetscCall(PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords, MPI_REPLACE));
4878: PetscCall(PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords, MPI_REPLACE));
4879: #endif
4880: PetscCallMPI(MPI_Type_free(&coordtype));
4881: }
4882: PetscCall(VecRestoreArray(coordinates, &coords));
4883: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
4884: PetscCall(VecDestroy(&coordinates));
4885: PetscCall(PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0));
4886: PetscFunctionReturn(PETSC_SUCCESS);
4887: }
4889: /*@
4890: DMPlexCreateFromCellListParallelPetsc - Create distributed `DMPLEX` from a list of vertices for each cell (common mesh generator output)
4892: Collective
4894: Input Parameters:
4895: + comm - The communicator
4896: . dim - The topological dimension of the mesh
4897: . numCells - The number of cells owned by this process
4898: . numVertices - The number of vertices owned by this process, or `PETSC_DECIDE`
4899: . NVertices - The global number of vertices, or `PETSC_DECIDE`
4900: . numCorners - The number of vertices for each cell
4901: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
4902: . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
4903: . spaceDim - The spatial dimension used for coordinates
4904: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
4906: Output Parameters:
4907: + dm - The `DM`
4908: . vertexSF - (Optional) `PetscSF` describing complete vertex ownership
4909: - verticesAdj - (Optional) vertex adjacency array
4911: Level: intermediate
4913: Notes:
4914: This function is just a convenient sequence of `DMCreate()`, `DMSetType()`, `DMSetDimension()`,
4915: `DMPlexBuildFromCellListParallel()`, `DMPlexInterpolate()`, `DMPlexBuildCoordinatesFromCellListParallel()`
4917: See `DMPlexBuildFromCellListParallel()` for an example and details about the topology-related parameters.
4919: See `DMPlexBuildCoordinatesFromCellListParallel()` for details about the geometry-related parameters.
4921: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListPetsc()`, `DMPlexBuildFromCellListParallel()`, `DMPlexBuildCoordinatesFromCellListParallel()`, `DMPlexCreateFromDAG()`, `DMPlexCreate()`
4922: @*/
4923: PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, PetscInt **verticesAdj, DM *dm)
4924: {
4925: PetscSF sfVert;
4927: PetscFunctionBegin;
4928: PetscCall(DMCreate(comm, dm));
4929: PetscCall(DMSetType(*dm, DMPLEX));
4932: PetscCall(DMSetDimension(*dm, dim));
4933: PetscCall(DMPlexBuildFromCellListParallel(*dm, numCells, numVertices, NVertices, numCorners, cells, &sfVert, verticesAdj));
4934: if (interpolate) {
4935: DM idm;
4937: PetscCall(DMPlexInterpolate(*dm, &idm));
4938: PetscCall(DMDestroy(dm));
4939: *dm = idm;
4940: }
4941: PetscCall(DMPlexBuildCoordinatesFromCellListParallel(*dm, spaceDim, sfVert, vertexCoords));
4942: if (vertexSF) *vertexSF = sfVert;
4943: else PetscCall(PetscSFDestroy(&sfVert));
4944: PetscFunctionReturn(PETSC_SUCCESS);
4945: }
4947: /*@C
4948: DMPlexBuildFromCellList - Build `DMPLEX` topology from a list of vertices for each cell (common mesh generator output)
4950: Collective; No Fortran Support
4952: Input Parameters:
4953: + dm - The `DM`
4954: . numCells - The number of cells owned by this process
4955: . numVertices - The number of vertices owned by this process, or `PETSC_DETERMINE`
4956: . numCorners - The number of vertices for each cell
4957: - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell
4959: Level: advanced
4961: Notes:
4962: Two triangles sharing a face
4963: .vb
4965: 2
4966: / | \
4967: / | \
4968: / | \
4969: 0 0 | 1 3
4970: \ | /
4971: \ | /
4972: \ | /
4973: 1
4974: .ve
4975: would have input
4976: .vb
4977: numCells = 2, numVertices = 4
4978: cells = [0 1 2 1 3 2]
4979: .ve
4980: which would result in the `DMPLEX`
4981: .vb
4983: 4
4984: / | \
4985: / | \
4986: / | \
4987: 2 0 | 1 5
4988: \ | /
4989: \ | /
4990: \ | /
4991: 3
4992: .ve
4994: If numVertices is `PETSC_DETERMINE`, it is computed by PETSc as the maximum vertex index in cells + 1.
4996: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildFromCellListParallel()`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromCellListPetsc()`
4997: @*/
4998: PetscErrorCode DMPlexBuildFromCellList(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const PetscInt cells[])
4999: {
5000: PetscInt *cones, c, p, dim;
5002: PetscFunctionBegin;
5003: PetscCall(PetscLogEventBegin(DMPLEX_BuildFromCellList, dm, 0, 0, 0));
5004: PetscCall(DMGetDimension(dm, &dim));
5005: /* Get/check global number of vertices */
5006: {
5007: PetscInt NVerticesInCells, i;
5008: const PetscInt len = numCells * numCorners;
5010: /* NVerticesInCells = max(cells) + 1 */
5011: NVerticesInCells = PETSC_MIN_INT;
5012: for (i = 0; i < len; i++)
5013: if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i];
5014: ++NVerticesInCells;
5016: if (numVertices == PETSC_DECIDE) numVertices = NVerticesInCells;
5017: else
5018: PetscCheck(numVertices >= NVerticesInCells, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Specified number of vertices %" PetscInt_FMT " must be greater than or equal to the number of vertices in cells %" PetscInt_FMT, numVertices, NVerticesInCells);
5019: }
5020: PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices));
5021: for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(dm, c, numCorners));
5022: PetscCall(DMSetUp(dm));
5023: PetscCall(DMPlexGetCones(dm, &cones));
5024: for (c = 0; c < numCells; ++c) {
5025: for (p = 0; p < numCorners; ++p) cones[c * numCorners + p] = cells[c * numCorners + p] + numCells;
5026: }
5027: PetscCall(DMPlexSymmetrize(dm));
5028: PetscCall(DMPlexStratify(dm));
5029: PetscCall(PetscLogEventEnd(DMPLEX_BuildFromCellList, dm, 0, 0, 0));
5030: PetscFunctionReturn(PETSC_SUCCESS);
5031: }
5033: /*@C
5034: DMPlexBuildCoordinatesFromCellList - Build `DM` coordinates from a list of coordinates for each owned vertex (common mesh generator output)
5036: Collective; No Fortran Support
5038: Input Parameters:
5039: + dm - The `DM`
5040: . spaceDim - The spatial dimension used for coordinates
5041: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex
5043: Level: advanced
5045: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildCoordinatesFromCellListParallel()`, `DMPlexCreateFromCellListPetsc()`, `DMPlexBuildFromCellList()`
5046: @*/
5047: PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM dm, PetscInt spaceDim, const PetscReal vertexCoords[])
5048: {
5049: PetscSection coordSection;
5050: Vec coordinates;
5051: DM cdm;
5052: PetscScalar *coords;
5053: PetscInt v, vStart, vEnd, d;
5055: PetscFunctionBegin;
5056: PetscCall(PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0));
5057: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
5058: PetscCheck(vStart >= 0 && vEnd >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first.");
5059: PetscCall(DMSetCoordinateDim(dm, spaceDim));
5060: PetscCall(DMGetCoordinateSection(dm, &coordSection));
5061: PetscCall(PetscSectionSetNumFields(coordSection, 1));
5062: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim));
5063: PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd));
5064: for (v = vStart; v < vEnd; ++v) {
5065: PetscCall(PetscSectionSetDof(coordSection, v, spaceDim));
5066: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim));
5067: }
5068: PetscCall(PetscSectionSetUp(coordSection));
5070: PetscCall(DMGetCoordinateDM(dm, &cdm));
5071: PetscCall(DMCreateLocalVector(cdm, &coordinates));
5072: PetscCall(VecSetBlockSize(coordinates, spaceDim));
5073: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
5074: PetscCall(VecGetArrayWrite(coordinates, &coords));
5075: for (v = 0; v < vEnd - vStart; ++v) {
5076: for (d = 0; d < spaceDim; ++d) coords[v * spaceDim + d] = vertexCoords[v * spaceDim + d];
5077: }
5078: PetscCall(VecRestoreArrayWrite(coordinates, &coords));
5079: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
5080: PetscCall(VecDestroy(&coordinates));
5081: PetscCall(PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0));
5082: PetscFunctionReturn(PETSC_SUCCESS);
5083: }
5085: /*@
5086: DMPlexCreateFromCellListPetsc - Create `DMPLEX` from a list of vertices for each cell (common mesh generator output), but only process 0 takes in the input
5088: Collective
5090: Input Parameters:
5091: + comm - The communicator
5092: . dim - The topological dimension of the mesh
5093: . numCells - The number of cells, only on process 0
5094: . numVertices - The number of vertices owned by this process, or `PETSC_DECIDE`, only on process 0
5095: . numCorners - The number of vertices for each cell, only on process 0
5096: . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
5097: . cells - An array of numCells*numCorners numbers, the vertices for each cell, only on process 0
5098: . spaceDim - The spatial dimension used for coordinates
5099: - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex, only on process 0
5101: Output Parameter:
5102: . dm - The `DM`, which only has points on process 0
5104: Level: intermediate
5106: Notes:
5107: This function is just a convenient sequence of `DMCreate()`, `DMSetType()`, `DMSetDimension()`, `DMPlexBuildFromCellList()`,
5108: `DMPlexInterpolate()`, `DMPlexBuildCoordinatesFromCellList()`
5110: See `DMPlexBuildFromCellList()` for an example and details about the topology-related parameters.
5111: See `DMPlexBuildCoordinatesFromCellList()` for details about the geometry-related parameters.
5112: See `DMPlexCreateFromCellListParallelPetsc()` for parallel input
5114: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildFromCellList()`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromDAG()`, `DMPlexCreate()`
5115: @*/
5116: PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], DM *dm)
5117: {
5118: PetscMPIInt rank;
5120: PetscFunctionBegin;
5121: PetscCheck(dim, comm, PETSC_ERR_ARG_OUTOFRANGE, "This is not appropriate for 0-dimensional meshes. Consider either creating the DM using DMPlexCreateFromDAG(), by hand, or using DMSwarm.");
5122: PetscCallMPI(MPI_Comm_rank(comm, &rank));
5123: PetscCall(DMCreate(comm, dm));
5124: PetscCall(DMSetType(*dm, DMPLEX));
5125: PetscCall(DMSetDimension(*dm, dim));
5126: if (rank == 0) PetscCall(DMPlexBuildFromCellList(*dm, numCells, numVertices, numCorners, cells));
5127: else PetscCall(DMPlexBuildFromCellList(*dm, 0, 0, 0, NULL));
5128: if (interpolate) {
5129: DM idm;
5131: PetscCall(DMPlexInterpolate(*dm, &idm));
5132: PetscCall(DMDestroy(dm));
5133: *dm = idm;
5134: }
5135: if (rank == 0) PetscCall(DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, vertexCoords));
5136: else PetscCall(DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, NULL));
5137: PetscFunctionReturn(PETSC_SUCCESS);
5138: }
5140: /*@
5141: DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a `DM`
5143: Input Parameters:
5144: + dm - The empty `DM` object, usually from `DMCreate()` and `DMSetDimension()`
5145: . depth - The depth of the DAG
5146: . numPoints - Array of size depth + 1 containing the number of points at each `depth`
5147: . coneSize - The cone size of each point
5148: . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point
5149: . coneOrientations - The orientation of each cone point
5150: - vertexCoords - An array of `numPoints`[0]*spacedim numbers representing the coordinates of each vertex, with spacedim the value set via `DMSetCoordinateDim()`
5152: Output Parameter:
5153: . dm - The `DM`
5155: Level: advanced
5157: Note:
5158: Two triangles sharing a face would have input
5159: .vb
5160: depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0]
5161: cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0]
5162: vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0]
5163: .ve
5164: which would result in the DMPlex
5165: .vb
5166: 4
5167: / | \
5168: / | \
5169: / | \
5170: 2 0 | 1 5
5171: \ | /
5172: \ | /
5173: \ | /
5174: 3
5175: .ve
5176: Notice that all points are numbered consecutively, unlike `DMPlexCreateFromCellListPetsc()`
5178: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListPetsc()`, `DMPlexCreate()`
5179: @*/
5180: PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[])
5181: {
5182: Vec coordinates;
5183: PetscSection coordSection;
5184: PetscScalar *coords;
5185: PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off;
5187: PetscFunctionBegin;
5188: PetscCall(DMGetDimension(dm, &dim));
5189: PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
5190: PetscCheck(dimEmbed >= dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Embedding dimension %" PetscInt_FMT " cannot be less than intrinsic dimension %" PetscInt_FMT, dimEmbed, dim);
5191: for (d = 0; d <= depth; ++d) pEnd += numPoints[d];
5192: PetscCall(DMPlexSetChart(dm, pStart, pEnd));
5193: for (p = pStart; p < pEnd; ++p) {
5194: PetscCall(DMPlexSetConeSize(dm, p, coneSize[p - pStart]));
5195: if (firstVertex < 0 && !coneSize[p - pStart]) firstVertex = p - pStart;
5196: }
5197: PetscCheck(firstVertex >= 0 || !numPoints[0], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected %" PetscInt_FMT " vertices but could not find any", numPoints[0]);
5198: PetscCall(DMSetUp(dm)); /* Allocate space for cones */
5199: for (p = pStart, off = 0; p < pEnd; off += coneSize[p - pStart], ++p) {
5200: PetscCall(DMPlexSetCone(dm, p, &cones[off]));
5201: PetscCall(DMPlexSetConeOrientation(dm, p, &coneOrientations[off]));
5202: }
5203: PetscCall(DMPlexSymmetrize(dm));
5204: PetscCall(DMPlexStratify(dm));
5205: /* Build coordinates */
5206: PetscCall(DMGetCoordinateSection(dm, &coordSection));
5207: PetscCall(PetscSectionSetNumFields(coordSection, 1));
5208: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed));
5209: PetscCall(PetscSectionSetChart(coordSection, firstVertex, firstVertex + numPoints[0]));
5210: for (v = firstVertex; v < firstVertex + numPoints[0]; ++v) {
5211: PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed));
5212: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed));
5213: }
5214: PetscCall(PetscSectionSetUp(coordSection));
5215: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
5216: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
5217: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
5218: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
5219: PetscCall(VecSetBlockSize(coordinates, dimEmbed));
5220: PetscCall(VecSetType(coordinates, VECSTANDARD));
5221: if (vertexCoords) {
5222: PetscCall(VecGetArray(coordinates, &coords));
5223: for (v = 0; v < numPoints[0]; ++v) {
5224: PetscInt off;
5226: PetscCall(PetscSectionGetOffset(coordSection, v + firstVertex, &off));
5227: for (d = 0; d < dimEmbed; ++d) coords[off + d] = vertexCoords[v * dimEmbed + d];
5228: }
5229: }
5230: PetscCall(VecRestoreArray(coordinates, &coords));
5231: PetscCall(DMSetCoordinatesLocal(dm, coordinates));
5232: PetscCall(VecDestroy(&coordinates));
5233: PetscFunctionReturn(PETSC_SUCCESS);
5234: }
5236: /*
5237: DMPlexCreateCellVertexFromFile - Create a `DMPLEX` mesh from a simple cell-vertex file.
5239: Collective
5241: + comm - The MPI communicator
5242: . filename - Name of the .dat file
5243: - interpolate - Create faces and edges in the mesh
5245: Output Parameter:
5246: . dm - The `DM` object representing the mesh
5248: Level: beginner
5250: Note:
5251: The format is the simplest possible:
5252: .vb
5253: Ne
5254: v0 v1 ... vk
5255: Nv
5256: x y z marker
5257: .ve
5259: Developer Note:
5260: Should use a `PetscViewer` not a filename
5262: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromFile()`, `DMPlexCreateMedFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()`
5263: */
5264: static PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
5265: {
5266: DMLabel marker;
5267: PetscViewer viewer;
5268: Vec coordinates;
5269: PetscSection coordSection;
5270: PetscScalar *coords;
5271: char line[PETSC_MAX_PATH_LEN];
5272: PetscInt dim = 3, cdim = 3, coordSize, v, c, d;
5273: PetscMPIInt rank;
5274: int snum, Nv, Nc, Ncn, Nl;
5276: PetscFunctionBegin;
5277: PetscCallMPI(MPI_Comm_rank(comm, &rank));
5278: PetscCall(PetscViewerCreate(comm, &viewer));
5279: PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
5280: PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
5281: PetscCall(PetscViewerFileSetName(viewer, filename));
5282: if (rank == 0) {
5283: PetscCall(PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING));
5284: snum = sscanf(line, "%d %d %d %d", &Nc, &Nv, &Ncn, &Nl);
5285: PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
5286: } else {
5287: Nc = Nv = Ncn = Nl = 0;
5288: }
5289: PetscCall(DMCreate(comm, dm));
5290: PetscCall(DMSetType(*dm, DMPLEX));
5291: PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv));
5292: PetscCall(DMSetDimension(*dm, dim));
5293: PetscCall(DMSetCoordinateDim(*dm, cdim));
5294: /* Read topology */
5295: if (rank == 0) {
5296: char format[PETSC_MAX_PATH_LEN];
5297: PetscInt cone[8];
5298: int vbuf[8], v;
5300: for (c = 0; c < Ncn; ++c) {
5301: format[c * 3 + 0] = '%';
5302: format[c * 3 + 1] = 'd';
5303: format[c * 3 + 2] = ' ';
5304: }
5305: format[Ncn * 3 - 1] = '\0';
5306: for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, Ncn));
5307: PetscCall(DMSetUp(*dm));
5308: for (c = 0; c < Nc; ++c) {
5309: PetscCall(PetscViewerRead(viewer, line, Ncn, NULL, PETSC_STRING));
5310: switch (Ncn) {
5311: case 2:
5312: snum = sscanf(line, format, &vbuf[0], &vbuf[1]);
5313: break;
5314: case 3:
5315: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2]);
5316: break;
5317: case 4:
5318: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3]);
5319: break;
5320: case 6:
5321: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5]);
5322: break;
5323: case 8:
5324: snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5], &vbuf[6], &vbuf[7]);
5325: break;
5326: default:
5327: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No cell shape with %d vertices", Ncn);
5328: }
5329: PetscCheck(snum == Ncn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
5330: for (v = 0; v < Ncn; ++v) cone[v] = vbuf[v] + Nc;
5331: /* Hexahedra are inverted */
5332: if (Ncn == 8) {
5333: PetscInt tmp = cone[1];
5334: cone[1] = cone[3];
5335: cone[3] = tmp;
5336: }
5337: PetscCall(DMPlexSetCone(*dm, c, cone));
5338: }
5339: }
5340: PetscCall(DMPlexSymmetrize(*dm));
5341: PetscCall(DMPlexStratify(*dm));
5342: /* Read coordinates */
5343: PetscCall(DMGetCoordinateSection(*dm, &coordSection));
5344: PetscCall(PetscSectionSetNumFields(coordSection, 1));
5345: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim));
5346: PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv));
5347: for (v = Nc; v < Nc + Nv; ++v) {
5348: PetscCall(PetscSectionSetDof(coordSection, v, cdim));
5349: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim));
5350: }
5351: PetscCall(PetscSectionSetUp(coordSection));
5352: PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
5353: PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
5354: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
5355: PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
5356: PetscCall(VecSetBlockSize(coordinates, cdim));
5357: PetscCall(VecSetType(coordinates, VECSTANDARD));
5358: PetscCall(VecGetArray(coordinates, &coords));
5359: if (rank == 0) {
5360: char format[PETSC_MAX_PATH_LEN];
5361: double x[3];
5362: int l, val[3];
5364: if (Nl) {
5365: for (l = 0; l < Nl; ++l) {
5366: format[l * 3 + 0] = '%';
5367: format[l * 3 + 1] = 'd';
5368: format[l * 3 + 2] = ' ';
5369: }
5370: format[Nl * 3 - 1] = '\0';
5371: PetscCall(DMCreateLabel(*dm, "marker"));
5372: PetscCall(DMGetLabel(*dm, "marker", &marker));
5373: }
5374: for (v = 0; v < Nv; ++v) {
5375: PetscCall(PetscViewerRead(viewer, line, 3 + Nl, NULL, PETSC_STRING));
5376: snum = sscanf(line, "%lg %lg %lg", &x[0], &x[1], &x[2]);
5377: PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
5378: switch (Nl) {
5379: case 0:
5380: snum = 0;
5381: break;
5382: case 1:
5383: snum = sscanf(line, format, &val[0]);
5384: break;
5385: case 2:
5386: snum = sscanf(line, format, &val[0], &val[1]);
5387: break;
5388: case 3:
5389: snum = sscanf(line, format, &val[0], &val[1], &val[2]);
5390: break;
5391: default:
5392: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Request support for %d labels", Nl);
5393: }
5394: PetscCheck(snum == Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line);
5395: for (d = 0; d < cdim; ++d) coords[v * cdim + d] = x[d];
5396: for (l = 0; l < Nl; ++l) PetscCall(DMLabelSetValue(marker, v + Nc, val[l]));
5397: }
5398: }
5399: PetscCall(VecRestoreArray(coordinates, &coords));
5400: PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
5401: PetscCall(VecDestroy(&coordinates));
5402: PetscCall(PetscViewerDestroy(&viewer));
5403: if (interpolate) {
5404: DM idm;
5405: DMLabel bdlabel;
5407: PetscCall(DMPlexInterpolate(*dm, &idm));
5408: PetscCall(DMDestroy(dm));
5409: *dm = idm;
5411: if (!Nl) {
5412: PetscCall(DMCreateLabel(*dm, "marker"));
5413: PetscCall(DMGetLabel(*dm, "marker", &bdlabel));
5414: PetscCall(DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel));
5415: PetscCall(DMPlexLabelComplete(*dm, bdlabel));
5416: }
5417: }
5418: PetscFunctionReturn(PETSC_SUCCESS);
5419: }
5421: /*@C
5422: DMPlexCreateFromFile - This takes a filename and produces a `DM`
5424: Collective
5426: Input Parameters:
5427: + comm - The communicator
5428: . filename - A file name
5429: . plexname - The object name of the resulting `DM`, also used for intra-datafile lookup by some formats
5430: - interpolate - Flag to create intermediate mesh pieces (edges, faces)
5432: Output Parameter:
5433: . dm - The `DM`
5435: Options Database Key:
5436: . -dm_plex_create_from_hdf5_xdmf - use the `PETSC_VIEWER_HDF5_XDMF` format for reading HDF5
5438: Use `-dm_plex_create_ prefix` to pass options to the internal `PetscViewer`, e.g.
5439: $ -dm_plex_create_viewer_hdf5_collective
5441: Level: beginner
5443: Notes:
5444: Using `PETSCVIEWERHDF5` type with `PETSC_VIEWER_HDF5_PETSC` format, one can save multiple `DMPLEX`
5445: meshes in a single HDF5 file. This in turn requires one to name the `DMPLEX` object with `PetscObjectSetName()`
5446: before saving it with `DMView()` and before loading it with `DMLoad()` for identification of the mesh object.
5447: The input parameter name is thus used to name the `DMPLEX` object when `DMPlexCreateFromFile()` internally
5448: calls `DMLoad()`. Currently, name is ignored for other viewer types and/or formats.
5450: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromDAG()`, `DMPlexCreateFromCellListPetsc()`, `DMPlexCreate()`, `PetscObjectSetName()`, `DMView()`, `DMLoad()`
5451: @*/
5452: PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], const char plexname[], PetscBool interpolate, DM *dm)
5453: {
5454: const char extGmsh[] = ".msh";
5455: const char extGmsh2[] = ".msh2";
5456: const char extGmsh4[] = ".msh4";
5457: const char extCGNS[] = ".cgns";
5458: const char extExodus[] = ".exo";
5459: const char extExodus_e[] = ".e";
5460: const char extGenesis[] = ".gen";
5461: const char extFluent[] = ".cas";
5462: const char extHDF5[] = ".h5";
5463: const char extXDMFHDF5[] = ".xdmf.h5";
5464: const char extMed[] = ".med";
5465: const char extPLY[] = ".ply";
5466: const char extEGADSLite[] = ".egadslite";
5467: const char extEGADS[] = ".egads";
5468: const char extIGES[] = ".igs";
5469: const char extSTEP[] = ".stp";
5470: const char extCV[] = ".dat";
5471: size_t len;
5472: PetscBool isGmsh, isGmsh2, isGmsh4, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isEGADSLite, isEGADS, isIGES, isSTEP, isCV, isXDMFHDF5;
5473: PetscMPIInt rank;
5475: PetscFunctionBegin;
5476: PetscAssertPointer(filename, 2);
5477: if (plexname) PetscAssertPointer(plexname, 3);
5478: PetscAssertPointer(dm, 5);
5479: PetscCall(DMInitializePackage());
5480: PetscCall(PetscLogEventBegin(DMPLEX_CreateFromFile, 0, 0, 0, 0));
5481: PetscCallMPI(MPI_Comm_rank(comm, &rank));
5482: PetscCall(PetscStrlen(filename, &len));
5483: PetscCheck(len, comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path");
5485: #define CheckExtension(extension__, is_extension__) \
5486: do { \
5487: PetscAssert(sizeof(extension__), comm, PETSC_ERR_PLIB, "Zero-size extension: %s", extension__); \
5488: /* don't count the null-terminator at the end */ \
5489: const size_t ext_len = sizeof(extension__) - 1; \
5490: if (len < ext_len) { \
5491: is_extension__ = PETSC_FALSE; \
5492: } else { \
5493: PetscCall(PetscStrncmp(filename + len - ext_len, extension__, ext_len, &is_extension__)); \
5494: } \
5495: } while (0)
5497: CheckExtension(extGmsh, isGmsh);
5498: CheckExtension(extGmsh2, isGmsh2);
5499: CheckExtension(extGmsh4, isGmsh4);
5500: CheckExtension(extCGNS, isCGNS);
5501: CheckExtension(extExodus, isExodus);
5502: if (!isExodus) CheckExtension(extExodus_e, isExodus);
5503: CheckExtension(extGenesis, isGenesis);
5504: CheckExtension(extFluent, isFluent);
5505: CheckExtension(extHDF5, isHDF5);
5506: CheckExtension(extMed, isMed);
5507: CheckExtension(extPLY, isPLY);
5508: CheckExtension(extEGADSLite, isEGADSLite);
5509: CheckExtension(extEGADS, isEGADS);
5510: CheckExtension(extIGES, isIGES);
5511: CheckExtension(extSTEP, isSTEP);
5512: CheckExtension(extCV, isCV);
5513: CheckExtension(extXDMFHDF5, isXDMFHDF5);
5515: #undef CheckExtension
5517: if (isGmsh || isGmsh2 || isGmsh4) {
5518: PetscCall(DMPlexCreateGmshFromFile(comm, filename, interpolate, dm));
5519: } else if (isCGNS) {
5520: PetscCall(DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm));
5521: } else if (isExodus || isGenesis) {
5522: PetscCall(DMPlexCreateExodusFromFile(comm, filename, interpolate, dm));
5523: } else if (isFluent) {
5524: PetscCall(DMPlexCreateFluentFromFile(comm, filename, interpolate, dm));
5525: } else if (isHDF5) {
5526: PetscViewer viewer;
5528: /* PETSC_VIEWER_HDF5_XDMF is used if the filename ends with .xdmf.h5, or if -dm_plex_create_from_hdf5_xdmf option is present */
5529: PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &isXDMFHDF5, NULL));
5530: PetscCall(PetscViewerCreate(comm, &viewer));
5531: PetscCall(PetscViewerSetType(viewer, PETSCVIEWERHDF5));
5532: PetscCall(PetscViewerSetOptionsPrefix(viewer, "dm_plex_create_"));
5533: PetscCall(PetscViewerSetFromOptions(viewer));
5534: PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
5535: PetscCall(PetscViewerFileSetName(viewer, filename));
5537: PetscCall(DMCreate(comm, dm));
5538: PetscCall(PetscObjectSetName((PetscObject)(*dm), plexname));
5539: PetscCall(DMSetType(*dm, DMPLEX));
5540: if (isXDMFHDF5) PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF));
5541: PetscCall(DMLoad(*dm, viewer));
5542: if (isXDMFHDF5) PetscCall(PetscViewerPopFormat(viewer));
5543: PetscCall(PetscViewerDestroy(&viewer));
5545: if (interpolate) {
5546: DM idm;
5548: PetscCall(DMPlexInterpolate(*dm, &idm));
5549: PetscCall(DMDestroy(dm));
5550: *dm = idm;
5551: }
5552: } else if (isMed) {
5553: PetscCall(DMPlexCreateMedFromFile(comm, filename, interpolate, dm));
5554: } else if (isPLY) {
5555: PetscCall(DMPlexCreatePLYFromFile(comm, filename, interpolate, dm));
5556: } else if (isEGADSLite || isEGADS || isIGES || isSTEP) {
5557: if (isEGADSLite) PetscCall(DMPlexCreateEGADSLiteFromFile(comm, filename, dm));
5558: else PetscCall(DMPlexCreateEGADSFromFile(comm, filename, dm));
5559: if (!interpolate) {
5560: DM udm;
5562: PetscCall(DMPlexUninterpolate(*dm, &udm));
5563: PetscCall(DMDestroy(dm));
5564: *dm = udm;
5565: }
5566: } else if (isCV) {
5567: PetscCall(DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm));
5568: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename);
5569: PetscCall(PetscStrlen(plexname, &len));
5570: if (len) PetscCall(PetscObjectSetName((PetscObject)(*dm), plexname));
5571: PetscCall(PetscLogEventEnd(DMPLEX_CreateFromFile, 0, 0, 0, 0));
5572: PetscFunctionReturn(PETSC_SUCCESS);
5573: }
5575: /*@C
5576: DMPlexCreateEphemeral - This takes a `DMPlexTransform` and a base `DMPlex` and produces an ephemeral `DM`, meaning one that is created on the fly in response to queries.
5578: Input Parameters:
5579: + tr - The `DMPlexTransform`
5580: - prefix - An options prefix, or NULL
5582: Output Parameter:
5583: . dm - The `DM`
5585: Level: beginner
5587: Notes:
5588: An emphemeral mesh is one that is not stored concretely, as in the default `DMPLEX` implementation, but rather is produced on the fly in response to queries, using information from the transform and the base mesh.
5590: .seealso: `DMPlexCreateFromFile`, `DMPlexCreateFromDAG()`, `DMPlexCreateFromCellListPetsc()`, `DMPlexCreate()`
5591: @*/
5592: PetscErrorCode DMPlexCreateEphemeral(DMPlexTransform tr, const char prefix[], DM *dm)
5593: {
5594: DM bdm, bcdm, cdm;
5595: Vec coordinates, coordinatesNew;
5596: PetscSection cs;
5597: PetscInt dim, cdim, Nl;
5599: PetscFunctionBegin;
5600: PetscCall(DMCreate(PetscObjectComm((PetscObject)tr), dm));
5601: PetscCall(DMSetType(*dm, DMPLEX));
5602: ((DM_Plex *)(*dm)->data)->interpolated = DMPLEX_INTERPOLATED_FULL;
5603: // Handle coordinates
5604: PetscCall(DMPlexTransformGetDM(tr, &bdm));
5605: PetscCall(DMGetCoordinateDim(bdm, &cdim));
5606: PetscCall(DMSetCoordinateDim(*dm, cdim));
5607: PetscCall(DMGetDimension(bdm, &dim));
5608: PetscCall(DMSetDimension(*dm, dim));
5609: PetscCall(DMGetCoordinateDM(bdm, &bcdm));
5610: PetscCall(DMGetCoordinateDM(*dm, &cdm));
5611: PetscCall(DMCopyDisc(bcdm, cdm));
5612: PetscCall(DMGetLocalSection(cdm, &cs));
5613: PetscCall(PetscSectionSetNumFields(cs, 1));
5614: PetscCall(PetscSectionSetFieldComponents(cs, 0, cdim));
5615: PetscCall(DMGetCoordinatesLocal(bdm, &coordinates));
5616: PetscCall(VecDuplicate(coordinates, &coordinatesNew));
5617: PetscCall(VecCopy(coordinates, coordinatesNew));
5618: PetscCall(DMSetCoordinatesLocal(*dm, coordinatesNew));
5619: PetscCall(VecDestroy(&coordinatesNew));
5621: PetscCall(PetscObjectReference((PetscObject)tr));
5622: PetscCall(DMPlexTransformDestroy(&((DM_Plex *)(*dm)->data)->tr));
5623: ((DM_Plex *)(*dm)->data)->tr = tr;
5624: PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
5625: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, prefix));
5626: PetscCall(DMSetFromOptions(*dm));
5628: PetscCall(DMGetNumLabels(bdm, &Nl));
5629: for (PetscInt l = 0; l < Nl; ++l) {
5630: DMLabel label, labelNew;
5631: const char *lname;
5632: PetscBool isDepth, isCellType;
5634: PetscCall(DMGetLabelName(bdm, l, &lname));
5635: PetscCall(PetscStrcmp(lname, "depth", &isDepth));
5636: if (isDepth) continue;
5637: PetscCall(PetscStrcmp(lname, "celltype", &isCellType));
5638: if (isCellType) continue;
5639: PetscCall(DMCreateLabel(*dm, lname));
5640: PetscCall(DMGetLabel(bdm, lname, &label));
5641: PetscCall(DMGetLabel(*dm, lname, &labelNew));
5642: PetscCall(DMLabelSetType(labelNew, DMLABELEPHEMERAL));
5643: PetscCall(DMLabelEphemeralSetLabel(labelNew, label));
5644: PetscCall(DMLabelEphemeralSetTransform(labelNew, tr));
5645: PetscCall(DMLabelSetUp(labelNew));
5646: }
5647: PetscFunctionReturn(PETSC_SUCCESS);
5648: }