xref: /petsc/src/dm/impls/plex/plexmetric.c (revision fe902acea6de04c8b00888cb5efe15799b75dfac)
120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
220282da2SJoe Wallwork #include <petscblaslapack.h>
320282da2SJoe Wallwork 
4*fe902aceSJoe Wallwork PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection;
5*fe902aceSJoe Wallwork 
631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
731b70465SJoe Wallwork {
831b70465SJoe Wallwork   MPI_Comm       comm;
993520041SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
10cc2eee55SJoe Wallwork   PetscBool      noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE;
1131b70465SJoe Wallwork   PetscErrorCode ierr;
12cc2eee55SJoe Wallwork   PetscInt       verbosity = -1, numIter = 3;
13cc2eee55SJoe 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;
1431b70465SJoe Wallwork 
1531b70465SJoe Wallwork   PetscFunctionBegin;
1631b70465SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1731b70465SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr);
1831b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr);
19cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, isotropic);CHKERRQ(ierr);
2093520041SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL);CHKERRQ(ierr);
2193520041SJoe Wallwork   ierr = DMPlexMetricSetUniform(dm, uniform);CHKERRQ(ierr);
2231b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr);
23cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst);CHKERRQ(ierr);
24cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL);CHKERRQ(ierr);
25cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoInsertion(dm, noInsert);CHKERRQ(ierr);
26cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL);CHKERRQ(ierr);
27cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoSwapping(dm, noSwap);CHKERRQ(ierr);
28cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL);CHKERRQ(ierr);
29cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoMovement(dm, noMove);CHKERRQ(ierr);
30cc2eee55SJoe Wallwork   ierr = PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0);CHKERRQ(ierr);
31cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNumIterations(dm, numIter);CHKERRQ(ierr);
32cc2eee55SJoe 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);
33cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetVerbosity(dm, verbosity);CHKERRQ(ierr);
3431b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr);
3531b70465SJoe Wallwork   ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr);
3631b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr);
3731b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr);
3831b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr);
3931b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr);
4031b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr);
4131b70465SJoe Wallwork   ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr);
4231b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr);
4331b70465SJoe Wallwork   ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr);
44cc2eee55SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL);CHKERRQ(ierr);
45cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetGradationFactor(dm, beta);CHKERRQ(ierr);
4631b70465SJoe Wallwork   ierr = PetscOptionsEnd();CHKERRQ(ierr);
4731b70465SJoe Wallwork   PetscFunctionReturn(0);
4831b70465SJoe Wallwork }
4931b70465SJoe Wallwork 
50cb7bfe8cSJoe Wallwork /*@
5131b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
5231b70465SJoe Wallwork 
5331b70465SJoe Wallwork   Input parameters:
5431b70465SJoe Wallwork + dm        - The DM
5531b70465SJoe Wallwork - isotropic - Is the metric isotropic?
5631b70465SJoe Wallwork 
5731b70465SJoe Wallwork   Level: beginner
5831b70465SJoe Wallwork 
5993520041SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst()
60cb7bfe8cSJoe Wallwork @*/
6131b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
6231b70465SJoe Wallwork {
6331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6431b70465SJoe Wallwork   PetscErrorCode ierr;
6531b70465SJoe Wallwork 
6631b70465SJoe Wallwork   PetscFunctionBegin;
6731b70465SJoe Wallwork   if (!plex->metricCtx) {
6831b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
6931b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
7031b70465SJoe Wallwork   }
7131b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
7231b70465SJoe Wallwork   PetscFunctionReturn(0);
7331b70465SJoe Wallwork }
7431b70465SJoe Wallwork 
75cb7bfe8cSJoe Wallwork /*@
7693520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7731b70465SJoe Wallwork 
7831b70465SJoe Wallwork   Input parameters:
7931b70465SJoe Wallwork . dm        - The DM
8031b70465SJoe Wallwork 
8131b70465SJoe Wallwork   Output parameters:
8231b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8331b70465SJoe Wallwork 
8431b70465SJoe Wallwork   Level: beginner
8531b70465SJoe Wallwork 
8693520041SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst()
87cb7bfe8cSJoe Wallwork @*/
8831b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
8931b70465SJoe Wallwork {
9031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
9131b70465SJoe Wallwork   PetscErrorCode ierr;
9231b70465SJoe Wallwork 
9331b70465SJoe Wallwork   PetscFunctionBegin;
9431b70465SJoe Wallwork   if (!plex->metricCtx) {
9531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
9631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
9731b70465SJoe Wallwork   }
9831b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
9931b70465SJoe Wallwork   PetscFunctionReturn(0);
10031b70465SJoe Wallwork }
10131b70465SJoe Wallwork 
102cb7bfe8cSJoe Wallwork /*@
10393520041SJoe Wallwork   DMPlexMetricSetUniform - Record whether a metric is uniform
10493520041SJoe Wallwork 
10593520041SJoe Wallwork   Input parameters:
10693520041SJoe Wallwork + dm      - The DM
10793520041SJoe Wallwork - uniform - Is the metric uniform?
10893520041SJoe Wallwork 
10993520041SJoe Wallwork   Level: beginner
11093520041SJoe Wallwork 
11193520041SJoe Wallwork   Notes:
11293520041SJoe Wallwork 
11393520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11493520041SJoe Wallwork 
11593520041SJoe Wallwork .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
11693520041SJoe Wallwork @*/
11793520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
11893520041SJoe Wallwork {
11993520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
12093520041SJoe Wallwork   PetscErrorCode ierr;
12193520041SJoe Wallwork 
12293520041SJoe Wallwork   PetscFunctionBegin;
12393520041SJoe Wallwork   if (!plex->metricCtx) {
12493520041SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
12593520041SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
12693520041SJoe Wallwork   }
12793520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
12893520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
12993520041SJoe Wallwork   PetscFunctionReturn(0);
13093520041SJoe Wallwork }
13193520041SJoe Wallwork 
13293520041SJoe Wallwork /*@
13393520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
13493520041SJoe Wallwork 
13593520041SJoe Wallwork   Input parameters:
13693520041SJoe Wallwork . dm      - The DM
13793520041SJoe Wallwork 
13893520041SJoe Wallwork   Output parameters:
13993520041SJoe Wallwork . uniform - Is the metric uniform?
14093520041SJoe Wallwork 
14193520041SJoe Wallwork   Level: beginner
14293520041SJoe Wallwork 
14393520041SJoe Wallwork .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
14493520041SJoe Wallwork @*/
14593520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
14693520041SJoe Wallwork {
14793520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
14893520041SJoe Wallwork   PetscErrorCode ierr;
14993520041SJoe Wallwork 
15093520041SJoe Wallwork   PetscFunctionBegin;
15193520041SJoe Wallwork   if (!plex->metricCtx) {
15293520041SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
15393520041SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
15493520041SJoe Wallwork   }
15593520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
15693520041SJoe Wallwork   PetscFunctionReturn(0);
15793520041SJoe Wallwork }
15893520041SJoe Wallwork 
15993520041SJoe Wallwork /*@
16031b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
16131b70465SJoe Wallwork 
16231b70465SJoe Wallwork   Input parameters:
16331b70465SJoe Wallwork + dm                      - The DM
16431b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
16531b70465SJoe Wallwork 
16631b70465SJoe Wallwork   Level: beginner
16731b70465SJoe Wallwork 
16831b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
169cb7bfe8cSJoe Wallwork @*/
17031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
17131b70465SJoe Wallwork {
17231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
17331b70465SJoe Wallwork   PetscErrorCode ierr;
17431b70465SJoe Wallwork 
17531b70465SJoe Wallwork   PetscFunctionBegin;
17631b70465SJoe Wallwork   if (!plex->metricCtx) {
17731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
17831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
17931b70465SJoe Wallwork   }
18031b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
18131b70465SJoe Wallwork   PetscFunctionReturn(0);
18231b70465SJoe Wallwork }
18331b70465SJoe Wallwork 
184cb7bfe8cSJoe Wallwork /*@
18531b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
18631b70465SJoe Wallwork 
18731b70465SJoe Wallwork   Input parameters:
18831b70465SJoe Wallwork . dm                      - The DM
18931b70465SJoe Wallwork 
19031b70465SJoe Wallwork   Output parameters:
19131b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
19231b70465SJoe Wallwork 
19331b70465SJoe Wallwork   Level: beginner
19431b70465SJoe Wallwork 
19531b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
196cb7bfe8cSJoe Wallwork @*/
19731b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
19831b70465SJoe Wallwork {
19931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
20031b70465SJoe Wallwork   PetscErrorCode ierr;
20131b70465SJoe Wallwork 
20231b70465SJoe Wallwork   PetscFunctionBegin;
20331b70465SJoe Wallwork   if (!plex->metricCtx) {
20431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
20531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
20631b70465SJoe Wallwork   }
20731b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
20831b70465SJoe Wallwork   PetscFunctionReturn(0);
20931b70465SJoe Wallwork }
21031b70465SJoe Wallwork 
211cb7bfe8cSJoe Wallwork /*@
212cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
213cc2eee55SJoe Wallwork 
214cc2eee55SJoe Wallwork   Input parameters:
215cc2eee55SJoe Wallwork + dm       - The DM
216cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
217cc2eee55SJoe Wallwork 
218cc2eee55SJoe Wallwork   Level: beginner
219cc2eee55SJoe Wallwork 
220cb7bfe8cSJoe Wallwork   Notes:
221cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
222cb7bfe8cSJoe Wallwork 
223cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement()
224cb7bfe8cSJoe Wallwork @*/
225cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
226cc2eee55SJoe Wallwork {
227cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
228cc2eee55SJoe Wallwork   PetscErrorCode ierr;
229cc2eee55SJoe Wallwork 
230cc2eee55SJoe Wallwork   PetscFunctionBegin;
231cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
232cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
233cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
234cc2eee55SJoe Wallwork   }
235cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
236cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
237cc2eee55SJoe Wallwork }
238cc2eee55SJoe Wallwork 
239cb7bfe8cSJoe Wallwork /*@
240cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
241cc2eee55SJoe Wallwork 
242cc2eee55SJoe Wallwork   Input parameters:
243cc2eee55SJoe Wallwork . dm       - The DM
244cc2eee55SJoe Wallwork 
245cc2eee55SJoe Wallwork   Output parameters:
246cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
247cc2eee55SJoe Wallwork 
248cc2eee55SJoe Wallwork   Level: beginner
249cc2eee55SJoe Wallwork 
250cb7bfe8cSJoe Wallwork   Notes:
251cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
252cb7bfe8cSJoe Wallwork 
253cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement()
254cb7bfe8cSJoe Wallwork @*/
255cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
256cc2eee55SJoe Wallwork {
257cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
258cc2eee55SJoe Wallwork   PetscErrorCode ierr;
259cc2eee55SJoe Wallwork 
260cc2eee55SJoe Wallwork   PetscFunctionBegin;
261cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
262cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
263cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
264cc2eee55SJoe Wallwork   }
265cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
266cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
267cc2eee55SJoe Wallwork }
268cc2eee55SJoe Wallwork 
269cb7bfe8cSJoe Wallwork /*@
270cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
271cc2eee55SJoe Wallwork 
272cc2eee55SJoe Wallwork   Input parameters:
273cc2eee55SJoe Wallwork + dm     - The DM
274cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
275cc2eee55SJoe Wallwork 
276cc2eee55SJoe Wallwork   Level: beginner
277cc2eee55SJoe Wallwork 
278cb7bfe8cSJoe Wallwork   Notes:
279cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
280cb7bfe8cSJoe Wallwork 
281cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement()
282cb7bfe8cSJoe Wallwork @*/
283cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
284cc2eee55SJoe Wallwork {
285cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
286cc2eee55SJoe Wallwork   PetscErrorCode ierr;
287cc2eee55SJoe Wallwork 
288cc2eee55SJoe Wallwork   PetscFunctionBegin;
289cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
290cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
291cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
292cc2eee55SJoe Wallwork   }
293cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
294cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
295cc2eee55SJoe Wallwork }
296cc2eee55SJoe Wallwork 
297cb7bfe8cSJoe Wallwork /*@
298cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
299cc2eee55SJoe Wallwork 
300cc2eee55SJoe Wallwork   Input parameters:
301cc2eee55SJoe Wallwork . dm     - The DM
302cc2eee55SJoe Wallwork 
303cc2eee55SJoe Wallwork   Output parameters:
304cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
305cc2eee55SJoe Wallwork 
306cc2eee55SJoe Wallwork   Level: beginner
307cc2eee55SJoe Wallwork 
308cb7bfe8cSJoe Wallwork   Notes:
309cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
310cb7bfe8cSJoe Wallwork 
311cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement()
312cb7bfe8cSJoe Wallwork @*/
313cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
314cc2eee55SJoe Wallwork {
315cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
316cc2eee55SJoe Wallwork   PetscErrorCode ierr;
317cc2eee55SJoe Wallwork 
318cc2eee55SJoe Wallwork   PetscFunctionBegin;
319cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
320cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
321cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
322cc2eee55SJoe Wallwork   }
323cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
324cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
325cc2eee55SJoe Wallwork }
326cc2eee55SJoe Wallwork 
327cb7bfe8cSJoe Wallwork /*@
328cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
329cc2eee55SJoe Wallwork 
330cc2eee55SJoe Wallwork   Input parameters:
331cc2eee55SJoe Wallwork + dm     - The DM
332cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
333cc2eee55SJoe Wallwork 
334cc2eee55SJoe Wallwork   Level: beginner
335cc2eee55SJoe Wallwork 
336cb7bfe8cSJoe Wallwork   Notes:
337cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
338cb7bfe8cSJoe Wallwork 
339cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping()
340cb7bfe8cSJoe Wallwork @*/
341cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
342cc2eee55SJoe Wallwork {
343cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
344cc2eee55SJoe Wallwork   PetscErrorCode ierr;
345cc2eee55SJoe Wallwork 
346cc2eee55SJoe Wallwork   PetscFunctionBegin;
347cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
348cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
349cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
350cc2eee55SJoe Wallwork   }
351cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
352cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
353cc2eee55SJoe Wallwork }
354cc2eee55SJoe Wallwork 
355cb7bfe8cSJoe Wallwork /*@
356cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
357cc2eee55SJoe Wallwork 
358cc2eee55SJoe Wallwork   Input parameters:
359cc2eee55SJoe Wallwork . dm     - The DM
360cc2eee55SJoe Wallwork 
361cc2eee55SJoe Wallwork   Output parameters:
362cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
363cc2eee55SJoe Wallwork 
364cc2eee55SJoe Wallwork   Level: beginner
365cc2eee55SJoe Wallwork 
366cb7bfe8cSJoe Wallwork   Notes:
367cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
368cb7bfe8cSJoe Wallwork 
369cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping()
370cb7bfe8cSJoe Wallwork @*/
371cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
372cc2eee55SJoe Wallwork {
373cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
374cc2eee55SJoe Wallwork   PetscErrorCode ierr;
375cc2eee55SJoe Wallwork 
376cc2eee55SJoe Wallwork   PetscFunctionBegin;
377cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
378cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
379cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
380cc2eee55SJoe Wallwork   }
381cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
382cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
383cc2eee55SJoe Wallwork }
384cc2eee55SJoe Wallwork 
385cb7bfe8cSJoe Wallwork /*@
38631b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
38731b70465SJoe Wallwork 
38831b70465SJoe Wallwork   Input parameters:
38931b70465SJoe Wallwork + dm    - The DM
39031b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
39131b70465SJoe Wallwork 
39231b70465SJoe Wallwork   Level: beginner
39331b70465SJoe Wallwork 
39431b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
395cb7bfe8cSJoe Wallwork @*/
39631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
39731b70465SJoe Wallwork {
39831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
39931b70465SJoe Wallwork   PetscErrorCode ierr;
40031b70465SJoe Wallwork 
40131b70465SJoe Wallwork   PetscFunctionBegin;
40231b70465SJoe Wallwork   if (!plex->metricCtx) {
40331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
40431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
40531b70465SJoe Wallwork   }
40631b70465SJoe Wallwork   if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min);
40731b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
40831b70465SJoe Wallwork   PetscFunctionReturn(0);
40931b70465SJoe Wallwork }
41031b70465SJoe Wallwork 
411cb7bfe8cSJoe Wallwork /*@
41231b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
41331b70465SJoe Wallwork 
41431b70465SJoe Wallwork   Input parameters:
41531b70465SJoe Wallwork . dm    - The DM
41631b70465SJoe Wallwork 
417cc2eee55SJoe Wallwork   Output parameters:
41831b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
41931b70465SJoe Wallwork 
42031b70465SJoe Wallwork   Level: beginner
42131b70465SJoe Wallwork 
42231b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
423cb7bfe8cSJoe Wallwork @*/
42431b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
42531b70465SJoe Wallwork {
42631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
42731b70465SJoe Wallwork   PetscErrorCode ierr;
42831b70465SJoe Wallwork 
42931b70465SJoe Wallwork   PetscFunctionBegin;
43031b70465SJoe Wallwork   if (!plex->metricCtx) {
43131b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
43231b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
43331b70465SJoe Wallwork   }
43431b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
43531b70465SJoe Wallwork   PetscFunctionReturn(0);
43631b70465SJoe Wallwork }
43731b70465SJoe Wallwork 
438cb7bfe8cSJoe Wallwork /*@
43931b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
44031b70465SJoe Wallwork 
44131b70465SJoe Wallwork   Input parameters:
44231b70465SJoe Wallwork + dm    - The DM
44331b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44431b70465SJoe Wallwork 
44531b70465SJoe Wallwork   Level: beginner
44631b70465SJoe Wallwork 
44731b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
448cb7bfe8cSJoe Wallwork @*/
44931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
45031b70465SJoe Wallwork {
45131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
45231b70465SJoe Wallwork   PetscErrorCode ierr;
45331b70465SJoe Wallwork 
45431b70465SJoe Wallwork   PetscFunctionBegin;
45531b70465SJoe Wallwork   if (!plex->metricCtx) {
45631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
45731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
45831b70465SJoe Wallwork   }
45931b70465SJoe Wallwork   if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max);
46031b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
46131b70465SJoe Wallwork   PetscFunctionReturn(0);
46231b70465SJoe Wallwork }
46331b70465SJoe Wallwork 
464cb7bfe8cSJoe Wallwork /*@
46531b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
46631b70465SJoe Wallwork 
46731b70465SJoe Wallwork   Input parameters:
46831b70465SJoe Wallwork . dm    - The DM
46931b70465SJoe Wallwork 
470cc2eee55SJoe Wallwork   Output parameters:
47131b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
47231b70465SJoe Wallwork 
47331b70465SJoe Wallwork   Level: beginner
47431b70465SJoe Wallwork 
47531b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
476cb7bfe8cSJoe Wallwork @*/
47731b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
47831b70465SJoe Wallwork {
47931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
48031b70465SJoe Wallwork   PetscErrorCode ierr;
48131b70465SJoe Wallwork 
48231b70465SJoe Wallwork   PetscFunctionBegin;
48331b70465SJoe Wallwork   if (!plex->metricCtx) {
48431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
48531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
48631b70465SJoe Wallwork   }
48731b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
48831b70465SJoe Wallwork   PetscFunctionReturn(0);
48931b70465SJoe Wallwork }
49031b70465SJoe Wallwork 
491cb7bfe8cSJoe Wallwork /*@
49231b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
49331b70465SJoe Wallwork 
49431b70465SJoe Wallwork   Input parameters:
49531b70465SJoe Wallwork + dm    - The DM
49631b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
49731b70465SJoe Wallwork 
49831b70465SJoe Wallwork   Level: beginner
49931b70465SJoe Wallwork 
50031b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
50131b70465SJoe Wallwork 
50231b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
503cb7bfe8cSJoe Wallwork @*/
50431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
50531b70465SJoe Wallwork {
50631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
50731b70465SJoe Wallwork   PetscErrorCode ierr;
50831b70465SJoe Wallwork 
50931b70465SJoe Wallwork   PetscFunctionBegin;
51031b70465SJoe Wallwork   if (!plex->metricCtx) {
51131b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
51231b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
51331b70465SJoe Wallwork   }
51431b70465SJoe Wallwork   if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max);
51531b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
51631b70465SJoe Wallwork   PetscFunctionReturn(0);
51731b70465SJoe Wallwork }
51831b70465SJoe Wallwork 
519cb7bfe8cSJoe Wallwork /*@
52031b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
52131b70465SJoe Wallwork 
52231b70465SJoe Wallwork   Input parameters:
52331b70465SJoe Wallwork . dm    - The DM
52431b70465SJoe Wallwork 
525cc2eee55SJoe Wallwork   Output parameters:
52631b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
52731b70465SJoe Wallwork 
52831b70465SJoe Wallwork   Level: beginner
52931b70465SJoe Wallwork 
53031b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
531cb7bfe8cSJoe Wallwork @*/
53231b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
53331b70465SJoe Wallwork {
53431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
53531b70465SJoe Wallwork   PetscErrorCode ierr;
53631b70465SJoe Wallwork 
53731b70465SJoe Wallwork   PetscFunctionBegin;
53831b70465SJoe Wallwork   if (!plex->metricCtx) {
53931b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
54031b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
54131b70465SJoe Wallwork   }
54231b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
54331b70465SJoe Wallwork   PetscFunctionReturn(0);
54431b70465SJoe Wallwork }
54531b70465SJoe Wallwork 
546cb7bfe8cSJoe Wallwork /*@
54731b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
54831b70465SJoe Wallwork 
54931b70465SJoe Wallwork   Input parameters:
55031b70465SJoe Wallwork + dm               - The DM
55131b70465SJoe Wallwork - targetComplexity - The target metric complexity
55231b70465SJoe Wallwork 
55331b70465SJoe Wallwork   Level: beginner
55431b70465SJoe Wallwork 
55531b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
556cb7bfe8cSJoe Wallwork @*/
55731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
55831b70465SJoe Wallwork {
55931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
56031b70465SJoe Wallwork   PetscErrorCode ierr;
56131b70465SJoe Wallwork 
56231b70465SJoe Wallwork   PetscFunctionBegin;
56331b70465SJoe Wallwork   if (!plex->metricCtx) {
56431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
56531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
56631b70465SJoe Wallwork   }
56731b70465SJoe Wallwork   if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive");
56831b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
56931b70465SJoe Wallwork   PetscFunctionReturn(0);
57031b70465SJoe Wallwork }
57131b70465SJoe Wallwork 
572cb7bfe8cSJoe Wallwork /*@
57331b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
57431b70465SJoe Wallwork 
57531b70465SJoe Wallwork   Input parameters:
57631b70465SJoe Wallwork . dm               - The DM
57731b70465SJoe Wallwork 
578cc2eee55SJoe Wallwork   Output parameters:
57931b70465SJoe Wallwork . targetComplexity - The target metric complexity
58031b70465SJoe Wallwork 
58131b70465SJoe Wallwork   Level: beginner
58231b70465SJoe Wallwork 
58331b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
584cb7bfe8cSJoe Wallwork @*/
58531b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
58631b70465SJoe Wallwork {
58731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
58831b70465SJoe Wallwork   PetscErrorCode ierr;
58931b70465SJoe Wallwork 
59031b70465SJoe Wallwork   PetscFunctionBegin;
59131b70465SJoe Wallwork   if (!plex->metricCtx) {
59231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
59331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
59431b70465SJoe Wallwork   }
59531b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
59631b70465SJoe Wallwork   PetscFunctionReturn(0);
59731b70465SJoe Wallwork }
59831b70465SJoe Wallwork 
599cb7bfe8cSJoe Wallwork /*@
60031b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
60131b70465SJoe Wallwork 
60231b70465SJoe Wallwork   Input parameters:
60331b70465SJoe Wallwork + dm - The DM
60431b70465SJoe Wallwork - p  - The normalization order
60531b70465SJoe Wallwork 
60631b70465SJoe Wallwork   Level: beginner
60731b70465SJoe Wallwork 
60831b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
609cb7bfe8cSJoe Wallwork @*/
61031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
61131b70465SJoe Wallwork {
61231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
61331b70465SJoe Wallwork   PetscErrorCode ierr;
61431b70465SJoe Wallwork 
61531b70465SJoe Wallwork   PetscFunctionBegin;
61631b70465SJoe Wallwork   if (!plex->metricCtx) {
61731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
61831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
61931b70465SJoe Wallwork   }
62031b70465SJoe Wallwork   if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater");
62131b70465SJoe Wallwork   plex->metricCtx->p = p;
62231b70465SJoe Wallwork   PetscFunctionReturn(0);
62331b70465SJoe Wallwork }
62431b70465SJoe Wallwork 
625cb7bfe8cSJoe Wallwork /*@
62631b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
62731b70465SJoe Wallwork 
62831b70465SJoe Wallwork   Input parameters:
62931b70465SJoe Wallwork . dm - The DM
63031b70465SJoe Wallwork 
631cc2eee55SJoe Wallwork   Output parameters:
63231b70465SJoe Wallwork . p - The normalization order
63331b70465SJoe Wallwork 
63431b70465SJoe Wallwork   Level: beginner
63531b70465SJoe Wallwork 
63631b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
637cb7bfe8cSJoe Wallwork @*/
63831b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
63931b70465SJoe Wallwork {
64031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
64131b70465SJoe Wallwork   PetscErrorCode ierr;
64231b70465SJoe Wallwork 
64331b70465SJoe Wallwork   PetscFunctionBegin;
64431b70465SJoe Wallwork   if (!plex->metricCtx) {
64531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
64631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
64731b70465SJoe Wallwork   }
64831b70465SJoe Wallwork   *p = plex->metricCtx->p;
64931b70465SJoe Wallwork   PetscFunctionReturn(0);
65031b70465SJoe Wallwork }
65131b70465SJoe Wallwork 
652cb7bfe8cSJoe Wallwork /*@
653cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
654cc2eee55SJoe Wallwork 
655cc2eee55SJoe Wallwork   Input parameters:
656cc2eee55SJoe Wallwork + dm   - The DM
657cc2eee55SJoe Wallwork - beta - The metric gradation factor
658cc2eee55SJoe Wallwork 
659cc2eee55SJoe Wallwork   Level: beginner
660cc2eee55SJoe Wallwork 
661cc2eee55SJoe Wallwork   Notes:
662cc2eee55SJoe Wallwork 
663cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
664cc2eee55SJoe Wallwork 
665cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
666cc2eee55SJoe Wallwork 
667cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
668cb7bfe8cSJoe Wallwork 
669cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetGradationFactor()
670cb7bfe8cSJoe Wallwork @*/
671cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
672cc2eee55SJoe Wallwork {
673cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
674cc2eee55SJoe Wallwork   PetscErrorCode ierr;
675cc2eee55SJoe Wallwork 
676cc2eee55SJoe Wallwork   PetscFunctionBegin;
677cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
678cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
679cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
680cc2eee55SJoe Wallwork   }
681cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
682cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
683cc2eee55SJoe Wallwork }
684cc2eee55SJoe Wallwork 
685cb7bfe8cSJoe Wallwork /*@
686cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
687cc2eee55SJoe Wallwork 
688cc2eee55SJoe Wallwork   Input parameters:
689cc2eee55SJoe Wallwork . dm   - The DM
690cc2eee55SJoe Wallwork 
691cc2eee55SJoe Wallwork   Output parameters:
692cc2eee55SJoe Wallwork . beta - The metric gradation factor
693cc2eee55SJoe Wallwork 
694cc2eee55SJoe Wallwork   Level: beginner
695cc2eee55SJoe Wallwork 
696cb7bfe8cSJoe Wallwork   Notes:
697cb7bfe8cSJoe Wallwork 
698cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
699cb7bfe8cSJoe Wallwork 
700cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
701cb7bfe8cSJoe Wallwork 
702cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
703cc2eee55SJoe Wallwork 
704cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetGradationFactor()
705cb7bfe8cSJoe Wallwork @*/
706cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
707cc2eee55SJoe Wallwork {
708cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
709cc2eee55SJoe Wallwork   PetscErrorCode ierr;
710cc2eee55SJoe Wallwork 
711cc2eee55SJoe Wallwork   PetscFunctionBegin;
712cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
713cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
714cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
715cc2eee55SJoe Wallwork   }
716cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
717cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
718cc2eee55SJoe Wallwork }
719cc2eee55SJoe Wallwork 
720cb7bfe8cSJoe Wallwork /*@
721cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
722cc2eee55SJoe Wallwork 
723cc2eee55SJoe Wallwork   Input parameters:
724cc2eee55SJoe Wallwork + dm        - The DM
725cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
726cc2eee55SJoe Wallwork 
727cb7bfe8cSJoe Wallwork   Level: beginner
728cb7bfe8cSJoe Wallwork 
729cb7bfe8cSJoe Wallwork   Notes:
730cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
731cb7bfe8cSJoe Wallwork 
732cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations()
733cb7bfe8cSJoe Wallwork @*/
734cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
735cc2eee55SJoe Wallwork {
736cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
737cc2eee55SJoe Wallwork   PetscErrorCode ierr;
738cc2eee55SJoe Wallwork 
739cc2eee55SJoe Wallwork   PetscFunctionBegin;
740cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
741cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
742cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
743cc2eee55SJoe Wallwork   }
744cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
745cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
746cc2eee55SJoe Wallwork }
747cc2eee55SJoe Wallwork 
748cb7bfe8cSJoe Wallwork /*@
749cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
750cc2eee55SJoe Wallwork 
751cc2eee55SJoe Wallwork   Input parameters:
752cc2eee55SJoe Wallwork . dm        - The DM
753cc2eee55SJoe Wallwork 
754cc2eee55SJoe Wallwork   Output parameters:
755cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
756cc2eee55SJoe Wallwork 
757cb7bfe8cSJoe Wallwork   Level: beginner
758cb7bfe8cSJoe Wallwork 
759cb7bfe8cSJoe Wallwork   Notes:
760cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
761cb7bfe8cSJoe Wallwork 
762cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
763cb7bfe8cSJoe Wallwork @*/
764cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
765cc2eee55SJoe Wallwork {
766cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
767cc2eee55SJoe Wallwork   PetscErrorCode ierr;
768cc2eee55SJoe Wallwork 
769cc2eee55SJoe Wallwork   PetscFunctionBegin;
770cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
771cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
772cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
773cc2eee55SJoe Wallwork   }
774cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
775cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
776cc2eee55SJoe Wallwork }
777cc2eee55SJoe Wallwork 
778cb7bfe8cSJoe Wallwork /*@
779cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
780cc2eee55SJoe Wallwork 
781cc2eee55SJoe Wallwork   Input parameters:
782cc2eee55SJoe Wallwork + dm      - The DM
783cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
784cc2eee55SJoe Wallwork 
785cb7bfe8cSJoe Wallwork   Level: beginner
786cb7bfe8cSJoe Wallwork 
787cb7bfe8cSJoe Wallwork   Notes:
788cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
789cc2eee55SJoe Wallwork 
790cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
791cb7bfe8cSJoe Wallwork @*/
792cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
793cc2eee55SJoe Wallwork {
794cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
795cc2eee55SJoe Wallwork   PetscErrorCode ierr;
796cc2eee55SJoe Wallwork 
797cc2eee55SJoe Wallwork   PetscFunctionBegin;
798cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
799cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
800cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
801cc2eee55SJoe Wallwork   }
802cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
803cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
804cc2eee55SJoe Wallwork }
805cc2eee55SJoe Wallwork 
806cb7bfe8cSJoe Wallwork /*@
807cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
808cc2eee55SJoe Wallwork 
809cc2eee55SJoe Wallwork   Input parameters:
810cc2eee55SJoe Wallwork . dm      - The DM
811cc2eee55SJoe Wallwork 
812cc2eee55SJoe Wallwork   Output parameters:
813cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
814cc2eee55SJoe Wallwork 
815cb7bfe8cSJoe Wallwork   Level: beginner
816cb7bfe8cSJoe Wallwork 
817cb7bfe8cSJoe Wallwork   Notes:
818cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
819cc2eee55SJoe Wallwork 
820cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity()
821cb7bfe8cSJoe Wallwork @*/
822cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
823cc2eee55SJoe Wallwork {
824cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
825cc2eee55SJoe Wallwork   PetscErrorCode ierr;
826cc2eee55SJoe Wallwork 
827cc2eee55SJoe Wallwork   PetscFunctionBegin;
828cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
829cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
830cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
831cc2eee55SJoe Wallwork   }
832cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
833cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
834cc2eee55SJoe Wallwork }
835cc2eee55SJoe Wallwork 
83620282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
83720282da2SJoe Wallwork {
83820282da2SJoe Wallwork   MPI_Comm       comm;
83920282da2SJoe Wallwork   PetscErrorCode ierr;
84020282da2SJoe Wallwork   PetscFE        fe;
84120282da2SJoe Wallwork   PetscInt       dim;
84220282da2SJoe Wallwork 
84320282da2SJoe Wallwork   PetscFunctionBegin;
84420282da2SJoe Wallwork 
84520282da2SJoe Wallwork   /* Extract metadata from dm */
84620282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
84720282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
84820282da2SJoe Wallwork 
84920282da2SJoe Wallwork   /* Create a P1 field of the requested size */
85020282da2SJoe Wallwork   ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
85120282da2SJoe Wallwork   ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr);
85220282da2SJoe Wallwork   ierr = DMCreateDS(dm);CHKERRQ(ierr);
85320282da2SJoe Wallwork   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
85420282da2SJoe Wallwork   ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr);
85520282da2SJoe Wallwork 
85620282da2SJoe Wallwork   PetscFunctionReturn(0);
85720282da2SJoe Wallwork }
85820282da2SJoe Wallwork 
859cb7bfe8cSJoe Wallwork /*@
86020282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
86120282da2SJoe Wallwork 
86220282da2SJoe Wallwork   Input parameters:
86320282da2SJoe Wallwork + dm     - The DM
86420282da2SJoe Wallwork - f      - The field number to use
86520282da2SJoe Wallwork 
86620282da2SJoe Wallwork   Output parameter:
86720282da2SJoe Wallwork . metric - The metric
86820282da2SJoe Wallwork 
86920282da2SJoe Wallwork   Level: beginner
87020282da2SJoe Wallwork 
87131b70465SJoe Wallwork   Notes:
87231b70465SJoe Wallwork 
87331b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
87431b70465SJoe Wallwork 
87531b70465SJoe Wallwork   Command line options for Riemannian metrics:
87631b70465SJoe Wallwork 
877cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
87893520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
879cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
880cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
881cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
882cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
883cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
884cb7bfe8cSJoe Wallwork . -dm_plex_metric_target_complexity         - Target metric complexity
885cb7bfe8cSJoe Wallwork 
886cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
887cb7bfe8cSJoe Wallwork 
888cb7bfe8cSJoe Wallwork . -dm_adaptor <pragmatic/mmg/parmmg>
889cb7bfe8cSJoe Wallwork 
890cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
891cb7bfe8cSJoe Wallwork 
892cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
893cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
894cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
895cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
896cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
897cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
89820282da2SJoe Wallwork 
89920282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
900cb7bfe8cSJoe Wallwork @*/
90120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
90220282da2SJoe Wallwork {
90331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
90493520041SJoe Wallwork   PetscBool      isotropic, uniform;
90520282da2SJoe Wallwork   PetscErrorCode ierr;
90620282da2SJoe Wallwork   PetscInt       coordDim, Nd;
90720282da2SJoe Wallwork 
90820282da2SJoe Wallwork   PetscFunctionBegin;
90931b70465SJoe Wallwork   if (!plex->metricCtx) {
91031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
91131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
91231b70465SJoe Wallwork   }
91320282da2SJoe Wallwork   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
91420282da2SJoe Wallwork   Nd = coordDim*coordDim;
91593520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
91693520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
91793520041SJoe Wallwork   if (uniform) {
91893520041SJoe Wallwork     MPI_Comm comm;
91993520041SJoe Wallwork 
92093520041SJoe Wallwork     ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
92193520041SJoe Wallwork     if (!isotropic) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
92293520041SJoe Wallwork     ierr = VecCreate(comm, metric);CHKERRQ(ierr);
92393520041SJoe Wallwork     ierr = VecSetSizes(*metric, 1, PETSC_DECIDE);CHKERRQ(ierr);
92493520041SJoe Wallwork     ierr = VecSetFromOptions(*metric);CHKERRQ(ierr);
92593520041SJoe Wallwork   } else if (isotropic) {
92693520041SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dm, f, 1, metric);CHKERRQ(ierr);
92793520041SJoe Wallwork   } else {
92820282da2SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr);
92993520041SJoe Wallwork   }
93020282da2SJoe Wallwork   PetscFunctionReturn(0);
93120282da2SJoe Wallwork }
93220282da2SJoe Wallwork 
933cb7bfe8cSJoe Wallwork /*@
93420282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
93520282da2SJoe Wallwork 
93620282da2SJoe Wallwork   Input parameters:
93720282da2SJoe Wallwork + dm     - The DM
93820282da2SJoe Wallwork . f      - The field number to use
93920282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
94020282da2SJoe Wallwork 
94120282da2SJoe Wallwork   Output parameter:
94220282da2SJoe Wallwork . metric - The uniform metric
94320282da2SJoe Wallwork 
94420282da2SJoe Wallwork   Level: beginner
94520282da2SJoe Wallwork 
94693520041SJoe Wallwork   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
94720282da2SJoe Wallwork 
94820282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
949cb7bfe8cSJoe Wallwork @*/
95020282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
95120282da2SJoe Wallwork {
95220282da2SJoe Wallwork   PetscErrorCode ierr;
95320282da2SJoe Wallwork 
95420282da2SJoe Wallwork   PetscFunctionBegin;
95593520041SJoe Wallwork   ierr = DMPlexMetricSetUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
95620282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
95720282da2SJoe Wallwork   if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
95820282da2SJoe Wallwork   if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
95993520041SJoe Wallwork   ierr = VecSet(*metric, alpha);CHKERRQ(ierr);
96093520041SJoe Wallwork   ierr = VecAssemblyBegin(*metric);CHKERRQ(ierr);
96193520041SJoe Wallwork   ierr = VecAssemblyEnd(*metric);CHKERRQ(ierr);
96220282da2SJoe Wallwork   PetscFunctionReturn(0);
96320282da2SJoe Wallwork }
96420282da2SJoe Wallwork 
96593520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
96693520041SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96793520041SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96893520041SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
96993520041SJoe Wallwork {
97093520041SJoe Wallwork   f0[0] = u[0];
97193520041SJoe Wallwork }
97293520041SJoe Wallwork 
973cb7bfe8cSJoe Wallwork /*@
97420282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
97520282da2SJoe Wallwork 
97620282da2SJoe Wallwork   Input parameters:
97720282da2SJoe Wallwork + dm        - The DM
97820282da2SJoe Wallwork . f         - The field number to use
97920282da2SJoe Wallwork - indicator - The error indicator
98020282da2SJoe Wallwork 
98120282da2SJoe Wallwork   Output parameter:
98220282da2SJoe Wallwork . metric    - The isotropic metric
98320282da2SJoe Wallwork 
98420282da2SJoe Wallwork   Level: beginner
98520282da2SJoe Wallwork 
98620282da2SJoe Wallwork   Notes:
98720282da2SJoe Wallwork 
98820282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
98920282da2SJoe Wallwork 
99093520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
99120282da2SJoe Wallwork 
99220282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
993cb7bfe8cSJoe Wallwork @*/
99420282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
99520282da2SJoe Wallwork {
99620282da2SJoe Wallwork   PetscErrorCode ierr;
99793520041SJoe Wallwork   PetscInt       m, n;
99820282da2SJoe Wallwork 
99920282da2SJoe Wallwork   PetscFunctionBegin;
100093520041SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, PETSC_TRUE);CHKERRQ(ierr);
100120282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
100293520041SJoe Wallwork   ierr = VecGetSize(indicator, &m);CHKERRQ(ierr);
100393520041SJoe Wallwork   ierr = VecGetSize(*metric, &n);CHKERRQ(ierr);
100493520041SJoe Wallwork   if (m == n) { ierr = VecCopy(indicator, *metric);CHKERRQ(ierr); }
100593520041SJoe Wallwork   else {
100693520041SJoe Wallwork     DM     dmIndi;
100793520041SJoe Wallwork     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
100893520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
100993520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101093520041SJoe Wallwork                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
101193520041SJoe Wallwork 
101220282da2SJoe Wallwork     ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr);
101393520041SJoe Wallwork     funcs[0] = identity;
101493520041SJoe Wallwork     ierr = DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric);CHKERRQ(ierr);
101520282da2SJoe Wallwork   }
101620282da2SJoe Wallwork   PetscFunctionReturn(0);
101720282da2SJoe Wallwork }
101820282da2SJoe Wallwork 
101982490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
102082490253SJoe Wallwork {
102182490253SJoe Wallwork   PetscInt i, j;
102282490253SJoe Wallwork 
102382490253SJoe Wallwork   PetscFunctionBegin;
102482490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
102582490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
102682490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
102782490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
102882490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
102982490253SJoe Wallwork       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]);
103082490253SJoe Wallwork       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]);
103182490253SJoe Wallwork     }
103282490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
103382490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
103482490253SJoe Wallwork   }
103582490253SJoe Wallwork   PetscFunctionReturn(0);
103682490253SJoe Wallwork }
103782490253SJoe Wallwork 
10383f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
103920282da2SJoe Wallwork {
104020282da2SJoe Wallwork   PetscErrorCode ierr;
104120282da2SJoe Wallwork   PetscInt       i, j, k;
104220282da2SJoe 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);
104320282da2SJoe Wallwork   PetscScalar   *Mpos;
104420282da2SJoe Wallwork 
104520282da2SJoe Wallwork   PetscFunctionBegin;
104620282da2SJoe Wallwork   ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr);
104720282da2SJoe Wallwork 
104820282da2SJoe Wallwork   /* Symmetrize */
104920282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
105020282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
105120282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
105220282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
105320282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
105420282da2SJoe Wallwork     }
105520282da2SJoe Wallwork   }
105620282da2SJoe Wallwork 
105720282da2SJoe Wallwork   /* Compute eigendecomposition */
105893520041SJoe Wallwork   if (dim == 1) {
105993520041SJoe Wallwork 
106093520041SJoe Wallwork     /* Isotropic case */
106193520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
106293520041SJoe Wallwork     Mpos[0] = 1.0;
106393520041SJoe Wallwork   } else {
106493520041SJoe Wallwork 
106593520041SJoe Wallwork     /* Anisotropic case */
106620282da2SJoe Wallwork     PetscScalar  *work;
106720282da2SJoe Wallwork     PetscBLASInt lwork;
106820282da2SJoe Wallwork 
106920282da2SJoe Wallwork     lwork = 5*dim;
107020282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
107120282da2SJoe Wallwork     {
107220282da2SJoe Wallwork       PetscBLASInt lierr;
107320282da2SJoe Wallwork       PetscBLASInt nb;
107420282da2SJoe Wallwork 
107520282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
107620282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
107720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
107820282da2SJoe Wallwork       {
107920282da2SJoe Wallwork         PetscReal *rwork;
108020282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
108120282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
108220282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
108320282da2SJoe Wallwork       }
108420282da2SJoe Wallwork #else
108520282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
108620282da2SJoe Wallwork #endif
108782490253SJoe Wallwork       if (lierr) {
108882490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
108982490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
109082490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
109182490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
109282490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
109382490253SJoe Wallwork           }
109482490253SJoe Wallwork         }
109582490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
109682490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
109782490253SJoe Wallwork       }
109820282da2SJoe Wallwork       ierr = PetscFPTrapPop();CHKERRQ(ierr);
109920282da2SJoe Wallwork     }
110020282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
110120282da2SJoe Wallwork   }
110220282da2SJoe Wallwork 
110320282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
110420282da2SJoe Wallwork   max_eig = 0.0;
110520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
110620282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
110720282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
110820282da2SJoe Wallwork   }
110920282da2SJoe Wallwork 
11103f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11113f07679eSJoe Wallwork   *detMp = 1.0;
111220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111320282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
11143f07679eSJoe Wallwork     *detMp *= eigs[i];
111520282da2SJoe Wallwork   }
111620282da2SJoe Wallwork 
111720282da2SJoe Wallwork   /* Reconstruct Hessian */
111820282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111920282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
112020282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
112120282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
112220282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
112320282da2SJoe Wallwork       }
112420282da2SJoe Wallwork     }
112520282da2SJoe Wallwork   }
112620282da2SJoe Wallwork   ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr);
112720282da2SJoe Wallwork 
112820282da2SJoe Wallwork   PetscFunctionReturn(0);
112920282da2SJoe Wallwork }
113020282da2SJoe Wallwork 
1131cb7bfe8cSJoe Wallwork /*@
113220282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
113320282da2SJoe Wallwork 
113420282da2SJoe Wallwork   Input parameters:
113520282da2SJoe Wallwork + dm                 - The DM
11363f07679eSJoe Wallwork . metricIn           - The metric
113799abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11383f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
113920282da2SJoe Wallwork 
114020282da2SJoe Wallwork   Output parameter:
11413f07679eSJoe Wallwork + metricOut          - The metric
11423f07679eSJoe Wallwork - determinant        - Its determinant
114320282da2SJoe Wallwork 
114420282da2SJoe Wallwork   Level: beginner
114520282da2SJoe Wallwork 
1146cb7bfe8cSJoe Wallwork   Notes:
1147cb7bfe8cSJoe Wallwork 
1148cb7bfe8cSJoe Wallwork   Relevant command line options:
1149cb7bfe8cSJoe Wallwork 
115093520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
115193520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
115293520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1153cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1154cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1155cb7bfe8cSJoe Wallwork 
115620282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
1157cb7bfe8cSJoe Wallwork @*/
11583f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant)
115920282da2SJoe Wallwork {
11603f07679eSJoe Wallwork   DM             dmDet;
116193520041SJoe Wallwork   PetscBool      isotropic, uniform;
116220282da2SJoe Wallwork   PetscErrorCode ierr;
116320282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
11643f07679eSJoe Wallwork   PetscScalar   *met, *det;
116520282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
116620282da2SJoe Wallwork 
116720282da2SJoe Wallwork   PetscFunctionBegin;
1168*fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0);CHKERRQ(ierr);
116920282da2SJoe Wallwork 
117020282da2SJoe Wallwork   /* Extract metadata from dm */
117120282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
117220282da2SJoe Wallwork   if (restrictSizes) {
117399abec2bSJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
117499abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
117599abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
117699abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
117799abec2bSJoe 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);
117899abec2bSJoe Wallwork   }
117999abec2bSJoe Wallwork   if (restrictAnisotropy) {
118099abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
118199abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
118220282da2SJoe Wallwork   }
118320282da2SJoe Wallwork 
118493520041SJoe Wallwork   /* Setup output metric */
11853f07679eSJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr);
11863f07679eSJoe Wallwork   ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr);
118793520041SJoe Wallwork 
118893520041SJoe Wallwork   /* Enforce SPD and extract determinant */
118993520041SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
119093520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
119193520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
119293520041SJoe Wallwork   if (uniform) {
119393520041SJoe Wallwork     if (!isotropic) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
119493520041SJoe Wallwork 
119593520041SJoe Wallwork     /* Uniform case */
119693520041SJoe Wallwork     ierr = VecDuplicate(metricIn, determinant);CHKERRQ(ierr);
119793520041SJoe Wallwork     ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
119893520041SJoe Wallwork     ierr = DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det);CHKERRQ(ierr);
119993520041SJoe Wallwork     ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
120093520041SJoe Wallwork   } else {
120193520041SJoe Wallwork 
120293520041SJoe Wallwork     /* Spatially varying case */
120393520041SJoe Wallwork     PetscInt nrow;
120493520041SJoe Wallwork 
120593520041SJoe Wallwork     if (isotropic) nrow = 1;
120693520041SJoe Wallwork     else nrow = dim;
12073f07679eSJoe Wallwork     ierr = DMClone(dm, &dmDet);CHKERRQ(ierr);
12083f07679eSJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant);CHKERRQ(ierr);
120920282da2SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12103f07679eSJoe Wallwork     ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
121120282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12123f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
121320282da2SJoe Wallwork       ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
12143f07679eSJoe Wallwork       ierr = DMPlexPointLocalRef(dmDet, v, det, &vdet);CHKERRQ(ierr);
121593520041SJoe Wallwork       ierr = DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet);CHKERRQ(ierr);
121620282da2SJoe Wallwork     }
12173f07679eSJoe Wallwork     ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
121893520041SJoe Wallwork   }
12193f07679eSJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
1220*fe902aceSJoe Wallwork 
1221*fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0);CHKERRQ(ierr);
122220282da2SJoe Wallwork   PetscFunctionReturn(0);
122320282da2SJoe Wallwork }
122420282da2SJoe Wallwork 
122520282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
122620282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
122720282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122820282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
122920282da2SJoe Wallwork {
123020282da2SJoe Wallwork   const PetscScalar p = constants[0];
123120282da2SJoe Wallwork 
12323f07679eSJoe Wallwork   f0[0] = PetscPowReal(u[0], p/(2.0*p + dim));
123320282da2SJoe Wallwork }
123420282da2SJoe Wallwork 
1235cb7bfe8cSJoe Wallwork /*@
123620282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
123720282da2SJoe Wallwork 
123820282da2SJoe Wallwork   Input parameters:
123920282da2SJoe Wallwork + dm                 - The DM
124020282da2SJoe Wallwork . metricIn           - The unnormalized metric
124116522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
124216522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
124320282da2SJoe Wallwork 
124420282da2SJoe Wallwork   Output parameter:
124520282da2SJoe Wallwork . metricOut          - The normalized metric
124620282da2SJoe Wallwork 
124720282da2SJoe Wallwork   Level: beginner
124820282da2SJoe Wallwork 
1249cb7bfe8cSJoe Wallwork   Notes:
1250cb7bfe8cSJoe Wallwork 
1251cb7bfe8cSJoe Wallwork   Relevant command line options:
1252cb7bfe8cSJoe Wallwork 
125393520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
125493520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
125593520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1256cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1257cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1258cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1259cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1260cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1261cb7bfe8cSJoe Wallwork 
126220282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
1263cb7bfe8cSJoe Wallwork @*/
126416522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut)
126520282da2SJoe Wallwork {
12663f07679eSJoe Wallwork   DM               dmDet;
126720282da2SJoe Wallwork   MPI_Comm         comm;
126893520041SJoe Wallwork   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
126920282da2SJoe Wallwork   PetscDS          ds;
127020282da2SJoe Wallwork   PetscErrorCode   ierr;
127120282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
12723f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
127393520041SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
12743f07679eSJoe Wallwork   Vec              determinant;
127520282da2SJoe Wallwork 
127620282da2SJoe Wallwork   PetscFunctionBegin;
1277*fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0);CHKERRQ(ierr);
127820282da2SJoe Wallwork 
127920282da2SJoe Wallwork   /* Extract metadata from dm */
128020282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
128120282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
128293520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
128393520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
128493520041SJoe Wallwork   if (isotropic) Nd = 1;
128593520041SJoe Wallwork   else Nd = dim*dim;
128620282da2SJoe Wallwork 
128720282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
128816522936SJoe Wallwork   ierr = DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst);CHKERRQ(ierr);
12893f07679eSJoe Wallwork   ierr = DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant);CHKERRQ(ierr);
129020282da2SJoe Wallwork 
129120282da2SJoe Wallwork   /* Compute global normalization factor */
129216522936SJoe Wallwork   ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
129316522936SJoe Wallwork   ierr = DMPlexMetricGetNormalizationOrder(dm, &p);CHKERRQ(ierr);
129416522936SJoe Wallwork   constants[0] = p;
129593520041SJoe Wallwork   if (uniform) {
129693520041SJoe Wallwork     if (!isotropic) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
129793520041SJoe Wallwork     DM  dmTmp;
129893520041SJoe Wallwork     Vec tmp;
129993520041SJoe Wallwork 
130093520041SJoe Wallwork     ierr = DMClone(dm, &dmTmp);CHKERRQ(ierr);
130193520041SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp);CHKERRQ(ierr);
130293520041SJoe Wallwork     ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
130393520041SJoe Wallwork     ierr = VecSet(tmp, det[0]);CHKERRQ(ierr);
130493520041SJoe Wallwork     ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
130593520041SJoe Wallwork     ierr = DMGetDS(dmTmp, &ds);CHKERRQ(ierr);
130693520041SJoe Wallwork     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
130793520041SJoe Wallwork     ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
130893520041SJoe Wallwork     ierr = DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL);CHKERRQ(ierr);
130993520041SJoe Wallwork     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
131093520041SJoe Wallwork     ierr = DMDestroy(&dmTmp);CHKERRQ(ierr);
131193520041SJoe Wallwork   } else {
13123f07679eSJoe Wallwork     ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr);
13133f07679eSJoe Wallwork     ierr = DMGetDS(dmDet, &ds);CHKERRQ(ierr);
131420282da2SJoe Wallwork     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
131520282da2SJoe Wallwork     ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
13163f07679eSJoe Wallwork     ierr = DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL);CHKERRQ(ierr);
131793520041SJoe Wallwork   }
13183f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
13193f07679eSJoe 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);
13203f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
132120282da2SJoe Wallwork 
132220282da2SJoe Wallwork   /* Apply local scaling */
132320282da2SJoe Wallwork   if (restrictSizes) {
132416522936SJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
132516522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
132616522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
132716522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
132816522936SJoe 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);
132916522936SJoe Wallwork   }
133016522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
133116522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
133216522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
133320282da2SJoe Wallwork   }
133420282da2SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
13353f07679eSJoe Wallwork   ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
133693520041SJoe Wallwork   if (uniform) {
133793520041SJoe Wallwork 
133893520041SJoe Wallwork     /* Uniform case */
133993520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
134093520041SJoe Wallwork     if (restrictSizes) { ierr = DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det);CHKERRQ(ierr); }
134193520041SJoe Wallwork   } else {
134293520041SJoe Wallwork 
134393520041SJoe Wallwork     /* Spatially varying case */
134493520041SJoe Wallwork     PetscInt nrow;
134593520041SJoe Wallwork 
134693520041SJoe Wallwork     if (isotropic) nrow = 1;
134793520041SJoe Wallwork     else nrow = dim;
134816522936SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
134920282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13503f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
135120282da2SJoe Wallwork 
135220282da2SJoe Wallwork       ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr);
13533f07679eSJoe Wallwork       ierr = DMPlexPointLocalRef(dmDet, v, det, &detM);CHKERRQ(ierr);
13543f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
135520282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
135693520041SJoe Wallwork       if (restrictSizes) { ierr = DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM);CHKERRQ(ierr); }
135793520041SJoe Wallwork     }
135820282da2SJoe Wallwork   }
13593f07679eSJoe Wallwork   ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
136020282da2SJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
13613f07679eSJoe Wallwork   ierr = VecDestroy(&determinant);CHKERRQ(ierr);
136293520041SJoe Wallwork   if (!uniform) { ierr = DMDestroy(&dmDet);CHKERRQ(ierr); }
136320282da2SJoe Wallwork 
1364*fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0);CHKERRQ(ierr);
136520282da2SJoe Wallwork   PetscFunctionReturn(0);
136620282da2SJoe Wallwork }
136720282da2SJoe Wallwork 
1368cb7bfe8cSJoe Wallwork /*@
136920282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
137020282da2SJoe Wallwork 
137120282da2SJoe Wallwork   Input Parameter:
137220282da2SJoe Wallwork + dm         - The DM
137320282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
137420282da2SJoe Wallwork . weights    - Weights for the average
137520282da2SJoe Wallwork - metrics    - The metrics to be averaged
137620282da2SJoe Wallwork 
137720282da2SJoe Wallwork   Output Parameter:
137820282da2SJoe Wallwork . metricAvg  - The averaged metric
137920282da2SJoe Wallwork 
138020282da2SJoe Wallwork   Level: beginner
138120282da2SJoe Wallwork 
138220282da2SJoe Wallwork   Notes:
138320282da2SJoe Wallwork   The weights should sum to unity.
138420282da2SJoe Wallwork 
138520282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
138620282da2SJoe Wallwork 
138720282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
1388cb7bfe8cSJoe Wallwork @*/
138920282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
139020282da2SJoe Wallwork {
139120282da2SJoe Wallwork   PetscBool      haveWeights = PETSC_TRUE;
139220282da2SJoe Wallwork   PetscErrorCode ierr;
139393520041SJoe Wallwork   PetscInt       i, m, n;
139420282da2SJoe Wallwork   PetscReal      sum = 0.0, tol = 1.0e-10;
139520282da2SJoe Wallwork 
139620282da2SJoe Wallwork   PetscFunctionBegin;
1397*fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0);CHKERRQ(ierr);
139893520041SJoe Wallwork   if (numMetrics < 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics);
139920282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr);
140020282da2SJoe Wallwork   ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr);
140193520041SJoe Wallwork   ierr = VecGetSize(*metricAvg, &m);CHKERRQ(ierr);
140293520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
140393520041SJoe Wallwork     ierr = VecGetSize(metrics[i], &n);CHKERRQ(ierr);
140493520041SJoe Wallwork     if (m != n) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
140593520041SJoe Wallwork   }
140620282da2SJoe Wallwork 
140720282da2SJoe Wallwork   /* Default to the unweighted case */
140820282da2SJoe Wallwork   if (!weights) {
140920282da2SJoe Wallwork     ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr);
141020282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
141120282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
141220282da2SJoe Wallwork   }
141320282da2SJoe Wallwork 
141420282da2SJoe Wallwork   /* Check weights sum to unity */
141593520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
141693520041SJoe Wallwork   if (PetscAbsReal(sum - 1) > tol) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
141720282da2SJoe Wallwork 
141820282da2SJoe Wallwork   /* Compute metric average */
141920282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); }
142020282da2SJoe Wallwork   if (!haveWeights) { ierr = PetscFree(weights); }
1421*fe902aceSJoe Wallwork 
1422*fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0);CHKERRQ(ierr);
142320282da2SJoe Wallwork   PetscFunctionReturn(0);
142420282da2SJoe Wallwork }
142520282da2SJoe Wallwork 
1426cb7bfe8cSJoe Wallwork /*@
142720282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
142820282da2SJoe Wallwork 
142920282da2SJoe Wallwork   Input Parameter:
143020282da2SJoe Wallwork + dm         - The DM
143120282da2SJoe Wallwork . metric1    - The first metric to be averaged
143220282da2SJoe Wallwork - metric2    - The second metric to be averaged
143320282da2SJoe Wallwork 
143420282da2SJoe Wallwork   Output Parameter:
143520282da2SJoe Wallwork . metricAvg  - The averaged metric
143620282da2SJoe Wallwork 
143720282da2SJoe Wallwork   Level: beginner
143820282da2SJoe Wallwork 
143920282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
1440cb7bfe8cSJoe Wallwork @*/
144120282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
144220282da2SJoe Wallwork {
144320282da2SJoe Wallwork   PetscErrorCode ierr;
144420282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
144520282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
144620282da2SJoe Wallwork 
144720282da2SJoe Wallwork   PetscFunctionBegin;
144820282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr);
144920282da2SJoe Wallwork   PetscFunctionReturn(0);
145020282da2SJoe Wallwork }
145120282da2SJoe Wallwork 
1452cb7bfe8cSJoe Wallwork /*@
145320282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
145420282da2SJoe Wallwork 
145520282da2SJoe Wallwork   Input Parameter:
145620282da2SJoe Wallwork + dm         - The DM
145720282da2SJoe Wallwork . metric1    - The first metric to be averaged
145820282da2SJoe Wallwork . metric2    - The second metric to be averaged
145920282da2SJoe Wallwork - metric3    - The third metric to be averaged
146020282da2SJoe Wallwork 
146120282da2SJoe Wallwork   Output Parameter:
146220282da2SJoe Wallwork . metricAvg  - The averaged metric
146320282da2SJoe Wallwork 
146420282da2SJoe Wallwork   Level: beginner
146520282da2SJoe Wallwork 
146620282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
1467cb7bfe8cSJoe Wallwork @*/
146820282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
146920282da2SJoe Wallwork {
147020282da2SJoe Wallwork   PetscErrorCode ierr;
147120282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
147220282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
147320282da2SJoe Wallwork 
147420282da2SJoe Wallwork   PetscFunctionBegin;
147520282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr);
147620282da2SJoe Wallwork   PetscFunctionReturn(0);
147720282da2SJoe Wallwork }
147820282da2SJoe Wallwork 
147920282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
148020282da2SJoe Wallwork {
148120282da2SJoe Wallwork   PetscErrorCode ierr;
148220282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
148320282da2SJoe Wallwork   PetscReal     *evals, *evals1;
148420282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
148520282da2SJoe Wallwork 
148620282da2SJoe Wallwork   PetscFunctionBegin;
148793520041SJoe Wallwork 
148893520041SJoe Wallwork   /* Isotropic case */
148993520041SJoe Wallwork   if (dim == 1) {
149093520041SJoe Wallwork     M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
149193520041SJoe Wallwork     PetscFunctionReturn(0);
149293520041SJoe Wallwork   }
149393520041SJoe Wallwork 
149493520041SJoe Wallwork   /* Anisotropic case */
149520282da2SJoe Wallwork   ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr);
149620282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
149720282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
149820282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
149920282da2SJoe Wallwork     }
150020282da2SJoe Wallwork   }
150120282da2SJoe Wallwork   {
150220282da2SJoe Wallwork     PetscScalar *work;
150320282da2SJoe Wallwork     PetscBLASInt lwork;
150420282da2SJoe Wallwork 
150520282da2SJoe Wallwork     lwork = 5*dim;
150620282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
150720282da2SJoe Wallwork     {
150820282da2SJoe Wallwork       PetscBLASInt lierr, nb;
150920282da2SJoe Wallwork       PetscReal    sqrtk;
151020282da2SJoe Wallwork 
151120282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
151220282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
151320282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
151420282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
151520282da2SJoe Wallwork       {
151620282da2SJoe Wallwork         PetscReal *rwork;
151720282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
151820282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
151920282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
152020282da2SJoe Wallwork       }
152120282da2SJoe Wallwork #else
152220282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
152320282da2SJoe Wallwork #endif
152482490253SJoe Wallwork       if (lierr) {
152582490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
152682490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
152782490253SJoe Wallwork       }
152820282da2SJoe Wallwork       ierr = PetscFPTrapPop();
152920282da2SJoe Wallwork 
153020282da2SJoe Wallwork       /* Compute square root and reciprocal */
153120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
153220282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
153320282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
153420282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
153520282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
153620282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
153720282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
153820282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
153920282da2SJoe Wallwork           }
154020282da2SJoe Wallwork         }
154120282da2SJoe Wallwork       }
154220282da2SJoe Wallwork 
154320282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
154420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
154520282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
154620282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
154720282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
154820282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
154920282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
155020282da2SJoe Wallwork             }
155120282da2SJoe Wallwork           }
155220282da2SJoe Wallwork         }
155320282da2SJoe Wallwork       }
155420282da2SJoe Wallwork 
155520282da2SJoe Wallwork       /* Compute eigendecomposition */
155620282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
155720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
155820282da2SJoe Wallwork       {
155920282da2SJoe Wallwork         PetscReal *rwork;
156020282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
156120282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
156220282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
156320282da2SJoe Wallwork       }
156420282da2SJoe Wallwork #else
156520282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
156620282da2SJoe Wallwork #endif
156782490253SJoe Wallwork       if (lierr) {
156882490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
156982490253SJoe Wallwork           for (j = 0; j < dim; ++j) {
157082490253SJoe Wallwork             evecs[i*dim+j] = 0.0;
157182490253SJoe Wallwork             for (k = 0; k < dim; ++k) {
157282490253SJoe Wallwork               for (l = 0; l < dim; ++l) {
157382490253SJoe Wallwork                 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
157482490253SJoe Wallwork               }
157582490253SJoe Wallwork             }
157682490253SJoe Wallwork           }
157782490253SJoe Wallwork         }
157882490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
157982490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
158082490253SJoe Wallwork       }
158120282da2SJoe Wallwork       ierr = PetscFPTrapPop();
158220282da2SJoe Wallwork 
158320282da2SJoe Wallwork       /* Modify eigenvalues */
158420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
158520282da2SJoe Wallwork 
158620282da2SJoe Wallwork       /* Map back to get the intersection */
158720282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
158820282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
158920282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
159020282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
159120282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
159220282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
159320282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
159420282da2SJoe Wallwork               }
159520282da2SJoe Wallwork             }
159620282da2SJoe Wallwork           }
159720282da2SJoe Wallwork         }
159820282da2SJoe Wallwork       }
159920282da2SJoe Wallwork     }
160020282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
160120282da2SJoe Wallwork   }
160220282da2SJoe Wallwork   ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr);
160320282da2SJoe Wallwork   PetscFunctionReturn(0);
160420282da2SJoe Wallwork }
160520282da2SJoe Wallwork 
1606cb7bfe8cSJoe Wallwork /*@
160720282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
160820282da2SJoe Wallwork 
160920282da2SJoe Wallwork   Input Parameter:
161020282da2SJoe Wallwork + dm         - The DM
161120282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
161220282da2SJoe Wallwork - metrics    - The metrics to be intersected
161320282da2SJoe Wallwork 
161420282da2SJoe Wallwork   Output Parameter:
161520282da2SJoe Wallwork . metricInt  - The intersected metric
161620282da2SJoe Wallwork 
161720282da2SJoe Wallwork   Level: beginner
161820282da2SJoe Wallwork 
161920282da2SJoe Wallwork   Notes:
162020282da2SJoe Wallwork 
162120282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
162220282da2SJoe Wallwork 
162320282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
162420282da2SJoe Wallwork 
162520282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
1626cb7bfe8cSJoe Wallwork @*/
162720282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
162820282da2SJoe Wallwork {
162993520041SJoe Wallwork   PetscBool      isotropic, uniform;
163020282da2SJoe Wallwork   PetscErrorCode ierr;
163193520041SJoe Wallwork   PetscInt       v, i, m, n;
163293520041SJoe Wallwork   PetscScalar   *met, *meti;
163320282da2SJoe Wallwork 
163420282da2SJoe Wallwork   PetscFunctionBegin;
1635*fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0);CHKERRQ(ierr);
163693520041SJoe Wallwork   if (numMetrics < 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics);
163720282da2SJoe Wallwork 
163820282da2SJoe Wallwork   /* Copy over the first metric */
163993520041SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr);
164020282da2SJoe Wallwork   ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr);
164193520041SJoe Wallwork   if (numMetrics == 1) PetscFunctionReturn(0);
164293520041SJoe Wallwork   ierr = VecGetSize(*metricInt, &m);CHKERRQ(ierr);
164393520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
164493520041SJoe Wallwork     ierr = VecGetSize(metrics[i], &n);CHKERRQ(ierr);
164593520041SJoe Wallwork     if (m != n) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
164693520041SJoe Wallwork   }
164720282da2SJoe Wallwork 
164820282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
164993520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
165093520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
165193520041SJoe Wallwork   if (uniform) {
165293520041SJoe Wallwork 
165393520041SJoe Wallwork     /* Uniform case */
165493520041SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
165593520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
165693520041SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
165793520041SJoe Wallwork       ierr = DMPlexMetricIntersection_Private(1, meti, met);CHKERRQ(ierr);
165893520041SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
165993520041SJoe Wallwork     }
166093520041SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
166193520041SJoe Wallwork   } else {
166293520041SJoe Wallwork 
166393520041SJoe Wallwork     /* Spatially varying case */
166493520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
166593520041SJoe Wallwork     PetscScalar *M, *Mi;
166693520041SJoe Wallwork 
166793520041SJoe Wallwork     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
166893520041SJoe Wallwork     if (isotropic) nrow = 1;
166993520041SJoe Wallwork     else nrow = dim;
167093520041SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
167120282da2SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
167220282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
167320282da2SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
167420282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
167520282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr);
167620282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr);
167793520041SJoe Wallwork         ierr = DMPlexMetricIntersection_Private(nrow, Mi, M);CHKERRQ(ierr);
167820282da2SJoe Wallwork       }
167920282da2SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
168020282da2SJoe Wallwork     }
168120282da2SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
168220282da2SJoe Wallwork   }
1683*fe902aceSJoe Wallwork 
1684*fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0);CHKERRQ(ierr);
168520282da2SJoe Wallwork   PetscFunctionReturn(0);
168620282da2SJoe Wallwork }
168720282da2SJoe Wallwork 
1688cb7bfe8cSJoe Wallwork /*@
168920282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
169020282da2SJoe Wallwork 
169120282da2SJoe Wallwork   Input Parameter:
169220282da2SJoe Wallwork + dm        - The DM
169320282da2SJoe Wallwork . metric1   - The first metric to be intersected
169420282da2SJoe Wallwork - metric2   - The second metric to be intersected
169520282da2SJoe Wallwork 
169620282da2SJoe Wallwork   Output Parameter:
169720282da2SJoe Wallwork . metricInt - The intersected metric
169820282da2SJoe Wallwork 
169920282da2SJoe Wallwork   Level: beginner
170020282da2SJoe Wallwork 
170120282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
1702cb7bfe8cSJoe Wallwork @*/
170320282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
170420282da2SJoe Wallwork {
170520282da2SJoe Wallwork   PetscErrorCode ierr;
170620282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
170720282da2SJoe Wallwork 
170820282da2SJoe Wallwork   PetscFunctionBegin;
170920282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr);
171020282da2SJoe Wallwork   PetscFunctionReturn(0);
171120282da2SJoe Wallwork }
171220282da2SJoe Wallwork 
1713cb7bfe8cSJoe Wallwork /*@
171420282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
171520282da2SJoe Wallwork 
171620282da2SJoe Wallwork   Input Parameter:
171720282da2SJoe Wallwork + dm        - The DM
171820282da2SJoe Wallwork . metric1   - The first metric to be intersected
171920282da2SJoe Wallwork . metric2   - The second metric to be intersected
172020282da2SJoe Wallwork - metric3   - The third metric to be intersected
172120282da2SJoe Wallwork 
172220282da2SJoe Wallwork   Output Parameter:
172320282da2SJoe Wallwork . metricInt - The intersected metric
172420282da2SJoe Wallwork 
172520282da2SJoe Wallwork   Level: beginner
172620282da2SJoe Wallwork 
172720282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
1728cb7bfe8cSJoe Wallwork @*/
172920282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
173020282da2SJoe Wallwork {
173120282da2SJoe Wallwork   PetscErrorCode ierr;
173220282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
173320282da2SJoe Wallwork 
173420282da2SJoe Wallwork   PetscFunctionBegin;
173520282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr);
173620282da2SJoe Wallwork   PetscFunctionReturn(0);
173720282da2SJoe Wallwork }
1738