Actual source code: plexrefine.c
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/petscfeimpl.h>
4: #include <petscdmplextransform.h>
5: #include <petscsf.h>
7: /*@
8: DMPlexCreateProcessSF - Create an `PetscSF` which just has process connectivity
10: Collective
12: Input Parameters:
13: + dm - The `DM`
14: - sfPoint - The `PetscSF` which encodes point connectivity
16: Output Parameters:
17: + processRanks - A list of process neighbors, or `NULL`
18: - sfProcess - An `PetscSF` encoding the process connectivity, or `NULL`
20: Level: developer
22: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSFCreate()`, `DMPlexCreateTwoSidedProcessSF()`
23: @*/
24: PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess)
25: {
26: PetscInt numRoots, numLeaves, l;
27: const PetscInt *localPoints;
28: const PetscSFNode *remotePoints;
29: PetscInt *localPointsNew;
30: PetscSFNode *remotePointsNew;
31: PetscMPIInt *ranks;
32: PetscInt *ranksNew;
33: PetscMPIInt size;
35: PetscFunctionBegin;
38: if (processRanks) PetscAssertPointer(processRanks, 3);
39: if (sfProcess) PetscAssertPointer(sfProcess, 4);
40: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
41: PetscCall(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints));
42: PetscCall(PetscMalloc1(numLeaves, &ranks));
43: for (l = 0; l < numLeaves; ++l) ranks[l] = (PetscMPIInt)remotePoints[l].rank;
44: PetscCall(PetscSortRemoveDupsMPIInt(&numLeaves, ranks));
45: PetscCall(PetscMalloc1(numLeaves, &ranksNew));
46: PetscCall(PetscMalloc1(numLeaves, &localPointsNew));
47: PetscCall(PetscMalloc1(numLeaves, &remotePointsNew));
48: for (l = 0; l < numLeaves; ++l) {
49: ranksNew[l] = ranks[l];
50: localPointsNew[l] = l;
51: remotePointsNew[l].index = 0;
52: remotePointsNew[l].rank = ranksNew[l];
53: }
54: PetscCall(PetscFree(ranks));
55: if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks));
56: else PetscCall(PetscFree(ranksNew));
57: if (sfProcess) {
58: PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess));
59: PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Process SF"));
60: PetscCall(PetscSFSetFromOptions(*sfProcess));
61: PetscCall(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
62: }
63: PetscFunctionReturn(PETSC_SUCCESS);
64: }
66: /*@
67: DMPlexCreateCoarsePointIS - Creates an `IS` covering the coarse `DM` chart with the fine points as data
69: Collective
71: Input Parameter:
72: . dm - The coarse `DM`
74: Output Parameter:
75: . fpointIS - The `IS` of all the fine points which exist in the original coarse mesh
77: Level: developer
79: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `IS`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetSubpointIS()`
80: @*/
81: PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS)
82: {
83: DMPlexTransform tr;
84: PetscInt *fpoints;
85: PetscInt pStart, pEnd, p, vStart, vEnd, v;
87: PetscFunctionBegin;
88: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
89: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
90: PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
91: PetscCall(DMPlexTransformSetUp(tr));
92: PetscCall(PetscMalloc1(pEnd - pStart, &fpoints));
93: for (p = 0; p < pEnd - pStart; ++p) fpoints[p] = -1;
94: for (v = vStart; v < vEnd; ++v) {
95: PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */
97: PetscCall(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew));
98: fpoints[v - pStart] = vNew;
99: }
100: PetscCall(DMPlexTransformDestroy(&tr));
101: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, fpoints, PETSC_OWN_POINTER, fpointIS));
102: PetscFunctionReturn(PETSC_SUCCESS);
103: }
105: /*@
106: DMPlexSetTransformType - Set the transform type for uniform refinement
108: Input Parameters:
109: + dm - The `DM`
110: - type - The transform type for uniform refinement
112: Level: developer
114: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexGetTransformType()`, `DMPlexSetRefinementUniform()`
115: @*/
116: PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type)
117: {
118: DM_Plex *mesh = (DM_Plex *)dm->data;
120: PetscFunctionBegin;
122: if (type) PetscAssertPointer(type, 2);
123: PetscCall(PetscFree(mesh->transformType));
124: PetscCall(PetscStrallocpy(type, &mesh->transformType));
125: PetscFunctionReturn(PETSC_SUCCESS);
126: }
128: /*@C
129: DMPlexGetTransformType - Retrieve the transform type for uniform refinement
131: Input Parameter:
132: . dm - The `DM`
134: Output Parameter:
135: . type - The transform type for uniform refinement
137: Level: developer
139: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexSetTransformType()`, `DMPlexGetRefinementUniform()`
140: @*/
141: PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type)
142: {
143: DM_Plex *mesh = (DM_Plex *)dm->data;
145: PetscFunctionBegin;
147: PetscAssertPointer(type, 2);
148: *type = mesh->transformType;
149: PetscFunctionReturn(PETSC_SUCCESS);
150: }
152: #include <petsc/private/dmplextransformimpl.h>
154: PetscErrorCode DMPlexSetTransform(DM dm, DMPlexTransform tr)
155: {
156: DM_Plex *mesh = (DM_Plex *)dm->data;
158: PetscFunctionBegin;
161: PetscCall(PetscObjectReference((PetscObject)tr));
162: PetscCall(DMPlexTransformDestroy(&mesh->transform));
163: // We need to remove the DM because we replace that exact DM with the transformed one in plexcreate.c
164: if (tr) PetscCall(DMPlexTransformSetDM(tr, NULL));
165: mesh->transform = tr;
166: PetscFunctionReturn(PETSC_SUCCESS);
167: }
169: PetscErrorCode DMPlexGetTransform(DM dm, DMPlexTransform *tr)
170: {
171: DM_Plex *mesh = (DM_Plex *)dm->data;
173: PetscFunctionBegin;
175: PetscAssertPointer(tr, 2);
176: *tr = mesh->transform;
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: PetscErrorCode DMPlexSetSaveTransform(DM dm, PetscBool save)
181: {
182: DM_Plex *mesh = (DM_Plex *)dm->data;
184: PetscFunctionBegin;
186: mesh->saveTransform = save;
187: PetscFunctionReturn(PETSC_SUCCESS);
188: }
190: PetscErrorCode DMPlexGetSaveTransform(DM dm, PetscBool *save)
191: {
192: DM_Plex *mesh = (DM_Plex *)dm->data;
194: PetscFunctionBegin;
196: PetscAssertPointer(save, 2);
197: *save = mesh->saveTransform;
198: PetscFunctionReturn(PETSC_SUCCESS);
199: }
201: /*@
202: DMPlexSetRefinementUniform - Set the flag for uniform refinement
204: Input Parameters:
205: + dm - The `DM`
206: - refinementUniform - The flag for uniform refinement
208: Level: developer
210: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
211: @*/
212: PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform)
213: {
214: DM_Plex *mesh = (DM_Plex *)dm->data;
216: PetscFunctionBegin;
218: mesh->refinementUniform = refinementUniform;
219: PetscFunctionReturn(PETSC_SUCCESS);
220: }
222: /*@
223: DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement
225: Input Parameter:
226: . dm - The `DM`
228: Output Parameter:
229: . refinementUniform - The flag for uniform refinement
231: Level: developer
233: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
234: @*/
235: PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform)
236: {
237: DM_Plex *mesh = (DM_Plex *)dm->data;
239: PetscFunctionBegin;
241: PetscAssertPointer(refinementUniform, 2);
242: *refinementUniform = mesh->refinementUniform;
243: PetscFunctionReturn(PETSC_SUCCESS);
244: }
246: /*@
247: DMPlexSetRefinementLimit - Set the maximum cell volume for refinement
249: Input Parameters:
250: + dm - The `DM`
251: - refinementLimit - The maximum cell volume in the refined mesh
253: Level: developer
255: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
256: @*/
257: PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit)
258: {
259: DM_Plex *mesh = (DM_Plex *)dm->data;
261: PetscFunctionBegin;
263: mesh->refinementLimit = refinementLimit;
264: PetscFunctionReturn(PETSC_SUCCESS);
265: }
267: /*@
268: DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement
270: Input Parameter:
271: . dm - The `DM`
273: Output Parameter:
274: . refinementLimit - The maximum cell volume in the refined mesh
276: Level: developer
278: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
279: @*/
280: PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit)
281: {
282: DM_Plex *mesh = (DM_Plex *)dm->data;
284: PetscFunctionBegin;
286: PetscAssertPointer(refinementLimit, 2);
287: /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */
288: *refinementLimit = mesh->refinementLimit;
289: PetscFunctionReturn(PETSC_SUCCESS);
290: }
292: /*@C
293: DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement
295: Input Parameters:
296: + dm - The `DM`
297: - refinementFunc - Function giving the maximum cell volume in the refined mesh
299: Calling Sequence of `refinementFunc`:
300: + coords - Coordinates of the current point, usually a cell centroid
301: - limit - The maximum cell volume for a cell containing this point
303: Level: developer
305: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
306: @*/
307: PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal coords[], PetscReal *limit))
308: {
309: DM_Plex *mesh = (DM_Plex *)dm->data;
311: PetscFunctionBegin;
313: mesh->refinementFunc = refinementFunc;
314: PetscFunctionReturn(PETSC_SUCCESS);
315: }
317: /*@C
318: DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement
320: Input Parameter:
321: . dm - The `DM`
323: Output Parameter:
324: . refinementFunc - Function giving the maximum cell volume in the refined mesh
326: Calling Sequence of `refinementFunc`:
327: + coords - Coordinates of the current point, usually a cell centroid
328: - limit - The maximum cell volume for a cell containing this point
330: Level: developer
332: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
333: @*/
334: PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal coords[], PetscReal *limit))
335: {
336: DM_Plex *mesh = (DM_Plex *)dm->data;
338: PetscFunctionBegin;
340: PetscAssertPointer(refinementFunc, 2);
341: *refinementFunc = mesh->refinementFunc;
342: PetscFunctionReturn(PETSC_SUCCESS);
343: }
345: PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm)
346: {
347: PetscBool isUniform;
349: PetscFunctionBegin;
350: PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
351: PetscCall(DMViewFromOptions(dm, NULL, "-initref_dm_view"));
352: if (isUniform) {
353: DMPlexTransform tr;
354: DM cdm, rcdm;
355: DMPlexTransformType trType;
356: const char *prefix;
357: PetscOptions options;
358: PetscInt cDegree;
359: PetscBool useCeed, save;
361: PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
362: PetscCall(DMPlexTransformSetDM(tr, dm));
363: PetscCall(DMPlexGetTransformType(dm, &trType));
364: if (trType) PetscCall(DMPlexTransformSetType(tr, trType));
365: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
366: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
367: PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
368: PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
369: PetscCall(DMPlexTransformSetFromOptions(tr));
370: PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
371: PetscCall(DMPlexTransformSetUp(tr));
372: PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
373: PetscCall(DMPlexTransformApply(tr, dm, rdm));
374: PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE));
375: PetscCall(DMPlexGetUseCeed(dm, &useCeed));
376: PetscCall(DMPlexSetUseCeed(*rdm, useCeed));
377: PetscCall(DMSetMatType(*rdm, dm->mattype));
378: PetscCall(DMCopyDisc(dm, *rdm));
379: PetscCall(DMGetCoordinateDM(dm, &cdm));
380: PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
381: PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree));
382: {
383: PetscDS cds, rcds;
385: PetscCall(DMPlexCreateCoordinateSpace(*rdm, cDegree, PETSC_FALSE, PETSC_TRUE));
386: PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
387: PetscCall(DMGetDS(cdm, &cds));
388: PetscCall(DMGetDS(rcdm, &rcds));
389: PetscCall(PetscDSCopyConstants(cds, rcds));
390: }
391: PetscCall(DMPlexGetUseCeed(cdm, &useCeed));
392: PetscCall(DMPlexSetUseCeed(rcdm, useCeed));
393: if (useCeed) {
394: PetscCall(DMPlexSetUseMatClosurePermutation(rcdm, PETSC_FALSE));
395: PetscCall(DMUseTensorOrder(rcdm, PETSC_TRUE));
396: }
397: PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm));
398: PetscCall(DMPlexGetSaveTransform(dm, &save));
399: if (save) PetscCall(DMPlexSetTransform(*rdm, tr));
400: PetscCall(DMPlexTransformDestroy(&tr));
401: } else {
402: PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm));
403: }
404: if (*rdm) {
405: ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
406: ((DM_Plex *)(*rdm)->data)->printL2 = ((DM_Plex *)dm->data)->printL2;
407: }
408: PetscCall(DMViewFromOptions(*rdm, NULL, "-postref_dm_view"));
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[])
413: {
414: DM cdm = dm;
415: PetscInt r;
416: PetscBool isUniform, localized, useCeed;
418: PetscFunctionBegin;
419: PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
420: PetscCall(DMGetCoordinatesLocalized(dm, &localized));
421: if (isUniform) {
422: for (r = 0; r < nlevels; ++r) {
423: DMPlexTransform tr;
424: DM codm, rcodm;
425: const char *prefix;
427: PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr));
428: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix));
429: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
430: PetscCall(DMPlexTransformSetDM(tr, cdm));
431: PetscCall(DMPlexTransformSetFromOptions(tr));
432: PetscCall(DMPlexTransformSetUp(tr));
433: PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r]));
434: PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown));
435: PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1));
436: PetscCall(DMSetMatType(rdm[r], dm->mattype));
437: PetscCall(DMPlexGetUseCeed(dm, &useCeed));
438: PetscCall(DMPlexSetUseCeed(rdm[r], useCeed));
439: PetscCall(DMCopyDisc(cdm, rdm[r]));
440: PetscCall(DMGetCoordinateDM(dm, &codm));
441: PetscCall(DMGetCoordinateDM(rdm[r], &rcodm));
442: PetscCall(DMCopyDisc(codm, rcodm));
443: PetscCall(DMPlexGetUseCeed(codm, &useCeed));
444: PetscCall(DMPlexSetUseCeed(rcodm, useCeed));
445: if (useCeed) {
446: PetscCall(DMPlexSetUseMatClosurePermutation(rcodm, PETSC_FALSE));
447: PetscCall(DMUseTensorOrder(rcodm, PETSC_TRUE));
448: }
449: PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r]));
450: PetscCall(DMSetCoarseDM(rdm[r], cdm));
451: PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE));
452: if (rdm[r]) {
453: ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
454: ((DM_Plex *)rdm[r]->data)->printL2 = ((DM_Plex *)dm->data)->printL2;
455: }
456: cdm = rdm[r];
457: PetscCall(DMPlexTransformDestroy(&tr));
458: }
459: } else {
460: for (r = 0; r < nlevels; ++r) {
461: PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r]));
462: PetscCall(DMPlexGetUseCeed(dm, &useCeed));
463: PetscCall(DMPlexSetUseCeed(rdm[r], useCeed));
464: PetscCall(DMCopyDisc(cdm, rdm[r]));
465: if (localized) PetscCall(DMLocalizeCoordinates(rdm[r]));
466: PetscCall(DMSetCoarseDM(rdm[r], cdm));
467: if (rdm[r]) {
468: ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
469: ((DM_Plex *)rdm[r]->data)->printL2 = ((DM_Plex *)dm->data)->printL2;
470: }
471: cdm = rdm[r];
472: }
473: }
474: PetscFunctionReturn(PETSC_SUCCESS);
475: }