Actual source code: dmperiodicity.c
1: #include <petsc/private/dmimpl.h>
3: #include <petscdmplex.h>
5: /*@C
6: DMGetPeriodicity - Get the description of mesh periodicity
8: Input Parameter:
9: . dm - The DM object
11: Output Parameters:
12: + maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates
13: . Lstart - If we assume the mesh is a torus, this is the start of each coordinate, or NULL for 0.0
14: - L - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0
16: Level: developer
18: .seealso: `DMGetPeriodicity()`
19: @*/
20: PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal **maxCell, const PetscReal **Lstart, const PetscReal **L)
21: {
22: PetscFunctionBegin;
24: if (maxCell) *maxCell = dm->maxCell;
25: if (Lstart) *Lstart = dm->Lstart;
26: if (L) *L = dm->L;
27: PetscFunctionReturn(PETSC_SUCCESS);
28: }
30: /*@C
31: DMSetPeriodicity - Set the description of mesh periodicity
33: Input Parameters:
34: + dm - The DM object
35: . maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates. Pass NULL to remove such information.
36: . Lstart - If we assume the mesh is a torus, this is the start of each coordinate, or NULL for 0.0
37: - L - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0
39: Level: developer
41: .seealso: `DMGetPeriodicity()`
42: @*/
43: PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal Lstart[], const PetscReal L[])
44: {
45: PetscInt dim, d;
47: PetscFunctionBegin;
52: PetscCall(DMGetDimension(dm, &dim));
53: if (maxCell) {
54: if (!dm->maxCell) PetscCall(PetscMalloc1(dim, &dm->maxCell));
55: for (d = 0; d < dim; ++d) dm->maxCell[d] = maxCell[d];
56: } else { /* remove maxCell information to disable automatic computation of localized vertices */
57: PetscCall(PetscFree(dm->maxCell));
58: dm->maxCell = NULL;
59: }
60: if (Lstart) {
61: if (!dm->Lstart) PetscCall(PetscMalloc1(dim, &dm->Lstart));
62: for (d = 0; d < dim; ++d) dm->Lstart[d] = Lstart[d];
63: } else { /* remove L information to disable automatic computation of localized vertices */
64: PetscCall(PetscFree(dm->Lstart));
65: dm->Lstart = NULL;
66: }
67: if (L) {
68: if (!dm->L) PetscCall(PetscMalloc1(dim, &dm->L));
69: for (d = 0; d < dim; ++d) dm->L[d] = L[d];
70: } else { /* remove L information to disable automatic computation of localized vertices */
71: PetscCall(PetscFree(dm->L));
72: dm->L = NULL;
73: }
74: PetscCheck((dm->maxCell && dm->L) || (!dm->maxCell && !dm->L), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot set only one of maxCell/L");
75: PetscFunctionReturn(PETSC_SUCCESS);
76: }
78: /*@
79: DMLocalizeCoordinate - If a mesh is periodic (a torus with lengths L_i, some of which can be infinite), project the coordinate onto [0, L_i) in each dimension.
81: Input Parameters:
82: + dm - The DM
83: . in - The input coordinate point (dim numbers)
84: - endpoint - Include the endpoint L_i
86: Output Parameter:
87: . out - The localized coordinate point
89: Level: developer
91: .seealso: `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()`
92: @*/
93: PetscErrorCode DMLocalizeCoordinate(DM dm, const PetscScalar in[], PetscBool endpoint, PetscScalar out[])
94: {
95: PetscInt dim, d;
97: PetscFunctionBegin;
98: PetscCall(DMGetCoordinateDim(dm, &dim));
99: if (!dm->maxCell) {
100: for (d = 0; d < dim; ++d) out[d] = in[d];
101: } else {
102: if (endpoint) {
103: for (d = 0; d < dim; ++d) {
104: if ((PetscAbsReal(PetscRealPart(in[d]) / dm->L[d] - PetscFloorReal(PetscRealPart(in[d]) / dm->L[d])) < PETSC_SMALL) && (PetscRealPart(in[d]) / dm->L[d] > PETSC_SMALL)) {
105: out[d] = in[d] - dm->L[d] * (PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]) - 1);
106: } else {
107: out[d] = in[d] - dm->L[d] * PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]);
108: }
109: }
110: } else {
111: for (d = 0; d < dim; ++d) out[d] = in[d] - dm->L[d] * PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]);
112: }
113: }
114: PetscFunctionReturn(PETSC_SUCCESS);
115: }
117: /*
118: DMLocalizeCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer.
120: Input Parameters:
121: + dm - The DM
122: . dim - The spatial dimension
123: . anchor - The anchor point, the input point can be no more than maxCell away from it
124: - in - The input coordinate point (dim numbers)
126: Output Parameter:
127: . out - The localized coordinate point
129: Level: developer
131: Note: This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell
133: .seealso: `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()`
134: */
135: PetscErrorCode DMLocalizeCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[])
136: {
137: PetscInt d;
139: PetscFunctionBegin;
140: if (!dm->maxCell) {
141: for (d = 0; d < dim; ++d) out[d] = in[d];
142: } else {
143: for (d = 0; d < dim; ++d) {
144: if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) {
145: out[d] = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d];
146: } else {
147: out[d] = in[d];
148: }
149: }
150: }
151: PetscFunctionReturn(PETSC_SUCCESS);
152: }
154: PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim, const PetscReal anchor[], const PetscReal in[], PetscReal out[])
155: {
156: PetscInt d;
158: PetscFunctionBegin;
159: if (!dm->maxCell) {
160: for (d = 0; d < dim; ++d) out[d] = in[d];
161: } else {
162: for (d = 0; d < dim; ++d) {
163: if ((dm->L[d] > 0.0) && (PetscAbsReal(anchor[d] - in[d]) > dm->maxCell[d])) {
164: out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d];
165: } else {
166: out[d] = in[d];
167: }
168: }
169: }
170: PetscFunctionReturn(PETSC_SUCCESS);
171: }
173: /*
174: DMLocalizeAddCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer.
176: Input Parameters:
177: + dm - The DM
178: . dim - The spatial dimension
179: . anchor - The anchor point, the input point can be no more than maxCell away from it
180: . in - The input coordinate delta (dim numbers)
181: - out - The input coordinate point (dim numbers)
183: Output Parameter:
184: . out - The localized coordinate in + out
186: Level: developer
188: Note: This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell
190: .seealso: `DMLocalizeCoordinates()`, `DMLocalizeCoordinate()`
191: */
192: PetscErrorCode DMLocalizeAddCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[])
193: {
194: PetscInt d;
196: PetscFunctionBegin;
197: if (!dm->maxCell) {
198: for (d = 0; d < dim; ++d) out[d] += in[d];
199: } else {
200: for (d = 0; d < dim; ++d) {
201: const PetscReal maxC = dm->maxCell[d];
203: if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > maxC)) {
204: const PetscScalar newCoord = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d];
206: if (PetscAbsScalar(newCoord - anchor[d]) > maxC)
207: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT "-Coordinate %g more than %g away from anchor %g", d, (double)PetscRealPart(in[d]), (double)maxC, (double)PetscRealPart(anchor[d]));
208: out[d] += newCoord;
209: } else {
210: out[d] += in[d];
211: }
212: }
213: }
214: PetscFunctionReturn(PETSC_SUCCESS);
215: }
217: /*@
218: DMGetCoordinatesLocalizedLocal - Check if the DM coordinates have been localized for cells on this process
220: Not collective
222: Input Parameter:
223: . dm - The DM
225: Output Parameter:
226: areLocalized - True if localized
228: Level: developer
230: .seealso: `DMLocalizeCoordinates()`, `DMGetCoordinatesLocalized()`, `DMSetPeriodicity()`
231: @*/
232: PetscErrorCode DMGetCoordinatesLocalizedLocal(DM dm, PetscBool *areLocalized)
233: {
234: PetscFunctionBegin;
237: *areLocalized = dm->coordinates[1].dim < 0 ? PETSC_FALSE : PETSC_TRUE;
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: /*@
242: DMGetCoordinatesLocalized - Check if the DM coordinates have been localized for cells
244: Collective on dm
246: Input Parameter:
247: . dm - The DM
249: Output Parameter:
250: areLocalized - True if localized
252: Level: developer
254: .seealso: `DMLocalizeCoordinates()`, `DMSetPeriodicity()`, `DMGetCoordinatesLocalizedLocal()`
255: @*/
256: PetscErrorCode DMGetCoordinatesLocalized(DM dm, PetscBool *areLocalized)
257: {
258: PetscBool localized;
260: PetscFunctionBegin;
263: PetscCall(DMGetCoordinatesLocalizedLocal(dm, &localized));
264: PetscCall(MPIU_Allreduce(&localized, areLocalized, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
265: PetscFunctionReturn(PETSC_SUCCESS);
266: }
268: /*@
269: DMLocalizeCoordinates - If a mesh is periodic, create local coordinates for cells having periodic faces
271: Collective on dm
273: Input Parameter:
274: . dm - The DM
276: Level: developer
278: .seealso: `DMSetPeriodicity()`, `DMLocalizeCoordinate()`, `DMLocalizeAddCoordinate()`
279: @*/
280: PetscErrorCode DMLocalizeCoordinates(DM dm)
281: {
282: DM cdm, cdgdm, cplex, plex;
283: PetscSection cs, csDG;
284: Vec coordinates, cVec;
285: PetscScalar *coordsDG, *anchor, *localized;
286: const PetscReal *Lstart, *L;
287: PetscInt Nc, vStart, vEnd, sStart, sEnd, newStart = PETSC_MAX_INT, newEnd = PETSC_MIN_INT, bs, coordSize;
288: PetscBool isLocalized, sparseLocalize = dm->sparseLocalize, useDG = PETSC_FALSE, useDGGlobal;
289: PetscInt maxHeight = 0, h;
290: PetscInt *pStart = NULL, *pEnd = NULL;
291: MPI_Comm comm;
293: PetscFunctionBegin;
295: PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
296: /* Cannot automatically localize without L and maxCell right now */
297: if (!L) PetscFunctionReturn(PETSC_SUCCESS);
298: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
299: PetscCall(DMGetCoordinatesLocalized(dm, &isLocalized));
300: if (isLocalized) PetscFunctionReturn(PETSC_SUCCESS);
302: PetscCall(DMGetCoordinateDM(dm, &cdm));
303: PetscCall(DMConvert(dm, DMPLEX, &plex));
304: PetscCall(DMConvert(cdm, DMPLEX, &cplex));
305: if (cplex) {
306: PetscCall(DMPlexGetDepthStratum(cplex, 0, &vStart, &vEnd));
307: PetscCall(DMPlexGetMaxProjectionHeight(cplex, &maxHeight));
308: PetscCall(DMGetWorkArray(dm, 2 * (maxHeight + 1), MPIU_INT, &pStart));
309: pEnd = &pStart[maxHeight + 1];
310: newStart = vStart;
311: newEnd = vEnd;
312: for (h = 0; h <= maxHeight; h++) {
313: PetscCall(DMPlexGetHeightStratum(cplex, h, &pStart[h], &pEnd[h]));
314: newStart = PetscMin(newStart, pStart[h]);
315: newEnd = PetscMax(newEnd, pEnd[h]);
316: }
317: } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM");
318: PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
319: PetscCheck(coordinates, comm, PETSC_ERR_SUP, "Missing local coordinate vector");
320: PetscCall(DMGetCoordinateSection(dm, &cs));
321: PetscCall(VecGetBlockSize(coordinates, &bs));
322: PetscCall(PetscSectionGetChart(cs, &sStart, &sEnd));
324: PetscCall(PetscSectionCreate(comm, &csDG));
325: PetscCall(PetscSectionSetNumFields(csDG, 1));
326: PetscCall(PetscSectionGetFieldComponents(cs, 0, &Nc));
327: PetscCall(PetscSectionSetFieldComponents(csDG, 0, Nc));
328: PetscCall(PetscSectionSetChart(csDG, newStart, newEnd));
329: PetscCheck(bs == Nc, comm, PETSC_ERR_ARG_INCOMP, "Coordinate block size %" PetscInt_FMT " != %" PetscInt_FMT " number of components", bs, Nc);
331: PetscCall(DMGetWorkArray(dm, 2 * Nc, MPIU_SCALAR, &anchor));
332: localized = &anchor[Nc];
333: for (h = 0; h <= maxHeight; h++) {
334: PetscInt cStart = pStart[h], cEnd = pEnd[h], c;
336: for (c = cStart; c < cEnd; ++c) {
337: PetscScalar *cellCoords = NULL;
338: DMPolytopeType ct;
339: PetscInt dof, d, p;
341: PetscCall(DMPlexGetCellType(plex, c, &ct));
342: if (ct == DM_POLYTOPE_FV_GHOST) continue;
343: PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
344: PetscCheck(!(dof % Nc), comm, PETSC_ERR_ARG_INCOMP, "Coordinate size on cell %" PetscInt_FMT " closure %" PetscInt_FMT " not divisible by %" PetscInt_FMT " number of components", c, dof, Nc);
345: for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[d];
346: for (p = 0; p < dof / Nc; ++p) {
347: PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p * Nc], localized));
348: for (d = 0; d < Nc; ++d)
349: if (cellCoords[p * Nc + d] != localized[d]) break;
350: if (d < Nc) break;
351: }
352: if (p < dof / Nc) useDG = PETSC_TRUE;
353: if (p < dof / Nc || !sparseLocalize) {
354: PetscCall(PetscSectionSetDof(csDG, c, dof));
355: PetscCall(PetscSectionSetFieldDof(csDG, c, 0, dof));
356: }
357: PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
358: }
359: }
360: PetscCallMPI(MPI_Allreduce(&useDG, &useDGGlobal, 1, MPIU_BOOL, MPI_LOR, comm));
361: if (!useDGGlobal) goto end;
363: PetscCall(PetscSectionSetUp(csDG));
364: PetscCall(PetscSectionGetStorageSize(csDG, &coordSize));
365: PetscCall(VecCreate(PETSC_COMM_SELF, &cVec));
366: PetscCall(PetscObjectSetName((PetscObject)cVec, "coordinates"));
367: PetscCall(VecSetBlockSize(cVec, bs));
368: PetscCall(VecSetSizes(cVec, coordSize, PETSC_DETERMINE));
369: PetscCall(VecSetType(cVec, VECSTANDARD));
370: PetscCall(VecGetArray(cVec, &coordsDG));
371: for (h = 0; h <= maxHeight; h++) {
372: PetscInt cStart = pStart[h], cEnd = pEnd[h], c;
374: for (c = cStart; c < cEnd; ++c) {
375: PetscScalar *cellCoords = NULL;
376: PetscInt p = 0, q, dof, cdof, d, offDG;
378: PetscCall(PetscSectionGetDof(csDG, c, &cdof));
379: if (!cdof) continue;
380: PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
381: PetscCall(PetscSectionGetOffset(csDG, c, &offDG));
382: // TODO The coordinates are set in closure order, which might not be the tensor order
383: for (q = 0; q < dof / Nc; ++q) {
384: // Select a trial anchor
385: for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[q * Nc + d];
386: for (p = 0; p < dof / Nc; ++p) {
387: PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p * Nc], &coordsDG[offDG + p * Nc]));
388: // We need the cell to fit into the torus [lower, lower+L)
389: for (d = 0; d < Nc; ++d)
390: if (L[d] > 0. && ((PetscRealPart(coordsDG[offDG + p * Nc + d]) < (Lstart ? Lstart[d] : 0.)) || (PetscRealPart(coordsDG[offDG + p * Nc + d]) > (Lstart ? Lstart[d] : 0.) + L[d]))) break;
391: if (d < Nc) break;
392: }
393: if (p == dof / Nc) break;
394: }
395: PetscCheck(p == dof / Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " does not fit into the torus %s[0, L]", c, Lstart ? "Lstart + " : "");
396: PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
397: }
398: }
399: PetscCall(VecRestoreArray(cVec, &coordsDG));
400: PetscCall(DMClone(cdm, &cdgdm));
401: PetscCall(DMSetCellCoordinateDM(dm, cdgdm));
402: PetscCall(DMSetCellCoordinateSection(dm, PETSC_DETERMINE, csDG));
403: PetscCall(DMSetCellCoordinatesLocal(dm, cVec));
404: PetscCall(VecDestroy(&cVec));
405: // Convert the discretization
406: {
407: PetscFE fe, dgfe;
408: PetscSpace P;
409: PetscDualSpace Q, dgQ;
410: PetscQuadrature q, fq;
411: PetscClassId id;
413: PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
414: PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
415: if (id == PETSCFE_CLASSID) {
416: PetscCall(PetscFEGetBasisSpace(fe, &P));
417: PetscCall(PetscObjectReference((PetscObject)P));
418: PetscCall(PetscFEGetDualSpace(fe, &Q));
419: PetscCall(PetscDualSpaceDuplicate(Q, &dgQ));
420: PetscCall(PetscDualSpaceLagrangeSetContinuity(dgQ, PETSC_FALSE));
421: PetscCall(PetscDualSpaceSetUp(dgQ));
422: PetscCall(PetscFEGetQuadrature(fe, &q));
423: PetscCall(PetscObjectReference((PetscObject)q));
424: PetscCall(PetscFEGetFaceQuadrature(fe, &fq));
425: PetscCall(PetscObjectReference((PetscObject)fq));
426: PetscCall(PetscFECreateFromSpaces(P, dgQ, q, fq, &dgfe));
427: PetscCall(DMSetField(cdgdm, 0, NULL, (PetscObject)dgfe));
428: PetscCall(PetscFEDestroy(&dgfe));
429: PetscCall(DMCreateDS(cdgdm));
430: }
431: }
432: PetscCall(DMDestroy(&cdgdm));
434: end:
435: PetscCall(DMRestoreWorkArray(dm, 2 * bs, MPIU_SCALAR, &anchor));
436: PetscCall(DMRestoreWorkArray(dm, 2 * (maxHeight + 1), MPIU_INT, &pStart));
437: PetscCall(PetscSectionDestroy(&csDG));
438: PetscCall(DMDestroy(&plex));
439: PetscCall(DMDestroy(&cplex));
440: PetscFunctionReturn(PETSC_SUCCESS);
441: }