1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2 #include <petscblaslapack.h>
3
4 PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection;
5
DMPlexMetricSetFromOptions(DM dm)6 PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
7 {
8 DM_Plex *plex = (DM_Plex *)dm->data;
9 MPI_Comm comm;
10 PetscBool isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
11 PetscBool noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE;
12 PetscInt verbosity = -1, numIter = 3;
13 PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0, beta = 1.3, hausd = 0.01;
14
15 PetscFunctionBegin;
16 if (!plex->metricCtx) PetscCall(PetscNew(&plex->metricCtx));
17 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
18 PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");
19 PetscCall(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL));
20 PetscCall(DMPlexMetricSetIsotropic(dm, isotropic));
21 PetscCall(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL));
22 PetscCall(DMPlexMetricSetUniform(dm, uniform));
23 PetscCall(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL));
24 PetscCall(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst));
25 PetscCall(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL));
26 PetscCall(DMPlexMetricSetNoInsertion(dm, noInsert));
27 PetscCall(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL));
28 PetscCall(DMPlexMetricSetNoSwapping(dm, noSwap));
29 PetscCall(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL));
30 PetscCall(DMPlexMetricSetNoMovement(dm, noMove));
31 PetscCall(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL));
32 PetscCall(DMPlexMetricSetNoSurf(dm, noSurf));
33 PetscCall(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0));
34 PetscCall(DMPlexMetricSetNumIterations(dm, numIter));
35 PetscCall(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10));
36 PetscCall(DMPlexMetricSetVerbosity(dm, verbosity));
37 PetscCall(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL));
38 PetscCall(DMPlexMetricSetMinimumMagnitude(dm, h_min));
39 PetscCall(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL));
40 PetscCall(DMPlexMetricSetMaximumMagnitude(dm, h_max));
41 PetscCall(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL));
42 PetscCall(DMPlexMetricSetMaximumAnisotropy(dm, a_max));
43 PetscCall(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL));
44 PetscCall(DMPlexMetricSetNormalizationOrder(dm, p));
45 PetscCall(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL));
46 PetscCall(DMPlexMetricSetTargetComplexity(dm, target));
47 PetscCall(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL));
48 PetscCall(DMPlexMetricSetGradationFactor(dm, beta));
49 PetscCall(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL));
50 PetscCall(DMPlexMetricSetHausdorffNumber(dm, hausd));
51 PetscOptionsEnd();
52 PetscFunctionReturn(PETSC_SUCCESS);
53 }
54
55 /*@
56 DMPlexMetricSetIsotropic - Record whether a metric is isotropic
57
58 Input Parameters:
59 + dm - The `DM`
60 - isotropic - Is the metric isotropic?
61
62 Level: beginner
63
64 .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetUniform()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
65 @*/
DMPlexMetricSetIsotropic(DM dm,PetscBool isotropic)66 PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
67 {
68 DM_Plex *plex = (DM_Plex *)dm->data;
69
70 PetscFunctionBegin;
71 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
72 plex->metricCtx->isotropic = isotropic;
73 PetscFunctionReturn(PETSC_SUCCESS);
74 }
75
76 /*@
77 DMPlexMetricIsIsotropic - Is a metric isotropic?
78
79 Input Parameters:
80 . dm - The `DM`
81
82 Output Parameters:
83 . isotropic - Is the metric isotropic?
84
85 Level: beginner
86
87 .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricIsUniform()`, `DMPlexMetricRestrictAnisotropyFirst()`
88 @*/
DMPlexMetricIsIsotropic(DM dm,PetscBool * isotropic)89 PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
90 {
91 DM_Plex *plex = (DM_Plex *)dm->data;
92
93 PetscFunctionBegin;
94 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
95 *isotropic = plex->metricCtx->isotropic;
96 PetscFunctionReturn(PETSC_SUCCESS);
97 }
98
99 /*@
100 DMPlexMetricSetUniform - Record whether a metric is uniform
101
102 Input Parameters:
103 + dm - The `DM`
104 - uniform - Is the metric uniform?
105
106 Level: beginner
107
108 Note:
109 If the metric is specified as uniform then it is assumed to be isotropic, too.
110
111 .seealso: `DMPLEX`, `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
112 @*/
DMPlexMetricSetUniform(DM dm,PetscBool uniform)113 PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
114 {
115 DM_Plex *plex = (DM_Plex *)dm->data;
116
117 PetscFunctionBegin;
118 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
119 plex->metricCtx->uniform = uniform;
120 if (uniform) plex->metricCtx->isotropic = uniform;
121 PetscFunctionReturn(PETSC_SUCCESS);
122 }
123
124 /*@
125 DMPlexMetricIsUniform - Is a metric uniform?
126
127 Input Parameters:
128 . dm - The `DM`
129
130 Output Parameters:
131 . uniform - Is the metric uniform?
132
133 Level: beginner
134
135 .seealso: `DMPLEX`, `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
136 @*/
DMPlexMetricIsUniform(DM dm,PetscBool * uniform)137 PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
138 {
139 DM_Plex *plex = (DM_Plex *)dm->data;
140
141 PetscFunctionBegin;
142 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
143 *uniform = plex->metricCtx->uniform;
144 PetscFunctionReturn(PETSC_SUCCESS);
145 }
146
147 /*@
148 DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
149
150 Input Parameters:
151 + dm - The `DM`
152 - restrictAnisotropyFirst - Should anisotropy be normalized first?
153
154 Level: beginner
155
156 .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
157 @*/
DMPlexMetricSetRestrictAnisotropyFirst(DM dm,PetscBool restrictAnisotropyFirst)158 PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
159 {
160 DM_Plex *plex = (DM_Plex *)dm->data;
161
162 PetscFunctionBegin;
163 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
164 plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
165 PetscFunctionReturn(PETSC_SUCCESS);
166 }
167
168 /*@
169 DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
170
171 Input Parameters:
172 . dm - The `DM`
173
174 Output Parameters:
175 . restrictAnisotropyFirst - Is anisotropy be normalized first?
176
177 Level: beginner
178
179 .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
180 @*/
DMPlexMetricRestrictAnisotropyFirst(DM dm,PetscBool * restrictAnisotropyFirst)181 PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
182 {
183 DM_Plex *plex = (DM_Plex *)dm->data;
184
185 PetscFunctionBegin;
186 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
187 *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
188 PetscFunctionReturn(PETSC_SUCCESS);
189 }
190
191 /*@
192 DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
193
194 Input Parameters:
195 + dm - The `DM`
196 - noInsert - Should node insertion and deletion be turned off?
197
198 Level: beginner
199
200 Note:
201 This is only used by Mmg and ParMmg (not Pragmatic).
202
203 .seealso: `DMPLEX`, `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
204 @*/
DMPlexMetricSetNoInsertion(DM dm,PetscBool noInsert)205 PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
206 {
207 DM_Plex *plex = (DM_Plex *)dm->data;
208
209 PetscFunctionBegin;
210 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
211 plex->metricCtx->noInsert = noInsert;
212 PetscFunctionReturn(PETSC_SUCCESS);
213 }
214
215 /*@
216 DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
217
218 Input Parameters:
219 . dm - The `DM`
220
221 Output Parameters:
222 . noInsert - Are node insertion and deletion turned off?
223
224 Level: beginner
225
226 Note:
227 This is only used by Mmg and ParMmg (not Pragmatic).
228
229 .seealso: `DMPLEX`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
230 @*/
DMPlexMetricNoInsertion(DM dm,PetscBool * noInsert)231 PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
232 {
233 DM_Plex *plex = (DM_Plex *)dm->data;
234
235 PetscFunctionBegin;
236 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
237 *noInsert = plex->metricCtx->noInsert;
238 PetscFunctionReturn(PETSC_SUCCESS);
239 }
240
241 /*@
242 DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
243
244 Input Parameters:
245 + dm - The `DM`
246 - noSwap - Should facet swapping be turned off?
247
248 Level: beginner
249
250 Note:
251 This is only used by Mmg and ParMmg (not Pragmatic).
252
253 .seealso: `DMPLEX`, `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
254 @*/
DMPlexMetricSetNoSwapping(DM dm,PetscBool noSwap)255 PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
256 {
257 DM_Plex *plex = (DM_Plex *)dm->data;
258
259 PetscFunctionBegin;
260 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
261 plex->metricCtx->noSwap = noSwap;
262 PetscFunctionReturn(PETSC_SUCCESS);
263 }
264
265 /*@
266 DMPlexMetricNoSwapping - Is facet swapping turned off?
267
268 Input Parameters:
269 . dm - The `DM`
270
271 Output Parameters:
272 . noSwap - Is facet swapping turned off?
273
274 Level: beginner
275
276 Note:
277 This is only used by Mmg and ParMmg (not Pragmatic).
278
279 .seealso: `DMPLEX`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
280 @*/
DMPlexMetricNoSwapping(DM dm,PetscBool * noSwap)281 PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
282 {
283 DM_Plex *plex = (DM_Plex *)dm->data;
284
285 PetscFunctionBegin;
286 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
287 *noSwap = plex->metricCtx->noSwap;
288 PetscFunctionReturn(PETSC_SUCCESS);
289 }
290
291 /*@
292 DMPlexMetricSetNoMovement - Should node movement be turned off?
293
294 Input Parameters:
295 + dm - The `DM`
296 - noMove - Should node movement be turned off?
297
298 Level: beginner
299
300 Note:
301 This is only used by Mmg and ParMmg (not Pragmatic).
302
303 .seealso: `DMPLEX`, `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()`
304 @*/
DMPlexMetricSetNoMovement(DM dm,PetscBool noMove)305 PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
306 {
307 DM_Plex *plex = (DM_Plex *)dm->data;
308
309 PetscFunctionBegin;
310 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
311 plex->metricCtx->noMove = noMove;
312 PetscFunctionReturn(PETSC_SUCCESS);
313 }
314
315 /*@
316 DMPlexMetricNoMovement - Is node movement turned off?
317
318 Input Parameters:
319 . dm - The `DM`
320
321 Output Parameters:
322 . noMove - Is node movement turned off?
323
324 Level: beginner
325
326 Note:
327 This is only used by Mmg and ParMmg (not Pragmatic).
328
329 .seealso: `DMPLEX`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()`
330 @*/
DMPlexMetricNoMovement(DM dm,PetscBool * noMove)331 PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
332 {
333 DM_Plex *plex = (DM_Plex *)dm->data;
334
335 PetscFunctionBegin;
336 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
337 *noMove = plex->metricCtx->noMove;
338 PetscFunctionReturn(PETSC_SUCCESS);
339 }
340
341 /*@
342 DMPlexMetricSetNoSurf - Should surface modification be turned off?
343
344 Input Parameters:
345 + dm - The `DM`
346 - noSurf - Should surface modification be turned off?
347
348 Level: beginner
349
350 Note:
351 This is only used by Mmg and ParMmg (not Pragmatic).
352
353 .seealso: `DMPLEX`, `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`
354 @*/
DMPlexMetricSetNoSurf(DM dm,PetscBool noSurf)355 PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
356 {
357 DM_Plex *plex = (DM_Plex *)dm->data;
358
359 PetscFunctionBegin;
360 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
361 plex->metricCtx->noSurf = noSurf;
362 PetscFunctionReturn(PETSC_SUCCESS);
363 }
364
365 /*@
366 DMPlexMetricNoSurf - Is surface modification turned off?
367
368 Input Parameters:
369 . dm - The `DM`
370
371 Output Parameters:
372 . noSurf - Is surface modification turned off?
373
374 Level: beginner
375
376 Note:
377 This is only used by Mmg and ParMmg (not Pragmatic).
378
379 .seealso: `DMPLEX`, `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`
380 @*/
DMPlexMetricNoSurf(DM dm,PetscBool * noSurf)381 PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
382 {
383 DM_Plex *plex = (DM_Plex *)dm->data;
384
385 PetscFunctionBegin;
386 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
387 *noSurf = plex->metricCtx->noSurf;
388 PetscFunctionReturn(PETSC_SUCCESS);
389 }
390
391 /*@
392 DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
393
394 Input Parameters:
395 + dm - The `DM`
396 - h_min - The minimum tolerated metric magnitude
397
398 Level: beginner
399
400 .seealso: `DMPLEX`, `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()`
401 @*/
DMPlexMetricSetMinimumMagnitude(DM dm,PetscReal h_min)402 PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
403 {
404 DM_Plex *plex = (DM_Plex *)dm->data;
405
406 PetscFunctionBegin;
407 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
408 PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
409 plex->metricCtx->h_min = h_min;
410 PetscFunctionReturn(PETSC_SUCCESS);
411 }
412
413 /*@
414 DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
415
416 Input Parameters:
417 . dm - The `DM`
418
419 Output Parameters:
420 . h_min - The minimum tolerated metric magnitude
421
422 Level: beginner
423
424 .seealso: `DMPLEX`, `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()`
425 @*/
DMPlexMetricGetMinimumMagnitude(DM dm,PetscReal * h_min)426 PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
427 {
428 DM_Plex *plex = (DM_Plex *)dm->data;
429
430 PetscFunctionBegin;
431 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
432 *h_min = plex->metricCtx->h_min;
433 PetscFunctionReturn(PETSC_SUCCESS);
434 }
435
436 /*@
437 DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
438
439 Input Parameters:
440 + dm - The `DM`
441 - h_max - The maximum tolerated metric magnitude
442
443 Level: beginner
444
445 .seealso: `DMPLEX`, `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()`
446 @*/
DMPlexMetricSetMaximumMagnitude(DM dm,PetscReal h_max)447 PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
448 {
449 DM_Plex *plex = (DM_Plex *)dm->data;
450
451 PetscFunctionBegin;
452 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
453 PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
454 plex->metricCtx->h_max = h_max;
455 PetscFunctionReturn(PETSC_SUCCESS);
456 }
457
458 /*@
459 DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
460
461 Input Parameters:
462 . dm - The `DM`
463
464 Output Parameters:
465 . h_max - The maximum tolerated metric magnitude
466
467 Level: beginner
468
469 .seealso: `DMPLEX`, `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()`
470 @*/
DMPlexMetricGetMaximumMagnitude(DM dm,PetscReal * h_max)471 PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
472 {
473 DM_Plex *plex = (DM_Plex *)dm->data;
474
475 PetscFunctionBegin;
476 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
477 *h_max = plex->metricCtx->h_max;
478 PetscFunctionReturn(PETSC_SUCCESS);
479 }
480
481 /*@
482 DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
483
484 Input Parameters:
485 + dm - The `DM`
486 - a_max - The maximum tolerated metric anisotropy
487
488 Level: beginner
489
490 Note:
491 If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
492
493 .seealso: `DMPLEX`, `DMPlexMetricGetMaximumAnisotropy()`, `DMPlexMetricSetMaximumMagnitude()`
494 @*/
DMPlexMetricSetMaximumAnisotropy(DM dm,PetscReal a_max)495 PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
496 {
497 DM_Plex *plex = (DM_Plex *)dm->data;
498
499 PetscFunctionBegin;
500 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
501 PetscCheck(a_max >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Anisotropy must be in [1, inf)");
502 plex->metricCtx->a_max = a_max;
503 PetscFunctionReturn(PETSC_SUCCESS);
504 }
505
506 /*@
507 DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
508
509 Input Parameters:
510 . dm - The `DM`
511
512 Output Parameters:
513 . a_max - The maximum tolerated metric anisotropy
514
515 Level: beginner
516
517 .seealso: `DMPLEX`, `DMPlexMetricSetMaximumAnisotropy()`, `DMPlexMetricGetMaximumMagnitude()`
518 @*/
DMPlexMetricGetMaximumAnisotropy(DM dm,PetscReal * a_max)519 PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
520 {
521 DM_Plex *plex = (DM_Plex *)dm->data;
522
523 PetscFunctionBegin;
524 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
525 *a_max = plex->metricCtx->a_max;
526 PetscFunctionReturn(PETSC_SUCCESS);
527 }
528
529 /*@
530 DMPlexMetricSetTargetComplexity - Set the target metric complexity
531
532 Input Parameters:
533 + dm - The `DM`
534 - targetComplexity - The target metric complexity
535
536 Level: beginner
537
538 .seealso: `DMPLEX`, `DMPlexMetricGetTargetComplexity()`, `DMPlexMetricSetNormalizationOrder()`
539 @*/
DMPlexMetricSetTargetComplexity(DM dm,PetscReal targetComplexity)540 PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
541 {
542 DM_Plex *plex = (DM_Plex *)dm->data;
543
544 PetscFunctionBegin;
545 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
546 PetscCheck(targetComplexity > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric complexity must be in (0, inf)");
547 plex->metricCtx->targetComplexity = targetComplexity;
548 PetscFunctionReturn(PETSC_SUCCESS);
549 }
550
551 /*@
552 DMPlexMetricGetTargetComplexity - Get the target metric complexity
553
554 Input Parameters:
555 . dm - The `DM`
556
557 Output Parameters:
558 . targetComplexity - The target metric complexity
559
560 Level: beginner
561
562 .seealso: `DMPLEX`, `DMPlexMetricSetTargetComplexity()`, `DMPlexMetricGetNormalizationOrder()`
563 @*/
DMPlexMetricGetTargetComplexity(DM dm,PetscReal * targetComplexity)564 PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
565 {
566 DM_Plex *plex = (DM_Plex *)dm->data;
567
568 PetscFunctionBegin;
569 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
570 *targetComplexity = plex->metricCtx->targetComplexity;
571 PetscFunctionReturn(PETSC_SUCCESS);
572 }
573
574 /*@
575 DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
576
577 Input Parameters:
578 + dm - The `DM`
579 - p - The normalization order
580
581 Level: beginner
582
583 .seealso: `DMPLEX`, `DMPlexMetricGetNormalizationOrder()`, `DMPlexMetricSetTargetComplexity()`
584 @*/
DMPlexMetricSetNormalizationOrder(DM dm,PetscReal p)585 PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
586 {
587 DM_Plex *plex = (DM_Plex *)dm->data;
588
589 PetscFunctionBegin;
590 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
591 PetscCheck(p >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Normalization order must be in [1, inf)");
592 plex->metricCtx->p = p;
593 PetscFunctionReturn(PETSC_SUCCESS);
594 }
595
596 /*@
597 DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
598
599 Input Parameters:
600 . dm - The `DM`
601
602 Output Parameters:
603 . p - The normalization order
604
605 Level: beginner
606
607 .seealso: `DMPLEX`, `DMPlexMetricSetNormalizationOrder()`, `DMPlexMetricGetTargetComplexity()`
608 @*/
DMPlexMetricGetNormalizationOrder(DM dm,PetscReal * p)609 PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
610 {
611 DM_Plex *plex = (DM_Plex *)dm->data;
612
613 PetscFunctionBegin;
614 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
615 *p = plex->metricCtx->p;
616 PetscFunctionReturn(PETSC_SUCCESS);
617 }
618
619 /*@
620 DMPlexMetricSetGradationFactor - Set the metric gradation factor
621
622 Input Parameters:
623 + dm - The `DM`
624 - beta - The metric gradation factor
625
626 Level: beginner
627
628 Notes:
629 The gradation factor is the maximum tolerated length ratio between adjacent edges.
630
631 Turn off gradation by passing the value -1. Otherwise, pass a positive value.
632
633 This is only used by Mmg and ParMmg (not Pragmatic).
634
635 .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
636 @*/
DMPlexMetricSetGradationFactor(DM dm,PetscReal beta)637 PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
638 {
639 DM_Plex *plex = (DM_Plex *)dm->data;
640
641 PetscFunctionBegin;
642 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
643 PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)");
644 plex->metricCtx->gradationFactor = beta;
645 PetscFunctionReturn(PETSC_SUCCESS);
646 }
647
648 /*@
649 DMPlexMetricGetGradationFactor - Get the metric gradation factor
650
651 Input Parameters:
652 . dm - The `DM`
653
654 Output Parameters:
655 . beta - The metric gradation factor
656
657 Level: beginner
658
659 Notes:
660 The gradation factor is the maximum tolerated length ratio between adjacent edges.
661
662 The value -1 implies that gradation is turned off.
663
664 This is only used by Mmg and ParMmg (not Pragmatic).
665
666 .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
667 @*/
DMPlexMetricGetGradationFactor(DM dm,PetscReal * beta)668 PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
669 {
670 DM_Plex *plex = (DM_Plex *)dm->data;
671
672 PetscFunctionBegin;
673 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
674 *beta = plex->metricCtx->gradationFactor;
675 PetscFunctionReturn(PETSC_SUCCESS);
676 }
677
678 /*@
679 DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
680
681 Input Parameters:
682 + dm - The `DM`
683 - hausd - The metric Hausdorff number
684
685 Level: beginner
686
687 Notes:
688 The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
689 boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
690 high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
691 an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
692 (resp. increase) the Hausdorff parameter. (Taken from
693 https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
694
695 This is only used by Mmg and ParMmg (not Pragmatic).
696
697 .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
698 @*/
DMPlexMetricSetHausdorffNumber(DM dm,PetscReal hausd)699 PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
700 {
701 DM_Plex *plex = (DM_Plex *)dm->data;
702
703 PetscFunctionBegin;
704 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
705 PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)");
706 plex->metricCtx->hausdorffNumber = hausd;
707 PetscFunctionReturn(PETSC_SUCCESS);
708 }
709
710 /*@
711 DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
712
713 Input Parameters:
714 . dm - The `DM`
715
716 Output Parameters:
717 . hausd - The metric Hausdorff number
718
719 Level: beginner
720
721 Notes:
722 The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
723 boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
724 high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
725 an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
726 (resp. increase) the Hausdorff parameter. (Taken from
727 https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
728
729 This is only used by Mmg and ParMmg (not Pragmatic).
730
731 .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
732 @*/
DMPlexMetricGetHausdorffNumber(DM dm,PetscReal * hausd)733 PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
734 {
735 DM_Plex *plex = (DM_Plex *)dm->data;
736
737 PetscFunctionBegin;
738 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
739 *hausd = plex->metricCtx->hausdorffNumber;
740 PetscFunctionReturn(PETSC_SUCCESS);
741 }
742
743 /*@
744 DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
745
746 Input Parameters:
747 + dm - The `DM`
748 - verbosity - The verbosity, where -1 is silent and 10 is maximum
749
750 Level: beginner
751
752 Note:
753 This is only used by Mmg and ParMmg (not Pragmatic).
754
755 .seealso: `DMPLEX`, `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()`
756 @*/
DMPlexMetricSetVerbosity(DM dm,PetscInt verbosity)757 PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
758 {
759 DM_Plex *plex = (DM_Plex *)dm->data;
760
761 PetscFunctionBegin;
762 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
763 plex->metricCtx->verbosity = verbosity;
764 PetscFunctionReturn(PETSC_SUCCESS);
765 }
766
767 /*@
768 DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
769
770 Input Parameters:
771 . dm - The `DM`
772
773 Output Parameters:
774 . verbosity - The verbosity, where -1 is silent and 10 is maximum
775
776 Level: beginner
777
778 Note:
779 This is only used by Mmg and ParMmg (not Pragmatic).
780
781 .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
782 @*/
DMPlexMetricGetVerbosity(DM dm,PetscInt * verbosity)783 PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
784 {
785 DM_Plex *plex = (DM_Plex *)dm->data;
786
787 PetscFunctionBegin;
788 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
789 *verbosity = plex->metricCtx->verbosity;
790 PetscFunctionReturn(PETSC_SUCCESS);
791 }
792
793 /*@
794 DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
795
796 Input Parameters:
797 + dm - The `DM`
798 - numIter - the number of parallel adaptation iterations
799
800 Level: beginner
801
802 Note:
803 This is only used by ParMmg (not Pragmatic or Mmg).
804
805 .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
806 @*/
DMPlexMetricSetNumIterations(DM dm,PetscInt numIter)807 PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
808 {
809 DM_Plex *plex = (DM_Plex *)dm->data;
810
811 PetscFunctionBegin;
812 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
813 plex->metricCtx->numIter = numIter;
814 PetscFunctionReturn(PETSC_SUCCESS);
815 }
816
817 /*@
818 DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
819
820 Input Parameters:
821 . dm - The `DM`
822
823 Output Parameters:
824 . numIter - the number of parallel adaptation iterations
825
826 Level: beginner
827
828 Note:
829 This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
830
831 .seealso: `DMPLEX`, `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()`
832 @*/
DMPlexMetricGetNumIterations(DM dm,PetscInt * numIter)833 PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
834 {
835 DM_Plex *plex = (DM_Plex *)dm->data;
836
837 PetscFunctionBegin;
838 if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
839 *numIter = plex->metricCtx->numIter;
840 PetscFunctionReturn(PETSC_SUCCESS);
841 }
842
DMPlexP1FieldCreate_Private(DM dm,PetscInt f,PetscInt size,Vec * metric)843 static PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
844 {
845 MPI_Comm comm;
846 PetscFE fe;
847 PetscInt dim;
848
849 PetscFunctionBegin;
850 /* Extract metadata from dm */
851 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
852 PetscCall(DMGetDimension(dm, &dim));
853
854 /* Create a P1 field of the requested size */
855 PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
856 PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
857 PetscCall(DMCreateDS(dm));
858 PetscCall(PetscFEDestroy(&fe));
859 PetscCall(DMCreateLocalVector(dm, metric));
860 PetscFunctionReturn(PETSC_SUCCESS);
861 }
862
863 /*@
864 DMPlexMetricCreate - Create a Riemannian metric field
865
866 Input Parameters:
867 + dm - The `DM`
868 - f - The field number to use
869
870 Output Parameter:
871 . metric - The metric
872
873 Options Database Key:
874 . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use
875
876 Options Database Keys for Mmg and ParMmg:
877 + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation
878 . -dm_plex_metric_num_iterations - Number of parallel mesh adaptation iterations for ParMmg
879 . -dm_plex_metric_no_insert - Should node insertion/deletion be turned off?
880 . -dm_plex_metric_no_swap - Should facet swapping be turned off?
881 . -dm_plex_metric_no_move - Should node movement be turned off?
882 - -dm_plex_metric_verbosity - Choose a verbosity level from -1 (silent) to 10 (maximum).
883
884 Options Database Keys for Riemannian metrics:
885 + -dm_plex_metric_isotropic - Is the metric isotropic?
886 . -dm_plex_metric_uniform - Is the metric uniform?
887 . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
888 . -dm_plex_metric_h_min - Minimum tolerated metric magnitude
889 . -dm_plex_metric_h_max - Maximum tolerated metric magnitude
890 . -dm_plex_metric_a_max - Maximum tolerated anisotropy
891 . -dm_plex_metric_p - L-p normalization order
892 - -dm_plex_metric_target_complexity - Target metric complexity
893
894 Level: beginner
895
896 Note:
897 It is assumed that the `DM` is comprised of simplices.
898
899 .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`
900 @*/
DMPlexMetricCreate(DM dm,PetscInt f,Vec * metric)901 PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
902 {
903 PetscBool isotropic, uniform;
904 PetscInt coordDim, Nd;
905
906 PetscFunctionBegin;
907 PetscCall(DMGetCoordinateDim(dm, &coordDim));
908 Nd = coordDim * coordDim;
909 PetscCall(DMPlexMetricIsUniform(dm, &uniform));
910 PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
911 if (uniform) {
912 MPI_Comm comm;
913
914 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
915 PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
916 PetscCall(VecCreate(comm, metric));
917 PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
918 PetscCall(VecSetFromOptions(*metric));
919 } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
920 else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
921 PetscFunctionReturn(PETSC_SUCCESS);
922 }
923
924 /*@
925 DMPlexMetricCreateUniform - Construct a uniform isotropic metric
926
927 Input Parameters:
928 + dm - The `DM`
929 . f - The field number to use
930 - alpha - Scaling parameter for the diagonal
931
932 Output Parameter:
933 . metric - The uniform metric
934
935 Level: beginner
936
937 Note:
938 In this case, the metric is constant in space and so is only specified for a single vertex.
939
940 .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
941 @*/
DMPlexMetricCreateUniform(DM dm,PetscInt f,PetscReal alpha,Vec * metric)942 PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
943 {
944 PetscFunctionBegin;
945 PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
946 PetscCall(DMPlexMetricCreate(dm, f, metric));
947 PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
948 PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
949 PetscCall(VecSet(*metric, alpha));
950 PetscCall(VecAssemblyBegin(*metric));
951 PetscCall(VecAssemblyEnd(*metric));
952 PetscFunctionReturn(PETSC_SUCCESS);
953 }
954
identity(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[])955 static void identity(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[])
956 {
957 f0[0] = u[0];
958 }
959
960 /*@
961 DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
962
963 Input Parameters:
964 + dm - The `DM`
965 . f - The field number to use
966 - indicator - The error indicator
967
968 Output Parameter:
969 . metric - The isotropic metric
970
971 Level: beginner
972
973 Notes:
974 It is assumed that the `DM` is comprised of simplices.
975
976 The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
977
978 .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
979 @*/
DMPlexMetricCreateIsotropic(DM dm,PetscInt f,Vec indicator,Vec * metric)980 PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
981 {
982 PetscInt m, n;
983
984 PetscFunctionBegin;
985 PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
986 PetscCall(DMPlexMetricCreate(dm, f, metric));
987 PetscCall(VecGetSize(indicator, &m));
988 PetscCall(VecGetSize(*metric, &n));
989 if (m == n) PetscCall(VecCopy(indicator, *metric));
990 else {
991 DM dmIndi;
992 void (*funcs[1])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
993
994 PetscCall(VecGetDM(indicator, &dmIndi));
995 funcs[0] = identity;
996 PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
997 }
998 PetscFunctionReturn(PETSC_SUCCESS);
999 }
1000
1001 /*@
1002 DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1003
1004 Input Parameters:
1005 + dm - The `DM` of the metric field
1006 - f - The field number to use
1007
1008 Output Parameters:
1009 + determinant - The determinant field
1010 - dmDet - The corresponding `DM`
1011
1012 Level: beginner
1013
1014 .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`, `DMPlexMetricCreate()`
1015 @*/
DMPlexMetricDeterminantCreate(DM dm,PetscInt f,Vec * determinant,DM * dmDet)1016 PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1017 {
1018 PetscBool uniform;
1019
1020 PetscFunctionBegin;
1021 PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1022 PetscCall(DMClone(dm, dmDet));
1023 if (uniform) {
1024 MPI_Comm comm;
1025
1026 PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm));
1027 PetscCall(VecCreate(comm, determinant));
1028 PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1029 PetscCall(VecSetFromOptions(*determinant));
1030 } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
1031 PetscFunctionReturn(PETSC_SUCCESS);
1032 }
1033
LAPACKsyevFail(PetscInt dim,PetscScalar Mpos[])1034 static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1035 {
1036 PetscInt i, j;
1037
1038 PetscFunctionBegin;
1039 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"));
1040 for (i = 0; i < dim; ++i) {
1041 if (i == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " [["));
1042 else PetscCall(PetscPrintf(PETSC_COMM_SELF, " ["));
1043 for (j = 0; j < dim; ++j) {
1044 if (j < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])));
1045 else PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])));
1046 }
1047 if (i < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
1048 else PetscCall(PetscPrintf(PETSC_COMM_SELF, "]]\n"));
1049 }
1050 PetscFunctionReturn(PETSC_SUCCESS);
1051 }
1052
DMPlexMetricModify_Private(PetscInt dim,PetscReal h_min,PetscReal h_max,PetscReal a_max,PetscScalar Mp[],PetscScalar * detMp)1053 static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1054 {
1055 PetscInt i, j, k;
1056 PetscReal *eigs, max_eig, l_min = 1.0 / (h_max * h_max), l_max = 1.0 / (h_min * h_min), la_min = 1.0 / (a_max * a_max);
1057 PetscScalar *Mpos;
1058
1059 PetscFunctionBegin;
1060 PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs));
1061
1062 /* Symmetrize */
1063 for (i = 0; i < dim; ++i) {
1064 Mpos[i * dim + i] = Mp[i * dim + i];
1065 for (j = i + 1; j < dim; ++j) {
1066 Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
1067 Mpos[j * dim + i] = Mpos[i * dim + j];
1068 }
1069 }
1070
1071 /* Compute eigendecomposition */
1072 if (dim == 1) {
1073 /* Isotropic case */
1074 eigs[0] = PetscRealPart(Mpos[0]);
1075 Mpos[0] = 1.0;
1076 } else {
1077 /* Anisotropic case */
1078 PetscScalar *work;
1079 PetscBLASInt lwork;
1080
1081 PetscCall(PetscBLASIntCast(5 * dim, &lwork));
1082 PetscCall(PetscMalloc1(5 * dim, &work));
1083 {
1084 PetscBLASInt lierr;
1085 PetscBLASInt nb;
1086
1087 PetscCall(PetscBLASIntCast(dim, &nb));
1088 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1089 #if defined(PETSC_USE_COMPLEX)
1090 {
1091 PetscReal *rwork;
1092 PetscCall(PetscMalloc1(3 * dim, &rwork));
1093 PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr));
1094 PetscCall(PetscFree(rwork));
1095 }
1096 #else
1097 PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr));
1098 #endif
1099 if (lierr) {
1100 for (i = 0; i < dim; ++i) {
1101 Mpos[i * dim + i] = Mp[i * dim + i];
1102 for (j = i + 1; j < dim; ++j) {
1103 Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
1104 Mpos[j * dim + i] = Mpos[i * dim + j];
1105 }
1106 }
1107 PetscCall(LAPACKsyevFail(dim, Mpos));
1108 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %" PetscBLASInt_FMT, lierr);
1109 }
1110 PetscCall(PetscFPTrapPop());
1111 }
1112 PetscCall(PetscFree(work));
1113 }
1114
1115 /* Reflect to positive orthant and enforce maximum and minimum size */
1116 max_eig = 0.0;
1117 for (i = 0; i < dim; ++i) {
1118 eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
1119 max_eig = PetscMax(eigs[i], max_eig);
1120 }
1121
1122 /* Enforce maximum anisotropy and compute determinant */
1123 *detMp = 1.0;
1124 for (i = 0; i < dim; ++i) {
1125 if (a_max >= 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min);
1126 *detMp *= eigs[i];
1127 }
1128
1129 /* Reconstruct Hessian */
1130 for (i = 0; i < dim; ++i) {
1131 for (j = 0; j < dim; ++j) {
1132 Mp[i * dim + j] = 0.0;
1133 for (k = 0; k < dim; ++k) Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j];
1134 }
1135 }
1136 PetscCall(PetscFree2(Mpos, eigs));
1137 PetscFunctionReturn(PETSC_SUCCESS);
1138 }
1139
1140 /*@
1141 DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
1142
1143 Input Parameters:
1144 + dm - The `DM`
1145 . metricIn - The metric
1146 . restrictSizes - Should maximum/minimum metric magnitudes be enforced?
1147 - restrictAnisotropy - Should maximum anisotropy be enforced?
1148
1149 Output Parameters:
1150 + metricOut - The metric
1151 - determinant - Its determinant
1152
1153 Options Database Keys:
1154 + -dm_plex_metric_isotropic - Is the metric isotropic?
1155 . -dm_plex_metric_uniform - Is the metric uniform?
1156 . -dm_plex_metric_h_min - Minimum tolerated metric magnitude
1157 . -dm_plex_metric_h_max - Maximum tolerated metric magnitude
1158 - -dm_plex_metric_a_max - Maximum tolerated anisotropy
1159
1160 Level: beginner
1161
1162 .seealso: `DMPLEX`, `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1163 @*/
DMPlexMetricEnforceSPD(DM dm,Vec metricIn,PetscBool restrictSizes,PetscBool restrictAnisotropy,Vec metricOut,Vec determinant)1164 PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1165 {
1166 DM dmDet;
1167 PetscBool isotropic, uniform;
1168 PetscInt dim, vStart, vEnd, v;
1169 PetscScalar *met, *det;
1170 PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+15;
1171
1172 PetscFunctionBegin;
1173 PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1174
1175 /* Extract metadata from dm */
1176 PetscCall(DMGetDimension(dm, &dim));
1177 if (restrictSizes) {
1178 PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1179 PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1180 h_min = PetscMax(h_min, 1.0e-30);
1181 h_max = PetscMin(h_max, 1.0e+30);
1182 PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1183 }
1184 if (restrictAnisotropy) {
1185 PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1186 a_max = PetscMin(a_max, 1.0e+30);
1187 }
1188
1189 /* Setup output metric */
1190 PetscCall(VecCopy(metricIn, metricOut));
1191
1192 /* Enforce SPD and extract determinant */
1193 PetscCall(VecGetArray(metricOut, &met));
1194 PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1195 PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1196 if (uniform) {
1197 PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
1198
1199 /* Uniform case */
1200 PetscCall(VecGetArray(determinant, &det));
1201 PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1202 PetscCall(VecRestoreArray(determinant, &det));
1203 } else {
1204 /* Spatially varying case */
1205 PetscInt nrow;
1206
1207 if (isotropic) nrow = 1;
1208 else nrow = dim;
1209 PetscCall(VecGetDM(determinant, &dmDet));
1210 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1211 PetscCall(VecGetArray(determinant, &det));
1212 for (v = vStart; v < vEnd; ++v) {
1213 PetscScalar *vmet, *vdet;
1214 PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
1215 PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
1216 PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
1217 }
1218 PetscCall(VecRestoreArray(determinant, &det));
1219 }
1220 PetscCall(VecRestoreArray(metricOut, &met));
1221
1222 PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1223 PetscFunctionReturn(PETSC_SUCCESS);
1224 }
1225
detMFunc(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[])1226 static void detMFunc(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[])
1227 {
1228 const PetscScalar p = constants[0];
1229
1230 f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim));
1231 }
1232
1233 /*@
1234 DMPlexMetricNormalize - Apply L-p normalization to a metric
1235
1236 Input Parameters:
1237 + dm - The `DM`
1238 . metricIn - The unnormalized metric
1239 . restrictSizes - Should maximum/minimum metric magnitudes be enforced?
1240 - restrictAnisotropy - Should maximum metric anisotropy be enforced?
1241
1242 Output Parameters:
1243 + metricOut - The normalized metric
1244 - determinant - computed determinant
1245
1246 Options Database Keys:
1247 + -dm_plex_metric_isotropic - Is the metric isotropic?
1248 . -dm_plex_metric_uniform - Is the metric uniform?
1249 . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1250 . -dm_plex_metric_h_min - Minimum tolerated metric magnitude
1251 . -dm_plex_metric_h_max - Maximum tolerated metric magnitude
1252 . -dm_plex_metric_a_max - Maximum tolerated anisotropy
1253 . -dm_plex_metric_p - L-p normalization order
1254 - -dm_plex_metric_target_complexity - Target metric complexity
1255
1256 Level: beginner
1257
1258 .seealso: `DMPLEX`, `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1259 @*/
DMPlexMetricNormalize(DM dm,Vec metricIn,PetscBool restrictSizes,PetscBool restrictAnisotropy,Vec metricOut,Vec determinant)1260 PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1261 {
1262 DM dmDet;
1263 MPI_Comm comm;
1264 PetscBool restrictAnisotropyFirst, isotropic, uniform;
1265 PetscDS ds;
1266 PetscInt dim, Nd, vStart, vEnd, v, i;
1267 PetscScalar *met, *det, integral, constants[1];
1268 PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
1269
1270 PetscFunctionBegin;
1271 PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1272
1273 /* Extract metadata from dm */
1274 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1275 PetscCall(DMGetDimension(dm, &dim));
1276 PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1277 PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1278 if (isotropic) Nd = 1;
1279 else Nd = dim * dim;
1280
1281 /* Set up metric and ensure it is SPD */
1282 PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
1283 PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
1284
1285 /* Compute global normalization factor */
1286 PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
1287 PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
1288 constants[0] = p;
1289 if (uniform) {
1290 PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
1291 DM dmTmp;
1292 Vec tmp;
1293
1294 PetscCall(DMClone(dm, &dmTmp));
1295 PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
1296 PetscCall(VecGetArray(determinant, &det));
1297 PetscCall(VecSet(tmp, det[0]));
1298 PetscCall(VecRestoreArray(determinant, &det));
1299 PetscCall(DMGetDS(dmTmp, &ds));
1300 PetscCall(PetscDSSetConstants(ds, 1, constants));
1301 PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1302 PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
1303 PetscCall(VecDestroy(&tmp));
1304 PetscCall(DMDestroy(&dmTmp));
1305 } else {
1306 PetscCall(VecGetDM(determinant, &dmDet));
1307 PetscCall(DMGetDS(dmDet, &ds));
1308 PetscCall(PetscDSSetConstants(ds, 1, constants));
1309 PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1310 PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
1311 }
1312 realIntegral = PetscRealPart(integral);
1313 PetscCheck(realIntegral > 1.0e-30, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Global metric normalization factor must be in (0, inf). Is the input metric positive-definite?");
1314 factGlob = PetscPowReal(target / realIntegral, 2.0 / dim);
1315
1316 /* Apply local scaling */
1317 if (restrictSizes) {
1318 PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1319 PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1320 h_min = PetscMax(h_min, 1.0e-30);
1321 h_max = PetscMin(h_max, 1.0e+30);
1322 PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1323 }
1324 if (restrictAnisotropy && !restrictAnisotropyFirst) {
1325 PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1326 a_max = PetscMin(a_max, 1.0e+30);
1327 }
1328 PetscCall(VecGetArray(metricOut, &met));
1329 PetscCall(VecGetArray(determinant, &det));
1330 if (uniform) {
1331 /* Uniform case */
1332 met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim));
1333 if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1334 } else {
1335 /* Spatially varying case */
1336 PetscInt nrow;
1337
1338 if (isotropic) nrow = 1;
1339 else nrow = dim;
1340 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1341 PetscCall(VecGetDM(determinant, &dmDet));
1342 for (v = vStart; v < vEnd; ++v) {
1343 PetscScalar *Mp, *detM;
1344
1345 PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
1346 PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
1347 fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim));
1348 for (i = 0; i < Nd; ++i) Mp[i] *= fact;
1349 if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
1350 }
1351 }
1352 PetscCall(VecRestoreArray(determinant, &det));
1353 PetscCall(VecRestoreArray(metricOut, &met));
1354
1355 PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1356 PetscFunctionReturn(PETSC_SUCCESS);
1357 }
1358
1359 /*@
1360 DMPlexMetricAverage - Compute the average of a list of metrics
1361
1362 Input Parameters:
1363 + dm - The `DM`
1364 . numMetrics - The number of metrics to be averaged
1365 . weights - Weights for the average
1366 - metrics - The metrics to be averaged
1367
1368 Output Parameter:
1369 . metricAvg - The averaged metric
1370
1371 Level: beginner
1372
1373 Notes:
1374 The weights should sum to unity.
1375
1376 If weights are not provided then an unweighted average is used.
1377
1378 .seealso: `DMPLEX`, `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1379 @*/
DMPlexMetricAverage(DM dm,PetscInt numMetrics,PetscReal weights[],Vec metrics[],Vec metricAvg)1380 PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1381 {
1382 PetscBool haveWeights = PETSC_TRUE;
1383 PetscInt i, m, n;
1384 PetscReal sum = 0.0, tol = 1.0e-10;
1385
1386 PetscFunctionBegin;
1387 PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0));
1388 PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
1389 PetscCall(VecSet(metricAvg, 0.0));
1390 PetscCall(VecGetSize(metricAvg, &m));
1391 for (i = 0; i < numMetrics; ++i) {
1392 PetscCall(VecGetSize(metrics[i], &n));
1393 PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
1394 }
1395
1396 /* Default to the unweighted case */
1397 if (!weights) {
1398 PetscCall(PetscMalloc1(numMetrics, &weights));
1399 haveWeights = PETSC_FALSE;
1400 for (i = 0; i < numMetrics; ++i) weights[i] = 1.0 / numMetrics;
1401 }
1402
1403 /* Check weights sum to unity */
1404 for (i = 0; i < numMetrics; ++i) sum += weights[i];
1405 PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
1406
1407 /* Compute metric average */
1408 for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
1409 if (!haveWeights) PetscCall(PetscFree(weights));
1410
1411 PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0));
1412 PetscFunctionReturn(PETSC_SUCCESS);
1413 }
1414
1415 /*@
1416 DMPlexMetricAverage2 - Compute the unweighted average of two metrics
1417
1418 Input Parameters:
1419 + dm - The `DM`
1420 . metric1 - The first metric to be averaged
1421 - metric2 - The second metric to be averaged
1422
1423 Output Parameter:
1424 . metricAvg - The averaged metric
1425
1426 Level: beginner
1427
1428 .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1429 @*/
DMPlexMetricAverage2(DM dm,Vec metric1,Vec metric2,Vec metricAvg)1430 PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1431 {
1432 PetscReal weights[2] = {0.5, 0.5};
1433 Vec metrics[2] = {metric1, metric2};
1434
1435 PetscFunctionBegin;
1436 PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
1437 PetscFunctionReturn(PETSC_SUCCESS);
1438 }
1439
1440 /*@
1441 DMPlexMetricAverage3 - Compute the unweighted average of three metrics
1442
1443 Input Parameters:
1444 + dm - The `DM`
1445 . metric1 - The first metric to be averaged
1446 . metric2 - The second metric to be averaged
1447 - metric3 - The third metric to be averaged
1448
1449 Output Parameter:
1450 . metricAvg - The averaged metric
1451
1452 Level: beginner
1453
1454 .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1455 @*/
DMPlexMetricAverage3(DM dm,Vec metric1,Vec metric2,Vec metric3,Vec metricAvg)1456 PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1457 {
1458 PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0};
1459 Vec metrics[3] = {metric1, metric2, metric3};
1460
1461 PetscFunctionBegin;
1462 PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
1463 PetscFunctionReturn(PETSC_SUCCESS);
1464 }
1465
DMPlexMetricIntersection_Private(PetscInt dim,PetscScalar M1[],PetscScalar M2[])1466 static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1467 {
1468 PetscInt i, j, k, l, m;
1469 PetscReal *evals;
1470 PetscScalar *evecs, *sqrtM1, *isqrtM1;
1471
1472 PetscFunctionBegin;
1473 /* Isotropic case */
1474 if (dim == 1) {
1475 M2[0] = PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
1476 PetscFunctionReturn(PETSC_SUCCESS);
1477 }
1478
1479 /* Anisotropic case */
1480 PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals));
1481 for (i = 0; i < dim; ++i) {
1482 for (j = 0; j < dim; ++j) evecs[i * dim + j] = M1[i * dim + j];
1483 }
1484 {
1485 PetscScalar *work;
1486 PetscBLASInt lwork;
1487
1488 PetscCall(PetscBLASIntCast(5 * dim, &lwork));
1489 PetscCall(PetscMalloc1(5 * dim, &work));
1490 {
1491 PetscBLASInt lierr, nb;
1492 PetscReal sqrtj;
1493
1494 /* Compute eigendecomposition of M1 */
1495 PetscCall(PetscBLASIntCast(dim, &nb));
1496 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1497 #if defined(PETSC_USE_COMPLEX)
1498 {
1499 PetscReal *rwork;
1500 PetscCall(PetscMalloc1(3 * dim, &rwork));
1501 PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1502 PetscCall(PetscFree(rwork));
1503 }
1504 #else
1505 PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1506 #endif
1507 if (lierr) {
1508 PetscCall(LAPACKsyevFail(dim, M1));
1509 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %" PetscBLASInt_FMT, lierr);
1510 }
1511 PetscCall(PetscFPTrapPop());
1512
1513 /* Compute square root and the reciprocal thereof */
1514 for (i = 0; i < dim; ++i) {
1515 for (k = 0; k < dim; ++k) {
1516 sqrtM1[i * dim + k] = 0.0;
1517 isqrtM1[i * dim + k] = 0.0;
1518 for (j = 0; j < dim; ++j) {
1519 sqrtj = PetscSqrtReal(evals[j]);
1520 sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k];
1521 isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k];
1522 }
1523 }
1524 }
1525
1526 /* Map M2 into the space spanned by the eigenvectors of M1 */
1527 for (i = 0; i < dim; ++i) {
1528 for (l = 0; l < dim; ++l) {
1529 evecs[i * dim + l] = 0.0;
1530 for (j = 0; j < dim; ++j) {
1531 for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
1532 }
1533 }
1534 }
1535
1536 /* Compute eigendecomposition */
1537 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1538 #if defined(PETSC_USE_COMPLEX)
1539 {
1540 PetscReal *rwork;
1541 PetscCall(PetscMalloc1(3 * dim, &rwork));
1542 PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1543 PetscCall(PetscFree(rwork));
1544 }
1545 #else
1546 PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1547 #endif
1548 if (lierr) {
1549 for (i = 0; i < dim; ++i) {
1550 for (l = 0; l < dim; ++l) {
1551 evecs[i * dim + l] = 0.0;
1552 for (j = 0; j < dim; ++j) {
1553 for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
1554 }
1555 }
1556 }
1557 PetscCall(LAPACKsyevFail(dim, evecs));
1558 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %" PetscBLASInt_FMT, lierr);
1559 }
1560 PetscCall(PetscFPTrapPop());
1561
1562 /* Modify eigenvalues */
1563 for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
1564
1565 /* Map back to get the intersection */
1566 for (i = 0; i < dim; ++i) {
1567 for (m = 0; m < dim; ++m) {
1568 M2[i * dim + m] = 0.0;
1569 for (j = 0; j < dim; ++j) {
1570 for (k = 0; k < dim; ++k) {
1571 for (l = 0; l < dim; ++l) M2[i * dim + m] += sqrtM1[j * dim + i] * evecs[j * dim + k] * evals[k] * evecs[l * dim + k] * sqrtM1[l * dim + m];
1572 }
1573 }
1574 }
1575 }
1576 }
1577 PetscCall(PetscFree(work));
1578 }
1579 PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
1580 PetscFunctionReturn(PETSC_SUCCESS);
1581 }
1582
1583 /*@
1584 DMPlexMetricIntersection - Compute the intersection of a list of metrics
1585
1586 Input Parameters:
1587 + dm - The `DM`
1588 . numMetrics - The number of metrics to be intersected
1589 - metrics - The metrics to be intersected
1590
1591 Output Parameter:
1592 . metricInt - The intersected metric
1593
1594 Level: beginner
1595
1596 Notes:
1597 The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
1598
1599 The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
1600
1601 .seealso: `DMPLEX`, `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1602 @*/
DMPlexMetricIntersection(DM dm,PetscInt numMetrics,Vec metrics[],Vec metricInt)1603 PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1604 {
1605 PetscBool isotropic, uniform;
1606 PetscInt v, i, m, n;
1607 PetscScalar *met, *meti;
1608
1609 PetscFunctionBegin;
1610 PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1611 PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
1612
1613 /* Copy over the first metric */
1614 PetscCall(VecCopy(metrics[0], metricInt));
1615 if (numMetrics == 1) PetscFunctionReturn(PETSC_SUCCESS);
1616 PetscCall(VecGetSize(metricInt, &m));
1617 for (i = 0; i < numMetrics; ++i) {
1618 PetscCall(VecGetSize(metrics[i], &n));
1619 PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
1620 }
1621
1622 /* Intersect subsequent metrics in turn */
1623 PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1624 PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1625 if (uniform) {
1626 /* Uniform case */
1627 PetscCall(VecGetArray(metricInt, &met));
1628 for (i = 1; i < numMetrics; ++i) {
1629 PetscCall(VecGetArray(metrics[i], &meti));
1630 PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
1631 PetscCall(VecRestoreArray(metrics[i], &meti));
1632 }
1633 PetscCall(VecRestoreArray(metricInt, &met));
1634 } else {
1635 /* Spatially varying case */
1636 PetscInt dim, vStart, vEnd, nrow;
1637 PetscScalar *M, *Mi;
1638
1639 PetscCall(DMGetDimension(dm, &dim));
1640 if (isotropic) nrow = 1;
1641 else nrow = dim;
1642 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1643 PetscCall(VecGetArray(metricInt, &met));
1644 for (i = 1; i < numMetrics; ++i) {
1645 PetscCall(VecGetArray(metrics[i], &meti));
1646 for (v = vStart; v < vEnd; ++v) {
1647 PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
1648 PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
1649 PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
1650 }
1651 PetscCall(VecRestoreArray(metrics[i], &meti));
1652 }
1653 PetscCall(VecRestoreArray(metricInt, &met));
1654 }
1655
1656 PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1657 PetscFunctionReturn(PETSC_SUCCESS);
1658 }
1659
1660 /*@
1661 DMPlexMetricIntersection2 - Compute the intersection of two metrics
1662
1663 Input Parameters:
1664 + dm - The `DM`
1665 . metric1 - The first metric to be intersected
1666 - metric2 - The second metric to be intersected
1667
1668 Output Parameter:
1669 . metricInt - The intersected metric
1670
1671 Level: beginner
1672
1673 .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1674 @*/
DMPlexMetricIntersection2(DM dm,Vec metric1,Vec metric2,Vec metricInt)1675 PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1676 {
1677 Vec metrics[2] = {metric1, metric2};
1678
1679 PetscFunctionBegin;
1680 PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
1681 PetscFunctionReturn(PETSC_SUCCESS);
1682 }
1683
1684 /*@
1685 DMPlexMetricIntersection3 - Compute the intersection of three metrics
1686
1687 Input Parameters:
1688 + dm - The `DM`
1689 . metric1 - The first metric to be intersected
1690 . metric2 - The second metric to be intersected
1691 - metric3 - The third metric to be intersected
1692
1693 Output Parameter:
1694 . metricInt - The intersected metric
1695
1696 Level: beginner
1697
1698 .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1699 @*/
DMPlexMetricIntersection3(DM dm,Vec metric1,Vec metric2,Vec metric3,Vec metricInt)1700 PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1701 {
1702 Vec metrics[3] = {metric1, metric2, metric3};
1703
1704 PetscFunctionBegin;
1705 PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
1706 PetscFunctionReturn(PETSC_SUCCESS);
1707 }
1708