xref: /petsc/src/dm/impls/plex/plexmetric.c (revision ae8b063e4f7e8ac6d2de6909aa6d539b126b38d9)
120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
220282da2SJoe Wallwork #include <petscblaslapack.h>
320282da2SJoe Wallwork 
4fe902aceSJoe Wallwork PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection;
5fe902aceSJoe 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;
13*ae8b063eSJoe 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, hausd = 0.01;
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);
46*ae8b063eSJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL);CHKERRQ(ierr);
47*ae8b063eSJoe Wallwork   ierr = DMPlexMetricSetHausdorffNumber(dm, hausd);CHKERRQ(ierr);
4831b70465SJoe Wallwork   ierr = PetscOptionsEnd();CHKERRQ(ierr);
4931b70465SJoe Wallwork   PetscFunctionReturn(0);
5031b70465SJoe Wallwork }
5131b70465SJoe Wallwork 
52cb7bfe8cSJoe Wallwork /*@
5331b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
5431b70465SJoe Wallwork 
5531b70465SJoe Wallwork   Input parameters:
5631b70465SJoe Wallwork + dm        - The DM
5731b70465SJoe Wallwork - isotropic - Is the metric isotropic?
5831b70465SJoe Wallwork 
5931b70465SJoe Wallwork   Level: beginner
6031b70465SJoe Wallwork 
6193520041SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst()
62cb7bfe8cSJoe Wallwork @*/
6331b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
6431b70465SJoe Wallwork {
6531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6631b70465SJoe Wallwork   PetscErrorCode ierr;
6731b70465SJoe Wallwork 
6831b70465SJoe Wallwork   PetscFunctionBegin;
6931b70465SJoe Wallwork   if (!plex->metricCtx) {
7031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
7131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
7231b70465SJoe Wallwork   }
7331b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
7431b70465SJoe Wallwork   PetscFunctionReturn(0);
7531b70465SJoe Wallwork }
7631b70465SJoe Wallwork 
77cb7bfe8cSJoe Wallwork /*@
7893520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7931b70465SJoe Wallwork 
8031b70465SJoe Wallwork   Input parameters:
8131b70465SJoe Wallwork . dm        - The DM
8231b70465SJoe Wallwork 
8331b70465SJoe Wallwork   Output parameters:
8431b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8531b70465SJoe Wallwork 
8631b70465SJoe Wallwork   Level: beginner
8731b70465SJoe Wallwork 
8893520041SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst()
89cb7bfe8cSJoe Wallwork @*/
9031b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
9131b70465SJoe Wallwork {
9231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
9331b70465SJoe Wallwork   PetscErrorCode ierr;
9431b70465SJoe Wallwork 
9531b70465SJoe Wallwork   PetscFunctionBegin;
9631b70465SJoe Wallwork   if (!plex->metricCtx) {
9731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
9831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
9931b70465SJoe Wallwork   }
10031b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
10131b70465SJoe Wallwork   PetscFunctionReturn(0);
10231b70465SJoe Wallwork }
10331b70465SJoe Wallwork 
104cb7bfe8cSJoe Wallwork /*@
10593520041SJoe Wallwork   DMPlexMetricSetUniform - Record whether a metric is uniform
10693520041SJoe Wallwork 
10793520041SJoe Wallwork   Input parameters:
10893520041SJoe Wallwork + dm      - The DM
10993520041SJoe Wallwork - uniform - Is the metric uniform?
11093520041SJoe Wallwork 
11193520041SJoe Wallwork   Level: beginner
11293520041SJoe Wallwork 
11393520041SJoe Wallwork   Notes:
11493520041SJoe Wallwork 
11593520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11693520041SJoe Wallwork 
11793520041SJoe Wallwork .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
11893520041SJoe Wallwork @*/
11993520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
12093520041SJoe Wallwork {
12193520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
12293520041SJoe Wallwork   PetscErrorCode ierr;
12393520041SJoe Wallwork 
12493520041SJoe Wallwork   PetscFunctionBegin;
12593520041SJoe Wallwork   if (!plex->metricCtx) {
12693520041SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
12793520041SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
12893520041SJoe Wallwork   }
12993520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
13093520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
13193520041SJoe Wallwork   PetscFunctionReturn(0);
13293520041SJoe Wallwork }
13393520041SJoe Wallwork 
13493520041SJoe Wallwork /*@
13593520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
13693520041SJoe Wallwork 
13793520041SJoe Wallwork   Input parameters:
13893520041SJoe Wallwork . dm      - The DM
13993520041SJoe Wallwork 
14093520041SJoe Wallwork   Output parameters:
14193520041SJoe Wallwork . uniform - Is the metric uniform?
14293520041SJoe Wallwork 
14393520041SJoe Wallwork   Level: beginner
14493520041SJoe Wallwork 
14593520041SJoe Wallwork .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
14693520041SJoe Wallwork @*/
14793520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
14893520041SJoe Wallwork {
14993520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
15093520041SJoe Wallwork   PetscErrorCode ierr;
15193520041SJoe Wallwork 
15293520041SJoe Wallwork   PetscFunctionBegin;
15393520041SJoe Wallwork   if (!plex->metricCtx) {
15493520041SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
15593520041SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
15693520041SJoe Wallwork   }
15793520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
15893520041SJoe Wallwork   PetscFunctionReturn(0);
15993520041SJoe Wallwork }
16093520041SJoe Wallwork 
16193520041SJoe Wallwork /*@
16231b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
16331b70465SJoe Wallwork 
16431b70465SJoe Wallwork   Input parameters:
16531b70465SJoe Wallwork + dm                      - The DM
16631b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
16731b70465SJoe Wallwork 
16831b70465SJoe Wallwork   Level: beginner
16931b70465SJoe Wallwork 
17031b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
171cb7bfe8cSJoe Wallwork @*/
17231b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
17331b70465SJoe Wallwork {
17431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
17531b70465SJoe Wallwork   PetscErrorCode ierr;
17631b70465SJoe Wallwork 
17731b70465SJoe Wallwork   PetscFunctionBegin;
17831b70465SJoe Wallwork   if (!plex->metricCtx) {
17931b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
18031b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
18131b70465SJoe Wallwork   }
18231b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
18331b70465SJoe Wallwork   PetscFunctionReturn(0);
18431b70465SJoe Wallwork }
18531b70465SJoe Wallwork 
186cb7bfe8cSJoe Wallwork /*@
18731b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
18831b70465SJoe Wallwork 
18931b70465SJoe Wallwork   Input parameters:
19031b70465SJoe Wallwork . dm                      - The DM
19131b70465SJoe Wallwork 
19231b70465SJoe Wallwork   Output parameters:
19331b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
19431b70465SJoe Wallwork 
19531b70465SJoe Wallwork   Level: beginner
19631b70465SJoe Wallwork 
19731b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
198cb7bfe8cSJoe Wallwork @*/
19931b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
20031b70465SJoe Wallwork {
20131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
20231b70465SJoe Wallwork   PetscErrorCode ierr;
20331b70465SJoe Wallwork 
20431b70465SJoe Wallwork   PetscFunctionBegin;
20531b70465SJoe Wallwork   if (!plex->metricCtx) {
20631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
20731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
20831b70465SJoe Wallwork   }
20931b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
21031b70465SJoe Wallwork   PetscFunctionReturn(0);
21131b70465SJoe Wallwork }
21231b70465SJoe Wallwork 
213cb7bfe8cSJoe Wallwork /*@
214cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
215cc2eee55SJoe Wallwork 
216cc2eee55SJoe Wallwork   Input parameters:
217cc2eee55SJoe Wallwork + dm       - The DM
218cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
219cc2eee55SJoe Wallwork 
220cc2eee55SJoe Wallwork   Level: beginner
221cc2eee55SJoe Wallwork 
222cb7bfe8cSJoe Wallwork   Notes:
223cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
224cb7bfe8cSJoe Wallwork 
225cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement()
226cb7bfe8cSJoe Wallwork @*/
227cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
228cc2eee55SJoe Wallwork {
229cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
230cc2eee55SJoe Wallwork   PetscErrorCode ierr;
231cc2eee55SJoe Wallwork 
232cc2eee55SJoe Wallwork   PetscFunctionBegin;
233cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
234cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
235cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
236cc2eee55SJoe Wallwork   }
237cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
238cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
239cc2eee55SJoe Wallwork }
240cc2eee55SJoe Wallwork 
241cb7bfe8cSJoe Wallwork /*@
242cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
243cc2eee55SJoe Wallwork 
244cc2eee55SJoe Wallwork   Input parameters:
245cc2eee55SJoe Wallwork . dm       - The DM
246cc2eee55SJoe Wallwork 
247cc2eee55SJoe Wallwork   Output parameters:
248cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
249cc2eee55SJoe Wallwork 
250cc2eee55SJoe Wallwork   Level: beginner
251cc2eee55SJoe Wallwork 
252cb7bfe8cSJoe Wallwork   Notes:
253cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
254cb7bfe8cSJoe Wallwork 
255cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement()
256cb7bfe8cSJoe Wallwork @*/
257cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
258cc2eee55SJoe Wallwork {
259cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
260cc2eee55SJoe Wallwork   PetscErrorCode ierr;
261cc2eee55SJoe Wallwork 
262cc2eee55SJoe Wallwork   PetscFunctionBegin;
263cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
264cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
265cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
266cc2eee55SJoe Wallwork   }
267cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
268cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
269cc2eee55SJoe Wallwork }
270cc2eee55SJoe Wallwork 
271cb7bfe8cSJoe Wallwork /*@
272cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
273cc2eee55SJoe Wallwork 
274cc2eee55SJoe Wallwork   Input parameters:
275cc2eee55SJoe Wallwork + dm     - The DM
276cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
277cc2eee55SJoe Wallwork 
278cc2eee55SJoe Wallwork   Level: beginner
279cc2eee55SJoe Wallwork 
280cb7bfe8cSJoe Wallwork   Notes:
281cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
282cb7bfe8cSJoe Wallwork 
283cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement()
284cb7bfe8cSJoe Wallwork @*/
285cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
286cc2eee55SJoe Wallwork {
287cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
288cc2eee55SJoe Wallwork   PetscErrorCode ierr;
289cc2eee55SJoe Wallwork 
290cc2eee55SJoe Wallwork   PetscFunctionBegin;
291cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
292cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
293cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
294cc2eee55SJoe Wallwork   }
295cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
296cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
297cc2eee55SJoe Wallwork }
298cc2eee55SJoe Wallwork 
299cb7bfe8cSJoe Wallwork /*@
300cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
301cc2eee55SJoe Wallwork 
302cc2eee55SJoe Wallwork   Input parameters:
303cc2eee55SJoe Wallwork . dm     - The DM
304cc2eee55SJoe Wallwork 
305cc2eee55SJoe Wallwork   Output parameters:
306cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
307cc2eee55SJoe Wallwork 
308cc2eee55SJoe Wallwork   Level: beginner
309cc2eee55SJoe Wallwork 
310cb7bfe8cSJoe Wallwork   Notes:
311cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
312cb7bfe8cSJoe Wallwork 
313cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement()
314cb7bfe8cSJoe Wallwork @*/
315cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
316cc2eee55SJoe Wallwork {
317cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
318cc2eee55SJoe Wallwork   PetscErrorCode ierr;
319cc2eee55SJoe Wallwork 
320cc2eee55SJoe Wallwork   PetscFunctionBegin;
321cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
322cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
323cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
324cc2eee55SJoe Wallwork   }
325cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
326cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
327cc2eee55SJoe Wallwork }
328cc2eee55SJoe Wallwork 
329cb7bfe8cSJoe Wallwork /*@
330cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
331cc2eee55SJoe Wallwork 
332cc2eee55SJoe Wallwork   Input parameters:
333cc2eee55SJoe Wallwork + dm     - The DM
334cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
335cc2eee55SJoe Wallwork 
336cc2eee55SJoe Wallwork   Level: beginner
337cc2eee55SJoe Wallwork 
338cb7bfe8cSJoe Wallwork   Notes:
339cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
340cb7bfe8cSJoe Wallwork 
341cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping()
342cb7bfe8cSJoe Wallwork @*/
343cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
344cc2eee55SJoe Wallwork {
345cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
346cc2eee55SJoe Wallwork   PetscErrorCode ierr;
347cc2eee55SJoe Wallwork 
348cc2eee55SJoe Wallwork   PetscFunctionBegin;
349cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
350cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
351cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
352cc2eee55SJoe Wallwork   }
353cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
354cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
355cc2eee55SJoe Wallwork }
356cc2eee55SJoe Wallwork 
357cb7bfe8cSJoe Wallwork /*@
358cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
359cc2eee55SJoe Wallwork 
360cc2eee55SJoe Wallwork   Input parameters:
361cc2eee55SJoe Wallwork . dm     - The DM
362cc2eee55SJoe Wallwork 
363cc2eee55SJoe Wallwork   Output parameters:
364cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
365cc2eee55SJoe Wallwork 
366cc2eee55SJoe Wallwork   Level: beginner
367cc2eee55SJoe Wallwork 
368cb7bfe8cSJoe Wallwork   Notes:
369cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
370cb7bfe8cSJoe Wallwork 
371cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping()
372cb7bfe8cSJoe Wallwork @*/
373cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
374cc2eee55SJoe Wallwork {
375cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
376cc2eee55SJoe Wallwork   PetscErrorCode ierr;
377cc2eee55SJoe Wallwork 
378cc2eee55SJoe Wallwork   PetscFunctionBegin;
379cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
380cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
381cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
382cc2eee55SJoe Wallwork   }
383cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
384cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
385cc2eee55SJoe Wallwork }
386cc2eee55SJoe Wallwork 
387cb7bfe8cSJoe Wallwork /*@
38831b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
38931b70465SJoe Wallwork 
39031b70465SJoe Wallwork   Input parameters:
39131b70465SJoe Wallwork + dm    - The DM
39231b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
39331b70465SJoe Wallwork 
39431b70465SJoe Wallwork   Level: beginner
39531b70465SJoe Wallwork 
39631b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
397cb7bfe8cSJoe Wallwork @*/
39831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
39931b70465SJoe Wallwork {
40031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
40131b70465SJoe Wallwork   PetscErrorCode ierr;
40231b70465SJoe Wallwork 
40331b70465SJoe Wallwork   PetscFunctionBegin;
40431b70465SJoe Wallwork   if (!plex->metricCtx) {
40531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
40631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
40731b70465SJoe Wallwork   }
4082c71b3e2SJacob Faibussowitsch   PetscCheckFalse(h_min <= 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min);
40931b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
41031b70465SJoe Wallwork   PetscFunctionReturn(0);
41131b70465SJoe Wallwork }
41231b70465SJoe Wallwork 
413cb7bfe8cSJoe Wallwork /*@
41431b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
41531b70465SJoe Wallwork 
41631b70465SJoe Wallwork   Input parameters:
41731b70465SJoe Wallwork . dm    - The DM
41831b70465SJoe Wallwork 
419cc2eee55SJoe Wallwork   Output parameters:
42031b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
42131b70465SJoe Wallwork 
42231b70465SJoe Wallwork   Level: beginner
42331b70465SJoe Wallwork 
42431b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
425cb7bfe8cSJoe Wallwork @*/
42631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
42731b70465SJoe Wallwork {
42831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
42931b70465SJoe Wallwork   PetscErrorCode ierr;
43031b70465SJoe Wallwork 
43131b70465SJoe Wallwork   PetscFunctionBegin;
43231b70465SJoe Wallwork   if (!plex->metricCtx) {
43331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
43431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
43531b70465SJoe Wallwork   }
43631b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
43731b70465SJoe Wallwork   PetscFunctionReturn(0);
43831b70465SJoe Wallwork }
43931b70465SJoe Wallwork 
440cb7bfe8cSJoe Wallwork /*@
44131b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
44231b70465SJoe Wallwork 
44331b70465SJoe Wallwork   Input parameters:
44431b70465SJoe Wallwork + dm    - The DM
44531b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44631b70465SJoe Wallwork 
44731b70465SJoe Wallwork   Level: beginner
44831b70465SJoe Wallwork 
44931b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
450cb7bfe8cSJoe Wallwork @*/
45131b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
45231b70465SJoe Wallwork {
45331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
45431b70465SJoe Wallwork   PetscErrorCode ierr;
45531b70465SJoe Wallwork 
45631b70465SJoe Wallwork   PetscFunctionBegin;
45731b70465SJoe Wallwork   if (!plex->metricCtx) {
45831b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
45931b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
46031b70465SJoe Wallwork   }
4612c71b3e2SJacob Faibussowitsch   PetscCheckFalse(h_max <= 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max);
46231b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
46331b70465SJoe Wallwork   PetscFunctionReturn(0);
46431b70465SJoe Wallwork }
46531b70465SJoe Wallwork 
466cb7bfe8cSJoe Wallwork /*@
46731b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
46831b70465SJoe Wallwork 
46931b70465SJoe Wallwork   Input parameters:
47031b70465SJoe Wallwork . dm    - The DM
47131b70465SJoe Wallwork 
472cc2eee55SJoe Wallwork   Output parameters:
47331b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
47431b70465SJoe Wallwork 
47531b70465SJoe Wallwork   Level: beginner
47631b70465SJoe Wallwork 
47731b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
478cb7bfe8cSJoe Wallwork @*/
47931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
48031b70465SJoe Wallwork {
48131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
48231b70465SJoe Wallwork   PetscErrorCode ierr;
48331b70465SJoe Wallwork 
48431b70465SJoe Wallwork   PetscFunctionBegin;
48531b70465SJoe Wallwork   if (!plex->metricCtx) {
48631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
48731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
48831b70465SJoe Wallwork   }
48931b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
49031b70465SJoe Wallwork   PetscFunctionReturn(0);
49131b70465SJoe Wallwork }
49231b70465SJoe Wallwork 
493cb7bfe8cSJoe Wallwork /*@
49431b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
49531b70465SJoe Wallwork 
49631b70465SJoe Wallwork   Input parameters:
49731b70465SJoe Wallwork + dm    - The DM
49831b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
49931b70465SJoe Wallwork 
50031b70465SJoe Wallwork   Level: beginner
50131b70465SJoe Wallwork 
50231b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
50331b70465SJoe Wallwork 
50431b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
505cb7bfe8cSJoe Wallwork @*/
50631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
50731b70465SJoe Wallwork {
50831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
50931b70465SJoe Wallwork   PetscErrorCode ierr;
51031b70465SJoe Wallwork 
51131b70465SJoe Wallwork   PetscFunctionBegin;
51231b70465SJoe Wallwork   if (!plex->metricCtx) {
51331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
51431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
51531b70465SJoe Wallwork   }
5162c71b3e2SJacob Faibussowitsch   PetscCheckFalse(a_max < 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max);
51731b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
51831b70465SJoe Wallwork   PetscFunctionReturn(0);
51931b70465SJoe Wallwork }
52031b70465SJoe Wallwork 
521cb7bfe8cSJoe Wallwork /*@
52231b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
52331b70465SJoe Wallwork 
52431b70465SJoe Wallwork   Input parameters:
52531b70465SJoe Wallwork . dm    - The DM
52631b70465SJoe Wallwork 
527cc2eee55SJoe Wallwork   Output parameters:
52831b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
52931b70465SJoe Wallwork 
53031b70465SJoe Wallwork   Level: beginner
53131b70465SJoe Wallwork 
53231b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
533cb7bfe8cSJoe Wallwork @*/
53431b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
53531b70465SJoe Wallwork {
53631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
53731b70465SJoe Wallwork   PetscErrorCode ierr;
53831b70465SJoe Wallwork 
53931b70465SJoe Wallwork   PetscFunctionBegin;
54031b70465SJoe Wallwork   if (!plex->metricCtx) {
54131b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
54231b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
54331b70465SJoe Wallwork   }
54431b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
54531b70465SJoe Wallwork   PetscFunctionReturn(0);
54631b70465SJoe Wallwork }
54731b70465SJoe Wallwork 
548cb7bfe8cSJoe Wallwork /*@
54931b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
55031b70465SJoe Wallwork 
55131b70465SJoe Wallwork   Input parameters:
55231b70465SJoe Wallwork + dm               - The DM
55331b70465SJoe Wallwork - targetComplexity - The target metric complexity
55431b70465SJoe Wallwork 
55531b70465SJoe Wallwork   Level: beginner
55631b70465SJoe Wallwork 
55731b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
558cb7bfe8cSJoe Wallwork @*/
55931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
56031b70465SJoe Wallwork {
56131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
56231b70465SJoe Wallwork   PetscErrorCode ierr;
56331b70465SJoe Wallwork 
56431b70465SJoe Wallwork   PetscFunctionBegin;
56531b70465SJoe Wallwork   if (!plex->metricCtx) {
56631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
56731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
56831b70465SJoe Wallwork   }
5692c71b3e2SJacob Faibussowitsch   PetscCheckFalse(targetComplexity <= 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive");
57031b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
57131b70465SJoe Wallwork   PetscFunctionReturn(0);
57231b70465SJoe Wallwork }
57331b70465SJoe Wallwork 
574cb7bfe8cSJoe Wallwork /*@
57531b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
57631b70465SJoe Wallwork 
57731b70465SJoe Wallwork   Input parameters:
57831b70465SJoe Wallwork . dm               - The DM
57931b70465SJoe Wallwork 
580cc2eee55SJoe Wallwork   Output parameters:
58131b70465SJoe Wallwork . targetComplexity - The target metric complexity
58231b70465SJoe Wallwork 
58331b70465SJoe Wallwork   Level: beginner
58431b70465SJoe Wallwork 
58531b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
586cb7bfe8cSJoe Wallwork @*/
58731b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
58831b70465SJoe Wallwork {
58931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
59031b70465SJoe Wallwork   PetscErrorCode ierr;
59131b70465SJoe Wallwork 
59231b70465SJoe Wallwork   PetscFunctionBegin;
59331b70465SJoe Wallwork   if (!plex->metricCtx) {
59431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
59531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
59631b70465SJoe Wallwork   }
59731b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
59831b70465SJoe Wallwork   PetscFunctionReturn(0);
59931b70465SJoe Wallwork }
60031b70465SJoe Wallwork 
601cb7bfe8cSJoe Wallwork /*@
60231b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
60331b70465SJoe Wallwork 
60431b70465SJoe Wallwork   Input parameters:
60531b70465SJoe Wallwork + dm - The DM
60631b70465SJoe Wallwork - p  - The normalization order
60731b70465SJoe Wallwork 
60831b70465SJoe Wallwork   Level: beginner
60931b70465SJoe Wallwork 
61031b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
611cb7bfe8cSJoe Wallwork @*/
61231b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
61331b70465SJoe Wallwork {
61431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
61531b70465SJoe Wallwork   PetscErrorCode ierr;
61631b70465SJoe Wallwork 
61731b70465SJoe Wallwork   PetscFunctionBegin;
61831b70465SJoe Wallwork   if (!plex->metricCtx) {
61931b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
62031b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
62131b70465SJoe Wallwork   }
6222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(p < 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater");
62331b70465SJoe Wallwork   plex->metricCtx->p = p;
62431b70465SJoe Wallwork   PetscFunctionReturn(0);
62531b70465SJoe Wallwork }
62631b70465SJoe Wallwork 
627cb7bfe8cSJoe Wallwork /*@
62831b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
62931b70465SJoe Wallwork 
63031b70465SJoe Wallwork   Input parameters:
63131b70465SJoe Wallwork . dm - The DM
63231b70465SJoe Wallwork 
633cc2eee55SJoe Wallwork   Output parameters:
63431b70465SJoe Wallwork . p - The normalization order
63531b70465SJoe Wallwork 
63631b70465SJoe Wallwork   Level: beginner
63731b70465SJoe Wallwork 
63831b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
639cb7bfe8cSJoe Wallwork @*/
64031b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
64131b70465SJoe Wallwork {
64231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
64331b70465SJoe Wallwork   PetscErrorCode ierr;
64431b70465SJoe Wallwork 
64531b70465SJoe Wallwork   PetscFunctionBegin;
64631b70465SJoe Wallwork   if (!plex->metricCtx) {
64731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
64831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
64931b70465SJoe Wallwork   }
65031b70465SJoe Wallwork   *p = plex->metricCtx->p;
65131b70465SJoe Wallwork   PetscFunctionReturn(0);
65231b70465SJoe Wallwork }
65331b70465SJoe Wallwork 
654cb7bfe8cSJoe Wallwork /*@
655cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
656cc2eee55SJoe Wallwork 
657cc2eee55SJoe Wallwork   Input parameters:
658cc2eee55SJoe Wallwork + dm   - The DM
659cc2eee55SJoe Wallwork - beta - The metric gradation factor
660cc2eee55SJoe Wallwork 
661cc2eee55SJoe Wallwork   Level: beginner
662cc2eee55SJoe Wallwork 
663cc2eee55SJoe Wallwork   Notes:
664cc2eee55SJoe Wallwork 
665cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
666cc2eee55SJoe Wallwork 
667cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
668cc2eee55SJoe Wallwork 
669cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
670cb7bfe8cSJoe Wallwork 
671*ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber()
672cb7bfe8cSJoe Wallwork @*/
673cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
674cc2eee55SJoe Wallwork {
675cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
676cc2eee55SJoe Wallwork   PetscErrorCode ierr;
677cc2eee55SJoe Wallwork 
678cc2eee55SJoe Wallwork   PetscFunctionBegin;
679cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
680cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
681cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
682cc2eee55SJoe Wallwork   }
683cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
684cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
685cc2eee55SJoe Wallwork }
686cc2eee55SJoe Wallwork 
687cb7bfe8cSJoe Wallwork /*@
688cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
689cc2eee55SJoe Wallwork 
690cc2eee55SJoe Wallwork   Input parameters:
691cc2eee55SJoe Wallwork . dm   - The DM
692cc2eee55SJoe Wallwork 
693cc2eee55SJoe Wallwork   Output parameters:
694cc2eee55SJoe Wallwork . beta - The metric gradation factor
695cc2eee55SJoe Wallwork 
696cc2eee55SJoe Wallwork   Level: beginner
697cc2eee55SJoe Wallwork 
698cb7bfe8cSJoe Wallwork   Notes:
699cb7bfe8cSJoe Wallwork 
700cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
701cb7bfe8cSJoe Wallwork 
702cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
703cb7bfe8cSJoe Wallwork 
704cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
705cc2eee55SJoe Wallwork 
706*ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber()
707cb7bfe8cSJoe Wallwork @*/
708cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
709cc2eee55SJoe Wallwork {
710cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
711cc2eee55SJoe Wallwork   PetscErrorCode ierr;
712cc2eee55SJoe Wallwork 
713cc2eee55SJoe Wallwork   PetscFunctionBegin;
714cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
715cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
716cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
717cc2eee55SJoe Wallwork   }
718cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
719cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
720cc2eee55SJoe Wallwork }
721cc2eee55SJoe Wallwork 
722cb7bfe8cSJoe Wallwork /*@
723*ae8b063eSJoe Wallwork   DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
724*ae8b063eSJoe Wallwork 
725*ae8b063eSJoe Wallwork   Input parameters:
726*ae8b063eSJoe Wallwork + dm    - The DM
727*ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number
728*ae8b063eSJoe Wallwork 
729*ae8b063eSJoe Wallwork   Level: beginner
730*ae8b063eSJoe Wallwork 
731*ae8b063eSJoe Wallwork   Notes:
732*ae8b063eSJoe Wallwork 
733*ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
734*ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
735*ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
736*ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
737*ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
738*ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
739*ae8b063eSJoe Wallwork 
740*ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
741*ae8b063eSJoe Wallwork 
742*ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber()
743*ae8b063eSJoe Wallwork @*/
744*ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
745*ae8b063eSJoe Wallwork {
746*ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
747*ae8b063eSJoe Wallwork   PetscErrorCode ierr;
748*ae8b063eSJoe Wallwork 
749*ae8b063eSJoe Wallwork   PetscFunctionBegin;
750*ae8b063eSJoe Wallwork   if (!plex->metricCtx) {
751*ae8b063eSJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
752*ae8b063eSJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
753*ae8b063eSJoe Wallwork   }
754*ae8b063eSJoe Wallwork   plex->metricCtx->hausdorffNumber = hausd;
755*ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
756*ae8b063eSJoe Wallwork }
757*ae8b063eSJoe Wallwork 
758*ae8b063eSJoe Wallwork /*@
759*ae8b063eSJoe Wallwork   DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
760*ae8b063eSJoe Wallwork 
761*ae8b063eSJoe Wallwork   Input parameters:
762*ae8b063eSJoe Wallwork . dm    - The DM
763*ae8b063eSJoe Wallwork 
764*ae8b063eSJoe Wallwork   Output parameters:
765*ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number
766*ae8b063eSJoe Wallwork 
767*ae8b063eSJoe Wallwork   Level: beginner
768*ae8b063eSJoe Wallwork 
769*ae8b063eSJoe Wallwork   Notes:
770*ae8b063eSJoe Wallwork 
771*ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
772*ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
773*ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
774*ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
775*ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
776*ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
777*ae8b063eSJoe Wallwork 
778*ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
779*ae8b063eSJoe Wallwork 
780*ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber()
781*ae8b063eSJoe Wallwork @*/
782*ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
783*ae8b063eSJoe Wallwork {
784*ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
785*ae8b063eSJoe Wallwork   PetscErrorCode ierr;
786*ae8b063eSJoe Wallwork 
787*ae8b063eSJoe Wallwork   PetscFunctionBegin;
788*ae8b063eSJoe Wallwork   if (!plex->metricCtx) {
789*ae8b063eSJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
790*ae8b063eSJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
791*ae8b063eSJoe Wallwork   }
792*ae8b063eSJoe Wallwork   *hausd = plex->metricCtx->hausdorffNumber;
793*ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
794*ae8b063eSJoe Wallwork }
795*ae8b063eSJoe Wallwork 
796*ae8b063eSJoe Wallwork /*@
797cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
798cc2eee55SJoe Wallwork 
799cc2eee55SJoe Wallwork   Input parameters:
800cc2eee55SJoe Wallwork + dm        - The DM
801cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
802cc2eee55SJoe Wallwork 
803cb7bfe8cSJoe Wallwork   Level: beginner
804cb7bfe8cSJoe Wallwork 
805cb7bfe8cSJoe Wallwork   Notes:
806cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
807cb7bfe8cSJoe Wallwork 
808cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations()
809cb7bfe8cSJoe Wallwork @*/
810cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
811cc2eee55SJoe Wallwork {
812cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
813cc2eee55SJoe Wallwork   PetscErrorCode ierr;
814cc2eee55SJoe Wallwork 
815cc2eee55SJoe Wallwork   PetscFunctionBegin;
816cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
817cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
818cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
819cc2eee55SJoe Wallwork   }
820cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
821cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
822cc2eee55SJoe Wallwork }
823cc2eee55SJoe Wallwork 
824cb7bfe8cSJoe Wallwork /*@
825cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
826cc2eee55SJoe Wallwork 
827cc2eee55SJoe Wallwork   Input parameters:
828cc2eee55SJoe Wallwork . dm        - The DM
829cc2eee55SJoe Wallwork 
830cc2eee55SJoe Wallwork   Output parameters:
831cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
832cc2eee55SJoe Wallwork 
833cb7bfe8cSJoe Wallwork   Level: beginner
834cb7bfe8cSJoe Wallwork 
835cb7bfe8cSJoe Wallwork   Notes:
836cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
837cb7bfe8cSJoe Wallwork 
838cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
839cb7bfe8cSJoe Wallwork @*/
840cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
841cc2eee55SJoe Wallwork {
842cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
843cc2eee55SJoe Wallwork   PetscErrorCode ierr;
844cc2eee55SJoe Wallwork 
845cc2eee55SJoe Wallwork   PetscFunctionBegin;
846cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
847cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
848cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
849cc2eee55SJoe Wallwork   }
850cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
851cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
852cc2eee55SJoe Wallwork }
853cc2eee55SJoe Wallwork 
854cb7bfe8cSJoe Wallwork /*@
855cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
856cc2eee55SJoe Wallwork 
857cc2eee55SJoe Wallwork   Input parameters:
858cc2eee55SJoe Wallwork + dm      - The DM
859cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
860cc2eee55SJoe Wallwork 
861cb7bfe8cSJoe Wallwork   Level: beginner
862cb7bfe8cSJoe Wallwork 
863cb7bfe8cSJoe Wallwork   Notes:
864cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
865cc2eee55SJoe Wallwork 
866cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
867cb7bfe8cSJoe Wallwork @*/
868cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
869cc2eee55SJoe Wallwork {
870cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
871cc2eee55SJoe Wallwork   PetscErrorCode ierr;
872cc2eee55SJoe Wallwork 
873cc2eee55SJoe Wallwork   PetscFunctionBegin;
874cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
875cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
876cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
877cc2eee55SJoe Wallwork   }
878cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
879cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
880cc2eee55SJoe Wallwork }
881cc2eee55SJoe Wallwork 
882cb7bfe8cSJoe Wallwork /*@
883cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
884cc2eee55SJoe Wallwork 
885cc2eee55SJoe Wallwork   Input parameters:
886cc2eee55SJoe Wallwork . dm      - The DM
887cc2eee55SJoe Wallwork 
888cc2eee55SJoe Wallwork   Output parameters:
889cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
890cc2eee55SJoe Wallwork 
891cb7bfe8cSJoe Wallwork   Level: beginner
892cb7bfe8cSJoe Wallwork 
893cb7bfe8cSJoe Wallwork   Notes:
894cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
895cc2eee55SJoe Wallwork 
896cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity()
897cb7bfe8cSJoe Wallwork @*/
898cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
899cc2eee55SJoe Wallwork {
900cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
901cc2eee55SJoe Wallwork   PetscErrorCode ierr;
902cc2eee55SJoe Wallwork 
903cc2eee55SJoe Wallwork   PetscFunctionBegin;
904cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
905cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
906cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
907cc2eee55SJoe Wallwork   }
908cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
909cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
910cc2eee55SJoe Wallwork }
911cc2eee55SJoe Wallwork 
91220282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
91320282da2SJoe Wallwork {
91420282da2SJoe Wallwork   MPI_Comm       comm;
91520282da2SJoe Wallwork   PetscErrorCode ierr;
91620282da2SJoe Wallwork   PetscFE        fe;
91720282da2SJoe Wallwork   PetscInt       dim;
91820282da2SJoe Wallwork 
91920282da2SJoe Wallwork   PetscFunctionBegin;
92020282da2SJoe Wallwork 
92120282da2SJoe Wallwork   /* Extract metadata from dm */
92220282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
92320282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
92420282da2SJoe Wallwork 
92520282da2SJoe Wallwork   /* Create a P1 field of the requested size */
92620282da2SJoe Wallwork   ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
92720282da2SJoe Wallwork   ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr);
92820282da2SJoe Wallwork   ierr = DMCreateDS(dm);CHKERRQ(ierr);
92920282da2SJoe Wallwork   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
93020282da2SJoe Wallwork   ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr);
93120282da2SJoe Wallwork 
93220282da2SJoe Wallwork   PetscFunctionReturn(0);
93320282da2SJoe Wallwork }
93420282da2SJoe Wallwork 
935cb7bfe8cSJoe Wallwork /*@
93620282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
93720282da2SJoe Wallwork 
93820282da2SJoe Wallwork   Input parameters:
93920282da2SJoe Wallwork + dm     - The DM
94020282da2SJoe Wallwork - f      - The field number to use
94120282da2SJoe Wallwork 
94220282da2SJoe Wallwork   Output parameter:
94320282da2SJoe Wallwork . metric - The metric
94420282da2SJoe Wallwork 
94520282da2SJoe Wallwork   Level: beginner
94620282da2SJoe Wallwork 
94731b70465SJoe Wallwork   Notes:
94831b70465SJoe Wallwork 
94931b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
95031b70465SJoe Wallwork 
95131b70465SJoe Wallwork   Command line options for Riemannian metrics:
95231b70465SJoe Wallwork 
953cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
95493520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
955cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
956cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
957cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
958cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
959cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
96067b8a455SSatish Balay - -dm_plex_metric_target_complexity         - Target metric complexity
961cb7bfe8cSJoe Wallwork 
962cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
963cb7bfe8cSJoe Wallwork 
96467b8a455SSatish Balay . -dm_adaptor <pragmatic/mmg/parmmg>        - specify dm adaptor to use
965cb7bfe8cSJoe Wallwork 
966cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
967cb7bfe8cSJoe Wallwork 
968cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
969cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
970cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
971cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
972cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
973cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
97420282da2SJoe Wallwork 
97520282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
976cb7bfe8cSJoe Wallwork @*/
97720282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
97820282da2SJoe Wallwork {
97931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
98093520041SJoe Wallwork   PetscBool      isotropic, uniform;
98120282da2SJoe Wallwork   PetscErrorCode ierr;
98220282da2SJoe Wallwork   PetscInt       coordDim, Nd;
98320282da2SJoe Wallwork 
98420282da2SJoe Wallwork   PetscFunctionBegin;
98531b70465SJoe Wallwork   if (!plex->metricCtx) {
98631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
98731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
98831b70465SJoe Wallwork   }
98920282da2SJoe Wallwork   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
99020282da2SJoe Wallwork   Nd = coordDim*coordDim;
99193520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
99293520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
99393520041SJoe Wallwork   if (uniform) {
99493520041SJoe Wallwork     MPI_Comm comm;
99593520041SJoe Wallwork 
99693520041SJoe Wallwork     ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
9972c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
99893520041SJoe Wallwork     ierr = VecCreate(comm, metric);CHKERRQ(ierr);
99993520041SJoe Wallwork     ierr = VecSetSizes(*metric, 1, PETSC_DECIDE);CHKERRQ(ierr);
100093520041SJoe Wallwork     ierr = VecSetFromOptions(*metric);CHKERRQ(ierr);
100193520041SJoe Wallwork   } else if (isotropic) {
100293520041SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dm, f, 1, metric);CHKERRQ(ierr);
100393520041SJoe Wallwork   } else {
100420282da2SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr);
100593520041SJoe Wallwork   }
100620282da2SJoe Wallwork   PetscFunctionReturn(0);
100720282da2SJoe Wallwork }
100820282da2SJoe Wallwork 
1009cb7bfe8cSJoe Wallwork /*@
101020282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
101120282da2SJoe Wallwork 
101220282da2SJoe Wallwork   Input parameters:
101320282da2SJoe Wallwork + dm     - The DM
101420282da2SJoe Wallwork . f      - The field number to use
101520282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
101620282da2SJoe Wallwork 
101720282da2SJoe Wallwork   Output parameter:
101820282da2SJoe Wallwork . metric - The uniform metric
101920282da2SJoe Wallwork 
102020282da2SJoe Wallwork   Level: beginner
102120282da2SJoe Wallwork 
102293520041SJoe Wallwork   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
102320282da2SJoe Wallwork 
102420282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
1025cb7bfe8cSJoe Wallwork @*/
102620282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
102720282da2SJoe Wallwork {
102820282da2SJoe Wallwork   PetscErrorCode ierr;
102920282da2SJoe Wallwork 
103020282da2SJoe Wallwork   PetscFunctionBegin;
103193520041SJoe Wallwork   ierr = DMPlexMetricSetUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
103220282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
10332c71b3e2SJacob Faibussowitsch   PetscCheck(alpha,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
10342c71b3e2SJacob Faibussowitsch   PetscCheck(alpha >= 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
103593520041SJoe Wallwork   ierr = VecSet(*metric, alpha);CHKERRQ(ierr);
103693520041SJoe Wallwork   ierr = VecAssemblyBegin(*metric);CHKERRQ(ierr);
103793520041SJoe Wallwork   ierr = VecAssemblyEnd(*metric);CHKERRQ(ierr);
103820282da2SJoe Wallwork   PetscFunctionReturn(0);
103920282da2SJoe Wallwork }
104020282da2SJoe Wallwork 
104193520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
104293520041SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
104393520041SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
104493520041SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
104593520041SJoe Wallwork {
104693520041SJoe Wallwork   f0[0] = u[0];
104793520041SJoe Wallwork }
104893520041SJoe Wallwork 
1049cb7bfe8cSJoe Wallwork /*@
105020282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
105120282da2SJoe Wallwork 
105220282da2SJoe Wallwork   Input parameters:
105320282da2SJoe Wallwork + dm        - The DM
105420282da2SJoe Wallwork . f         - The field number to use
105520282da2SJoe Wallwork - indicator - The error indicator
105620282da2SJoe Wallwork 
105720282da2SJoe Wallwork   Output parameter:
105820282da2SJoe Wallwork . metric    - The isotropic metric
105920282da2SJoe Wallwork 
106020282da2SJoe Wallwork   Level: beginner
106120282da2SJoe Wallwork 
106220282da2SJoe Wallwork   Notes:
106320282da2SJoe Wallwork 
106420282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
106520282da2SJoe Wallwork 
106693520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
106720282da2SJoe Wallwork 
106820282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
1069cb7bfe8cSJoe Wallwork @*/
107020282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
107120282da2SJoe Wallwork {
107220282da2SJoe Wallwork   PetscErrorCode ierr;
107393520041SJoe Wallwork   PetscInt       m, n;
107420282da2SJoe Wallwork 
107520282da2SJoe Wallwork   PetscFunctionBegin;
107693520041SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, PETSC_TRUE);CHKERRQ(ierr);
107720282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
107893520041SJoe Wallwork   ierr = VecGetSize(indicator, &m);CHKERRQ(ierr);
107993520041SJoe Wallwork   ierr = VecGetSize(*metric, &n);CHKERRQ(ierr);
108093520041SJoe Wallwork   if (m == n) { ierr = VecCopy(indicator, *metric);CHKERRQ(ierr); }
108193520041SJoe Wallwork   else {
108293520041SJoe Wallwork     DM     dmIndi;
108393520041SJoe Wallwork     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
108493520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
108593520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
108693520041SJoe Wallwork                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
108793520041SJoe Wallwork 
108820282da2SJoe Wallwork     ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr);
108993520041SJoe Wallwork     funcs[0] = identity;
109093520041SJoe Wallwork     ierr = DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric);CHKERRQ(ierr);
109120282da2SJoe Wallwork   }
109220282da2SJoe Wallwork   PetscFunctionReturn(0);
109320282da2SJoe Wallwork }
109420282da2SJoe Wallwork 
109582490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
109682490253SJoe Wallwork {
109782490253SJoe Wallwork   PetscInt i, j;
109882490253SJoe Wallwork 
109982490253SJoe Wallwork   PetscFunctionBegin;
110082490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
110182490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
110282490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
110382490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
110482490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
110582490253SJoe Wallwork       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]);
110682490253SJoe Wallwork       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]);
110782490253SJoe Wallwork     }
110882490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
110982490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
111082490253SJoe Wallwork   }
111182490253SJoe Wallwork   PetscFunctionReturn(0);
111282490253SJoe Wallwork }
111382490253SJoe Wallwork 
11143f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
111520282da2SJoe Wallwork {
111620282da2SJoe Wallwork   PetscErrorCode ierr;
111720282da2SJoe Wallwork   PetscInt       i, j, k;
111820282da2SJoe 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);
111920282da2SJoe Wallwork   PetscScalar   *Mpos;
112020282da2SJoe Wallwork 
112120282da2SJoe Wallwork   PetscFunctionBegin;
112220282da2SJoe Wallwork   ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr);
112320282da2SJoe Wallwork 
112420282da2SJoe Wallwork   /* Symmetrize */
112520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
112620282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
112720282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
112820282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
112920282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
113020282da2SJoe Wallwork     }
113120282da2SJoe Wallwork   }
113220282da2SJoe Wallwork 
113320282da2SJoe Wallwork   /* Compute eigendecomposition */
113493520041SJoe Wallwork   if (dim == 1) {
113593520041SJoe Wallwork 
113693520041SJoe Wallwork     /* Isotropic case */
113793520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
113893520041SJoe Wallwork     Mpos[0] = 1.0;
113993520041SJoe Wallwork   } else {
114093520041SJoe Wallwork 
114193520041SJoe Wallwork     /* Anisotropic case */
114220282da2SJoe Wallwork     PetscScalar  *work;
114320282da2SJoe Wallwork     PetscBLASInt lwork;
114420282da2SJoe Wallwork 
114520282da2SJoe Wallwork     lwork = 5*dim;
114620282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
114720282da2SJoe Wallwork     {
114820282da2SJoe Wallwork       PetscBLASInt lierr;
114920282da2SJoe Wallwork       PetscBLASInt nb;
115020282da2SJoe Wallwork 
115120282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
115220282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
115320282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
115420282da2SJoe Wallwork       {
115520282da2SJoe Wallwork         PetscReal *rwork;
115620282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
115720282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
115820282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
115920282da2SJoe Wallwork       }
116020282da2SJoe Wallwork #else
116120282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
116220282da2SJoe Wallwork #endif
116382490253SJoe Wallwork       if (lierr) {
116482490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
116582490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
116682490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
116782490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
116882490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
116982490253SJoe Wallwork           }
117082490253SJoe Wallwork         }
117182490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
117298921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
117382490253SJoe Wallwork       }
117420282da2SJoe Wallwork       ierr = PetscFPTrapPop();CHKERRQ(ierr);
117520282da2SJoe Wallwork     }
117620282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
117720282da2SJoe Wallwork   }
117820282da2SJoe Wallwork 
117920282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
118020282da2SJoe Wallwork   max_eig = 0.0;
118120282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
118220282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
118320282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
118420282da2SJoe Wallwork   }
118520282da2SJoe Wallwork 
11863f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11873f07679eSJoe Wallwork   *detMp = 1.0;
118820282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
118920282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
11903f07679eSJoe Wallwork     *detMp *= eigs[i];
119120282da2SJoe Wallwork   }
119220282da2SJoe Wallwork 
119320282da2SJoe Wallwork   /* Reconstruct Hessian */
119420282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
119520282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
119620282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
119720282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
119820282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
119920282da2SJoe Wallwork       }
120020282da2SJoe Wallwork     }
120120282da2SJoe Wallwork   }
120220282da2SJoe Wallwork   ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr);
120320282da2SJoe Wallwork 
120420282da2SJoe Wallwork   PetscFunctionReturn(0);
120520282da2SJoe Wallwork }
120620282da2SJoe Wallwork 
1207cb7bfe8cSJoe Wallwork /*@
120820282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
120920282da2SJoe Wallwork 
121020282da2SJoe Wallwork   Input parameters:
121120282da2SJoe Wallwork + dm                 - The DM
12123f07679eSJoe Wallwork . metricIn           - The metric
121399abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
12143f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
121520282da2SJoe Wallwork 
121620282da2SJoe Wallwork   Output parameter:
12173f07679eSJoe Wallwork + metricOut          - The metric
12183f07679eSJoe Wallwork - determinant        - Its determinant
121920282da2SJoe Wallwork 
122020282da2SJoe Wallwork   Level: beginner
122120282da2SJoe Wallwork 
1222cb7bfe8cSJoe Wallwork   Notes:
1223cb7bfe8cSJoe Wallwork 
1224cb7bfe8cSJoe Wallwork   Relevant command line options:
1225cb7bfe8cSJoe Wallwork 
122693520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
122793520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
122893520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1229cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1230cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1231cb7bfe8cSJoe Wallwork 
123220282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
1233cb7bfe8cSJoe Wallwork @*/
12343f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant)
123520282da2SJoe Wallwork {
12363f07679eSJoe Wallwork   DM             dmDet;
123793520041SJoe Wallwork   PetscBool      isotropic, uniform;
123820282da2SJoe Wallwork   PetscErrorCode ierr;
123920282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
12403f07679eSJoe Wallwork   PetscScalar   *met, *det;
124120282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
124220282da2SJoe Wallwork 
124320282da2SJoe Wallwork   PetscFunctionBegin;
1244fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0);CHKERRQ(ierr);
124520282da2SJoe Wallwork 
124620282da2SJoe Wallwork   /* Extract metadata from dm */
124720282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
124820282da2SJoe Wallwork   if (restrictSizes) {
124999abec2bSJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
125099abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
125199abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
125299abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
125354c59aa7SJacob Faibussowitsch     PetscCheck(h_min < h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
125499abec2bSJoe Wallwork   }
125599abec2bSJoe Wallwork   if (restrictAnisotropy) {
125699abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
125799abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
125820282da2SJoe Wallwork   }
125920282da2SJoe Wallwork 
126093520041SJoe Wallwork   /* Setup output metric */
12613f07679eSJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr);
12623f07679eSJoe Wallwork   ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr);
126393520041SJoe Wallwork 
126493520041SJoe Wallwork   /* Enforce SPD and extract determinant */
126593520041SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
126693520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
126793520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
126893520041SJoe Wallwork   if (uniform) {
1269d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
127093520041SJoe Wallwork 
127193520041SJoe Wallwork     /* Uniform case */
127293520041SJoe Wallwork     ierr = VecDuplicate(metricIn, determinant);CHKERRQ(ierr);
127393520041SJoe Wallwork     ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
127493520041SJoe Wallwork     ierr = DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det);CHKERRQ(ierr);
127593520041SJoe Wallwork     ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
127693520041SJoe Wallwork   } else {
127793520041SJoe Wallwork 
127893520041SJoe Wallwork     /* Spatially varying case */
127993520041SJoe Wallwork     PetscInt nrow;
128093520041SJoe Wallwork 
128193520041SJoe Wallwork     if (isotropic) nrow = 1;
128293520041SJoe Wallwork     else nrow = dim;
12833f07679eSJoe Wallwork     ierr = DMClone(dm, &dmDet);CHKERRQ(ierr);
12843f07679eSJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant);CHKERRQ(ierr);
128520282da2SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12863f07679eSJoe Wallwork     ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
128720282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12883f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
128920282da2SJoe Wallwork       ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
12903f07679eSJoe Wallwork       ierr = DMPlexPointLocalRef(dmDet, v, det, &vdet);CHKERRQ(ierr);
129193520041SJoe Wallwork       ierr = DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet);CHKERRQ(ierr);
129220282da2SJoe Wallwork     }
12933f07679eSJoe Wallwork     ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
129493520041SJoe Wallwork   }
12953f07679eSJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
1296fe902aceSJoe Wallwork 
1297fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0);CHKERRQ(ierr);
129820282da2SJoe Wallwork   PetscFunctionReturn(0);
129920282da2SJoe Wallwork }
130020282da2SJoe Wallwork 
130120282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
130220282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
130320282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130420282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
130520282da2SJoe Wallwork {
130620282da2SJoe Wallwork   const PetscScalar p = constants[0];
130720282da2SJoe Wallwork 
13083f07679eSJoe Wallwork   f0[0] = PetscPowReal(u[0], p/(2.0*p + dim));
130920282da2SJoe Wallwork }
131020282da2SJoe Wallwork 
1311cb7bfe8cSJoe Wallwork /*@
131220282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
131320282da2SJoe Wallwork 
131420282da2SJoe Wallwork   Input parameters:
131520282da2SJoe Wallwork + dm                 - The DM
131620282da2SJoe Wallwork . metricIn           - The unnormalized metric
131716522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
131816522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
131920282da2SJoe Wallwork 
132020282da2SJoe Wallwork   Output parameter:
132120282da2SJoe Wallwork . metricOut          - The normalized metric
132220282da2SJoe Wallwork 
132320282da2SJoe Wallwork   Level: beginner
132420282da2SJoe Wallwork 
1325cb7bfe8cSJoe Wallwork   Notes:
1326cb7bfe8cSJoe Wallwork 
1327cb7bfe8cSJoe Wallwork   Relevant command line options:
1328cb7bfe8cSJoe Wallwork 
132993520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
133093520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
133193520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1332cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1333cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1334cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1335cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1336cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1337cb7bfe8cSJoe Wallwork 
133820282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
1339cb7bfe8cSJoe Wallwork @*/
134016522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut)
134120282da2SJoe Wallwork {
13423f07679eSJoe Wallwork   DM               dmDet;
134320282da2SJoe Wallwork   MPI_Comm         comm;
134493520041SJoe Wallwork   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
134520282da2SJoe Wallwork   PetscDS          ds;
134620282da2SJoe Wallwork   PetscErrorCode   ierr;
134720282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
13483f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
134993520041SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
13503f07679eSJoe Wallwork   Vec              determinant;
135120282da2SJoe Wallwork 
135220282da2SJoe Wallwork   PetscFunctionBegin;
1353fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0);CHKERRQ(ierr);
135420282da2SJoe Wallwork 
135520282da2SJoe Wallwork   /* Extract metadata from dm */
135620282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
135720282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
135893520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
135993520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
136093520041SJoe Wallwork   if (isotropic) Nd = 1;
136193520041SJoe Wallwork   else Nd = dim*dim;
136220282da2SJoe Wallwork 
136320282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
136416522936SJoe Wallwork   ierr = DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst);CHKERRQ(ierr);
13653f07679eSJoe Wallwork   ierr = DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant);CHKERRQ(ierr);
136620282da2SJoe Wallwork 
136720282da2SJoe Wallwork   /* Compute global normalization factor */
136816522936SJoe Wallwork   ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
136916522936SJoe Wallwork   ierr = DMPlexMetricGetNormalizationOrder(dm, &p);CHKERRQ(ierr);
137016522936SJoe Wallwork   constants[0] = p;
137193520041SJoe Wallwork   if (uniform) {
13722c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
137393520041SJoe Wallwork     DM  dmTmp;
137493520041SJoe Wallwork     Vec tmp;
137593520041SJoe Wallwork 
137693520041SJoe Wallwork     ierr = DMClone(dm, &dmTmp);CHKERRQ(ierr);
137793520041SJoe Wallwork     ierr = DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp);CHKERRQ(ierr);
137893520041SJoe Wallwork     ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
137993520041SJoe Wallwork     ierr = VecSet(tmp, det[0]);CHKERRQ(ierr);
138093520041SJoe Wallwork     ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
138193520041SJoe Wallwork     ierr = DMGetDS(dmTmp, &ds);CHKERRQ(ierr);
138293520041SJoe Wallwork     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
138393520041SJoe Wallwork     ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
138493520041SJoe Wallwork     ierr = DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL);CHKERRQ(ierr);
138593520041SJoe Wallwork     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
138693520041SJoe Wallwork     ierr = DMDestroy(&dmTmp);CHKERRQ(ierr);
138793520041SJoe Wallwork   } else {
13883f07679eSJoe Wallwork     ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr);
13893f07679eSJoe Wallwork     ierr = DMGetDS(dmDet, &ds);CHKERRQ(ierr);
139020282da2SJoe Wallwork     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
139120282da2SJoe Wallwork     ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
13923f07679eSJoe Wallwork     ierr = DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL);CHKERRQ(ierr);
139393520041SJoe Wallwork   }
13943f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
13952c71b3e2SJacob Faibussowitsch   PetscCheckFalse(realIntegral < 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global metric normalization factor should be strictly positive, not %.4e Is the input metric positive-definite?", realIntegral);
13963f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
139720282da2SJoe Wallwork 
139820282da2SJoe Wallwork   /* Apply local scaling */
139920282da2SJoe Wallwork   if (restrictSizes) {
140016522936SJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
140116522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
140216522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
140316522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
14042c71b3e2SJacob Faibussowitsch     PetscCheckFalse(h_min >= h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
140516522936SJoe Wallwork   }
140616522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
140716522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
140816522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
140920282da2SJoe Wallwork   }
141020282da2SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
14113f07679eSJoe Wallwork   ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
141293520041SJoe Wallwork   if (uniform) {
141393520041SJoe Wallwork 
141493520041SJoe Wallwork     /* Uniform case */
141593520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
141693520041SJoe Wallwork     if (restrictSizes) { ierr = DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det);CHKERRQ(ierr); }
141793520041SJoe Wallwork   } else {
141893520041SJoe Wallwork 
141993520041SJoe Wallwork     /* Spatially varying case */
142093520041SJoe Wallwork     PetscInt nrow;
142193520041SJoe Wallwork 
142293520041SJoe Wallwork     if (isotropic) nrow = 1;
142393520041SJoe Wallwork     else nrow = dim;
142416522936SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
142520282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
14263f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
142720282da2SJoe Wallwork 
142820282da2SJoe Wallwork       ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr);
14293f07679eSJoe Wallwork       ierr = DMPlexPointLocalRef(dmDet, v, det, &detM);CHKERRQ(ierr);
14303f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
143120282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
143293520041SJoe Wallwork       if (restrictSizes) { ierr = DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM);CHKERRQ(ierr); }
143393520041SJoe Wallwork     }
143420282da2SJoe Wallwork   }
14353f07679eSJoe Wallwork   ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
143620282da2SJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
14373f07679eSJoe Wallwork   ierr = VecDestroy(&determinant);CHKERRQ(ierr);
143893520041SJoe Wallwork   if (!uniform) { ierr = DMDestroy(&dmDet);CHKERRQ(ierr); }
143920282da2SJoe Wallwork 
1440fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0);CHKERRQ(ierr);
144120282da2SJoe Wallwork   PetscFunctionReturn(0);
144220282da2SJoe Wallwork }
144320282da2SJoe Wallwork 
1444cb7bfe8cSJoe Wallwork /*@
144520282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
144620282da2SJoe Wallwork 
1447f1a722f8SMatthew G. Knepley   Input Parameters:
144820282da2SJoe Wallwork + dm         - The DM
144920282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
145020282da2SJoe Wallwork . weights    - Weights for the average
145120282da2SJoe Wallwork - metrics    - The metrics to be averaged
145220282da2SJoe Wallwork 
145320282da2SJoe Wallwork   Output Parameter:
145420282da2SJoe Wallwork . metricAvg  - The averaged metric
145520282da2SJoe Wallwork 
145620282da2SJoe Wallwork   Level: beginner
145720282da2SJoe Wallwork 
145820282da2SJoe Wallwork   Notes:
145920282da2SJoe Wallwork   The weights should sum to unity.
146020282da2SJoe Wallwork 
146120282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
146220282da2SJoe Wallwork 
146320282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
1464cb7bfe8cSJoe Wallwork @*/
146520282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
146620282da2SJoe Wallwork {
146720282da2SJoe Wallwork   PetscBool      haveWeights = PETSC_TRUE;
146820282da2SJoe Wallwork   PetscErrorCode ierr;
146993520041SJoe Wallwork   PetscInt       i, m, n;
147020282da2SJoe Wallwork   PetscReal      sum = 0.0, tol = 1.0e-10;
147120282da2SJoe Wallwork 
147220282da2SJoe Wallwork   PetscFunctionBegin;
1473fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0);CHKERRQ(ierr);
14742c71b3e2SJacob Faibussowitsch   PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics);
147520282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr);
147620282da2SJoe Wallwork   ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr);
147793520041SJoe Wallwork   ierr = VecGetSize(*metricAvg, &m);CHKERRQ(ierr);
147893520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
147993520041SJoe Wallwork     ierr = VecGetSize(metrics[i], &n);CHKERRQ(ierr);
14802c71b3e2SJacob Faibussowitsch     PetscCheckFalse(m != n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
148193520041SJoe Wallwork   }
148220282da2SJoe Wallwork 
148320282da2SJoe Wallwork   /* Default to the unweighted case */
148420282da2SJoe Wallwork   if (!weights) {
148520282da2SJoe Wallwork     ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr);
148620282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
148720282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
148820282da2SJoe Wallwork   }
148920282da2SJoe Wallwork 
149020282da2SJoe Wallwork   /* Check weights sum to unity */
149193520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
14922c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsReal(sum - 1) > tol,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
149320282da2SJoe Wallwork 
149420282da2SJoe Wallwork   /* Compute metric average */
149520282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); }
149620282da2SJoe Wallwork   if (!haveWeights) { ierr = PetscFree(weights); }
1497fe902aceSJoe Wallwork 
1498fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0);CHKERRQ(ierr);
149920282da2SJoe Wallwork   PetscFunctionReturn(0);
150020282da2SJoe Wallwork }
150120282da2SJoe Wallwork 
1502cb7bfe8cSJoe Wallwork /*@
150320282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
150420282da2SJoe Wallwork 
1505f1a722f8SMatthew G. Knepley   Input Parameters:
150620282da2SJoe Wallwork + dm         - The DM
150720282da2SJoe Wallwork . metric1    - The first metric to be averaged
150820282da2SJoe Wallwork - metric2    - The second metric to be averaged
150920282da2SJoe Wallwork 
151020282da2SJoe Wallwork   Output Parameter:
151120282da2SJoe Wallwork . metricAvg  - The averaged metric
151220282da2SJoe Wallwork 
151320282da2SJoe Wallwork   Level: beginner
151420282da2SJoe Wallwork 
151520282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
1516cb7bfe8cSJoe Wallwork @*/
151720282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
151820282da2SJoe Wallwork {
151920282da2SJoe Wallwork   PetscErrorCode ierr;
152020282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
152120282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
152220282da2SJoe Wallwork 
152320282da2SJoe Wallwork   PetscFunctionBegin;
152420282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr);
152520282da2SJoe Wallwork   PetscFunctionReturn(0);
152620282da2SJoe Wallwork }
152720282da2SJoe Wallwork 
1528cb7bfe8cSJoe Wallwork /*@
152920282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
153020282da2SJoe Wallwork 
1531f1a722f8SMatthew G. Knepley   Input Parameters:
153220282da2SJoe Wallwork + dm         - The DM
153320282da2SJoe Wallwork . metric1    - The first metric to be averaged
153420282da2SJoe Wallwork . metric2    - The second metric to be averaged
153520282da2SJoe Wallwork - metric3    - The third metric to be averaged
153620282da2SJoe Wallwork 
153720282da2SJoe Wallwork   Output Parameter:
153820282da2SJoe Wallwork . metricAvg  - The averaged metric
153920282da2SJoe Wallwork 
154020282da2SJoe Wallwork   Level: beginner
154120282da2SJoe Wallwork 
154220282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
1543cb7bfe8cSJoe Wallwork @*/
154420282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
154520282da2SJoe Wallwork {
154620282da2SJoe Wallwork   PetscErrorCode ierr;
154720282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
154820282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
154920282da2SJoe Wallwork 
155020282da2SJoe Wallwork   PetscFunctionBegin;
155120282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr);
155220282da2SJoe Wallwork   PetscFunctionReturn(0);
155320282da2SJoe Wallwork }
155420282da2SJoe Wallwork 
155520282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
155620282da2SJoe Wallwork {
155720282da2SJoe Wallwork   PetscErrorCode ierr;
155820282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
155920282da2SJoe Wallwork   PetscReal     *evals, *evals1;
156020282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
156120282da2SJoe Wallwork 
156220282da2SJoe Wallwork   PetscFunctionBegin;
156393520041SJoe Wallwork 
156493520041SJoe Wallwork   /* Isotropic case */
156593520041SJoe Wallwork   if (dim == 1) {
156693520041SJoe Wallwork     M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
156793520041SJoe Wallwork     PetscFunctionReturn(0);
156893520041SJoe Wallwork   }
156993520041SJoe Wallwork 
157093520041SJoe Wallwork   /* Anisotropic case */
157120282da2SJoe Wallwork   ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr);
157220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
157320282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
157420282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
157520282da2SJoe Wallwork     }
157620282da2SJoe Wallwork   }
157720282da2SJoe Wallwork   {
157820282da2SJoe Wallwork     PetscScalar *work;
157920282da2SJoe Wallwork     PetscBLASInt lwork;
158020282da2SJoe Wallwork 
158120282da2SJoe Wallwork     lwork = 5*dim;
158220282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
158320282da2SJoe Wallwork     {
158420282da2SJoe Wallwork       PetscBLASInt lierr, nb;
158520282da2SJoe Wallwork       PetscReal    sqrtk;
158620282da2SJoe Wallwork 
158720282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
158820282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
158920282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
159020282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
159120282da2SJoe Wallwork       {
159220282da2SJoe Wallwork         PetscReal *rwork;
159320282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
159420282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
159520282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
159620282da2SJoe Wallwork       }
159720282da2SJoe Wallwork #else
159820282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
159920282da2SJoe Wallwork #endif
160082490253SJoe Wallwork       if (lierr) {
160182490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
160298921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
160382490253SJoe Wallwork       }
160420282da2SJoe Wallwork       ierr = PetscFPTrapPop();
160520282da2SJoe Wallwork 
160620282da2SJoe Wallwork       /* Compute square root and reciprocal */
160720282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
160820282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
160920282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
161020282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
161120282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
161220282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
161320282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
161420282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
161520282da2SJoe Wallwork           }
161620282da2SJoe Wallwork         }
161720282da2SJoe Wallwork       }
161820282da2SJoe Wallwork 
161920282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
162020282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
162120282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
162220282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
162320282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
162420282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
162520282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
162620282da2SJoe Wallwork             }
162720282da2SJoe Wallwork           }
162820282da2SJoe Wallwork         }
162920282da2SJoe Wallwork       }
163020282da2SJoe Wallwork 
163120282da2SJoe Wallwork       /* Compute eigendecomposition */
163220282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
163320282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
163420282da2SJoe Wallwork       {
163520282da2SJoe Wallwork         PetscReal *rwork;
163620282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
163720282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
163820282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
163920282da2SJoe Wallwork       }
164020282da2SJoe Wallwork #else
164120282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
164220282da2SJoe Wallwork #endif
164382490253SJoe Wallwork       if (lierr) {
164482490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
164582490253SJoe Wallwork           for (j = 0; j < dim; ++j) {
164682490253SJoe Wallwork             evecs[i*dim+j] = 0.0;
164782490253SJoe Wallwork             for (k = 0; k < dim; ++k) {
164882490253SJoe Wallwork               for (l = 0; l < dim; ++l) {
164982490253SJoe Wallwork                 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
165082490253SJoe Wallwork               }
165182490253SJoe Wallwork             }
165282490253SJoe Wallwork           }
165382490253SJoe Wallwork         }
165482490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
165598921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
165682490253SJoe Wallwork       }
165720282da2SJoe Wallwork       ierr = PetscFPTrapPop();
165820282da2SJoe Wallwork 
165920282da2SJoe Wallwork       /* Modify eigenvalues */
166020282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
166120282da2SJoe Wallwork 
166220282da2SJoe Wallwork       /* Map back to get the intersection */
166320282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
166420282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
166520282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
166620282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
166720282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
166820282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
166920282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
167020282da2SJoe Wallwork               }
167120282da2SJoe Wallwork             }
167220282da2SJoe Wallwork           }
167320282da2SJoe Wallwork         }
167420282da2SJoe Wallwork       }
167520282da2SJoe Wallwork     }
167620282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
167720282da2SJoe Wallwork   }
167820282da2SJoe Wallwork   ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr);
167920282da2SJoe Wallwork   PetscFunctionReturn(0);
168020282da2SJoe Wallwork }
168120282da2SJoe Wallwork 
1682cb7bfe8cSJoe Wallwork /*@
168320282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
168420282da2SJoe Wallwork 
1685f1a722f8SMatthew G. Knepley   Input Parameters:
168620282da2SJoe Wallwork + dm         - The DM
168720282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
168820282da2SJoe Wallwork - metrics    - The metrics to be intersected
168920282da2SJoe Wallwork 
169020282da2SJoe Wallwork   Output Parameter:
169120282da2SJoe Wallwork . metricInt  - The intersected metric
169220282da2SJoe Wallwork 
169320282da2SJoe Wallwork   Level: beginner
169420282da2SJoe Wallwork 
169520282da2SJoe Wallwork   Notes:
169620282da2SJoe Wallwork 
169720282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
169820282da2SJoe Wallwork 
169920282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
170020282da2SJoe Wallwork 
170120282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
1702cb7bfe8cSJoe Wallwork @*/
170320282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
170420282da2SJoe Wallwork {
170593520041SJoe Wallwork   PetscBool      isotropic, uniform;
170620282da2SJoe Wallwork   PetscErrorCode ierr;
170793520041SJoe Wallwork   PetscInt       v, i, m, n;
170893520041SJoe Wallwork   PetscScalar   *met, *meti;
170920282da2SJoe Wallwork 
171020282da2SJoe Wallwork   PetscFunctionBegin;
1711fe902aceSJoe Wallwork   ierr = PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0);CHKERRQ(ierr);
17122c71b3e2SJacob Faibussowitsch   PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics);
171320282da2SJoe Wallwork 
171420282da2SJoe Wallwork   /* Copy over the first metric */
171593520041SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr);
171620282da2SJoe Wallwork   ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr);
171793520041SJoe Wallwork   if (numMetrics == 1) PetscFunctionReturn(0);
171893520041SJoe Wallwork   ierr = VecGetSize(*metricInt, &m);CHKERRQ(ierr);
171993520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
172093520041SJoe Wallwork     ierr = VecGetSize(metrics[i], &n);CHKERRQ(ierr);
17212c71b3e2SJacob Faibussowitsch     PetscCheckFalse(m != n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
172293520041SJoe Wallwork   }
172320282da2SJoe Wallwork 
172420282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
172593520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
172693520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
172793520041SJoe Wallwork   if (uniform) {
172893520041SJoe Wallwork 
172993520041SJoe Wallwork     /* Uniform case */
173093520041SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
173193520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
173293520041SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
173393520041SJoe Wallwork       ierr = DMPlexMetricIntersection_Private(1, meti, met);CHKERRQ(ierr);
173493520041SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
173593520041SJoe Wallwork     }
173693520041SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
173793520041SJoe Wallwork   } else {
173893520041SJoe Wallwork 
173993520041SJoe Wallwork     /* Spatially varying case */
174093520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
174193520041SJoe Wallwork     PetscScalar *M, *Mi;
174293520041SJoe Wallwork 
174393520041SJoe Wallwork     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
174493520041SJoe Wallwork     if (isotropic) nrow = 1;
174593520041SJoe Wallwork     else nrow = dim;
174693520041SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
174720282da2SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
174820282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
174920282da2SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
175020282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
175120282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr);
175220282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr);
175393520041SJoe Wallwork         ierr = DMPlexMetricIntersection_Private(nrow, Mi, M);CHKERRQ(ierr);
175420282da2SJoe Wallwork       }
175520282da2SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
175620282da2SJoe Wallwork     }
175720282da2SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
175820282da2SJoe Wallwork   }
1759fe902aceSJoe Wallwork 
1760fe902aceSJoe Wallwork   ierr = PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0);CHKERRQ(ierr);
176120282da2SJoe Wallwork   PetscFunctionReturn(0);
176220282da2SJoe Wallwork }
176320282da2SJoe Wallwork 
1764cb7bfe8cSJoe Wallwork /*@
176520282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
176620282da2SJoe Wallwork 
1767f1a722f8SMatthew G. Knepley   Input Parameters:
176820282da2SJoe Wallwork + dm        - The DM
176920282da2SJoe Wallwork . metric1   - The first metric to be intersected
177020282da2SJoe Wallwork - metric2   - The second metric to be intersected
177120282da2SJoe Wallwork 
177220282da2SJoe Wallwork   Output Parameter:
177320282da2SJoe Wallwork . metricInt - The intersected metric
177420282da2SJoe Wallwork 
177520282da2SJoe Wallwork   Level: beginner
177620282da2SJoe Wallwork 
177720282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
1778cb7bfe8cSJoe Wallwork @*/
177920282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
178020282da2SJoe Wallwork {
178120282da2SJoe Wallwork   PetscErrorCode ierr;
178220282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
178320282da2SJoe Wallwork 
178420282da2SJoe Wallwork   PetscFunctionBegin;
178520282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr);
178620282da2SJoe Wallwork   PetscFunctionReturn(0);
178720282da2SJoe Wallwork }
178820282da2SJoe Wallwork 
1789cb7bfe8cSJoe Wallwork /*@
179020282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
179120282da2SJoe Wallwork 
1792f1a722f8SMatthew G. Knepley   Input Parameters:
179320282da2SJoe Wallwork + dm        - The DM
179420282da2SJoe Wallwork . metric1   - The first metric to be intersected
179520282da2SJoe Wallwork . metric2   - The second metric to be intersected
179620282da2SJoe Wallwork - metric3   - The third metric to be intersected
179720282da2SJoe Wallwork 
179820282da2SJoe Wallwork   Output Parameter:
179920282da2SJoe Wallwork . metricInt - The intersected metric
180020282da2SJoe Wallwork 
180120282da2SJoe Wallwork   Level: beginner
180220282da2SJoe Wallwork 
180320282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
1804cb7bfe8cSJoe Wallwork @*/
180520282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
180620282da2SJoe Wallwork {
180720282da2SJoe Wallwork   PetscErrorCode ierr;
180820282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
180920282da2SJoe Wallwork 
181020282da2SJoe Wallwork   PetscFunctionBegin;
181120282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr);
181220282da2SJoe Wallwork   PetscFunctionReturn(0);
181320282da2SJoe Wallwork }
1814