xref: /petsc/src/dm/impls/plex/plexrefine.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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 SF which just has process connectivity
9 
10   Collective on dm
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 SF encoding the process connectivity, or NULL
19 
20   Level: developer
21 
22 .seealso: `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) PetscValidPointer(processRanks, 3);
38   if (sfProcess) PetscValidPointer(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(0);
63 }
64 
65 /*@
66   DMPlexCreateCoarsePointIS - Creates an IS covering the coarse DM chart with the fine points as data
67 
68   Collective on dm
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: `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(0);
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: `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) PetscValidCharPointer(type, 2);
122   PetscCall(PetscFree(mesh->transformType));
123   PetscCall(PetscStrallocpy(type, &mesh->transformType));
124   PetscFunctionReturn(0);
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: `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   PetscValidPointer(type, 2);
147   *type = mesh->transformType;
148   PetscFunctionReturn(0);
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: `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(0);
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: `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   PetscValidBoolPointer(refinementUniform, 2);
192   *refinementUniform = mesh->refinementUniform;
193   PetscFunctionReturn(0);
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: `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(0);
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: `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   PetscValidRealPointer(refinementLimit, 2);
237   /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */
238   *refinementLimit = mesh->refinementLimit;
239   PetscFunctionReturn(0);
240 }
241 
242 /*@
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   Note: The calling sequence is refinementFunc(coords, limit)
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: `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
256 @*/
257 PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal[], PetscReal *))
258 {
259   DM_Plex *mesh = (DM_Plex *)dm->data;
260 
261   PetscFunctionBegin;
262   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
263   mesh->refinementFunc = refinementFunc;
264   PetscFunctionReturn(0);
265 }
266 
267 /*@
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   Note: The calling sequence is refinementFunc(coords, limit)
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: `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
283 @*/
284 PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal[], PetscReal *))
285 {
286   DM_Plex *mesh = (DM_Plex *)dm->data;
287 
288   PetscFunctionBegin;
289   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
290   PetscValidPointer(refinementFunc, 2);
291   *refinementFunc = mesh->refinementFunc;
292   PetscFunctionReturn(0);
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 
309     PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
310     PetscCall(DMPlexTransformSetDM(tr, dm));
311     PetscCall(DMPlexGetTransformType(dm, &trType));
312     if (trType) PetscCall(DMPlexTransformSetType(tr, trType));
313     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
314     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
315     PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
316     PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
317     PetscCall(DMPlexTransformSetFromOptions(tr));
318     PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
319     PetscCall(DMPlexTransformSetUp(tr));
320     PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
321     PetscCall(DMPlexTransformApply(tr, dm, rdm));
322     PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE));
323     PetscCall(DMCopyDisc(dm, *rdm));
324     PetscCall(DMGetCoordinateDM(dm, &cdm));
325     PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
326     PetscCall(DMCopyDisc(cdm, rcdm));
327     PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm));
328     PetscCall(DMPlexTransformDestroy(&tr));
329   } else {
330     PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm));
331   }
332   if (*rdm) {
333     ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
334     ((DM_Plex *)(*rdm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
335   }
336   PetscCall(DMViewFromOptions(dm, NULL, "-postref_dm_view"));
337   PetscFunctionReturn(0);
338 }
339 
340 PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[])
341 {
342   DM        cdm = dm;
343   PetscInt  r;
344   PetscBool isUniform, localized;
345 
346   PetscFunctionBegin;
347   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
348   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
349   if (isUniform) {
350     for (r = 0; r < nlevels; ++r) {
351       DMPlexTransform tr;
352       DM              codm, rcodm;
353       const char     *prefix;
354 
355       PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr));
356       PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix));
357       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
358       PetscCall(DMPlexTransformSetDM(tr, cdm));
359       PetscCall(DMPlexTransformSetFromOptions(tr));
360       PetscCall(DMPlexTransformSetUp(tr));
361       PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r]));
362       PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown));
363       PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1));
364       PetscCall(DMCopyDisc(cdm, rdm[r]));
365       PetscCall(DMGetCoordinateDM(dm, &codm));
366       PetscCall(DMGetCoordinateDM(rdm[r], &rcodm));
367       PetscCall(DMCopyDisc(codm, rcodm));
368       PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r]));
369       PetscCall(DMSetCoarseDM(rdm[r], cdm));
370       PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE));
371       if (rdm[r]) {
372         ((DM_Plex *)(rdm[r])->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
373         ((DM_Plex *)(rdm[r])->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
374       }
375       cdm = rdm[r];
376       PetscCall(DMPlexTransformDestroy(&tr));
377     }
378   } else {
379     for (r = 0; r < nlevels; ++r) {
380       PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r]));
381       PetscCall(DMCopyDisc(cdm, rdm[r]));
382       if (localized) PetscCall(DMLocalizeCoordinates(rdm[r]));
383       PetscCall(DMSetCoarseDM(rdm[r], cdm));
384       if (rdm[r]) {
385         ((DM_Plex *)(rdm[r])->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
386         ((DM_Plex *)(rdm[r])->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
387       }
388       cdm = rdm[r];
389     }
390   }
391   PetscFunctionReturn(0);
392 }
393