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