xref: /petsc/src/dm/impls/plex/plexmetric.c (revision 93520041193d7b76aece9b436de8c2aa14914128)
120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
220282da2SJoe Wallwork #include <petscblaslapack.h>
320282da2SJoe Wallwork 
431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
531b70465SJoe Wallwork {
631b70465SJoe Wallwork   MPI_Comm       comm;
7*93520041SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
8cc2eee55SJoe Wallwork   PetscBool      noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE;
931b70465SJoe Wallwork   PetscErrorCode ierr;
10cc2eee55SJoe Wallwork   PetscInt       verbosity = -1, numIter = 3;
11cc2eee55SJoe Wallwork   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;
1231b70465SJoe Wallwork 
1331b70465SJoe Wallwork   PetscFunctionBegin;
1431b70465SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1531b70465SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr);
1631b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr);
17cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, isotropic);CHKERRQ(ierr);
18*93520041SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL);CHKERRQ(ierr);
19*93520041SJoe Wallwork   ierr = DMPlexMetricSetUniform(dm, uniform);CHKERRQ(ierr);
2031b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr);
21cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst);CHKERRQ(ierr);
22cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL);CHKERRQ(ierr);
23cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoInsertion(dm, noInsert);CHKERRQ(ierr);
24cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL);CHKERRQ(ierr);
25cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoSwapping(dm, noSwap);CHKERRQ(ierr);
26cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL);CHKERRQ(ierr);
27cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoMovement(dm, noMove);CHKERRQ(ierr);
28cc2eee55SJoe Wallwork   ierr = PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0);CHKERRQ(ierr);
29cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNumIterations(dm, numIter);CHKERRQ(ierr);
30cc2eee55SJoe Wallwork   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);
31cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetVerbosity(dm, verbosity);CHKERRQ(ierr);
3231b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr);
3331b70465SJoe Wallwork   ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr);
3431b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr);
3531b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr);
3631b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr);
3731b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr);
3831b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr);
3931b70465SJoe Wallwork   ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr);
4031b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr);
4131b70465SJoe Wallwork   ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr);
42cc2eee55SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL);CHKERRQ(ierr);
43cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetGradationFactor(dm, beta);CHKERRQ(ierr);
4431b70465SJoe Wallwork   ierr = PetscOptionsEnd();CHKERRQ(ierr);
4531b70465SJoe Wallwork   PetscFunctionReturn(0);
4631b70465SJoe Wallwork }
4731b70465SJoe Wallwork 
48cb7bfe8cSJoe Wallwork /*@
4931b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
5031b70465SJoe Wallwork 
5131b70465SJoe Wallwork   Input parameters:
5231b70465SJoe Wallwork + dm        - The DM
5331b70465SJoe Wallwork - isotropic - Is the metric isotropic?
5431b70465SJoe Wallwork 
5531b70465SJoe Wallwork   Level: beginner
5631b70465SJoe Wallwork 
57*93520041SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst()
58cb7bfe8cSJoe Wallwork @*/
5931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
6031b70465SJoe Wallwork {
6131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6231b70465SJoe Wallwork   PetscErrorCode ierr;
6331b70465SJoe Wallwork 
6431b70465SJoe Wallwork   PetscFunctionBegin;
6531b70465SJoe Wallwork   if (!plex->metricCtx) {
6631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
6731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
6831b70465SJoe Wallwork   }
6931b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
7031b70465SJoe Wallwork   PetscFunctionReturn(0);
7131b70465SJoe Wallwork }
7231b70465SJoe Wallwork 
73cb7bfe8cSJoe Wallwork /*@
74*93520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7531b70465SJoe Wallwork 
7631b70465SJoe Wallwork   Input parameters:
7731b70465SJoe Wallwork . dm        - The DM
7831b70465SJoe Wallwork 
7931b70465SJoe Wallwork   Output parameters:
8031b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8131b70465SJoe Wallwork 
8231b70465SJoe Wallwork   Level: beginner
8331b70465SJoe Wallwork 
84*93520041SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst()
85cb7bfe8cSJoe Wallwork @*/
8631b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
8731b70465SJoe Wallwork {
8831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
8931b70465SJoe Wallwork   PetscErrorCode ierr;
9031b70465SJoe Wallwork 
9131b70465SJoe Wallwork   PetscFunctionBegin;
9231b70465SJoe Wallwork   if (!plex->metricCtx) {
9331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
9431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
9531b70465SJoe Wallwork   }
9631b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
9731b70465SJoe Wallwork   PetscFunctionReturn(0);
9831b70465SJoe Wallwork }
9931b70465SJoe Wallwork 
100cb7bfe8cSJoe Wallwork /*@
101*93520041SJoe Wallwork   DMPlexMetricSetUniform - Record whether a metric is uniform
102*93520041SJoe Wallwork 
103*93520041SJoe Wallwork   Input parameters:
104*93520041SJoe Wallwork + dm      - The DM
105*93520041SJoe Wallwork - uniform - Is the metric uniform?
106*93520041SJoe Wallwork 
107*93520041SJoe Wallwork   Level: beginner
108*93520041SJoe Wallwork 
109*93520041SJoe Wallwork   Notes:
110*93520041SJoe Wallwork 
111*93520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
112*93520041SJoe Wallwork 
113*93520041SJoe Wallwork .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
114*93520041SJoe Wallwork @*/
115*93520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
116*93520041SJoe Wallwork {
117*93520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
118*93520041SJoe Wallwork   PetscErrorCode ierr;
119*93520041SJoe Wallwork 
120*93520041SJoe Wallwork   PetscFunctionBegin;
121*93520041SJoe Wallwork   if (!plex->metricCtx) {
122*93520041SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
123*93520041SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
124*93520041SJoe Wallwork   }
125*93520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
126*93520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
127*93520041SJoe Wallwork   PetscFunctionReturn(0);
128*93520041SJoe Wallwork }
129*93520041SJoe Wallwork 
130*93520041SJoe Wallwork /*@
131*93520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
132*93520041SJoe Wallwork 
133*93520041SJoe Wallwork   Input parameters:
134*93520041SJoe Wallwork . dm      - The DM
135*93520041SJoe Wallwork 
136*93520041SJoe Wallwork   Output parameters:
137*93520041SJoe Wallwork . uniform - Is the metric uniform?
138*93520041SJoe Wallwork 
139*93520041SJoe Wallwork   Level: beginner
140*93520041SJoe Wallwork 
141*93520041SJoe Wallwork .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
142*93520041SJoe Wallwork @*/
143*93520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
144*93520041SJoe Wallwork {
145*93520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
146*93520041SJoe Wallwork   PetscErrorCode ierr;
147*93520041SJoe Wallwork 
148*93520041SJoe Wallwork   PetscFunctionBegin;
149*93520041SJoe Wallwork   if (!plex->metricCtx) {
150*93520041SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
151*93520041SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
152*93520041SJoe Wallwork   }
153*93520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
154*93520041SJoe Wallwork   PetscFunctionReturn(0);
155*93520041SJoe Wallwork }
156*93520041SJoe Wallwork 
157*93520041SJoe Wallwork /*@
15831b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
15931b70465SJoe Wallwork 
16031b70465SJoe Wallwork   Input parameters:
16131b70465SJoe Wallwork + dm                      - The DM
16231b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
16331b70465SJoe Wallwork 
16431b70465SJoe Wallwork   Level: beginner
16531b70465SJoe Wallwork 
16631b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
167cb7bfe8cSJoe Wallwork @*/
16831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
16931b70465SJoe Wallwork {
17031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
17131b70465SJoe Wallwork   PetscErrorCode ierr;
17231b70465SJoe Wallwork 
17331b70465SJoe Wallwork   PetscFunctionBegin;
17431b70465SJoe Wallwork   if (!plex->metricCtx) {
17531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
17631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
17731b70465SJoe Wallwork   }
17831b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
17931b70465SJoe Wallwork   PetscFunctionReturn(0);
18031b70465SJoe Wallwork }
18131b70465SJoe Wallwork 
182cb7bfe8cSJoe Wallwork /*@
18331b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
18431b70465SJoe Wallwork 
18531b70465SJoe Wallwork   Input parameters:
18631b70465SJoe Wallwork . dm                      - The DM
18731b70465SJoe Wallwork 
18831b70465SJoe Wallwork   Output parameters:
18931b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
19031b70465SJoe Wallwork 
19131b70465SJoe Wallwork   Level: beginner
19231b70465SJoe Wallwork 
19331b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
194cb7bfe8cSJoe Wallwork @*/
19531b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
19631b70465SJoe Wallwork {
19731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
19831b70465SJoe Wallwork   PetscErrorCode ierr;
19931b70465SJoe Wallwork 
20031b70465SJoe Wallwork   PetscFunctionBegin;
20131b70465SJoe Wallwork   if (!plex->metricCtx) {
20231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
20331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
20431b70465SJoe Wallwork   }
20531b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
20631b70465SJoe Wallwork   PetscFunctionReturn(0);
20731b70465SJoe Wallwork }
20831b70465SJoe Wallwork 
209cb7bfe8cSJoe Wallwork /*@
210cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
211cc2eee55SJoe Wallwork 
212cc2eee55SJoe Wallwork   Input parameters:
213cc2eee55SJoe Wallwork + dm       - The DM
214cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
215cc2eee55SJoe Wallwork 
216cc2eee55SJoe Wallwork   Level: beginner
217cc2eee55SJoe Wallwork 
218cb7bfe8cSJoe Wallwork   Notes:
219cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
220cb7bfe8cSJoe Wallwork 
221cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement()
222cb7bfe8cSJoe Wallwork @*/
223cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
224cc2eee55SJoe Wallwork {
225cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
226cc2eee55SJoe Wallwork   PetscErrorCode ierr;
227cc2eee55SJoe Wallwork 
228cc2eee55SJoe Wallwork   PetscFunctionBegin;
229cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
230cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
231cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
232cc2eee55SJoe Wallwork   }
233cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
234cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
235cc2eee55SJoe Wallwork }
236cc2eee55SJoe Wallwork 
237cb7bfe8cSJoe Wallwork /*@
238cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
239cc2eee55SJoe Wallwork 
240cc2eee55SJoe Wallwork   Input parameters:
241cc2eee55SJoe Wallwork . dm       - The DM
242cc2eee55SJoe Wallwork 
243cc2eee55SJoe Wallwork   Output parameters:
244cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
245cc2eee55SJoe Wallwork 
246cc2eee55SJoe Wallwork   Level: beginner
247cc2eee55SJoe Wallwork 
248cb7bfe8cSJoe Wallwork   Notes:
249cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
250cb7bfe8cSJoe Wallwork 
251cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement()
252cb7bfe8cSJoe Wallwork @*/
253cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
254cc2eee55SJoe Wallwork {
255cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
256cc2eee55SJoe Wallwork   PetscErrorCode ierr;
257cc2eee55SJoe Wallwork 
258cc2eee55SJoe Wallwork   PetscFunctionBegin;
259cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
260cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
261cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
262cc2eee55SJoe Wallwork   }
263cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
264cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
265cc2eee55SJoe Wallwork }
266cc2eee55SJoe Wallwork 
267cb7bfe8cSJoe Wallwork /*@
268cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
269cc2eee55SJoe Wallwork 
270cc2eee55SJoe Wallwork   Input parameters:
271cc2eee55SJoe Wallwork + dm     - The DM
272cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
273cc2eee55SJoe Wallwork 
274cc2eee55SJoe Wallwork   Level: beginner
275cc2eee55SJoe Wallwork 
276cb7bfe8cSJoe Wallwork   Notes:
277cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
278cb7bfe8cSJoe Wallwork 
279cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement()
280cb7bfe8cSJoe Wallwork @*/
281cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
282cc2eee55SJoe Wallwork {
283cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
284cc2eee55SJoe Wallwork   PetscErrorCode ierr;
285cc2eee55SJoe Wallwork 
286cc2eee55SJoe Wallwork   PetscFunctionBegin;
287cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
288cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
289cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
290cc2eee55SJoe Wallwork   }
291cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
292cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
293cc2eee55SJoe Wallwork }
294cc2eee55SJoe Wallwork 
295cb7bfe8cSJoe Wallwork /*@
296cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
297cc2eee55SJoe Wallwork 
298cc2eee55SJoe Wallwork   Input parameters:
299cc2eee55SJoe Wallwork . dm     - The DM
300cc2eee55SJoe Wallwork 
301cc2eee55SJoe Wallwork   Output parameters:
302cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
303cc2eee55SJoe Wallwork 
304cc2eee55SJoe Wallwork   Level: beginner
305cc2eee55SJoe Wallwork 
306cb7bfe8cSJoe Wallwork   Notes:
307cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
308cb7bfe8cSJoe Wallwork 
309cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement()
310cb7bfe8cSJoe Wallwork @*/
311cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
312cc2eee55SJoe Wallwork {
313cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
314cc2eee55SJoe Wallwork   PetscErrorCode ierr;
315cc2eee55SJoe Wallwork 
316cc2eee55SJoe Wallwork   PetscFunctionBegin;
317cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
318cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
319cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
320cc2eee55SJoe Wallwork   }
321cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
322cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
323cc2eee55SJoe Wallwork }
324cc2eee55SJoe Wallwork 
325cb7bfe8cSJoe Wallwork /*@
326cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
327cc2eee55SJoe Wallwork 
328cc2eee55SJoe Wallwork   Input parameters:
329cc2eee55SJoe Wallwork + dm     - The DM
330cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
331cc2eee55SJoe Wallwork 
332cc2eee55SJoe Wallwork   Level: beginner
333cc2eee55SJoe Wallwork 
334cb7bfe8cSJoe Wallwork   Notes:
335cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
336cb7bfe8cSJoe Wallwork 
337cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping()
338cb7bfe8cSJoe Wallwork @*/
339cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
340cc2eee55SJoe Wallwork {
341cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
342cc2eee55SJoe Wallwork   PetscErrorCode ierr;
343cc2eee55SJoe Wallwork 
344cc2eee55SJoe Wallwork   PetscFunctionBegin;
345cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
346cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
347cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
348cc2eee55SJoe Wallwork   }
349cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
350cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
351cc2eee55SJoe Wallwork }
352cc2eee55SJoe Wallwork 
353cb7bfe8cSJoe Wallwork /*@
354cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
355cc2eee55SJoe Wallwork 
356cc2eee55SJoe Wallwork   Input parameters:
357cc2eee55SJoe Wallwork . dm     - The DM
358cc2eee55SJoe Wallwork 
359cc2eee55SJoe Wallwork   Output parameters:
360cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
361cc2eee55SJoe Wallwork 
362cc2eee55SJoe Wallwork   Level: beginner
363cc2eee55SJoe Wallwork 
364cb7bfe8cSJoe Wallwork   Notes:
365cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
366cb7bfe8cSJoe Wallwork 
367cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping()
368cb7bfe8cSJoe Wallwork @*/
369cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
370cc2eee55SJoe Wallwork {
371cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
372cc2eee55SJoe Wallwork   PetscErrorCode ierr;
373cc2eee55SJoe Wallwork 
374cc2eee55SJoe Wallwork   PetscFunctionBegin;
375cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
376cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
377cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
378cc2eee55SJoe Wallwork   }
379cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
380cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
381cc2eee55SJoe Wallwork }
382cc2eee55SJoe Wallwork 
383cb7bfe8cSJoe Wallwork /*@
38431b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
38531b70465SJoe Wallwork 
38631b70465SJoe Wallwork   Input parameters:
38731b70465SJoe Wallwork + dm    - The DM
38831b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
38931b70465SJoe Wallwork 
39031b70465SJoe Wallwork   Level: beginner
39131b70465SJoe Wallwork 
39231b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
393cb7bfe8cSJoe Wallwork @*/
39431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
39531b70465SJoe Wallwork {
39631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
39731b70465SJoe Wallwork   PetscErrorCode ierr;
39831b70465SJoe Wallwork 
39931b70465SJoe Wallwork   PetscFunctionBegin;
40031b70465SJoe Wallwork   if (!plex->metricCtx) {
40131b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
40231b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
40331b70465SJoe Wallwork   }
40431b70465SJoe Wallwork   if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min);
40531b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
40631b70465SJoe Wallwork   PetscFunctionReturn(0);
40731b70465SJoe Wallwork }
40831b70465SJoe Wallwork 
409cb7bfe8cSJoe Wallwork /*@
41031b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
41131b70465SJoe Wallwork 
41231b70465SJoe Wallwork   Input parameters:
41331b70465SJoe Wallwork . dm    - The DM
41431b70465SJoe Wallwork 
415cc2eee55SJoe Wallwork   Output parameters:
41631b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
41731b70465SJoe Wallwork 
41831b70465SJoe Wallwork   Level: beginner
41931b70465SJoe Wallwork 
42031b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
421cb7bfe8cSJoe Wallwork @*/
42231b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
42331b70465SJoe Wallwork {
42431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
42531b70465SJoe Wallwork   PetscErrorCode ierr;
42631b70465SJoe Wallwork 
42731b70465SJoe Wallwork   PetscFunctionBegin;
42831b70465SJoe Wallwork   if (!plex->metricCtx) {
42931b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
43031b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
43131b70465SJoe Wallwork   }
43231b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
43331b70465SJoe Wallwork   PetscFunctionReturn(0);
43431b70465SJoe Wallwork }
43531b70465SJoe Wallwork 
436cb7bfe8cSJoe Wallwork /*@
43731b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
43831b70465SJoe Wallwork 
43931b70465SJoe Wallwork   Input parameters:
44031b70465SJoe Wallwork + dm    - The DM
44131b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44231b70465SJoe Wallwork 
44331b70465SJoe Wallwork   Level: beginner
44431b70465SJoe Wallwork 
44531b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
446cb7bfe8cSJoe Wallwork @*/
44731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
44831b70465SJoe Wallwork {
44931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
45031b70465SJoe Wallwork   PetscErrorCode ierr;
45131b70465SJoe Wallwork 
45231b70465SJoe Wallwork   PetscFunctionBegin;
45331b70465SJoe Wallwork   if (!plex->metricCtx) {
45431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
45531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
45631b70465SJoe Wallwork   }
45731b70465SJoe Wallwork   if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max);
45831b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
45931b70465SJoe Wallwork   PetscFunctionReturn(0);
46031b70465SJoe Wallwork }
46131b70465SJoe Wallwork 
462cb7bfe8cSJoe Wallwork /*@
46331b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
46431b70465SJoe Wallwork 
46531b70465SJoe Wallwork   Input parameters:
46631b70465SJoe Wallwork . dm    - The DM
46731b70465SJoe Wallwork 
468cc2eee55SJoe Wallwork   Output parameters:
46931b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
47031b70465SJoe Wallwork 
47131b70465SJoe Wallwork   Level: beginner
47231b70465SJoe Wallwork 
47331b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
474cb7bfe8cSJoe Wallwork @*/
47531b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
47631b70465SJoe Wallwork {
47731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
47831b70465SJoe Wallwork   PetscErrorCode ierr;
47931b70465SJoe Wallwork 
48031b70465SJoe Wallwork   PetscFunctionBegin;
48131b70465SJoe Wallwork   if (!plex->metricCtx) {
48231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
48331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
48431b70465SJoe Wallwork   }
48531b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
48631b70465SJoe Wallwork   PetscFunctionReturn(0);
48731b70465SJoe Wallwork }
48831b70465SJoe Wallwork 
489cb7bfe8cSJoe Wallwork /*@
49031b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
49131b70465SJoe Wallwork 
49231b70465SJoe Wallwork   Input parameters:
49331b70465SJoe Wallwork + dm    - The DM
49431b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
49531b70465SJoe Wallwork 
49631b70465SJoe Wallwork   Level: beginner
49731b70465SJoe Wallwork 
49831b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
49931b70465SJoe Wallwork 
50031b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
501cb7bfe8cSJoe Wallwork @*/
50231b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
50331b70465SJoe Wallwork {
50431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
50531b70465SJoe Wallwork   PetscErrorCode ierr;
50631b70465SJoe Wallwork 
50731b70465SJoe Wallwork   PetscFunctionBegin;
50831b70465SJoe Wallwork   if (!plex->metricCtx) {
50931b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
51031b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
51131b70465SJoe Wallwork   }
51231b70465SJoe Wallwork   if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max);
51331b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
51431b70465SJoe Wallwork   PetscFunctionReturn(0);
51531b70465SJoe Wallwork }
51631b70465SJoe Wallwork 
517cb7bfe8cSJoe Wallwork /*@
51831b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
51931b70465SJoe Wallwork 
52031b70465SJoe Wallwork   Input parameters:
52131b70465SJoe Wallwork . dm    - The DM
52231b70465SJoe Wallwork 
523cc2eee55SJoe Wallwork   Output parameters:
52431b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
52531b70465SJoe Wallwork 
52631b70465SJoe Wallwork   Level: beginner
52731b70465SJoe Wallwork 
52831b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
529cb7bfe8cSJoe Wallwork @*/
53031b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
53131b70465SJoe Wallwork {
53231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
53331b70465SJoe Wallwork   PetscErrorCode ierr;
53431b70465SJoe Wallwork 
53531b70465SJoe Wallwork   PetscFunctionBegin;
53631b70465SJoe Wallwork   if (!plex->metricCtx) {
53731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
53831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
53931b70465SJoe Wallwork   }
54031b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
54131b70465SJoe Wallwork   PetscFunctionReturn(0);
54231b70465SJoe Wallwork }
54331b70465SJoe Wallwork 
544cb7bfe8cSJoe Wallwork /*@
54531b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
54631b70465SJoe Wallwork 
54731b70465SJoe Wallwork   Input parameters:
54831b70465SJoe Wallwork + dm               - The DM
54931b70465SJoe Wallwork - targetComplexity - The target metric complexity
55031b70465SJoe Wallwork 
55131b70465SJoe Wallwork   Level: beginner
55231b70465SJoe Wallwork 
55331b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
554cb7bfe8cSJoe Wallwork @*/
55531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
55631b70465SJoe Wallwork {
55731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
55831b70465SJoe Wallwork   PetscErrorCode ierr;
55931b70465SJoe Wallwork 
56031b70465SJoe Wallwork   PetscFunctionBegin;
56131b70465SJoe Wallwork   if (!plex->metricCtx) {
56231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
56331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
56431b70465SJoe Wallwork   }
56531b70465SJoe Wallwork   if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive");
56631b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
56731b70465SJoe Wallwork   PetscFunctionReturn(0);
56831b70465SJoe Wallwork }
56931b70465SJoe Wallwork 
570cb7bfe8cSJoe Wallwork /*@
57131b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
57231b70465SJoe Wallwork 
57331b70465SJoe Wallwork   Input parameters:
57431b70465SJoe Wallwork . dm               - The DM
57531b70465SJoe Wallwork 
576cc2eee55SJoe Wallwork   Output parameters:
57731b70465SJoe Wallwork . targetComplexity - The target metric complexity
57831b70465SJoe Wallwork 
57931b70465SJoe Wallwork   Level: beginner
58031b70465SJoe Wallwork 
58131b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
582cb7bfe8cSJoe Wallwork @*/
58331b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
58431b70465SJoe Wallwork {
58531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
58631b70465SJoe Wallwork   PetscErrorCode ierr;
58731b70465SJoe Wallwork 
58831b70465SJoe Wallwork   PetscFunctionBegin;
58931b70465SJoe Wallwork   if (!plex->metricCtx) {
59031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
59131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
59231b70465SJoe Wallwork   }
59331b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
59431b70465SJoe Wallwork   PetscFunctionReturn(0);
59531b70465SJoe Wallwork }
59631b70465SJoe Wallwork 
597cb7bfe8cSJoe Wallwork /*@
59831b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
59931b70465SJoe Wallwork 
60031b70465SJoe Wallwork   Input parameters:
60131b70465SJoe Wallwork + dm - The DM
60231b70465SJoe Wallwork - p  - The normalization order
60331b70465SJoe Wallwork 
60431b70465SJoe Wallwork   Level: beginner
60531b70465SJoe Wallwork 
60631b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
607cb7bfe8cSJoe Wallwork @*/
60831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
60931b70465SJoe Wallwork {
61031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
61131b70465SJoe Wallwork   PetscErrorCode ierr;
61231b70465SJoe Wallwork 
61331b70465SJoe Wallwork   PetscFunctionBegin;
61431b70465SJoe Wallwork   if (!plex->metricCtx) {
61531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
61631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
61731b70465SJoe Wallwork   }
61831b70465SJoe Wallwork   if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater");
61931b70465SJoe Wallwork   plex->metricCtx->p = p;
62031b70465SJoe Wallwork   PetscFunctionReturn(0);
62131b70465SJoe Wallwork }
62231b70465SJoe Wallwork 
623cb7bfe8cSJoe Wallwork /*@
62431b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
62531b70465SJoe Wallwork 
62631b70465SJoe Wallwork   Input parameters:
62731b70465SJoe Wallwork . dm - The DM
62831b70465SJoe Wallwork 
629cc2eee55SJoe Wallwork   Output parameters:
63031b70465SJoe Wallwork . p - The normalization order
63131b70465SJoe Wallwork 
63231b70465SJoe Wallwork   Level: beginner
63331b70465SJoe Wallwork 
63431b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
635cb7bfe8cSJoe Wallwork @*/
63631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
63731b70465SJoe Wallwork {
63831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
63931b70465SJoe Wallwork   PetscErrorCode ierr;
64031b70465SJoe Wallwork 
64131b70465SJoe Wallwork   PetscFunctionBegin;
64231b70465SJoe Wallwork   if (!plex->metricCtx) {
64331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
64431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
64531b70465SJoe Wallwork   }
64631b70465SJoe Wallwork   *p = plex->metricCtx->p;
64731b70465SJoe Wallwork   PetscFunctionReturn(0);
64831b70465SJoe Wallwork }
64931b70465SJoe Wallwork 
650cb7bfe8cSJoe Wallwork /*@
651cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
652cc2eee55SJoe Wallwork 
653cc2eee55SJoe Wallwork   Input parameters:
654cc2eee55SJoe Wallwork + dm   - The DM
655cc2eee55SJoe Wallwork - beta - The metric gradation factor
656cc2eee55SJoe Wallwork 
657cc2eee55SJoe Wallwork   Level: beginner
658cc2eee55SJoe Wallwork 
659cc2eee55SJoe Wallwork   Notes:
660cc2eee55SJoe Wallwork 
661cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
662cc2eee55SJoe Wallwork 
663cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
664cc2eee55SJoe Wallwork 
665cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
666cb7bfe8cSJoe Wallwork 
667cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetGradationFactor()
668cb7bfe8cSJoe Wallwork @*/
669cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
670cc2eee55SJoe Wallwork {
671cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
672cc2eee55SJoe Wallwork   PetscErrorCode ierr;
673cc2eee55SJoe Wallwork 
674cc2eee55SJoe Wallwork   PetscFunctionBegin;
675cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
676cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
677cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
678cc2eee55SJoe Wallwork   }
679cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
680cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
681cc2eee55SJoe Wallwork }
682cc2eee55SJoe Wallwork 
683cb7bfe8cSJoe Wallwork /*@
684cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
685cc2eee55SJoe Wallwork 
686cc2eee55SJoe Wallwork   Input parameters:
687cc2eee55SJoe Wallwork . dm   - The DM
688cc2eee55SJoe Wallwork 
689cc2eee55SJoe Wallwork   Output parameters:
690cc2eee55SJoe Wallwork . beta - The metric gradation factor
691cc2eee55SJoe Wallwork 
692cc2eee55SJoe Wallwork   Level: beginner
693cc2eee55SJoe Wallwork 
694cb7bfe8cSJoe Wallwork   Notes:
695cb7bfe8cSJoe Wallwork 
696cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
697cb7bfe8cSJoe Wallwork 
698cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
699cb7bfe8cSJoe Wallwork 
700cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
701cc2eee55SJoe Wallwork 
702cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetGradationFactor()
703cb7bfe8cSJoe Wallwork @*/
704cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
705cc2eee55SJoe Wallwork {
706cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
707cc2eee55SJoe Wallwork   PetscErrorCode ierr;
708cc2eee55SJoe Wallwork 
709cc2eee55SJoe Wallwork   PetscFunctionBegin;
710cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
711cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
712cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
713cc2eee55SJoe Wallwork   }
714cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
715cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
716cc2eee55SJoe Wallwork }
717cc2eee55SJoe Wallwork 
718cb7bfe8cSJoe Wallwork /*@
719cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
720cc2eee55SJoe Wallwork 
721cc2eee55SJoe Wallwork   Input parameters:
722cc2eee55SJoe Wallwork + dm        - The DM
723cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
724cc2eee55SJoe Wallwork 
725cb7bfe8cSJoe Wallwork   Level: beginner
726cb7bfe8cSJoe Wallwork 
727cb7bfe8cSJoe Wallwork   Notes:
728cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
729cb7bfe8cSJoe Wallwork 
730cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations()
731cb7bfe8cSJoe Wallwork @*/
732cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
733cc2eee55SJoe Wallwork {
734cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
735cc2eee55SJoe Wallwork   PetscErrorCode ierr;
736cc2eee55SJoe Wallwork 
737cc2eee55SJoe Wallwork   PetscFunctionBegin;
738cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
739cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
740cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
741cc2eee55SJoe Wallwork   }
742cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
743cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
744cc2eee55SJoe Wallwork }
745cc2eee55SJoe Wallwork 
746cb7bfe8cSJoe Wallwork /*@
747cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
748cc2eee55SJoe Wallwork 
749cc2eee55SJoe Wallwork   Input parameters:
750cc2eee55SJoe Wallwork . dm        - The DM
751cc2eee55SJoe Wallwork 
752cc2eee55SJoe Wallwork   Output parameters:
753cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
754cc2eee55SJoe Wallwork 
755cb7bfe8cSJoe Wallwork   Level: beginner
756cb7bfe8cSJoe Wallwork 
757cb7bfe8cSJoe Wallwork   Notes:
758cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
759cb7bfe8cSJoe Wallwork 
760cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
761cb7bfe8cSJoe Wallwork @*/
762cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
763cc2eee55SJoe Wallwork {
764cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
765cc2eee55SJoe Wallwork   PetscErrorCode ierr;
766cc2eee55SJoe Wallwork 
767cc2eee55SJoe Wallwork   PetscFunctionBegin;
768cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
769cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
770cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
771cc2eee55SJoe Wallwork   }
772cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
773cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
774cc2eee55SJoe Wallwork }
775cc2eee55SJoe Wallwork 
776cb7bfe8cSJoe Wallwork /*@
777cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
778cc2eee55SJoe Wallwork 
779cc2eee55SJoe Wallwork   Input parameters:
780cc2eee55SJoe Wallwork + dm      - The DM
781cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
782cc2eee55SJoe Wallwork 
783cb7bfe8cSJoe Wallwork   Level: beginner
784cb7bfe8cSJoe Wallwork 
785cb7bfe8cSJoe Wallwork   Notes:
786cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
787cc2eee55SJoe Wallwork 
788cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
789cb7bfe8cSJoe Wallwork @*/
790cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
791cc2eee55SJoe Wallwork {
792cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
793cc2eee55SJoe Wallwork   PetscErrorCode ierr;
794cc2eee55SJoe Wallwork 
795cc2eee55SJoe Wallwork   PetscFunctionBegin;
796cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
797cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
798cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
799cc2eee55SJoe Wallwork   }
800cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
801cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
802cc2eee55SJoe Wallwork }
803cc2eee55SJoe Wallwork 
804cb7bfe8cSJoe Wallwork /*@
805cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
806cc2eee55SJoe Wallwork 
807cc2eee55SJoe Wallwork   Input parameters:
808cc2eee55SJoe Wallwork . dm      - The DM
809cc2eee55SJoe Wallwork 
810cc2eee55SJoe Wallwork   Output parameters:
811cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
812cc2eee55SJoe Wallwork 
813cb7bfe8cSJoe Wallwork   Level: beginner
814cb7bfe8cSJoe Wallwork 
815cb7bfe8cSJoe Wallwork   Notes:
816cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
817cc2eee55SJoe Wallwork 
818cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity()
819cb7bfe8cSJoe Wallwork @*/
820cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
821cc2eee55SJoe Wallwork {
822cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
823cc2eee55SJoe Wallwork   PetscErrorCode ierr;
824cc2eee55SJoe Wallwork 
825cc2eee55SJoe Wallwork   PetscFunctionBegin;
826cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
827cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
828cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
829cc2eee55SJoe Wallwork   }
830cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
831cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
832cc2eee55SJoe Wallwork }
833cc2eee55SJoe Wallwork 
83420282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
83520282da2SJoe Wallwork {
83620282da2SJoe Wallwork   MPI_Comm       comm;
83720282da2SJoe Wallwork   PetscErrorCode ierr;
83820282da2SJoe Wallwork   PetscFE        fe;
83920282da2SJoe Wallwork   PetscInt       dim;
84020282da2SJoe Wallwork 
84120282da2SJoe Wallwork   PetscFunctionBegin;
84220282da2SJoe Wallwork 
84320282da2SJoe Wallwork   /* Extract metadata from dm */
84420282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
84520282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
84620282da2SJoe Wallwork 
84720282da2SJoe Wallwork   /* Create a P1 field of the requested size */
84820282da2SJoe Wallwork   ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
84920282da2SJoe Wallwork   ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr);
85020282da2SJoe Wallwork   ierr = DMCreateDS(dm);CHKERRQ(ierr);
85120282da2SJoe Wallwork   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
85220282da2SJoe Wallwork   ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr);
85320282da2SJoe Wallwork 
85420282da2SJoe Wallwork   PetscFunctionReturn(0);
85520282da2SJoe Wallwork }
85620282da2SJoe Wallwork 
857cb7bfe8cSJoe Wallwork /*@
85820282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
85920282da2SJoe Wallwork 
86020282da2SJoe Wallwork   Input parameters:
86120282da2SJoe Wallwork + dm     - The DM
86220282da2SJoe Wallwork - f      - The field number to use
86320282da2SJoe Wallwork 
86420282da2SJoe Wallwork   Output parameter:
86520282da2SJoe Wallwork . metric - The metric
86620282da2SJoe Wallwork 
86720282da2SJoe Wallwork   Level: beginner
86820282da2SJoe Wallwork 
86931b70465SJoe Wallwork   Notes:
87031b70465SJoe Wallwork 
87131b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
87231b70465SJoe Wallwork 
87331b70465SJoe Wallwork   Command line options for Riemannian metrics:
87431b70465SJoe Wallwork 
875cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
876*93520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
877cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
878cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
879cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
880cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
881cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
882cb7bfe8cSJoe Wallwork . -dm_plex_metric_target_complexity         - Target metric complexity
883cb7bfe8cSJoe Wallwork 
884cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
885cb7bfe8cSJoe Wallwork 
886cb7bfe8cSJoe Wallwork . -dm_adaptor <pragmatic/mmg/parmmg>
887cb7bfe8cSJoe Wallwork 
888cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
889cb7bfe8cSJoe Wallwork 
890cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
891cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
892cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
893cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
894cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
895cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
89620282da2SJoe Wallwork 
89720282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
898cb7bfe8cSJoe Wallwork @*/
89920282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
90020282da2SJoe Wallwork {
90131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
902*93520041SJoe Wallwork   PetscBool      isotropic, uniform;
90320282da2SJoe Wallwork   PetscErrorCode ierr;
90420282da2SJoe Wallwork   PetscInt       coordDim, Nd;
90520282da2SJoe Wallwork 
90620282da2SJoe Wallwork   PetscFunctionBegin;
90731b70465SJoe Wallwork   if (!plex->metricCtx) {
90831b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
90931b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
91031b70465SJoe Wallwork   }
91120282da2SJoe Wallwork   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
91220282da2SJoe Wallwork   Nd = coordDim*coordDim;
913*93520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
914*93520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
915*93520041SJoe Wallwork   if (uniform) {
916*93520041SJoe Wallwork     MPI_Comm comm;
917*93520041SJoe Wallwork 
918*93520041SJoe Wallwork     ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
919*93520041SJoe Wallwork     if (!isotropic) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
920*93520041SJoe Wallwork     ierr = VecCreate(comm, metric);CHKERRQ(ierr);
921*93520041SJoe Wallwork     ierr = VecSetSizes(*metric, 1, PETSC_DECIDE);CHKERRQ(ierr);
922*93520041SJoe Wallwork     ierr = VecSetFromOptions(*metric);CHKERRQ(ierr);
923*93520041SJoe Wallwork   } else if (isotropic) {
924*93520041SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dm, f, 1, metric);CHKERRQ(ierr);
925*93520041SJoe Wallwork   } else {
92620282da2SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr);
927*93520041SJoe Wallwork   }
92820282da2SJoe Wallwork   PetscFunctionReturn(0);
92920282da2SJoe Wallwork }
93020282da2SJoe Wallwork 
931cb7bfe8cSJoe Wallwork /*@
93220282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
93320282da2SJoe Wallwork 
93420282da2SJoe Wallwork   Input parameters:
93520282da2SJoe Wallwork + dm     - The DM
93620282da2SJoe Wallwork . f      - The field number to use
93720282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
93820282da2SJoe Wallwork 
93920282da2SJoe Wallwork   Output parameter:
94020282da2SJoe Wallwork . metric - The uniform metric
94120282da2SJoe Wallwork 
94220282da2SJoe Wallwork   Level: beginner
94320282da2SJoe Wallwork 
944*93520041SJoe Wallwork   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
94520282da2SJoe Wallwork 
94620282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
947cb7bfe8cSJoe Wallwork @*/
94820282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
94920282da2SJoe Wallwork {
95020282da2SJoe Wallwork   PetscErrorCode ierr;
95120282da2SJoe Wallwork 
95220282da2SJoe Wallwork   PetscFunctionBegin;
953*93520041SJoe Wallwork   ierr = DMPlexMetricSetUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
95420282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
95520282da2SJoe Wallwork   if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
95620282da2SJoe Wallwork   if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
957*93520041SJoe Wallwork   ierr = VecSet(*metric, alpha);CHKERRQ(ierr);
958*93520041SJoe Wallwork   ierr = VecAssemblyBegin(*metric);CHKERRQ(ierr);
959*93520041SJoe Wallwork   ierr = VecAssemblyEnd(*metric);CHKERRQ(ierr);
96020282da2SJoe Wallwork   PetscFunctionReturn(0);
96120282da2SJoe Wallwork }
96220282da2SJoe Wallwork 
963*93520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
964*93520041SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
965*93520041SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
966*93520041SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
967*93520041SJoe Wallwork {
968*93520041SJoe Wallwork   f0[0] = u[0];
969*93520041SJoe Wallwork }
970*93520041SJoe Wallwork 
971cb7bfe8cSJoe Wallwork /*@
97220282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
97320282da2SJoe Wallwork 
97420282da2SJoe Wallwork   Input parameters:
97520282da2SJoe Wallwork + dm        - The DM
97620282da2SJoe Wallwork . f         - The field number to use
97720282da2SJoe Wallwork - indicator - The error indicator
97820282da2SJoe Wallwork 
97920282da2SJoe Wallwork   Output parameter:
98020282da2SJoe Wallwork . metric    - The isotropic metric
98120282da2SJoe Wallwork 
98220282da2SJoe Wallwork   Level: beginner
98320282da2SJoe Wallwork 
98420282da2SJoe Wallwork   Notes:
98520282da2SJoe Wallwork 
98620282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
98720282da2SJoe Wallwork 
988*93520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
98920282da2SJoe Wallwork 
99020282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
991cb7bfe8cSJoe Wallwork @*/
99220282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
99320282da2SJoe Wallwork {
99420282da2SJoe Wallwork   PetscErrorCode ierr;
995*93520041SJoe Wallwork   PetscInt       m, n;
99620282da2SJoe Wallwork 
99720282da2SJoe Wallwork   PetscFunctionBegin;
998*93520041SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, PETSC_TRUE);CHKERRQ(ierr);
99920282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
1000*93520041SJoe Wallwork   ierr = VecGetSize(indicator, &m);CHKERRQ(ierr);
1001*93520041SJoe Wallwork   ierr = VecGetSize(*metric, &n);CHKERRQ(ierr);
1002*93520041SJoe Wallwork   if (m == n) { ierr = VecCopy(indicator, *metric);CHKERRQ(ierr); }
1003*93520041SJoe Wallwork   else {
1004*93520041SJoe Wallwork     DM     dmIndi;
1005*93520041SJoe Wallwork     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
1006*93520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1007*93520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1008*93520041SJoe Wallwork                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
1009*93520041SJoe Wallwork 
101020282da2SJoe Wallwork     ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr);
1011*93520041SJoe Wallwork     funcs[0] = identity;
1012*93520041SJoe Wallwork     ierr = DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric);CHKERRQ(ierr);
101320282da2SJoe Wallwork   }
101420282da2SJoe Wallwork   PetscFunctionReturn(0);
101520282da2SJoe Wallwork }
101620282da2SJoe Wallwork 
101782490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
101882490253SJoe Wallwork {
101982490253SJoe Wallwork   PetscInt i, j;
102082490253SJoe Wallwork 
102182490253SJoe Wallwork   PetscFunctionBegin;
102282490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
102382490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
102482490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
102582490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
102682490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
102782490253SJoe Wallwork       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]);
102882490253SJoe Wallwork       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]);
102982490253SJoe Wallwork     }
103082490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
103182490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
103282490253SJoe Wallwork   }
103382490253SJoe Wallwork   PetscFunctionReturn(0);
103482490253SJoe Wallwork }
103582490253SJoe Wallwork 
10363f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
103720282da2SJoe Wallwork {
103820282da2SJoe Wallwork   PetscErrorCode ierr;
103920282da2SJoe Wallwork   PetscInt       i, j, k;
104020282da2SJoe Wallwork   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);
104120282da2SJoe Wallwork   PetscScalar   *Mpos;
104220282da2SJoe Wallwork 
104320282da2SJoe Wallwork   PetscFunctionBegin;
104420282da2SJoe Wallwork   ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr);
104520282da2SJoe Wallwork 
104620282da2SJoe Wallwork   /* Symmetrize */
104720282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
104820282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
104920282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
105020282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
105120282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
105220282da2SJoe Wallwork     }
105320282da2SJoe Wallwork   }
105420282da2SJoe Wallwork 
105520282da2SJoe Wallwork   /* Compute eigendecomposition */
1056*93520041SJoe Wallwork   if (dim == 1) {
1057*93520041SJoe Wallwork 
1058*93520041SJoe Wallwork     /* Isotropic case */
1059*93520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
1060*93520041SJoe Wallwork     Mpos[0] = 1.0;
1061*93520041SJoe Wallwork   } else {
1062*93520041SJoe Wallwork 
1063*93520041SJoe Wallwork     /* Anisotropic case */
106420282da2SJoe Wallwork     PetscScalar  *work;
106520282da2SJoe Wallwork     PetscBLASInt lwork;
106620282da2SJoe Wallwork 
106720282da2SJoe Wallwork     lwork = 5*dim;
106820282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
106920282da2SJoe Wallwork     {
107020282da2SJoe Wallwork       PetscBLASInt lierr;
107120282da2SJoe Wallwork       PetscBLASInt nb;
107220282da2SJoe Wallwork 
107320282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
107420282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
107520282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
107620282da2SJoe Wallwork       {
107720282da2SJoe Wallwork         PetscReal *rwork;
107820282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
107920282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
108020282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
108120282da2SJoe Wallwork       }
108220282da2SJoe Wallwork #else
108320282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
108420282da2SJoe Wallwork #endif
108582490253SJoe Wallwork       if (lierr) {
108682490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
108782490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
108882490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
108982490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
109082490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
109182490253SJoe Wallwork           }
109282490253SJoe Wallwork         }
109382490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
109482490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
109582490253SJoe Wallwork       }
109620282da2SJoe Wallwork       ierr = PetscFPTrapPop();CHKERRQ(ierr);
109720282da2SJoe Wallwork     }
109820282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
109920282da2SJoe Wallwork   }
110020282da2SJoe Wallwork 
110120282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
110220282da2SJoe Wallwork   max_eig = 0.0;
110320282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
110420282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
110520282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
110620282da2SJoe Wallwork   }
110720282da2SJoe Wallwork 
11083f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11093f07679eSJoe Wallwork   *detMp = 1.0;
111020282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111120282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
11123f07679eSJoe Wallwork     *detMp *= eigs[i];
111320282da2SJoe Wallwork   }
111420282da2SJoe Wallwork 
111520282da2SJoe Wallwork   /* Reconstruct Hessian */
111620282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111720282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
111820282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
111920282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
112020282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
112120282da2SJoe Wallwork       }
112220282da2SJoe Wallwork     }
112320282da2SJoe Wallwork   }
112420282da2SJoe Wallwork   ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr);
112520282da2SJoe Wallwork 
112620282da2SJoe Wallwork   PetscFunctionReturn(0);
112720282da2SJoe Wallwork }
112820282da2SJoe Wallwork 
1129cb7bfe8cSJoe Wallwork /*@
113020282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
113120282da2SJoe Wallwork 
113220282da2SJoe Wallwork   Input parameters:
113320282da2SJoe Wallwork + dm                 - The DM
11343f07679eSJoe Wallwork . metricIn           - The metric
113599abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11363f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
113720282da2SJoe Wallwork 
113820282da2SJoe Wallwork   Output parameter:
11393f07679eSJoe Wallwork + metricOut          - The metric
11403f07679eSJoe Wallwork - determinant        - Its determinant
114120282da2SJoe Wallwork 
114220282da2SJoe Wallwork   Level: beginner
114320282da2SJoe Wallwork 
1144cb7bfe8cSJoe Wallwork   Notes:
1145cb7bfe8cSJoe Wallwork 
1146cb7bfe8cSJoe Wallwork   Relevant command line options:
1147cb7bfe8cSJoe Wallwork 
1148*93520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
1149*93520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
1150*93520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1151cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1152cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1153cb7bfe8cSJoe Wallwork 
115420282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
1155cb7bfe8cSJoe Wallwork @*/
11563f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant)
115720282da2SJoe Wallwork {
11583f07679eSJoe Wallwork   DM             dmDet;
1159*93520041SJoe Wallwork   PetscBool      isotropic, uniform;
116020282da2SJoe Wallwork   PetscErrorCode ierr;
116120282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
11623f07679eSJoe Wallwork   PetscScalar   *met, *det;
116320282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
116420282da2SJoe Wallwork 
116520282da2SJoe Wallwork   PetscFunctionBegin;
116620282da2SJoe Wallwork 
116720282da2SJoe Wallwork   /* Extract metadata from dm */
116820282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
116920282da2SJoe Wallwork   if (restrictSizes) {
117099abec2bSJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
117199abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
117299abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
117399abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
117499abec2bSJoe Wallwork     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);
117599abec2bSJoe Wallwork   }
117699abec2bSJoe Wallwork   if (restrictAnisotropy) {
117799abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
117899abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
117920282da2SJoe Wallwork   }
118020282da2SJoe Wallwork 
1181*93520041SJoe Wallwork   /* Setup output metric */
11823f07679eSJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr);
11833f07679eSJoe Wallwork   ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr);
1184*93520041SJoe Wallwork 
1185*93520041SJoe Wallwork   /* Enforce SPD and extract determinant */
1186*93520041SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
1187*93520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
1188*93520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
1189*93520041SJoe Wallwork   if (uniform) {
1190*93520041SJoe Wallwork     if (!isotropic) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
1191*93520041SJoe Wallwork 
1192*93520041SJoe Wallwork     /* Uniform case */
1193*93520041SJoe Wallwork     ierr = VecDuplicate(metricIn, determinant);CHKERRQ(ierr);
1194*93520041SJoe Wallwork     ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
1195*93520041SJoe Wallwork     ierr = DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det);CHKERRQ(ierr);
1196*93520041SJoe Wallwork     ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
1197*93520041SJoe Wallwork   } else {
1198*93520041SJoe Wallwork 
1199*93520041SJoe Wallwork     /* Spatially varying case */
1200*93520041SJoe Wallwork     PetscInt nrow;
1201*93520041SJoe Wallwork 
1202*93520041SJoe Wallwork     if (isotropic) nrow = 1;
1203*93520041SJoe Wallwork     else nrow = dim;
12043f07679eSJoe Wallwork     ierr = DMClone(dm, &dmDet);CHKERRQ(ierr);
12053f07679eSJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant);CHKERRQ(ierr);
120620282da2SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12073f07679eSJoe Wallwork     ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
120820282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12093f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
121020282da2SJoe Wallwork       ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
12113f07679eSJoe Wallwork       ierr = DMPlexPointLocalRef(dmDet, v, det, &vdet);CHKERRQ(ierr);
1212*93520041SJoe Wallwork       ierr = DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet);CHKERRQ(ierr);
121320282da2SJoe Wallwork     }
12143f07679eSJoe Wallwork     ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
1215*93520041SJoe Wallwork   }
12163f07679eSJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
121720282da2SJoe Wallwork   PetscFunctionReturn(0);
121820282da2SJoe Wallwork }
121920282da2SJoe Wallwork 
122020282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
122120282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
122220282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122320282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
122420282da2SJoe Wallwork {
122520282da2SJoe Wallwork   const PetscScalar p = constants[0];
122620282da2SJoe Wallwork 
12273f07679eSJoe Wallwork   f0[0] = PetscPowReal(u[0], p/(2.0*p + dim));
122820282da2SJoe Wallwork }
122920282da2SJoe Wallwork 
1230cb7bfe8cSJoe Wallwork /*@
123120282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
123220282da2SJoe Wallwork 
123320282da2SJoe Wallwork   Input parameters:
123420282da2SJoe Wallwork + dm                 - The DM
123520282da2SJoe Wallwork . metricIn           - The unnormalized metric
123616522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
123716522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
123820282da2SJoe Wallwork 
123920282da2SJoe Wallwork   Output parameter:
124020282da2SJoe Wallwork . metricOut          - The normalized metric
124120282da2SJoe Wallwork 
124220282da2SJoe Wallwork   Level: beginner
124320282da2SJoe Wallwork 
1244cb7bfe8cSJoe Wallwork   Notes:
1245cb7bfe8cSJoe Wallwork 
1246cb7bfe8cSJoe Wallwork   Relevant command line options:
1247cb7bfe8cSJoe Wallwork 
1248*93520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
1249*93520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
1250*93520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1251cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1252cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1253cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1254cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1255cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1256cb7bfe8cSJoe Wallwork 
125720282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
1258cb7bfe8cSJoe Wallwork @*/
125916522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut)
126020282da2SJoe Wallwork {
12613f07679eSJoe Wallwork   DM               dmDet;
126220282da2SJoe Wallwork   MPI_Comm         comm;
1263*93520041SJoe Wallwork   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
126420282da2SJoe Wallwork   PetscDS          ds;
126520282da2SJoe Wallwork   PetscErrorCode   ierr;
126620282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
12673f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
1268*93520041SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
12693f07679eSJoe Wallwork   Vec              determinant;
127020282da2SJoe Wallwork 
127120282da2SJoe Wallwork   PetscFunctionBegin;
127220282da2SJoe Wallwork 
127320282da2SJoe Wallwork   /* Extract metadata from dm */
127420282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
127520282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1276*93520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
1277*93520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
1278*93520041SJoe Wallwork   if (isotropic) Nd = 1;
1279*93520041SJoe Wallwork   else Nd = dim*dim;
128020282da2SJoe Wallwork 
128120282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
128216522936SJoe Wallwork   ierr = DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst);CHKERRQ(ierr);
12833f07679eSJoe Wallwork   ierr = DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant);CHKERRQ(ierr);
128420282da2SJoe Wallwork 
128520282da2SJoe Wallwork   /* Compute global normalization factor */
128616522936SJoe Wallwork   ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
128716522936SJoe Wallwork   ierr = DMPlexMetricGetNormalizationOrder(dm, &p);CHKERRQ(ierr);
128816522936SJoe Wallwork   constants[0] = p;
1289*93520041SJoe Wallwork   if (uniform) {
1290*93520041SJoe Wallwork     if (!isotropic) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
1291*93520041SJoe Wallwork     DM  dmTmp;
1292*93520041SJoe Wallwork     Vec tmp;
1293*93520041SJoe Wallwork 
1294*93520041SJoe Wallwork     ierr = DMClone(dm, &dmTmp);CHKERRQ(ierr);
1295*93520041SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp);CHKERRQ(ierr);
1296*93520041SJoe Wallwork     ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
1297*93520041SJoe Wallwork     ierr = VecSet(tmp, det[0]);CHKERRQ(ierr);
1298*93520041SJoe Wallwork     ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
1299*93520041SJoe Wallwork     ierr = DMGetDS(dmTmp, &ds);CHKERRQ(ierr);
1300*93520041SJoe Wallwork     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
1301*93520041SJoe Wallwork     ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
1302*93520041SJoe Wallwork     ierr = DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL);CHKERRQ(ierr);
1303*93520041SJoe Wallwork     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
1304*93520041SJoe Wallwork     ierr = DMDestroy(&dmTmp);CHKERRQ(ierr);
1305*93520041SJoe Wallwork   } else {
13063f07679eSJoe Wallwork     ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr);
13073f07679eSJoe Wallwork     ierr = DMGetDS(dmDet, &ds);CHKERRQ(ierr);
130820282da2SJoe Wallwork     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
130920282da2SJoe Wallwork     ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
13103f07679eSJoe Wallwork     ierr = DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL);CHKERRQ(ierr);
1311*93520041SJoe Wallwork   }
13123f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
13133f07679eSJoe Wallwork   if (realIntegral < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global metric normalization factor should be strictly positive, not %.4e Is the input metric positive-definite?", realIntegral);
13143f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
131520282da2SJoe Wallwork 
131620282da2SJoe Wallwork   /* Apply local scaling */
131720282da2SJoe Wallwork   if (restrictSizes) {
131816522936SJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
131916522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
132016522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
132116522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
132216522936SJoe Wallwork     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);
132316522936SJoe Wallwork   }
132416522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
132516522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
132616522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
132720282da2SJoe Wallwork   }
132820282da2SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
13293f07679eSJoe Wallwork   ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
1330*93520041SJoe Wallwork   if (uniform) {
1331*93520041SJoe Wallwork 
1332*93520041SJoe Wallwork     /* Uniform case */
1333*93520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
1334*93520041SJoe Wallwork     if (restrictSizes) { ierr = DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det);CHKERRQ(ierr); }
1335*93520041SJoe Wallwork   } else {
1336*93520041SJoe Wallwork 
1337*93520041SJoe Wallwork     /* Spatially varying case */
1338*93520041SJoe Wallwork     PetscInt nrow;
1339*93520041SJoe Wallwork 
1340*93520041SJoe Wallwork     if (isotropic) nrow = 1;
1341*93520041SJoe Wallwork     else nrow = dim;
134216522936SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
134320282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13443f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
134520282da2SJoe Wallwork 
134620282da2SJoe Wallwork       ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr);
13473f07679eSJoe Wallwork       ierr = DMPlexPointLocalRef(dmDet, v, det, &detM);CHKERRQ(ierr);
13483f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
134920282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
1350*93520041SJoe Wallwork       if (restrictSizes) { ierr = DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM);CHKERRQ(ierr); }
1351*93520041SJoe Wallwork     }
135220282da2SJoe Wallwork   }
13533f07679eSJoe Wallwork   ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
135420282da2SJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
13553f07679eSJoe Wallwork   ierr = VecDestroy(&determinant);CHKERRQ(ierr);
1356*93520041SJoe Wallwork   if (!uniform) { ierr = DMDestroy(&dmDet);CHKERRQ(ierr); }
135720282da2SJoe Wallwork 
135820282da2SJoe Wallwork   PetscFunctionReturn(0);
135920282da2SJoe Wallwork }
136020282da2SJoe Wallwork 
1361cb7bfe8cSJoe Wallwork /*@
136220282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
136320282da2SJoe Wallwork 
136420282da2SJoe Wallwork   Input Parameter:
136520282da2SJoe Wallwork + dm         - The DM
136620282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
136720282da2SJoe Wallwork . weights    - Weights for the average
136820282da2SJoe Wallwork - metrics    - The metrics to be averaged
136920282da2SJoe Wallwork 
137020282da2SJoe Wallwork   Output Parameter:
137120282da2SJoe Wallwork . metricAvg  - The averaged metric
137220282da2SJoe Wallwork 
137320282da2SJoe Wallwork   Level: beginner
137420282da2SJoe Wallwork 
137520282da2SJoe Wallwork   Notes:
137620282da2SJoe Wallwork   The weights should sum to unity.
137720282da2SJoe Wallwork 
137820282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
137920282da2SJoe Wallwork 
138020282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
1381cb7bfe8cSJoe Wallwork @*/
138220282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
138320282da2SJoe Wallwork {
138420282da2SJoe Wallwork   PetscBool      haveWeights = PETSC_TRUE;
138520282da2SJoe Wallwork   PetscErrorCode ierr;
1386*93520041SJoe Wallwork   PetscInt       i, m, n;
138720282da2SJoe Wallwork   PetscReal      sum = 0.0, tol = 1.0e-10;
138820282da2SJoe Wallwork 
138920282da2SJoe Wallwork   PetscFunctionBegin;
1390*93520041SJoe Wallwork   if (numMetrics < 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics);
139120282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr);
139220282da2SJoe Wallwork   ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr);
1393*93520041SJoe Wallwork   ierr = VecGetSize(*metricAvg, &m);CHKERRQ(ierr);
1394*93520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
1395*93520041SJoe Wallwork     ierr = VecGetSize(metrics[i], &n);CHKERRQ(ierr);
1396*93520041SJoe Wallwork     if (m != n) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
1397*93520041SJoe Wallwork   }
139820282da2SJoe Wallwork 
139920282da2SJoe Wallwork   /* Default to the unweighted case */
140020282da2SJoe Wallwork   if (!weights) {
140120282da2SJoe Wallwork     ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr);
140220282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
140320282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
140420282da2SJoe Wallwork   }
140520282da2SJoe Wallwork 
140620282da2SJoe Wallwork   /* Check weights sum to unity */
1407*93520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
1408*93520041SJoe Wallwork   if (PetscAbsReal(sum - 1) > tol) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
140920282da2SJoe Wallwork 
141020282da2SJoe Wallwork   /* Compute metric average */
141120282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); }
141220282da2SJoe Wallwork   if (!haveWeights) { ierr = PetscFree(weights); }
141320282da2SJoe Wallwork   PetscFunctionReturn(0);
141420282da2SJoe Wallwork }
141520282da2SJoe Wallwork 
1416cb7bfe8cSJoe Wallwork /*@
141720282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
141820282da2SJoe Wallwork 
141920282da2SJoe Wallwork   Input Parameter:
142020282da2SJoe Wallwork + dm         - The DM
142120282da2SJoe Wallwork . metric1    - The first metric to be averaged
142220282da2SJoe Wallwork - metric2    - The second metric to be averaged
142320282da2SJoe Wallwork 
142420282da2SJoe Wallwork   Output Parameter:
142520282da2SJoe Wallwork . metricAvg  - The averaged metric
142620282da2SJoe Wallwork 
142720282da2SJoe Wallwork   Level: beginner
142820282da2SJoe Wallwork 
142920282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
1430cb7bfe8cSJoe Wallwork @*/
143120282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
143220282da2SJoe Wallwork {
143320282da2SJoe Wallwork   PetscErrorCode ierr;
143420282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
143520282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
143620282da2SJoe Wallwork 
143720282da2SJoe Wallwork   PetscFunctionBegin;
143820282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr);
143920282da2SJoe Wallwork   PetscFunctionReturn(0);
144020282da2SJoe Wallwork }
144120282da2SJoe Wallwork 
1442cb7bfe8cSJoe Wallwork /*@
144320282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
144420282da2SJoe Wallwork 
144520282da2SJoe Wallwork   Input Parameter:
144620282da2SJoe Wallwork + dm         - The DM
144720282da2SJoe Wallwork . metric1    - The first metric to be averaged
144820282da2SJoe Wallwork . metric2    - The second metric to be averaged
144920282da2SJoe Wallwork - metric3    - The third metric to be averaged
145020282da2SJoe Wallwork 
145120282da2SJoe Wallwork   Output Parameter:
145220282da2SJoe Wallwork . metricAvg  - The averaged metric
145320282da2SJoe Wallwork 
145420282da2SJoe Wallwork   Level: beginner
145520282da2SJoe Wallwork 
145620282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
1457cb7bfe8cSJoe Wallwork @*/
145820282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
145920282da2SJoe Wallwork {
146020282da2SJoe Wallwork   PetscErrorCode ierr;
146120282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
146220282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
146320282da2SJoe Wallwork 
146420282da2SJoe Wallwork   PetscFunctionBegin;
146520282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr);
146620282da2SJoe Wallwork   PetscFunctionReturn(0);
146720282da2SJoe Wallwork }
146820282da2SJoe Wallwork 
146920282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
147020282da2SJoe Wallwork {
147120282da2SJoe Wallwork   PetscErrorCode ierr;
147220282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
147320282da2SJoe Wallwork   PetscReal     *evals, *evals1;
147420282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
147520282da2SJoe Wallwork 
147620282da2SJoe Wallwork   PetscFunctionBegin;
1477*93520041SJoe Wallwork 
1478*93520041SJoe Wallwork   /* Isotropic case */
1479*93520041SJoe Wallwork   if (dim == 1) {
1480*93520041SJoe Wallwork     M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
1481*93520041SJoe Wallwork     PetscFunctionReturn(0);
1482*93520041SJoe Wallwork   }
1483*93520041SJoe Wallwork 
1484*93520041SJoe Wallwork   /* Anisotropic case */
148520282da2SJoe Wallwork   ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr);
148620282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
148720282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
148820282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
148920282da2SJoe Wallwork     }
149020282da2SJoe Wallwork   }
149120282da2SJoe Wallwork   {
149220282da2SJoe Wallwork     PetscScalar *work;
149320282da2SJoe Wallwork     PetscBLASInt lwork;
149420282da2SJoe Wallwork 
149520282da2SJoe Wallwork     lwork = 5*dim;
149620282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
149720282da2SJoe Wallwork     {
149820282da2SJoe Wallwork       PetscBLASInt lierr, nb;
149920282da2SJoe Wallwork       PetscReal    sqrtk;
150020282da2SJoe Wallwork 
150120282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
150220282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
150320282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
150420282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
150520282da2SJoe Wallwork       {
150620282da2SJoe Wallwork         PetscReal *rwork;
150720282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
150820282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
150920282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
151020282da2SJoe Wallwork       }
151120282da2SJoe Wallwork #else
151220282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
151320282da2SJoe Wallwork #endif
151482490253SJoe Wallwork       if (lierr) {
151582490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
151682490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
151782490253SJoe Wallwork       }
151820282da2SJoe Wallwork       ierr = PetscFPTrapPop();
151920282da2SJoe Wallwork 
152020282da2SJoe Wallwork       /* Compute square root and reciprocal */
152120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
152220282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
152320282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
152420282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
152520282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
152620282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
152720282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
152820282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
152920282da2SJoe Wallwork           }
153020282da2SJoe Wallwork         }
153120282da2SJoe Wallwork       }
153220282da2SJoe Wallwork 
153320282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
153420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
153520282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
153620282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
153720282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
153820282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
153920282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
154020282da2SJoe Wallwork             }
154120282da2SJoe Wallwork           }
154220282da2SJoe Wallwork         }
154320282da2SJoe Wallwork       }
154420282da2SJoe Wallwork 
154520282da2SJoe Wallwork       /* Compute eigendecomposition */
154620282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
154720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
154820282da2SJoe Wallwork       {
154920282da2SJoe Wallwork         PetscReal *rwork;
155020282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
155120282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
155220282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
155320282da2SJoe Wallwork       }
155420282da2SJoe Wallwork #else
155520282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
155620282da2SJoe Wallwork #endif
155782490253SJoe Wallwork       if (lierr) {
155882490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
155982490253SJoe Wallwork           for (j = 0; j < dim; ++j) {
156082490253SJoe Wallwork             evecs[i*dim+j] = 0.0;
156182490253SJoe Wallwork             for (k = 0; k < dim; ++k) {
156282490253SJoe Wallwork               for (l = 0; l < dim; ++l) {
156382490253SJoe Wallwork                 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
156482490253SJoe Wallwork               }
156582490253SJoe Wallwork             }
156682490253SJoe Wallwork           }
156782490253SJoe Wallwork         }
156882490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
156982490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
157082490253SJoe Wallwork       }
157120282da2SJoe Wallwork       ierr = PetscFPTrapPop();
157220282da2SJoe Wallwork 
157320282da2SJoe Wallwork       /* Modify eigenvalues */
157420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
157520282da2SJoe Wallwork 
157620282da2SJoe Wallwork       /* Map back to get the intersection */
157720282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
157820282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
157920282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
158020282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
158120282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
158220282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
158320282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
158420282da2SJoe Wallwork               }
158520282da2SJoe Wallwork             }
158620282da2SJoe Wallwork           }
158720282da2SJoe Wallwork         }
158820282da2SJoe Wallwork       }
158920282da2SJoe Wallwork     }
159020282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
159120282da2SJoe Wallwork   }
159220282da2SJoe Wallwork   ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr);
159320282da2SJoe Wallwork   PetscFunctionReturn(0);
159420282da2SJoe Wallwork }
159520282da2SJoe Wallwork 
1596cb7bfe8cSJoe Wallwork /*@
159720282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
159820282da2SJoe Wallwork 
159920282da2SJoe Wallwork   Input Parameter:
160020282da2SJoe Wallwork + dm         - The DM
160120282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
160220282da2SJoe Wallwork - metrics    - The metrics to be intersected
160320282da2SJoe Wallwork 
160420282da2SJoe Wallwork   Output Parameter:
160520282da2SJoe Wallwork . metricInt  - The intersected metric
160620282da2SJoe Wallwork 
160720282da2SJoe Wallwork   Level: beginner
160820282da2SJoe Wallwork 
160920282da2SJoe Wallwork   Notes:
161020282da2SJoe Wallwork 
161120282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
161220282da2SJoe Wallwork 
161320282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
161420282da2SJoe Wallwork 
161520282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
1616cb7bfe8cSJoe Wallwork @*/
161720282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
161820282da2SJoe Wallwork {
1619*93520041SJoe Wallwork   PetscBool      isotropic, uniform;
162020282da2SJoe Wallwork   PetscErrorCode ierr;
1621*93520041SJoe Wallwork   PetscInt       v, i, m, n;
1622*93520041SJoe Wallwork   PetscScalar   *met, *meti;
162320282da2SJoe Wallwork 
162420282da2SJoe Wallwork   PetscFunctionBegin;
1625*93520041SJoe Wallwork   if (numMetrics < 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics);
162620282da2SJoe Wallwork 
162720282da2SJoe Wallwork   /* Copy over the first metric */
1628*93520041SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr);
162920282da2SJoe Wallwork   ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr);
1630*93520041SJoe Wallwork   if (numMetrics == 1) PetscFunctionReturn(0);
1631*93520041SJoe Wallwork   ierr = VecGetSize(*metricInt, &m);CHKERRQ(ierr);
1632*93520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
1633*93520041SJoe Wallwork     ierr = VecGetSize(metrics[i], &n);CHKERRQ(ierr);
1634*93520041SJoe Wallwork     if (m != n) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
1635*93520041SJoe Wallwork   }
163620282da2SJoe Wallwork 
163720282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
1638*93520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
1639*93520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
1640*93520041SJoe Wallwork   if (uniform) {
1641*93520041SJoe Wallwork 
1642*93520041SJoe Wallwork     /* Uniform case */
1643*93520041SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
1644*93520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
1645*93520041SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
1646*93520041SJoe Wallwork       ierr = DMPlexMetricIntersection_Private(1, meti, met);CHKERRQ(ierr);
1647*93520041SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
1648*93520041SJoe Wallwork     }
1649*93520041SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
1650*93520041SJoe Wallwork   } else {
1651*93520041SJoe Wallwork 
1652*93520041SJoe Wallwork     /* Spatially varying case */
1653*93520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
1654*93520041SJoe Wallwork     PetscScalar *M, *Mi;
1655*93520041SJoe Wallwork 
1656*93520041SJoe Wallwork     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1657*93520041SJoe Wallwork     if (isotropic) nrow = 1;
1658*93520041SJoe Wallwork     else nrow = dim;
1659*93520041SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
166020282da2SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
166120282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
166220282da2SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
166320282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
166420282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr);
166520282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr);
1666*93520041SJoe Wallwork         ierr = DMPlexMetricIntersection_Private(nrow, Mi, M);CHKERRQ(ierr);
166720282da2SJoe Wallwork       }
166820282da2SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
166920282da2SJoe Wallwork     }
167020282da2SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
167120282da2SJoe Wallwork   }
167220282da2SJoe Wallwork   PetscFunctionReturn(0);
167320282da2SJoe Wallwork }
167420282da2SJoe Wallwork 
1675cb7bfe8cSJoe Wallwork /*@
167620282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
167720282da2SJoe Wallwork 
167820282da2SJoe Wallwork   Input Parameter:
167920282da2SJoe Wallwork + dm        - The DM
168020282da2SJoe Wallwork . metric1   - The first metric to be intersected
168120282da2SJoe Wallwork - metric2   - The second metric to be intersected
168220282da2SJoe Wallwork 
168320282da2SJoe Wallwork   Output Parameter:
168420282da2SJoe Wallwork . metricInt - The intersected metric
168520282da2SJoe Wallwork 
168620282da2SJoe Wallwork   Level: beginner
168720282da2SJoe Wallwork 
168820282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
1689cb7bfe8cSJoe Wallwork @*/
169020282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
169120282da2SJoe Wallwork {
169220282da2SJoe Wallwork   PetscErrorCode ierr;
169320282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
169420282da2SJoe Wallwork 
169520282da2SJoe Wallwork   PetscFunctionBegin;
169620282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr);
169720282da2SJoe Wallwork   PetscFunctionReturn(0);
169820282da2SJoe Wallwork }
169920282da2SJoe Wallwork 
1700cb7bfe8cSJoe Wallwork /*@
170120282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
170220282da2SJoe Wallwork 
170320282da2SJoe Wallwork   Input Parameter:
170420282da2SJoe Wallwork + dm        - The DM
170520282da2SJoe Wallwork . metric1   - The first metric to be intersected
170620282da2SJoe Wallwork . metric2   - The second metric to be intersected
170720282da2SJoe Wallwork - metric3   - The third metric to be intersected
170820282da2SJoe Wallwork 
170920282da2SJoe Wallwork   Output Parameter:
171020282da2SJoe Wallwork . metricInt - The intersected metric
171120282da2SJoe Wallwork 
171220282da2SJoe Wallwork   Level: beginner
171320282da2SJoe Wallwork 
171420282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
1715cb7bfe8cSJoe Wallwork @*/
171620282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
171720282da2SJoe Wallwork {
171820282da2SJoe Wallwork   PetscErrorCode ierr;
171920282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
172020282da2SJoe Wallwork 
172120282da2SJoe Wallwork   PetscFunctionBegin;
172220282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr);
172320282da2SJoe Wallwork   PetscFunctionReturn(0);
172420282da2SJoe Wallwork }
1725