xref: /petsc/src/dm/impls/plex/plexmetric.c (revision 24ded41b4e3afbef0dd5eaa1b3d8dd0172f6dba2)
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   PetscBool isotropic, uniform;
914   PetscInt  coordDim, Nd;
915 
916   PetscFunctionBegin;
917   PetscCall(DMGetCoordinateDim(dm, &coordDim));
918   Nd = coordDim*coordDim;
919   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
920   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
921   if (uniform) {
922     MPI_Comm comm;
923 
924     PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
925     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
926     PetscCall(VecCreate(comm, metric));
927     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
928     PetscCall(VecSetFromOptions(*metric));
929   } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
930   else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
931   PetscFunctionReturn(0);
932 }
933 
934 /*@
935   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
936 
937   Input parameters:
938 + dm     - The DM
939 . f      - The field number to use
940 - alpha  - Scaling parameter for the diagonal
941 
942   Output parameter:
943 . metric - The uniform metric
944 
945   Level: beginner
946 
947   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
948 
949 .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
950 @*/
951 PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
952 {
953   PetscFunctionBegin;
954   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
955   PetscCall(DMPlexMetricCreate(dm, f, metric));
956   PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
957   PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
958   PetscCall(VecSet(*metric, alpha));
959   PetscCall(VecAssemblyBegin(*metric));
960   PetscCall(VecAssemblyEnd(*metric));
961   PetscFunctionReturn(0);
962 }
963 
964 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
965                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
966                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
967                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
968 {
969   f0[0] = u[0];
970 }
971 
972 /*@
973   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
974 
975   Input parameters:
976 + dm        - The DM
977 . f         - The field number to use
978 - indicator - The error indicator
979 
980   Output parameter:
981 . metric    - The isotropic metric
982 
983   Level: beginner
984 
985   Notes:
986 
987   It is assumed that the DM is comprised of simplices.
988 
989   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
990 
991 .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
992 @*/
993 PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
994 {
995   PetscInt       m, n;
996 
997   PetscFunctionBegin;
998   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
999   PetscCall(DMPlexMetricCreate(dm, f, metric));
1000   PetscCall(VecGetSize(indicator, &m));
1001   PetscCall(VecGetSize(*metric, &n));
1002   if (m == n) PetscCall(VecCopy(indicator, *metric));
1003   else {
1004     DM     dmIndi;
1005     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
1006                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1007                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1008                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
1009 
1010     PetscCall(VecGetDM(indicator, &dmIndi));
1011     funcs[0] = identity;
1012     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
1013   }
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 /*@
1018   DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1019 
1020   Input parameters:
1021 + dm          - The DM of the metric field
1022 - f           - The field number to use
1023 
1024   Output parameter:
1025 + determinant - The determinant field
1026 - dmDet       - The corresponding DM
1027 
1028   Level: beginner
1029 
1030 .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic(), DMPlexMetricCreate()
1031 @*/
1032 PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1033 {
1034   PetscBool uniform;
1035 
1036   PetscFunctionBegin;
1037   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1038   PetscCall(DMClone(dm, dmDet));
1039   if (uniform) {
1040     MPI_Comm comm;
1041 
1042     PetscCall(PetscObjectGetComm((PetscObject) *dmDet, &comm));
1043     PetscCall(VecCreate(comm, determinant));
1044     PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1045     PetscCall(VecSetFromOptions(*determinant));
1046   } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
1047   PetscFunctionReturn(0);
1048 }
1049 
1050 static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1051 {
1052   PetscInt i, j;
1053 
1054   PetscFunctionBegin;
1055   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
1056   for (i = 0; i < dim; ++i) {
1057     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
1058     else        PetscPrintf(PETSC_COMM_SELF, "     [");
1059     for (j = 0; j < dim; ++j) {
1060       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i*dim+j]));
1061       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i*dim+j]));
1062     }
1063     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
1064     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
1065   }
1066   PetscFunctionReturn(0);
1067 }
1068 
1069 static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1070 {
1071   PetscInt       i, j, k;
1072   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);
1073   PetscScalar   *Mpos;
1074 
1075   PetscFunctionBegin;
1076   PetscCall(PetscMalloc2(dim*dim, &Mpos, dim, &eigs));
1077 
1078   /* Symmetrize */
1079   for (i = 0; i < dim; ++i) {
1080     Mpos[i*dim+i] = Mp[i*dim+i];
1081     for (j = i+1; j < dim; ++j) {
1082       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
1083       Mpos[j*dim+i] = Mpos[i*dim+j];
1084     }
1085   }
1086 
1087   /* Compute eigendecomposition */
1088   if (dim == 1) {
1089 
1090     /* Isotropic case */
1091     eigs[0] = PetscRealPart(Mpos[0]);
1092     Mpos[0] = 1.0;
1093   } else {
1094 
1095     /* Anisotropic case */
1096     PetscScalar  *work;
1097     PetscBLASInt lwork;
1098 
1099     lwork = 5*dim;
1100     PetscCall(PetscMalloc1(5*dim, &work));
1101     {
1102       PetscBLASInt lierr;
1103       PetscBLASInt nb;
1104 
1105       PetscCall(PetscBLASIntCast(dim, &nb));
1106       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1107 #if defined(PETSC_USE_COMPLEX)
1108       {
1109         PetscReal *rwork;
1110         PetscCall(PetscMalloc1(3*dim, &rwork));
1111         PetscCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
1112         PetscCall(PetscFree(rwork));
1113       }
1114 #else
1115       PetscCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
1116 #endif
1117       if (lierr) {
1118         for (i = 0; i < dim; ++i) {
1119           Mpos[i*dim+i] = Mp[i*dim+i];
1120           for (j = i+1; j < dim; ++j) {
1121             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
1122             Mpos[j*dim+i] = Mpos[i*dim+j];
1123           }
1124         }
1125         LAPACKsyevFail(dim, Mpos);
1126         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
1127       }
1128       PetscCall(PetscFPTrapPop());
1129     }
1130     PetscCall(PetscFree(work));
1131   }
1132 
1133   /* Reflect to positive orthant and enforce maximum and minimum size */
1134   max_eig = 0.0;
1135   for (i = 0; i < dim; ++i) {
1136     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
1137     max_eig = PetscMax(eigs[i], max_eig);
1138   }
1139 
1140   /* Enforce maximum anisotropy and compute determinant */
1141   *detMp = 1.0;
1142   for (i = 0; i < dim; ++i) {
1143     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
1144     *detMp *= eigs[i];
1145   }
1146 
1147   /* Reconstruct Hessian */
1148   for (i = 0; i < dim; ++i) {
1149     for (j = 0; j < dim; ++j) {
1150       Mp[i*dim+j] = 0.0;
1151       for (k = 0; k < dim; ++k) {
1152         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
1153       }
1154     }
1155   }
1156   PetscCall(PetscFree2(Mpos, eigs));
1157 
1158   PetscFunctionReturn(0);
1159 }
1160 
1161 /*@
1162   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
1163 
1164   Input parameters:
1165 + dm                 - The DM
1166 . metricIn           - The metric
1167 . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
1168 - restrictAnisotropy - Should maximum anisotropy be enforced?
1169 
1170   Output parameter:
1171 + metricOut          - The metric
1172 - determinant        - Its determinant
1173 
1174   Level: beginner
1175 
1176   Notes:
1177 
1178   Relevant command line options:
1179 
1180 + -dm_plex_metric_isotropic - Is the metric isotropic?
1181 . -dm_plex_metric_uniform   - Is the metric uniform?
1182 . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1183 . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1184 - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1185 
1186 .seealso: `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1187 @*/
1188 PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1189 {
1190   DM             dmDet;
1191   PetscBool      isotropic, uniform;
1192   PetscInt       dim, vStart, vEnd, v;
1193   PetscScalar   *met, *det;
1194   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
1195 
1196   PetscFunctionBegin;
1197   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0));
1198 
1199   /* Extract metadata from dm */
1200   PetscCall(DMGetDimension(dm, &dim));
1201   if (restrictSizes) {
1202     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1203     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1204     h_min = PetscMax(h_min, 1.0e-30);
1205     h_max = PetscMin(h_max, 1.0e+30);
1206     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1207   }
1208   if (restrictAnisotropy) {
1209     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1210     a_max = PetscMin(a_max, 1.0e+30);
1211   }
1212 
1213   /* Setup output metric */
1214   PetscCall(VecCopy(metricIn, metricOut));
1215 
1216   /* Enforce SPD and extract determinant */
1217   PetscCall(VecGetArray(metricOut, &met));
1218   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1219   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1220   if (uniform) {
1221     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
1222 
1223     /* Uniform case */
1224     PetscCall(VecGetArray(determinant, &det));
1225     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1226     PetscCall(VecRestoreArray(determinant, &det));
1227   } else {
1228 
1229     /* Spatially varying case */
1230     PetscInt nrow;
1231 
1232     if (isotropic) nrow = 1;
1233     else nrow = dim;
1234     PetscCall(VecGetDM(determinant, &dmDet));
1235     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1236     PetscCall(VecGetArray(determinant, &det));
1237     for (v = vStart; v < vEnd; ++v) {
1238       PetscScalar *vmet, *vdet;
1239       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
1240       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
1241       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
1242     }
1243     PetscCall(VecRestoreArray(determinant, &det));
1244   }
1245   PetscCall(VecRestoreArray(metricOut, &met));
1246 
1247   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0));
1248   PetscFunctionReturn(0);
1249 }
1250 
1251 static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1252                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1253                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1254                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1255 {
1256   const PetscScalar p = constants[0];
1257 
1258   f0[0] = PetscPowScalar(u[0], p/(2.0*p + dim));
1259 }
1260 
1261 /*@
1262   DMPlexMetricNormalize - Apply L-p normalization to a metric
1263 
1264   Input parameters:
1265 + dm                 - The DM
1266 . metricIn           - The unnormalized metric
1267 . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
1268 - restrictAnisotropy - Should maximum metric anisotropy be enforced?
1269 
1270   Output parameter:
1271 . metricOut          - The normalized metric
1272 
1273   Level: beginner
1274 
1275   Notes:
1276 
1277   Relevant command line options:
1278 
1279 + -dm_plex_metric_isotropic                 - Is the metric isotropic?
1280 . -dm_plex_metric_uniform                   - Is the metric uniform?
1281 . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1282 . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1283 . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1284 . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1285 . -dm_plex_metric_p                         - L-p normalization order
1286 - -dm_plex_metric_target_complexity         - Target metric complexity
1287 
1288 .seealso: `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1289 @*/
1290 PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1291 {
1292   DM               dmDet;
1293   MPI_Comm         comm;
1294   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
1295   PetscDS          ds;
1296   PetscInt         dim, Nd, vStart, vEnd, v, i;
1297   PetscScalar     *met, *det, integral, constants[1];
1298   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
1299 
1300   PetscFunctionBegin;
1301   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0));
1302 
1303   /* Extract metadata from dm */
1304   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
1305   PetscCall(DMGetDimension(dm, &dim));
1306   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1307   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1308   if (isotropic) Nd = 1;
1309   else Nd = dim*dim;
1310 
1311   /* Set up metric and ensure it is SPD */
1312   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
1313   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
1314 
1315   /* Compute global normalization factor */
1316   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
1317   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
1318   constants[0] = p;
1319   if (uniform) {
1320     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
1321     DM  dmTmp;
1322     Vec tmp;
1323 
1324     PetscCall(DMClone(dm, &dmTmp));
1325     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
1326     PetscCall(VecGetArray(determinant, &det));
1327     PetscCall(VecSet(tmp, det[0]));
1328     PetscCall(VecRestoreArray(determinant, &det));
1329     PetscCall(DMGetDS(dmTmp, &ds));
1330     PetscCall(PetscDSSetConstants(ds, 1, constants));
1331     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1332     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
1333     PetscCall(VecDestroy(&tmp));
1334     PetscCall(DMDestroy(&dmTmp));
1335   } else {
1336     PetscCall(VecGetDM(determinant, &dmDet));
1337     PetscCall(DMGetDS(dmDet, &ds));
1338     PetscCall(PetscDSSetConstants(ds, 1, constants));
1339     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1340     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
1341   }
1342   realIntegral = PetscRealPart(integral);
1343   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?");
1344   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
1345 
1346   /* Apply local scaling */
1347   if (restrictSizes) {
1348     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1349     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1350     h_min = PetscMax(h_min, 1.0e-30);
1351     h_max = PetscMin(h_max, 1.0e+30);
1352     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1353   }
1354   if (restrictAnisotropy && !restrictAnisotropyFirst) {
1355     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1356     a_max = PetscMin(a_max, 1.0e+30);
1357   }
1358   PetscCall(VecGetArray(metricOut, &met));
1359   PetscCall(VecGetArray(determinant, &det));
1360   if (uniform) {
1361 
1362     /* Uniform case */
1363     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
1364     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1365   } else {
1366 
1367     /* Spatially varying case */
1368     PetscInt nrow;
1369 
1370     if (isotropic) nrow = 1;
1371     else nrow = dim;
1372     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1373     PetscCall(VecGetDM(determinant, &dmDet));
1374     for (v = vStart; v < vEnd; ++v) {
1375       PetscScalar *Mp, *detM;
1376 
1377       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
1378       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
1379       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
1380       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
1381       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
1382     }
1383   }
1384   PetscCall(VecRestoreArray(determinant, &det));
1385   PetscCall(VecRestoreArray(metricOut, &met));
1386 
1387   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0));
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 /*@
1392   DMPlexMetricAverage - Compute the average of a list of metrics
1393 
1394   Input Parameters:
1395 + dm         - The DM
1396 . numMetrics - The number of metrics to be averaged
1397 . weights    - Weights for the average
1398 - metrics    - The metrics to be averaged
1399 
1400   Output Parameter:
1401 . metricAvg  - The averaged metric
1402 
1403   Level: beginner
1404 
1405   Notes:
1406   The weights should sum to unity.
1407 
1408   If weights are not provided then an unweighted average is used.
1409 
1410 .seealso: `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1411 @*/
1412 PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1413 {
1414   PetscBool haveWeights = PETSC_TRUE;
1415   PetscInt  i, m, n;
1416   PetscReal sum = 0.0, tol = 1.0e-10;
1417 
1418   PetscFunctionBegin;
1419   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0));
1420   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
1421   PetscCall(VecSet(metricAvg, 0.0));
1422   PetscCall(VecGetSize(metricAvg, &m));
1423   for (i = 0; i < numMetrics; ++i) {
1424     PetscCall(VecGetSize(metrics[i], &n));
1425     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
1426   }
1427 
1428   /* Default to the unweighted case */
1429   if (!weights) {
1430     PetscCall(PetscMalloc1(numMetrics, &weights));
1431     haveWeights = PETSC_FALSE;
1432     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
1433   }
1434 
1435   /* Check weights sum to unity */
1436   for (i = 0; i < numMetrics; ++i) sum += weights[i];
1437   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
1438 
1439   /* Compute metric average */
1440   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
1441   if (!haveWeights) PetscCall(PetscFree(weights));
1442 
1443   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0));
1444   PetscFunctionReturn(0);
1445 }
1446 
1447 /*@
1448   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
1449 
1450   Input Parameters:
1451 + dm         - The DM
1452 . metric1    - The first metric to be averaged
1453 - metric2    - The second metric to be averaged
1454 
1455   Output Parameter:
1456 . metricAvg  - The averaged metric
1457 
1458   Level: beginner
1459 
1460 .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1461 @*/
1462 PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1463 {
1464   PetscReal      weights[2] = {0.5, 0.5};
1465   Vec            metrics[2] = {metric1, metric2};
1466 
1467   PetscFunctionBegin;
1468   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
1469   PetscFunctionReturn(0);
1470 }
1471 
1472 /*@
1473   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
1474 
1475   Input Parameters:
1476 + dm         - The DM
1477 . metric1    - The first metric to be averaged
1478 . metric2    - The second metric to be averaged
1479 - metric3    - The third metric to be averaged
1480 
1481   Output Parameter:
1482 . metricAvg  - The averaged metric
1483 
1484   Level: beginner
1485 
1486 .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1487 @*/
1488 PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1489 {
1490   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
1491   Vec            metrics[3] = {metric1, metric2, metric3};
1492 
1493   PetscFunctionBegin;
1494   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
1495   PetscFunctionReturn(0);
1496 }
1497 
1498 static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1499 {
1500   PetscInt       i, j, k, l, m;
1501   PetscReal     *evals;
1502   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
1503 
1504   PetscFunctionBegin;
1505 
1506   /* Isotropic case */
1507   if (dim == 1) {
1508     M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
1509     PetscFunctionReturn(0);
1510   }
1511 
1512   /* Anisotropic case */
1513   PetscCall(PetscMalloc4(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals));
1514   for (i = 0; i < dim; ++i) {
1515     for (j = 0; j < dim; ++j) {
1516       evecs[i*dim+j] = M1[i*dim+j];
1517     }
1518   }
1519   {
1520     PetscScalar *work;
1521     PetscBLASInt lwork;
1522 
1523     lwork = 5*dim;
1524     PetscCall(PetscMalloc1(5*dim, &work));
1525     {
1526       PetscBLASInt lierr, nb;
1527       PetscReal    sqrtj;
1528 
1529       /* Compute eigendecomposition of M1 */
1530       PetscCall(PetscBLASIntCast(dim, &nb));
1531       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1532 #if defined(PETSC_USE_COMPLEX)
1533       {
1534         PetscReal *rwork;
1535         PetscCall(PetscMalloc1(3*dim, &rwork));
1536         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1537         PetscCall(PetscFree(rwork));
1538       }
1539 #else
1540       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1541 #endif
1542       if (lierr) {
1543         LAPACKsyevFail(dim, M1);
1544         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
1545       }
1546       PetscCall(PetscFPTrapPop());
1547 
1548       /* Compute square root and the reciprocal thereof */
1549       for (i = 0; i < dim; ++i) {
1550         for (k = 0; k < dim; ++k) {
1551           sqrtM1[i*dim+k] = 0.0;
1552           isqrtM1[i*dim+k] = 0.0;
1553           for (j = 0; j < dim; ++j) {
1554             sqrtj = PetscSqrtReal(evals[j]);
1555             sqrtM1[i*dim+k] += evecs[j*dim+i] * sqrtj * evecs[j*dim+k];
1556             isqrtM1[i*dim+k] += evecs[j*dim+i] * (1.0/sqrtj) * evecs[j*dim+k];
1557           }
1558         }
1559       }
1560 
1561       /* Map M2 into the space spanned by the eigenvectors of M1 */
1562       for (i = 0; i < dim; ++i) {
1563         for (l = 0; l < dim; ++l) {
1564           evecs[i*dim+l] = 0.0;
1565           for (j = 0; j < dim; ++j) {
1566             for (k = 0; k < dim; ++k) {
1567               evecs[i*dim+l] += isqrtM1[j*dim+i] * M2[j*dim+k] * isqrtM1[k*dim+l];
1568             }
1569           }
1570         }
1571       }
1572 
1573       /* Compute eigendecomposition */
1574       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1575 #if defined(PETSC_USE_COMPLEX)
1576       {
1577         PetscReal *rwork;
1578         PetscCall(PetscMalloc1(3*dim, &rwork));
1579         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1580         PetscCall(PetscFree(rwork));
1581       }
1582 #else
1583       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1584 #endif
1585       if (lierr) {
1586         for (i = 0; i < dim; ++i) {
1587           for (l = 0; l < dim; ++l) {
1588             evecs[i*dim+l] = 0.0;
1589             for (j = 0; j < dim; ++j) {
1590               for (k = 0; k < dim; ++k) {
1591                 evecs[i*dim+l] += isqrtM1[j*dim+i] * M2[j*dim+k] * isqrtM1[k*dim+l];
1592               }
1593             }
1594           }
1595         }
1596         LAPACKsyevFail(dim, evecs);
1597         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
1598       }
1599       PetscCall(PetscFPTrapPop());
1600 
1601       /* Modify eigenvalues */
1602       for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
1603 
1604       /* Map back to get the intersection */
1605       for (i = 0; i < dim; ++i) {
1606         for (m = 0; m < dim; ++m) {
1607           M2[i*dim+m] = 0.0;
1608           for (j = 0; j < dim; ++j) {
1609             for (k = 0; k < dim; ++k) {
1610               for (l = 0; l < dim; ++l) {
1611                 M2[i*dim+m] += sqrtM1[j*dim+i] * evecs[j*dim+k] * evals[k] * evecs[l*dim+k] * sqrtM1[l*dim+m];
1612               }
1613             }
1614           }
1615         }
1616       }
1617     }
1618     PetscCall(PetscFree(work));
1619   }
1620   PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
1621   PetscFunctionReturn(0);
1622 }
1623 
1624 /*@
1625   DMPlexMetricIntersection - Compute the intersection of a list of metrics
1626 
1627   Input Parameters:
1628 + dm         - The DM
1629 . numMetrics - The number of metrics to be intersected
1630 - metrics    - The metrics to be intersected
1631 
1632   Output Parameter:
1633 . metricInt  - The intersected metric
1634 
1635   Level: beginner
1636 
1637   Notes:
1638 
1639   The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
1640 
1641   The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
1642 
1643 .seealso: `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1644 @*/
1645 PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1646 {
1647   PetscBool      isotropic, uniform;
1648   PetscInt       v, i, m, n;
1649   PetscScalar   *met, *meti;
1650 
1651   PetscFunctionBegin;
1652   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0));
1653   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
1654 
1655   /* Copy over the first metric */
1656   PetscCall(VecCopy(metrics[0], metricInt));
1657   if (numMetrics == 1) PetscFunctionReturn(0);
1658   PetscCall(VecGetSize(metricInt, &m));
1659   for (i = 0; i < numMetrics; ++i) {
1660     PetscCall(VecGetSize(metrics[i], &n));
1661     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
1662   }
1663 
1664   /* Intersect subsequent metrics in turn */
1665   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1666   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1667   if (uniform) {
1668 
1669     /* Uniform case */
1670     PetscCall(VecGetArray(metricInt, &met));
1671     for (i = 1; i < numMetrics; ++i) {
1672       PetscCall(VecGetArray(metrics[i], &meti));
1673       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
1674       PetscCall(VecRestoreArray(metrics[i], &meti));
1675     }
1676     PetscCall(VecRestoreArray(metricInt, &met));
1677   } else {
1678 
1679     /* Spatially varying case */
1680     PetscInt     dim, vStart, vEnd, nrow;
1681     PetscScalar *M, *Mi;
1682 
1683     PetscCall(DMGetDimension(dm, &dim));
1684     if (isotropic) nrow = 1;
1685     else nrow = dim;
1686     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1687     PetscCall(VecGetArray(metricInt, &met));
1688     for (i = 1; i < numMetrics; ++i) {
1689       PetscCall(VecGetArray(metrics[i], &meti));
1690       for (v = vStart; v < vEnd; ++v) {
1691         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
1692         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
1693         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
1694       }
1695       PetscCall(VecRestoreArray(metrics[i], &meti));
1696     }
1697     PetscCall(VecRestoreArray(metricInt, &met));
1698   }
1699 
1700   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0));
1701   PetscFunctionReturn(0);
1702 }
1703 
1704 /*@
1705   DMPlexMetricIntersection2 - Compute the intersection of two metrics
1706 
1707   Input Parameters:
1708 + dm        - The DM
1709 . metric1   - The first metric to be intersected
1710 - metric2   - The second metric to be intersected
1711 
1712   Output Parameter:
1713 . metricInt - The intersected metric
1714 
1715   Level: beginner
1716 
1717 .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1718 @*/
1719 PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1720 {
1721   Vec            metrics[2] = {metric1, metric2};
1722 
1723   PetscFunctionBegin;
1724   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
1725   PetscFunctionReturn(0);
1726 }
1727 
1728 /*@
1729   DMPlexMetricIntersection3 - Compute the intersection of three metrics
1730 
1731   Input Parameters:
1732 + dm        - The DM
1733 . metric1   - The first metric to be intersected
1734 . metric2   - The second metric to be intersected
1735 - metric3   - The third metric to be intersected
1736 
1737   Output Parameter:
1738 . metricInt - The intersected metric
1739 
1740   Level: beginner
1741 
1742 .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1743 @*/
1744 PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1745 {
1746   Vec            metrics[3] = {metric1, metric2, metric3};
1747 
1748   PetscFunctionBegin;
1749   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
1750   PetscFunctionReturn(0);
1751 }
1752