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