xref: /petsc/src/dm/impls/plex/plexmetric.c (revision 99abec2b429b33f09a86a719839caa30dbc4af3f)
120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
220282da2SJoe Wallwork #include <petscblaslapack.h>
320282da2SJoe Wallwork 
431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
531b70465SJoe Wallwork {
631b70465SJoe Wallwork   MPI_Comm       comm;
731b70465SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
831b70465SJoe Wallwork   PetscErrorCode ierr;
931b70465SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0;
1031b70465SJoe Wallwork 
1131b70465SJoe Wallwork   PetscFunctionBegin;
1231b70465SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1331b70465SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr);
1431b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr);
1531b70465SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, isotropic);
1631b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr);
1731b70465SJoe Wallwork   ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst);
1831b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr);
1931b70465SJoe Wallwork   ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr);
2031b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr);
2131b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr);
2231b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr);
2331b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr);
2431b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr);
2531b70465SJoe Wallwork   ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr);
2631b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr);
2731b70465SJoe Wallwork   ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr);
2831b70465SJoe Wallwork   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2931b70465SJoe Wallwork   PetscFunctionReturn(0);
3031b70465SJoe Wallwork }
3131b70465SJoe Wallwork 
3231b70465SJoe Wallwork /*
3331b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
3431b70465SJoe Wallwork 
3531b70465SJoe Wallwork   Input parameters:
3631b70465SJoe Wallwork + dm        - The DM
3731b70465SJoe Wallwork - isotropic - Is the metric isotropic?
3831b70465SJoe Wallwork 
3931b70465SJoe Wallwork   Level: beginner
4031b70465SJoe Wallwork 
4131b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
4231b70465SJoe Wallwork */
4331b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
4431b70465SJoe Wallwork {
4531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
4631b70465SJoe Wallwork   PetscErrorCode ierr;
4731b70465SJoe Wallwork 
4831b70465SJoe Wallwork   PetscFunctionBegin;
4931b70465SJoe Wallwork   if (!plex->metricCtx) {
5031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
5131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
5231b70465SJoe Wallwork   }
5331b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
5431b70465SJoe Wallwork   PetscFunctionReturn(0);
5531b70465SJoe Wallwork }
5631b70465SJoe Wallwork 
5731b70465SJoe Wallwork /*
5831b70465SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric is isotropic?
5931b70465SJoe Wallwork 
6031b70465SJoe Wallwork   Input parameters:
6131b70465SJoe Wallwork . dm        - The DM
6231b70465SJoe Wallwork 
6331b70465SJoe Wallwork   Output parameters:
6431b70465SJoe Wallwork . isotropic - Is the metric isotropic?
6531b70465SJoe Wallwork 
6631b70465SJoe Wallwork   Level: beginner
6731b70465SJoe Wallwork 
6831b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
6931b70465SJoe Wallwork */
7031b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
7131b70465SJoe Wallwork {
7231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
7331b70465SJoe Wallwork   PetscErrorCode ierr;
7431b70465SJoe Wallwork 
7531b70465SJoe Wallwork   PetscFunctionBegin;
7631b70465SJoe Wallwork   if (!plex->metricCtx) {
7731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
7831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
7931b70465SJoe Wallwork   }
8031b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
8131b70465SJoe Wallwork   PetscFunctionReturn(0);
8231b70465SJoe Wallwork }
8331b70465SJoe Wallwork 
8431b70465SJoe Wallwork /*
8531b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
8631b70465SJoe Wallwork 
8731b70465SJoe Wallwork   Input parameters:
8831b70465SJoe Wallwork + dm                      - The DM
8931b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
9031b70465SJoe Wallwork 
9131b70465SJoe Wallwork   Level: beginner
9231b70465SJoe Wallwork 
9331b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
9431b70465SJoe Wallwork */
9531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
9631b70465SJoe Wallwork {
9731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
9831b70465SJoe Wallwork   PetscErrorCode ierr;
9931b70465SJoe Wallwork 
10031b70465SJoe Wallwork   PetscFunctionBegin;
10131b70465SJoe Wallwork   if (!plex->metricCtx) {
10231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
10331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
10431b70465SJoe Wallwork   }
10531b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
10631b70465SJoe Wallwork   PetscFunctionReturn(0);
10731b70465SJoe Wallwork }
10831b70465SJoe Wallwork 
10931b70465SJoe Wallwork /*
11031b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
11131b70465SJoe Wallwork 
11231b70465SJoe Wallwork   Input parameters:
11331b70465SJoe Wallwork . dm                      - The DM
11431b70465SJoe Wallwork 
11531b70465SJoe Wallwork   Output parameters:
11631b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
11731b70465SJoe Wallwork 
11831b70465SJoe Wallwork   Level: beginner
11931b70465SJoe Wallwork 
12031b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
12131b70465SJoe Wallwork */
12231b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
12331b70465SJoe Wallwork {
12431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
12531b70465SJoe Wallwork   PetscErrorCode ierr;
12631b70465SJoe Wallwork 
12731b70465SJoe Wallwork   PetscFunctionBegin;
12831b70465SJoe Wallwork   if (!plex->metricCtx) {
12931b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
13031b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
13131b70465SJoe Wallwork   }
13231b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
13331b70465SJoe Wallwork   PetscFunctionReturn(0);
13431b70465SJoe Wallwork }
13531b70465SJoe Wallwork 
13631b70465SJoe Wallwork /*
13731b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
13831b70465SJoe Wallwork 
13931b70465SJoe Wallwork   Input parameters:
14031b70465SJoe Wallwork + dm    - The DM
14131b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
14231b70465SJoe Wallwork 
14331b70465SJoe Wallwork   Level: beginner
14431b70465SJoe Wallwork 
14531b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
14631b70465SJoe Wallwork */
14731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
14831b70465SJoe Wallwork {
14931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
15031b70465SJoe Wallwork   PetscErrorCode ierr;
15131b70465SJoe Wallwork 
15231b70465SJoe Wallwork   PetscFunctionBegin;
15331b70465SJoe Wallwork   if (!plex->metricCtx) {
15431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
15531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
15631b70465SJoe Wallwork   }
15731b70465SJoe Wallwork   if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min);
15831b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
15931b70465SJoe Wallwork   PetscFunctionReturn(0);
16031b70465SJoe Wallwork }
16131b70465SJoe Wallwork 
16231b70465SJoe Wallwork /*
16331b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
16431b70465SJoe Wallwork 
16531b70465SJoe Wallwork   Input parameters:
16631b70465SJoe Wallwork . dm    - The DM
16731b70465SJoe Wallwork 
16831b70465SJoe Wallwork   Input parameters:
16931b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
17031b70465SJoe Wallwork 
17131b70465SJoe Wallwork   Level: beginner
17231b70465SJoe Wallwork 
17331b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
17431b70465SJoe Wallwork */
17531b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
17631b70465SJoe Wallwork {
17731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
17831b70465SJoe Wallwork   PetscErrorCode ierr;
17931b70465SJoe Wallwork 
18031b70465SJoe Wallwork   PetscFunctionBegin;
18131b70465SJoe Wallwork   if (!plex->metricCtx) {
18231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
18331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
18431b70465SJoe Wallwork   }
18531b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
18631b70465SJoe Wallwork   PetscFunctionReturn(0);
18731b70465SJoe Wallwork }
18831b70465SJoe Wallwork 
18931b70465SJoe Wallwork /*
19031b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
19131b70465SJoe Wallwork 
19231b70465SJoe Wallwork   Input parameters:
19331b70465SJoe Wallwork + dm    - The DM
19431b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
19531b70465SJoe Wallwork 
19631b70465SJoe Wallwork   Level: beginner
19731b70465SJoe Wallwork 
19831b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
19931b70465SJoe Wallwork */
20031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
20131b70465SJoe Wallwork {
20231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
20331b70465SJoe Wallwork   PetscErrorCode ierr;
20431b70465SJoe Wallwork 
20531b70465SJoe Wallwork   PetscFunctionBegin;
20631b70465SJoe Wallwork   if (!plex->metricCtx) {
20731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
20831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
20931b70465SJoe Wallwork   }
21031b70465SJoe Wallwork   if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max);
21131b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
21231b70465SJoe Wallwork   PetscFunctionReturn(0);
21331b70465SJoe Wallwork }
21431b70465SJoe Wallwork 
21531b70465SJoe Wallwork /*
21631b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
21731b70465SJoe Wallwork 
21831b70465SJoe Wallwork   Input parameters:
21931b70465SJoe Wallwork . dm    - The DM
22031b70465SJoe Wallwork 
22131b70465SJoe Wallwork   Input parameters:
22231b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
22331b70465SJoe Wallwork 
22431b70465SJoe Wallwork   Level: beginner
22531b70465SJoe Wallwork 
22631b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
22731b70465SJoe Wallwork */
22831b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
22931b70465SJoe Wallwork {
23031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
23131b70465SJoe Wallwork   PetscErrorCode ierr;
23231b70465SJoe Wallwork 
23331b70465SJoe Wallwork   PetscFunctionBegin;
23431b70465SJoe Wallwork   if (!plex->metricCtx) {
23531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
23631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
23731b70465SJoe Wallwork   }
23831b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
23931b70465SJoe Wallwork   PetscFunctionReturn(0);
24031b70465SJoe Wallwork }
24131b70465SJoe Wallwork 
24231b70465SJoe Wallwork /*
24331b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
24431b70465SJoe Wallwork 
24531b70465SJoe Wallwork   Input parameters:
24631b70465SJoe Wallwork + dm    - The DM
24731b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
24831b70465SJoe Wallwork 
24931b70465SJoe Wallwork   Level: beginner
25031b70465SJoe Wallwork 
25131b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
25231b70465SJoe Wallwork 
25331b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
25431b70465SJoe Wallwork */
25531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
25631b70465SJoe Wallwork {
25731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
25831b70465SJoe Wallwork   PetscErrorCode ierr;
25931b70465SJoe Wallwork 
26031b70465SJoe Wallwork   PetscFunctionBegin;
26131b70465SJoe Wallwork   if (!plex->metricCtx) {
26231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
26331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
26431b70465SJoe Wallwork   }
26531b70465SJoe Wallwork   if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max);
26631b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
26731b70465SJoe Wallwork   PetscFunctionReturn(0);
26831b70465SJoe Wallwork }
26931b70465SJoe Wallwork 
27031b70465SJoe Wallwork /*
27131b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
27231b70465SJoe Wallwork 
27331b70465SJoe Wallwork   Input parameters:
27431b70465SJoe Wallwork . dm    - The DM
27531b70465SJoe Wallwork 
27631b70465SJoe Wallwork   Input parameters:
27731b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
27831b70465SJoe Wallwork 
27931b70465SJoe Wallwork   Level: beginner
28031b70465SJoe Wallwork 
28131b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
28231b70465SJoe Wallwork */
28331b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
28431b70465SJoe Wallwork {
28531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
28631b70465SJoe Wallwork   PetscErrorCode ierr;
28731b70465SJoe Wallwork 
28831b70465SJoe Wallwork   PetscFunctionBegin;
28931b70465SJoe Wallwork   if (!plex->metricCtx) {
29031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
29131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
29231b70465SJoe Wallwork   }
29331b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
29431b70465SJoe Wallwork   PetscFunctionReturn(0);
29531b70465SJoe Wallwork }
29631b70465SJoe Wallwork 
29731b70465SJoe Wallwork /*
29831b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
29931b70465SJoe Wallwork 
30031b70465SJoe Wallwork   Input parameters:
30131b70465SJoe Wallwork + dm               - The DM
30231b70465SJoe Wallwork - targetComplexity - The target metric complexity
30331b70465SJoe Wallwork 
30431b70465SJoe Wallwork   Level: beginner
30531b70465SJoe Wallwork 
30631b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
30731b70465SJoe Wallwork */
30831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
30931b70465SJoe Wallwork {
31031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
31131b70465SJoe Wallwork   PetscErrorCode ierr;
31231b70465SJoe Wallwork 
31331b70465SJoe Wallwork   PetscFunctionBegin;
31431b70465SJoe Wallwork   if (!plex->metricCtx) {
31531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
31631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
31731b70465SJoe Wallwork   }
31831b70465SJoe Wallwork   if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive");
31931b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
32031b70465SJoe Wallwork   PetscFunctionReturn(0);
32131b70465SJoe Wallwork }
32231b70465SJoe Wallwork 
32331b70465SJoe Wallwork /*
32431b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
32531b70465SJoe Wallwork 
32631b70465SJoe Wallwork   Input parameters:
32731b70465SJoe Wallwork . dm               - The DM
32831b70465SJoe Wallwork 
32931b70465SJoe Wallwork   Input parameters:
33031b70465SJoe Wallwork . targetComplexity - The target metric complexity
33131b70465SJoe Wallwork 
33231b70465SJoe Wallwork   Level: beginner
33331b70465SJoe Wallwork 
33431b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
33531b70465SJoe Wallwork */
33631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
33731b70465SJoe Wallwork {
33831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
33931b70465SJoe Wallwork   PetscErrorCode ierr;
34031b70465SJoe Wallwork 
34131b70465SJoe Wallwork   PetscFunctionBegin;
34231b70465SJoe Wallwork   if (!plex->metricCtx) {
34331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
34431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
34531b70465SJoe Wallwork   }
34631b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
34731b70465SJoe Wallwork   PetscFunctionReturn(0);
34831b70465SJoe Wallwork }
34931b70465SJoe Wallwork 
35031b70465SJoe Wallwork /*
35131b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
35231b70465SJoe Wallwork 
35331b70465SJoe Wallwork   Input parameters:
35431b70465SJoe Wallwork + dm - The DM
35531b70465SJoe Wallwork - p  - The normalization order
35631b70465SJoe Wallwork 
35731b70465SJoe Wallwork   Level: beginner
35831b70465SJoe Wallwork 
35931b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
36031b70465SJoe Wallwork */
36131b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
36231b70465SJoe Wallwork {
36331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
36431b70465SJoe Wallwork   PetscErrorCode ierr;
36531b70465SJoe Wallwork 
36631b70465SJoe Wallwork   PetscFunctionBegin;
36731b70465SJoe Wallwork   if (!plex->metricCtx) {
36831b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
36931b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
37031b70465SJoe Wallwork   }
37131b70465SJoe Wallwork   if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater");
37231b70465SJoe Wallwork   plex->metricCtx->p = p;
37331b70465SJoe Wallwork   PetscFunctionReturn(0);
37431b70465SJoe Wallwork }
37531b70465SJoe Wallwork 
37631b70465SJoe Wallwork /*
37731b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
37831b70465SJoe Wallwork 
37931b70465SJoe Wallwork   Input parameters:
38031b70465SJoe Wallwork . dm - The DM
38131b70465SJoe Wallwork 
38231b70465SJoe Wallwork   Input parameters:
38331b70465SJoe Wallwork . p - The normalization order
38431b70465SJoe Wallwork 
38531b70465SJoe Wallwork   Level: beginner
38631b70465SJoe Wallwork 
38731b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
38831b70465SJoe Wallwork */
38931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
39031b70465SJoe Wallwork {
39131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
39231b70465SJoe Wallwork   PetscErrorCode ierr;
39331b70465SJoe Wallwork 
39431b70465SJoe Wallwork   PetscFunctionBegin;
39531b70465SJoe Wallwork   if (!plex->metricCtx) {
39631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
39731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
39831b70465SJoe Wallwork   }
39931b70465SJoe Wallwork   *p = plex->metricCtx->p;
40031b70465SJoe Wallwork   PetscFunctionReturn(0);
40131b70465SJoe Wallwork }
40231b70465SJoe Wallwork 
40320282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
40420282da2SJoe Wallwork {
40520282da2SJoe Wallwork   MPI_Comm       comm;
40620282da2SJoe Wallwork   PetscErrorCode ierr;
40720282da2SJoe Wallwork   PetscFE        fe;
40820282da2SJoe Wallwork   PetscInt       dim;
40920282da2SJoe Wallwork 
41020282da2SJoe Wallwork   PetscFunctionBegin;
41120282da2SJoe Wallwork 
41220282da2SJoe Wallwork   /* Extract metadata from dm */
41320282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
41420282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
41520282da2SJoe Wallwork 
41620282da2SJoe Wallwork   /* Create a P1 field of the requested size */
41720282da2SJoe Wallwork   ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
41820282da2SJoe Wallwork   ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr);
41920282da2SJoe Wallwork   ierr = DMCreateDS(dm);CHKERRQ(ierr);
42020282da2SJoe Wallwork   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
42120282da2SJoe Wallwork   ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr);
42220282da2SJoe Wallwork 
42320282da2SJoe Wallwork   PetscFunctionReturn(0);
42420282da2SJoe Wallwork }
42520282da2SJoe Wallwork 
42620282da2SJoe Wallwork /*
42720282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
42820282da2SJoe Wallwork 
42920282da2SJoe Wallwork   Input parameters:
43020282da2SJoe Wallwork + dm     - The DM
43120282da2SJoe Wallwork - f      - The field number to use
43220282da2SJoe Wallwork 
43320282da2SJoe Wallwork   Output parameter:
43420282da2SJoe Wallwork . metric - The metric
43520282da2SJoe Wallwork 
43620282da2SJoe Wallwork   Level: beginner
43720282da2SJoe Wallwork 
43831b70465SJoe Wallwork   Notes:
43931b70465SJoe Wallwork 
44031b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
44131b70465SJoe Wallwork 
44231b70465SJoe Wallwork   Command line options for Riemannian metrics:
44331b70465SJoe Wallwork 
44431b70465SJoe Wallwork   -dm_plex_metric_isotropic                 - Is the metric isotropic?
44531b70465SJoe Wallwork   -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
44631b70465SJoe Wallwork   -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
44731b70465SJoe Wallwork   -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
44831b70465SJoe Wallwork   -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
44931b70465SJoe Wallwork   -dm_plex_metric_p                         - L-p normalization order
45031b70465SJoe Wallwork   -dm_plex_metric_target_complexity         - Target metric complexity
45120282da2SJoe Wallwork 
45220282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
45320282da2SJoe Wallwork */
45420282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
45520282da2SJoe Wallwork {
45631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
45720282da2SJoe Wallwork   PetscErrorCode ierr;
45820282da2SJoe Wallwork   PetscInt       coordDim, Nd;
45920282da2SJoe Wallwork 
46020282da2SJoe Wallwork   PetscFunctionBegin;
46131b70465SJoe Wallwork   if (!plex->metricCtx) {
46231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
46331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
46431b70465SJoe Wallwork   }
46520282da2SJoe Wallwork   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
46620282da2SJoe Wallwork   Nd = coordDim*coordDim;
46720282da2SJoe Wallwork   ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr);
46820282da2SJoe Wallwork   PetscFunctionReturn(0);
46920282da2SJoe Wallwork }
47020282da2SJoe Wallwork 
47120282da2SJoe Wallwork typedef struct {
47220282da2SJoe Wallwork   PetscReal scaling;  /* Scaling for uniform metric diagonal */
47320282da2SJoe Wallwork } DMPlexMetricUniformCtx;
47420282da2SJoe Wallwork 
47520282da2SJoe Wallwork static PetscErrorCode diagonal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
47620282da2SJoe Wallwork {
47720282da2SJoe Wallwork   DMPlexMetricUniformCtx *user = (DMPlexMetricUniformCtx*)ctx;
47820282da2SJoe Wallwork   PetscInt                i, j;
47920282da2SJoe Wallwork 
48020282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
48120282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
48220282da2SJoe Wallwork       if (i == j) u[i+dim*j] = user->scaling;
48320282da2SJoe Wallwork       else u[i+dim*j] = 0.0;
48420282da2SJoe Wallwork     }
48520282da2SJoe Wallwork   }
48620282da2SJoe Wallwork   return 0;
48720282da2SJoe Wallwork }
48820282da2SJoe Wallwork 
48920282da2SJoe Wallwork /*
49020282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
49120282da2SJoe Wallwork 
49220282da2SJoe Wallwork   Input parameters:
49320282da2SJoe Wallwork + dm     - The DM
49420282da2SJoe Wallwork . f      - The field number to use
49520282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
49620282da2SJoe Wallwork 
49720282da2SJoe Wallwork   Output parameter:
49820282da2SJoe Wallwork . metric - The uniform metric
49920282da2SJoe Wallwork 
50020282da2SJoe Wallwork   Level: beginner
50120282da2SJoe Wallwork 
50220282da2SJoe Wallwork   Note: It is assumed that the DM is comprised of simplices.
50320282da2SJoe Wallwork 
50420282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
50520282da2SJoe Wallwork */
50620282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
50720282da2SJoe Wallwork {
50820282da2SJoe Wallwork   DMPlexMetricUniformCtx user;
50920282da2SJoe Wallwork   PetscErrorCode       (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
51020282da2SJoe Wallwork   PetscErrorCode         ierr;
51120282da2SJoe Wallwork   void                  *ctxs[1];
51220282da2SJoe Wallwork 
51320282da2SJoe Wallwork   PetscFunctionBegin;
51420282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
51520282da2SJoe Wallwork   if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
51620282da2SJoe Wallwork   if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
51720282da2SJoe Wallwork   else user.scaling = alpha;
51820282da2SJoe Wallwork   funcs[0] = diagonal;
51920282da2SJoe Wallwork   ctxs[0] = &user;
52020282da2SJoe Wallwork   ierr = DMProjectFunctionLocal(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, *metric);CHKERRQ(ierr);
52120282da2SJoe Wallwork   PetscFunctionReturn(0);
52220282da2SJoe Wallwork }
52320282da2SJoe Wallwork 
52420282da2SJoe Wallwork /*
52520282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
52620282da2SJoe Wallwork 
52720282da2SJoe Wallwork   Input parameters:
52820282da2SJoe Wallwork + dm        - The DM
52920282da2SJoe Wallwork . f         - The field number to use
53020282da2SJoe Wallwork - indicator - The error indicator
53120282da2SJoe Wallwork 
53220282da2SJoe Wallwork   Output parameter:
53320282da2SJoe Wallwork . metric    - The isotropic metric
53420282da2SJoe Wallwork 
53520282da2SJoe Wallwork   Level: beginner
53620282da2SJoe Wallwork 
53720282da2SJoe Wallwork   Notes:
53820282da2SJoe Wallwork 
53920282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
54020282da2SJoe Wallwork 
54120282da2SJoe Wallwork   The indicator needs to be a scalar field defined at *vertices*.
54220282da2SJoe Wallwork 
54320282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
54420282da2SJoe Wallwork */
54520282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
54620282da2SJoe Wallwork {
54720282da2SJoe Wallwork   DM                 dmIndi;
54820282da2SJoe Wallwork   PetscErrorCode     ierr;
54920282da2SJoe Wallwork   PetscInt           dim, d, vStart, vEnd, v;
55020282da2SJoe Wallwork   const PetscScalar *indi;
55120282da2SJoe Wallwork   PetscScalar       *met;
55220282da2SJoe Wallwork 
55320282da2SJoe Wallwork   PetscFunctionBegin;
55420282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
55520282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
55620282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
55720282da2SJoe Wallwork   ierr = VecGetArrayRead(indicator, &indi);CHKERRQ(ierr);
55820282da2SJoe Wallwork   ierr = VecGetArrayWrite(*metric, &met);CHKERRQ(ierr);
55920282da2SJoe Wallwork   ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr);
56020282da2SJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
56120282da2SJoe Wallwork     PetscScalar *vindi, *vmet;
56220282da2SJoe Wallwork     ierr = DMPlexPointLocalRead(dmIndi, v, indi, &vindi);CHKERRQ(ierr);
56320282da2SJoe Wallwork     ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
56420282da2SJoe Wallwork     for (d = 0; d < dim; ++d) vmet[d*(dim+1)] = vindi[0];
56520282da2SJoe Wallwork   }
56620282da2SJoe Wallwork   ierr = VecRestoreArrayWrite(*metric, &met);CHKERRQ(ierr);
56720282da2SJoe Wallwork   ierr = VecRestoreArrayRead(indicator, &indi);CHKERRQ(ierr);
56820282da2SJoe Wallwork   PetscFunctionReturn(0);
56920282da2SJoe Wallwork }
57020282da2SJoe Wallwork 
57120282da2SJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[])
57220282da2SJoe Wallwork {
57320282da2SJoe Wallwork   PetscErrorCode ierr;
57420282da2SJoe Wallwork   PetscInt       i, j, k;
57520282da2SJoe 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);
57620282da2SJoe Wallwork   PetscScalar   *Mpos;
57720282da2SJoe Wallwork 
57820282da2SJoe Wallwork   PetscFunctionBegin;
57920282da2SJoe Wallwork   ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr);
58020282da2SJoe Wallwork 
58120282da2SJoe Wallwork   /* Symmetrize */
58220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
58320282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
58420282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
58520282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
58620282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
58720282da2SJoe Wallwork     }
58820282da2SJoe Wallwork   }
58920282da2SJoe Wallwork 
59020282da2SJoe Wallwork   /* Compute eigendecomposition */
59120282da2SJoe Wallwork   {
59220282da2SJoe Wallwork     PetscScalar  *work;
59320282da2SJoe Wallwork     PetscBLASInt lwork;
59420282da2SJoe Wallwork 
59520282da2SJoe Wallwork     lwork = 5*dim;
59620282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
59720282da2SJoe Wallwork     {
59820282da2SJoe Wallwork       PetscBLASInt lierr;
59920282da2SJoe Wallwork       PetscBLASInt nb;
60020282da2SJoe Wallwork 
60120282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
60220282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
60320282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
60420282da2SJoe Wallwork       {
60520282da2SJoe Wallwork         PetscReal *rwork;
60620282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
60720282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
60820282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
60920282da2SJoe Wallwork       }
61020282da2SJoe Wallwork #else
61120282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
61220282da2SJoe Wallwork #endif
61320282da2SJoe Wallwork       if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
61420282da2SJoe Wallwork       ierr = PetscFPTrapPop();CHKERRQ(ierr);
61520282da2SJoe Wallwork     }
61620282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
61720282da2SJoe Wallwork   }
61820282da2SJoe Wallwork 
61920282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
62020282da2SJoe Wallwork   max_eig = 0.0;
62120282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
62220282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
62320282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
62420282da2SJoe Wallwork   }
62520282da2SJoe Wallwork 
62620282da2SJoe Wallwork   /* Enforce maximum anisotropy */
62720282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
62820282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
62920282da2SJoe Wallwork   }
63020282da2SJoe Wallwork 
63120282da2SJoe Wallwork   /* Reconstruct Hessian */
63220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
63320282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
63420282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
63520282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
63620282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
63720282da2SJoe Wallwork       }
63820282da2SJoe Wallwork     }
63920282da2SJoe Wallwork   }
64020282da2SJoe Wallwork   ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr);
64120282da2SJoe Wallwork 
64220282da2SJoe Wallwork   PetscFunctionReturn(0);
64320282da2SJoe Wallwork }
64420282da2SJoe Wallwork 
64520282da2SJoe Wallwork /*
64620282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
64720282da2SJoe Wallwork 
64820282da2SJoe Wallwork   Input parameters:
64920282da2SJoe Wallwork + dm                 - The DM
650*99abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
651*99abec2bSJoe Wallwork . restrictAnisotropy - Should maximum anisotropy be enforced?
65220282da2SJoe Wallwork - metric             - The metric
65320282da2SJoe Wallwork 
65420282da2SJoe Wallwork   Output parameter:
65520282da2SJoe Wallwork . metric             - The metric
65620282da2SJoe Wallwork 
65720282da2SJoe Wallwork   Level: beginner
65820282da2SJoe Wallwork 
65920282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
66020282da2SJoe Wallwork */
661*99abec2bSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metric)
66220282da2SJoe Wallwork {
66320282da2SJoe Wallwork   PetscErrorCode ierr;
66420282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
66520282da2SJoe Wallwork   PetscScalar   *met;
66620282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
66720282da2SJoe Wallwork 
66820282da2SJoe Wallwork   PetscFunctionBegin;
66920282da2SJoe Wallwork 
67020282da2SJoe Wallwork   /* Extract metadata from dm */
67120282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
67220282da2SJoe Wallwork   if (restrictSizes) {
673*99abec2bSJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
674*99abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
675*99abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
676*99abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
677*99abec2bSJoe Wallwork     if (h_min >= h_max) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
678*99abec2bSJoe Wallwork   }
679*99abec2bSJoe Wallwork   if (restrictAnisotropy) {
680*99abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
681*99abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
68220282da2SJoe Wallwork   }
68320282da2SJoe Wallwork 
68420282da2SJoe Wallwork   /* Enforce SPD */
68520282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
68620282da2SJoe Wallwork   ierr = VecGetArray(metric, &met);CHKERRQ(ierr);
68720282da2SJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
68820282da2SJoe Wallwork     PetscScalar *vmet;
68920282da2SJoe Wallwork     ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
69020282da2SJoe Wallwork     ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, vmet);CHKERRQ(ierr);
69120282da2SJoe Wallwork   }
69220282da2SJoe Wallwork   ierr = VecRestoreArray(metric, &met);CHKERRQ(ierr);
69320282da2SJoe Wallwork 
69420282da2SJoe Wallwork   PetscFunctionReturn(0);
69520282da2SJoe Wallwork }
69620282da2SJoe Wallwork 
69720282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
69820282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69920282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70020282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
70120282da2SJoe Wallwork {
70220282da2SJoe Wallwork   const PetscScalar p = constants[0];
70320282da2SJoe Wallwork   PetscReal         detH = 0.0;
70420282da2SJoe Wallwork 
70520282da2SJoe Wallwork   if      (dim == 2) DMPlex_Det2D_Scalar_Internal(&detH, u);
70620282da2SJoe Wallwork   else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detH, u);
70720282da2SJoe Wallwork   f0[0] = PetscPowReal(detH, p/(2.0*p + dim));
70820282da2SJoe Wallwork }
70920282da2SJoe Wallwork 
71020282da2SJoe Wallwork /*
71120282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
71220282da2SJoe Wallwork 
71320282da2SJoe Wallwork   Input parameters:
71420282da2SJoe Wallwork + dm            - The DM
71520282da2SJoe Wallwork . metricIn      - The unnormalized metric
71620282da2SJoe Wallwork - restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced?
71720282da2SJoe Wallwork 
71820282da2SJoe Wallwork   Output parameter:
71920282da2SJoe Wallwork . metricOut     - The normalized metric
72020282da2SJoe Wallwork 
72120282da2SJoe Wallwork   Level: beginner
72220282da2SJoe Wallwork 
72320282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
72420282da2SJoe Wallwork */
72520282da2SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, Vec *metricOut)
72620282da2SJoe Wallwork {
72720282da2SJoe Wallwork   DMPlexMetricCtx *user;
72820282da2SJoe Wallwork   MPI_Comm         comm;
72920282da2SJoe Wallwork   PetscDS          ds;
73020282da2SJoe Wallwork   PetscErrorCode   ierr;
73120282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
73220282da2SJoe Wallwork   PetscScalar     *met, integral, constants[1];
73320282da2SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, target;
73420282da2SJoe Wallwork 
73520282da2SJoe Wallwork   PetscFunctionBegin;
73620282da2SJoe Wallwork 
73720282da2SJoe Wallwork   /* Extract metadata from dm */
73820282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
73920282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
74020282da2SJoe Wallwork   Nd = dim*dim;
74120282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
74220282da2SJoe Wallwork   ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr);
74320282da2SJoe Wallwork   if (restrictSizes && user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max;
74420282da2SJoe Wallwork   if (PetscAbsReal(user->p) >= 1.0) p = user->p;
74520282da2SJoe Wallwork   else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric normalization order %f should be greater than one.", user->p);
74620282da2SJoe Wallwork   constants[0] = p;
74720282da2SJoe Wallwork   if (user->targetComplexity > 0.0) target = user->targetComplexity;
74820282da2SJoe Wallwork   else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Target metric complexity %f should be positive.", user->targetComplexity);
74920282da2SJoe Wallwork 
75020282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
75120282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr);
75220282da2SJoe Wallwork   ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr);
75320282da2SJoe Wallwork   ierr = DMPlexMetricEnforceSPD(dm, PETSC_FALSE, *metricOut);CHKERRQ(ierr);
75420282da2SJoe Wallwork 
75520282da2SJoe Wallwork   /* Compute global normalization factor */
75620282da2SJoe Wallwork   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
75720282da2SJoe Wallwork   ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
75820282da2SJoe Wallwork   ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
75920282da2SJoe Wallwork   ierr = DMPlexComputeIntegralFEM(dm, *metricOut, &integral, NULL);CHKERRQ(ierr);
76020282da2SJoe Wallwork   factGlob = PetscPowReal(target/PetscRealPart(integral), 2.0/dim);
76120282da2SJoe Wallwork 
76220282da2SJoe Wallwork   /* Apply local scaling */
76320282da2SJoe Wallwork   a_max = 0.0;
76420282da2SJoe Wallwork   if (restrictSizes) {
76520282da2SJoe Wallwork     if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max);
76620282da2SJoe Wallwork     if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min);
76720282da2SJoe Wallwork     if (!user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max;
76820282da2SJoe Wallwork   }
76920282da2SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
77020282da2SJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
77120282da2SJoe Wallwork     PetscScalar       *Mp;
77220282da2SJoe Wallwork     PetscReal          detM, fact;
77320282da2SJoe Wallwork 
77420282da2SJoe Wallwork     ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr);
77520282da2SJoe Wallwork     if      (dim == 2) DMPlex_Det2D_Scalar_Internal(&detM, Mp);
77620282da2SJoe Wallwork     else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detM, Mp);
77720282da2SJoe Wallwork     else SETERRQ1(comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
77820282da2SJoe Wallwork     fact = factGlob * PetscPowReal(PetscAbsReal(detM), -1.0/(2*p+dim));
77920282da2SJoe Wallwork     for (i = 0; i < Nd; ++i) Mp[i] *= fact;
78020282da2SJoe Wallwork     if (restrictSizes) { ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, Mp);CHKERRQ(ierr); }
78120282da2SJoe Wallwork   }
78220282da2SJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
78320282da2SJoe Wallwork 
78420282da2SJoe Wallwork   PetscFunctionReturn(0);
78520282da2SJoe Wallwork }
78620282da2SJoe Wallwork 
78720282da2SJoe Wallwork /*
78820282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
78920282da2SJoe Wallwork 
79020282da2SJoe Wallwork   Input Parameter:
79120282da2SJoe Wallwork + dm         - The DM
79220282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
79320282da2SJoe Wallwork . weights    - Weights for the average
79420282da2SJoe Wallwork - metrics    - The metrics to be averaged
79520282da2SJoe Wallwork 
79620282da2SJoe Wallwork   Output Parameter:
79720282da2SJoe Wallwork . metricAvg  - The averaged metric
79820282da2SJoe Wallwork 
79920282da2SJoe Wallwork   Level: beginner
80020282da2SJoe Wallwork 
80120282da2SJoe Wallwork   Notes:
80220282da2SJoe Wallwork   The weights should sum to unity.
80320282da2SJoe Wallwork 
80420282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
80520282da2SJoe Wallwork 
80620282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
80720282da2SJoe Wallwork */
80820282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
80920282da2SJoe Wallwork {
81020282da2SJoe Wallwork   PetscBool      haveWeights = PETSC_TRUE;
81120282da2SJoe Wallwork   PetscErrorCode ierr;
81220282da2SJoe Wallwork   PetscInt       i;
81320282da2SJoe Wallwork   PetscReal      sum = 0.0, tol = 1.0e-10;
81420282da2SJoe Wallwork 
81520282da2SJoe Wallwork   PetscFunctionBegin;
81620282da2SJoe Wallwork   if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); }
81720282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr);
81820282da2SJoe Wallwork   ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr);
81920282da2SJoe Wallwork 
82020282da2SJoe Wallwork   /* Default to the unweighted case */
82120282da2SJoe Wallwork   if (!weights) {
82220282da2SJoe Wallwork     ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr);
82320282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
82420282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
82520282da2SJoe Wallwork   }
82620282da2SJoe Wallwork 
82720282da2SJoe Wallwork   /* Check weights sum to unity */
82820282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { sum += weights[i]; }
82920282da2SJoe Wallwork   if (PetscAbsReal(sum - 1) > tol) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); }
83020282da2SJoe Wallwork 
83120282da2SJoe Wallwork   /* Compute metric average */
83220282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); }
83320282da2SJoe Wallwork   if (!haveWeights) {ierr = PetscFree(weights); }
83420282da2SJoe Wallwork   PetscFunctionReturn(0);
83520282da2SJoe Wallwork }
83620282da2SJoe Wallwork 
83720282da2SJoe Wallwork /*
83820282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
83920282da2SJoe Wallwork 
84020282da2SJoe Wallwork   Input Parameter:
84120282da2SJoe Wallwork + dm         - The DM
84220282da2SJoe Wallwork . metric1    - The first metric to be averaged
84320282da2SJoe Wallwork - metric2    - The second metric to be averaged
84420282da2SJoe Wallwork 
84520282da2SJoe Wallwork   Output Parameter:
84620282da2SJoe Wallwork . metricAvg  - The averaged metric
84720282da2SJoe Wallwork 
84820282da2SJoe Wallwork   Level: beginner
84920282da2SJoe Wallwork 
85020282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
85120282da2SJoe Wallwork */
85220282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
85320282da2SJoe Wallwork {
85420282da2SJoe Wallwork   PetscErrorCode ierr;
85520282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
85620282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
85720282da2SJoe Wallwork 
85820282da2SJoe Wallwork   PetscFunctionBegin;
85920282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr);
86020282da2SJoe Wallwork   PetscFunctionReturn(0);
86120282da2SJoe Wallwork }
86220282da2SJoe Wallwork 
86320282da2SJoe Wallwork /*
86420282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
86520282da2SJoe Wallwork 
86620282da2SJoe Wallwork   Input Parameter:
86720282da2SJoe Wallwork + dm         - The DM
86820282da2SJoe Wallwork . metric1    - The first metric to be averaged
86920282da2SJoe Wallwork . metric2    - The second metric to be averaged
87020282da2SJoe Wallwork - metric3    - The third metric to be averaged
87120282da2SJoe Wallwork 
87220282da2SJoe Wallwork   Output Parameter:
87320282da2SJoe Wallwork . metricAvg  - The averaged metric
87420282da2SJoe Wallwork 
87520282da2SJoe Wallwork   Level: beginner
87620282da2SJoe Wallwork 
87720282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
87820282da2SJoe Wallwork */
87920282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
88020282da2SJoe Wallwork {
88120282da2SJoe Wallwork   PetscErrorCode ierr;
88220282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
88320282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
88420282da2SJoe Wallwork 
88520282da2SJoe Wallwork   PetscFunctionBegin;
88620282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr);
88720282da2SJoe Wallwork   PetscFunctionReturn(0);
88820282da2SJoe Wallwork }
88920282da2SJoe Wallwork 
89020282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
89120282da2SJoe Wallwork {
89220282da2SJoe Wallwork   PetscErrorCode ierr;
89320282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
89420282da2SJoe Wallwork   PetscReal     *evals, *evals1;
89520282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
89620282da2SJoe Wallwork 
89720282da2SJoe Wallwork   PetscFunctionBegin;
89820282da2SJoe Wallwork   ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr);
89920282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
90020282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
90120282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
90220282da2SJoe Wallwork     }
90320282da2SJoe Wallwork   }
90420282da2SJoe Wallwork   {
90520282da2SJoe Wallwork     PetscScalar *work;
90620282da2SJoe Wallwork     PetscBLASInt lwork;
90720282da2SJoe Wallwork 
90820282da2SJoe Wallwork     lwork = 5*dim;
90920282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
91020282da2SJoe Wallwork     {
91120282da2SJoe Wallwork       PetscBLASInt lierr, nb;
91220282da2SJoe Wallwork       PetscReal    sqrtk;
91320282da2SJoe Wallwork 
91420282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
91520282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
91620282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
91720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
91820282da2SJoe Wallwork       {
91920282da2SJoe Wallwork         PetscReal *rwork;
92020282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
92120282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
92220282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
92320282da2SJoe Wallwork       }
92420282da2SJoe Wallwork #else
92520282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
92620282da2SJoe Wallwork #endif
92720282da2SJoe Wallwork       if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
92820282da2SJoe Wallwork       ierr = PetscFPTrapPop();
92920282da2SJoe Wallwork 
93020282da2SJoe Wallwork       /* Compute square root and reciprocal */
93120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
93220282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
93320282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
93420282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
93520282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
93620282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
93720282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
93820282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
93920282da2SJoe Wallwork           }
94020282da2SJoe Wallwork         }
94120282da2SJoe Wallwork       }
94220282da2SJoe Wallwork 
94320282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
94420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
94520282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
94620282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
94720282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
94820282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
94920282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
95020282da2SJoe Wallwork             }
95120282da2SJoe Wallwork           }
95220282da2SJoe Wallwork         }
95320282da2SJoe Wallwork       }
95420282da2SJoe Wallwork 
95520282da2SJoe Wallwork       /* Compute eigendecomposition */
95620282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
95720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
95820282da2SJoe Wallwork       {
95920282da2SJoe Wallwork         PetscReal *rwork;
96020282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
96120282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
96220282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
96320282da2SJoe Wallwork       }
96420282da2SJoe Wallwork #else
96520282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
96620282da2SJoe Wallwork #endif
96720282da2SJoe Wallwork       if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
96820282da2SJoe Wallwork       ierr = PetscFPTrapPop();
96920282da2SJoe Wallwork 
97020282da2SJoe Wallwork       /* Modify eigenvalues */
97120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
97220282da2SJoe Wallwork 
97320282da2SJoe Wallwork       /* Map back to get the intersection */
97420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
97520282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
97620282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
97720282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
97820282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
97920282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
98020282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
98120282da2SJoe Wallwork               }
98220282da2SJoe Wallwork             }
98320282da2SJoe Wallwork           }
98420282da2SJoe Wallwork         }
98520282da2SJoe Wallwork       }
98620282da2SJoe Wallwork     }
98720282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
98820282da2SJoe Wallwork   }
98920282da2SJoe Wallwork   ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr);
99020282da2SJoe Wallwork   PetscFunctionReturn(0);
99120282da2SJoe Wallwork }
99220282da2SJoe Wallwork 
99320282da2SJoe Wallwork /*
99420282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
99520282da2SJoe Wallwork 
99620282da2SJoe Wallwork   Input Parameter:
99720282da2SJoe Wallwork + dm         - The DM
99820282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
99920282da2SJoe Wallwork - metrics    - The metrics to be intersected
100020282da2SJoe Wallwork 
100120282da2SJoe Wallwork   Output Parameter:
100220282da2SJoe Wallwork . metricInt  - The intersected metric
100320282da2SJoe Wallwork 
100420282da2SJoe Wallwork   Level: beginner
100520282da2SJoe Wallwork 
100620282da2SJoe Wallwork   Notes:
100720282da2SJoe Wallwork 
100820282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
100920282da2SJoe Wallwork 
101020282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
101120282da2SJoe Wallwork 
101220282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
101320282da2SJoe Wallwork */
101420282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
101520282da2SJoe Wallwork {
101620282da2SJoe Wallwork   PetscErrorCode ierr;
101720282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v, i;
101820282da2SJoe Wallwork   PetscScalar   *met, *meti, *M, *Mi;
101920282da2SJoe Wallwork 
102020282da2SJoe Wallwork   PetscFunctionBegin;
102120282da2SJoe Wallwork   if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); }
102220282da2SJoe Wallwork 
102320282da2SJoe Wallwork   /* Extract metadata from dm */
102420282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
102520282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr);
102620282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
102720282da2SJoe Wallwork 
102820282da2SJoe Wallwork   /* Copy over the first metric */
102920282da2SJoe Wallwork   ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr);
103020282da2SJoe Wallwork 
103120282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
103220282da2SJoe Wallwork   if (numMetrics > 1) {
103320282da2SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
103420282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
103520282da2SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
103620282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
103720282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr);
103820282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr);
103920282da2SJoe Wallwork         ierr = DMPlexMetricIntersection_Private(dim, Mi, M);CHKERRQ(ierr);
104020282da2SJoe Wallwork       }
104120282da2SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
104220282da2SJoe Wallwork     }
104320282da2SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
104420282da2SJoe Wallwork   }
104520282da2SJoe Wallwork 
104620282da2SJoe Wallwork   PetscFunctionReturn(0);
104720282da2SJoe Wallwork }
104820282da2SJoe Wallwork 
104920282da2SJoe Wallwork /*
105020282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
105120282da2SJoe Wallwork 
105220282da2SJoe Wallwork   Input Parameter:
105320282da2SJoe Wallwork + dm        - The DM
105420282da2SJoe Wallwork . metric1   - The first metric to be intersected
105520282da2SJoe Wallwork - metric2   - The second metric to be intersected
105620282da2SJoe Wallwork 
105720282da2SJoe Wallwork   Output Parameter:
105820282da2SJoe Wallwork . metricInt - The intersected metric
105920282da2SJoe Wallwork 
106020282da2SJoe Wallwork   Level: beginner
106120282da2SJoe Wallwork 
106220282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
106320282da2SJoe Wallwork */
106420282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
106520282da2SJoe Wallwork {
106620282da2SJoe Wallwork   PetscErrorCode ierr;
106720282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
106820282da2SJoe Wallwork 
106920282da2SJoe Wallwork   PetscFunctionBegin;
107020282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr);
107120282da2SJoe Wallwork   PetscFunctionReturn(0);
107220282da2SJoe Wallwork }
107320282da2SJoe Wallwork 
107420282da2SJoe Wallwork /*
107520282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
107620282da2SJoe Wallwork 
107720282da2SJoe Wallwork   Input Parameter:
107820282da2SJoe Wallwork + dm        - The DM
107920282da2SJoe Wallwork . metric1   - The first metric to be intersected
108020282da2SJoe Wallwork . metric2   - The second metric to be intersected
108120282da2SJoe Wallwork - metric3   - The third metric to be intersected
108220282da2SJoe Wallwork 
108320282da2SJoe Wallwork   Output Parameter:
108420282da2SJoe Wallwork . metricInt - The intersected metric
108520282da2SJoe Wallwork 
108620282da2SJoe Wallwork   Level: beginner
108720282da2SJoe Wallwork 
108820282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
108920282da2SJoe Wallwork */
109020282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
109120282da2SJoe Wallwork {
109220282da2SJoe Wallwork   PetscErrorCode ierr;
109320282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
109420282da2SJoe Wallwork 
109520282da2SJoe Wallwork   PetscFunctionBegin;
109620282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr);
109720282da2SJoe Wallwork   PetscFunctionReturn(0);
109820282da2SJoe Wallwork }
1099