xref: /petsc/src/dm/impls/plex/plexrefine.c (revision 9140fee14acbea959c6ed74f4836a1a89f708038)
1 #include <petsc/private/dmpleximpl.h>  /*I      "petscdmplex.h"   I*/
2 #include <petsc/private/petscfeimpl.h> /* For PetscFEInterpolate_Static() */
3 
4 #include <petscdmplextransform.h>
5 #include <petscsf.h>
6 
7 /*@
8   DMPlexCreateProcessSF - Create an `PetscSF` which just has process connectivity
9 
10   Collective
11 
12   Input Parameters:
13 + dm      - The `DM`
14 - sfPoint - The `PetscSF` which encodes point connectivity
15 
16   Output Parameters:
17 + processRanks - A list of process neighbors, or `NULL`
18 - sfProcess    - An `PetscSF` encoding the process connectivity, or `NULL`
19 
20   Level: developer
21 
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   PetscInt          *ranks, *ranksNew;
32   PetscMPIInt        size;
33 
34   PetscFunctionBegin;
35   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
37   if (processRanks) PetscAssertPointer(processRanks, 3);
38   if (sfProcess) PetscAssertPointer(sfProcess, 4);
39   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
40   PetscCall(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints));
41   PetscCall(PetscMalloc1(numLeaves, &ranks));
42   for (l = 0; l < numLeaves; ++l) ranks[l] = remotePoints[l].rank;
43   PetscCall(PetscSortRemoveDupsInt(&numLeaves, ranks));
44   PetscCall(PetscMalloc1(numLeaves, &ranksNew));
45   PetscCall(PetscMalloc1(numLeaves, &localPointsNew));
46   PetscCall(PetscMalloc1(numLeaves, &remotePointsNew));
47   for (l = 0; l < numLeaves; ++l) {
48     ranksNew[l]              = ranks[l];
49     localPointsNew[l]        = l;
50     remotePointsNew[l].index = 0;
51     remotePointsNew[l].rank  = ranksNew[l];
52   }
53   PetscCall(PetscFree(ranks));
54   if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks));
55   else PetscCall(PetscFree(ranksNew));
56   if (sfProcess) {
57     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess));
58     PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Process SF"));
59     PetscCall(PetscSFSetFromOptions(*sfProcess));
60     PetscCall(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
61   }
62   PetscFunctionReturn(PETSC_SUCCESS);
63 }
64 
65 /*@
66   DMPlexCreateCoarsePointIS - Creates an `IS` covering the coarse `DM` chart with the fine points as data
67 
68   Collective
69 
70   Input Parameter:
71 . dm - The coarse `DM`
72 
73   Output Parameter:
74 . fpointIS - The `IS` of all the fine points which exist in the original coarse mesh
75 
76   Level: developer
77 
78 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `IS`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetSubpointIS()`
79 @*/
80 PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS)
81 {
82   DMPlexTransform tr;
83   PetscInt       *fpoints;
84   PetscInt        pStart, pEnd, p, vStart, vEnd, v;
85 
86   PetscFunctionBegin;
87   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
88   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
89   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
90   PetscCall(DMPlexTransformSetUp(tr));
91   PetscCall(PetscMalloc1(pEnd - pStart, &fpoints));
92   for (p = 0; p < pEnd - pStart; ++p) fpoints[p] = -1;
93   for (v = vStart; v < vEnd; ++v) {
94     PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */
95 
96     PetscCall(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew));
97     fpoints[v - pStart] = vNew;
98   }
99   PetscCall(DMPlexTransformDestroy(&tr));
100   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, fpoints, PETSC_OWN_POINTER, fpointIS));
101   PetscFunctionReturn(PETSC_SUCCESS);
102 }
103 
104 /*@C
105   DMPlexSetTransformType - Set the transform type for uniform refinement
106 
107   Input Parameters:
108 + dm   - The `DM`
109 - type - The transform type for uniform refinement
110 
111   Level: developer
112 
113 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexGetTransformType()`, `DMPlexSetRefinementUniform()`
114 @*/
115 PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type)
116 {
117   DM_Plex *mesh = (DM_Plex *)dm->data;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
121   if (type) PetscAssertPointer(type, 2);
122   PetscCall(PetscFree(mesh->transformType));
123   PetscCall(PetscStrallocpy(type, &mesh->transformType));
124   PetscFunctionReturn(PETSC_SUCCESS);
125 }
126 
127 /*@C
128   DMPlexGetTransformType - Retrieve the transform type for uniform refinement
129 
130   Input Parameter:
131 . dm - The `DM`
132 
133   Output Parameter:
134 . type - The transform type for uniform refinement
135 
136   Level: developer
137 
138 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexSetTransformType()`, `DMPlexGetRefinementUniform()`
139 @*/
140 PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type)
141 {
142   DM_Plex *mesh = (DM_Plex *)dm->data;
143 
144   PetscFunctionBegin;
145   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
146   PetscAssertPointer(type, 2);
147   *type = mesh->transformType;
148   PetscFunctionReturn(PETSC_SUCCESS);
149 }
150 
151 /*@
152   DMPlexSetRefinementUniform - Set the flag for uniform refinement
153 
154   Input Parameters:
155 + dm                - The `DM`
156 - refinementUniform - The flag for uniform refinement
157 
158   Level: developer
159 
160 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
161 @*/
162 PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform)
163 {
164   DM_Plex *mesh = (DM_Plex *)dm->data;
165 
166   PetscFunctionBegin;
167   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
168   mesh->refinementUniform = refinementUniform;
169   PetscFunctionReturn(PETSC_SUCCESS);
170 }
171 
172 /*@
173   DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement
174 
175   Input Parameter:
176 . dm - The `DM`
177 
178   Output Parameter:
179 . refinementUniform - The flag for uniform refinement
180 
181   Level: developer
182 
183 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
184 @*/
185 PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform)
186 {
187   DM_Plex *mesh = (DM_Plex *)dm->data;
188 
189   PetscFunctionBegin;
190   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
191   PetscAssertPointer(refinementUniform, 2);
192   *refinementUniform = mesh->refinementUniform;
193   PetscFunctionReturn(PETSC_SUCCESS);
194 }
195 
196 /*@
197   DMPlexSetRefinementLimit - Set the maximum cell volume for refinement
198 
199   Input Parameters:
200 + dm              - The `DM`
201 - refinementLimit - The maximum cell volume in the refined mesh
202 
203   Level: developer
204 
205 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
206 @*/
207 PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit)
208 {
209   DM_Plex *mesh = (DM_Plex *)dm->data;
210 
211   PetscFunctionBegin;
212   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
213   mesh->refinementLimit = refinementLimit;
214   PetscFunctionReturn(PETSC_SUCCESS);
215 }
216 
217 /*@
218   DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement
219 
220   Input Parameter:
221 . dm - The `DM`
222 
223   Output Parameter:
224 . refinementLimit - The maximum cell volume in the refined mesh
225 
226   Level: developer
227 
228 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
229 @*/
230 PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit)
231 {
232   DM_Plex *mesh = (DM_Plex *)dm->data;
233 
234   PetscFunctionBegin;
235   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
236   PetscAssertPointer(refinementLimit, 2);
237   /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */
238   *refinementLimit = mesh->refinementLimit;
239   PetscFunctionReturn(PETSC_SUCCESS);
240 }
241 
242 /*@C
243   DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement
244 
245   Input Parameters:
246 + dm             - The `DM`
247 - refinementFunc - Function giving the maximum cell volume in the refined mesh
248 
249   Calling Sequence of `refinementFunc`:
250 + coords - Coordinates of the current point, usually a cell centroid
251 - limit  - The maximum cell volume for a cell containing this point
252 
253   Level: developer
254 
255 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
256 @*/
257 PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal coords[], PetscReal *limit))
258 {
259   DM_Plex *mesh = (DM_Plex *)dm->data;
260 
261   PetscFunctionBegin;
262   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
263   mesh->refinementFunc = refinementFunc;
264   PetscFunctionReturn(PETSC_SUCCESS);
265 }
266 
267 /*@C
268   DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement
269 
270   Input Parameter:
271 . dm - The `DM`
272 
273   Output Parameter:
274 . refinementFunc - Function giving the maximum cell volume in the refined mesh
275 
276   Calling Sequence of `refinementFunc`:
277 + coords - Coordinates of the current point, usually a cell centroid
278 - limit  - The maximum cell volume for a cell containing this point
279 
280   Level: developer
281 
282 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
283 @*/
284 PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal coords[], PetscReal *limit))
285 {
286   DM_Plex *mesh = (DM_Plex *)dm->data;
287 
288   PetscFunctionBegin;
289   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
290   PetscAssertPointer(refinementFunc, 2);
291   *refinementFunc = mesh->refinementFunc;
292   PetscFunctionReturn(PETSC_SUCCESS);
293 }
294 
295 PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm)
296 {
297   PetscBool isUniform;
298 
299   PetscFunctionBegin;
300   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
301   PetscCall(DMViewFromOptions(dm, NULL, "-initref_dm_view"));
302   if (isUniform) {
303     DMPlexTransform     tr;
304     DM                  cdm, rcdm;
305     DMPlexTransformType trType;
306     const char         *prefix;
307     PetscOptions        options;
308     PetscInt            cDegree;
309     PetscBool           useCeed;
310 
311     PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
312     PetscCall(DMPlexTransformSetDM(tr, dm));
313     PetscCall(DMPlexGetTransformType(dm, &trType));
314     if (trType) PetscCall(DMPlexTransformSetType(tr, trType));
315     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
316     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
317     PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
318     PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
319     PetscCall(DMPlexTransformSetFromOptions(tr));
320     PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
321     PetscCall(DMPlexTransformSetUp(tr));
322     PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
323     PetscCall(DMPlexTransformApply(tr, dm, rdm));
324     PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE));
325     PetscCall(DMPlexGetUseCeed(dm, &useCeed));
326     PetscCall(DMPlexSetUseCeed(*rdm, useCeed));
327     PetscCall(DMSetMatType((*rdm), dm->mattype));
328     PetscCall(DMCopyDisc(dm, *rdm));
329     PetscCall(DMGetCoordinateDM(dm, &cdm));
330     PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
331     PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree));
332     if (cDegree <= 1) {
333       PetscCall(DMCopyDisc(cdm, rcdm));
334     } else {
335       PetscCall(DMPlexCreateCoordinateSpace(*rdm, 1, PETSC_TRUE, NULL));
336       PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
337     }
338     PetscCall(DMPlexGetUseCeed(cdm, &useCeed));
339     PetscCall(DMPlexSetUseCeed(rcdm, useCeed));
340     if (useCeed) {
341       PetscCall(DMPlexSetUseMatClosurePermutation(rcdm, PETSC_FALSE));
342       PetscCall(DMUseTensorOrder(rcdm, PETSC_TRUE));
343     }
344     PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm));
345     PetscCall(DMPlexTransformDestroy(&tr));
346   } else {
347     PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm));
348   }
349   if (*rdm) {
350     ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
351     ((DM_Plex *)(*rdm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
352   }
353   PetscCall(DMViewFromOptions(dm, NULL, "-postref_dm_view"));
354   PetscFunctionReturn(PETSC_SUCCESS);
355 }
356 
357 PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[])
358 {
359   DM        cdm = dm;
360   PetscInt  r;
361   PetscBool isUniform, localized, useCeed;
362 
363   PetscFunctionBegin;
364   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
365   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
366   if (isUniform) {
367     for (r = 0; r < nlevels; ++r) {
368       DMPlexTransform tr;
369       DM              codm, rcodm;
370       const char     *prefix;
371 
372       PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr));
373       PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix));
374       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
375       PetscCall(DMPlexTransformSetDM(tr, cdm));
376       PetscCall(DMPlexTransformSetFromOptions(tr));
377       PetscCall(DMPlexTransformSetUp(tr));
378       PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r]));
379       PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown));
380       PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1));
381       PetscCall(DMSetMatType(rdm[r], dm->mattype));
382       PetscCall(DMPlexGetUseCeed(dm, &useCeed));
383       PetscCall(DMPlexSetUseCeed(rdm[r], useCeed));
384       PetscCall(DMCopyDisc(cdm, rdm[r]));
385       PetscCall(DMGetCoordinateDM(dm, &codm));
386       PetscCall(DMGetCoordinateDM(rdm[r], &rcodm));
387       PetscCall(DMCopyDisc(codm, rcodm));
388       PetscCall(DMPlexGetUseCeed(codm, &useCeed));
389       PetscCall(DMPlexSetUseCeed(rcodm, useCeed));
390       if (useCeed) {
391         PetscCall(DMPlexSetUseMatClosurePermutation(rcodm, PETSC_FALSE));
392         PetscCall(DMUseTensorOrder(rcodm, PETSC_TRUE));
393       }
394       PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r]));
395       PetscCall(DMSetCoarseDM(rdm[r], cdm));
396       PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE));
397       if (rdm[r]) {
398         ((DM_Plex *)(rdm[r])->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
399         ((DM_Plex *)(rdm[r])->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
400       }
401       cdm = rdm[r];
402       PetscCall(DMPlexTransformDestroy(&tr));
403     }
404   } else {
405     for (r = 0; r < nlevels; ++r) {
406       PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r]));
407       PetscCall(DMPlexGetUseCeed(dm, &useCeed));
408       PetscCall(DMPlexSetUseCeed(rdm[r], useCeed));
409       PetscCall(DMCopyDisc(cdm, rdm[r]));
410       if (localized) PetscCall(DMLocalizeCoordinates(rdm[r]));
411       PetscCall(DMSetCoarseDM(rdm[r], cdm));
412       if (rdm[r]) {
413         ((DM_Plex *)(rdm[r])->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
414         ((DM_Plex *)(rdm[r])->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
415       }
416       cdm = rdm[r];
417     }
418   }
419   PetscFunctionReturn(PETSC_SUCCESS);
420 }
421