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