xref: /petsc/src/dm/impls/plex/plexmetric.c (revision bc00d7af965a8f4b5f2e4adea421fc706cfaca52)
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 {
8*bc00d7afSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
931b70465SJoe Wallwork   MPI_Comm       comm;
1093520041SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
1176f360caSJoe Wallwork   PetscBool      noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE;
12cc2eee55SJoe Wallwork   PetscInt       verbosity = -1, numIter = 3;
13ae8b063eSJoe 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;
16*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(PetscNew(&plex->metricCtx));
179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
18d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL));
209566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, isotropic));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL));
229566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, uniform));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL));
249566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL));
269566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoInsertion(dm, noInsert));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL));
289566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoSwapping(dm, noSwap));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL));
309566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoMovement(dm, noMove));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL));
329566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoSurf(dm, noSurf));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0));
349566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNumIterations(dm, numIter));
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10));
369566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetVerbosity(dm, verbosity));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL));
389566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMinimumMagnitude(dm, h_min));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL));
409566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMaximumMagnitude(dm, h_max));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL));
429566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMaximumAnisotropy(dm, a_max));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL));
449566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNormalizationOrder(dm, p));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL));
469566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetTargetComplexity(dm, target));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL));
489566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetGradationFactor(dm, beta));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL));
509566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetHausdorffNumber(dm, hausd));
51d0609cedSBarry Smith   PetscOptionsEnd();
5231b70465SJoe Wallwork   PetscFunctionReturn(0);
5331b70465SJoe Wallwork }
5431b70465SJoe Wallwork 
55cb7bfe8cSJoe Wallwork /*@
5631b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
5731b70465SJoe Wallwork 
5831b70465SJoe Wallwork   Input parameters:
5931b70465SJoe Wallwork + dm        - The DM
6031b70465SJoe Wallwork - isotropic - Is the metric isotropic?
6131b70465SJoe Wallwork 
6231b70465SJoe Wallwork   Level: beginner
6331b70465SJoe Wallwork 
6493520041SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst()
65cb7bfe8cSJoe Wallwork @*/
6631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
6731b70465SJoe Wallwork {
6831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6931b70465SJoe Wallwork 
7031b70465SJoe Wallwork   PetscFunctionBegin;
71*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
7231b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
7331b70465SJoe Wallwork   PetscFunctionReturn(0);
7431b70465SJoe Wallwork }
7531b70465SJoe Wallwork 
76cb7bfe8cSJoe Wallwork /*@
7793520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7831b70465SJoe Wallwork 
7931b70465SJoe Wallwork   Input parameters:
8031b70465SJoe Wallwork . dm        - The DM
8131b70465SJoe Wallwork 
8231b70465SJoe Wallwork   Output parameters:
8331b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8431b70465SJoe Wallwork 
8531b70465SJoe Wallwork   Level: beginner
8631b70465SJoe Wallwork 
8793520041SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst()
88cb7bfe8cSJoe Wallwork @*/
8931b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
9031b70465SJoe Wallwork {
9131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
9231b70465SJoe Wallwork 
9331b70465SJoe Wallwork   PetscFunctionBegin;
94*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
9531b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
9631b70465SJoe Wallwork   PetscFunctionReturn(0);
9731b70465SJoe Wallwork }
9831b70465SJoe Wallwork 
99cb7bfe8cSJoe Wallwork /*@
10093520041SJoe Wallwork   DMPlexMetricSetUniform - Record whether a metric is uniform
10193520041SJoe Wallwork 
10293520041SJoe Wallwork   Input parameters:
10393520041SJoe Wallwork + dm      - The DM
10493520041SJoe Wallwork - uniform - Is the metric uniform?
10593520041SJoe Wallwork 
10693520041SJoe Wallwork   Level: beginner
10793520041SJoe Wallwork 
10893520041SJoe Wallwork   Notes:
10993520041SJoe Wallwork 
11093520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11193520041SJoe Wallwork 
11293520041SJoe Wallwork .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
11393520041SJoe Wallwork @*/
11493520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
11593520041SJoe Wallwork {
11693520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
11793520041SJoe Wallwork 
11893520041SJoe Wallwork   PetscFunctionBegin;
119*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
12093520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
12193520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
12293520041SJoe Wallwork   PetscFunctionReturn(0);
12393520041SJoe Wallwork }
12493520041SJoe Wallwork 
12593520041SJoe Wallwork /*@
12693520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
12793520041SJoe Wallwork 
12893520041SJoe Wallwork   Input parameters:
12993520041SJoe Wallwork . dm      - The DM
13093520041SJoe Wallwork 
13193520041SJoe Wallwork   Output parameters:
13293520041SJoe Wallwork . uniform - Is the metric uniform?
13393520041SJoe Wallwork 
13493520041SJoe Wallwork   Level: beginner
13593520041SJoe Wallwork 
13693520041SJoe Wallwork .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
13793520041SJoe Wallwork @*/
13893520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
13993520041SJoe Wallwork {
14093520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
14193520041SJoe Wallwork 
14293520041SJoe Wallwork   PetscFunctionBegin;
143*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
14493520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
14593520041SJoe Wallwork   PetscFunctionReturn(0);
14693520041SJoe Wallwork }
14793520041SJoe Wallwork 
14893520041SJoe Wallwork /*@
14931b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
15031b70465SJoe Wallwork 
15131b70465SJoe Wallwork   Input parameters:
15231b70465SJoe Wallwork + dm                      - The DM
15331b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
15431b70465SJoe Wallwork 
15531b70465SJoe Wallwork   Level: beginner
15631b70465SJoe Wallwork 
15731b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
158cb7bfe8cSJoe Wallwork @*/
15931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
16031b70465SJoe Wallwork {
16131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
16231b70465SJoe Wallwork 
16331b70465SJoe Wallwork   PetscFunctionBegin;
164*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
16531b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
16631b70465SJoe Wallwork   PetscFunctionReturn(0);
16731b70465SJoe Wallwork }
16831b70465SJoe Wallwork 
169cb7bfe8cSJoe Wallwork /*@
17031b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
17131b70465SJoe Wallwork 
17231b70465SJoe Wallwork   Input parameters:
17331b70465SJoe Wallwork . dm                      - The DM
17431b70465SJoe Wallwork 
17531b70465SJoe Wallwork   Output parameters:
17631b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
17731b70465SJoe Wallwork 
17831b70465SJoe Wallwork   Level: beginner
17931b70465SJoe Wallwork 
18031b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
181cb7bfe8cSJoe Wallwork @*/
18231b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
18331b70465SJoe Wallwork {
18431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
18531b70465SJoe Wallwork 
18631b70465SJoe Wallwork   PetscFunctionBegin;
187*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
18831b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
18931b70465SJoe Wallwork   PetscFunctionReturn(0);
19031b70465SJoe Wallwork }
19131b70465SJoe Wallwork 
192cb7bfe8cSJoe Wallwork /*@
193cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
194cc2eee55SJoe Wallwork 
195cc2eee55SJoe Wallwork   Input parameters:
196cc2eee55SJoe Wallwork + dm       - The DM
197cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
198cc2eee55SJoe Wallwork 
199cc2eee55SJoe Wallwork   Level: beginner
200cc2eee55SJoe Wallwork 
201cb7bfe8cSJoe Wallwork   Notes:
202cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
203cb7bfe8cSJoe Wallwork 
20476f360caSJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf()
205cb7bfe8cSJoe Wallwork @*/
206cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
207cc2eee55SJoe Wallwork {
208cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
209cc2eee55SJoe Wallwork 
210cc2eee55SJoe Wallwork   PetscFunctionBegin;
211*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
212cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
213cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
214cc2eee55SJoe Wallwork }
215cc2eee55SJoe Wallwork 
216cb7bfe8cSJoe Wallwork /*@
217cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
218cc2eee55SJoe Wallwork 
219cc2eee55SJoe Wallwork   Input parameters:
220cc2eee55SJoe Wallwork . dm       - The DM
221cc2eee55SJoe Wallwork 
222cc2eee55SJoe Wallwork   Output parameters:
223cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
224cc2eee55SJoe Wallwork 
225cc2eee55SJoe Wallwork   Level: beginner
226cc2eee55SJoe Wallwork 
227cb7bfe8cSJoe Wallwork   Notes:
228cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
229cb7bfe8cSJoe Wallwork 
23076f360caSJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf()
231cb7bfe8cSJoe Wallwork @*/
232cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
233cc2eee55SJoe Wallwork {
234cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
235cc2eee55SJoe Wallwork 
236cc2eee55SJoe Wallwork   PetscFunctionBegin;
237*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
238cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
239cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
240cc2eee55SJoe Wallwork }
241cc2eee55SJoe Wallwork 
242cb7bfe8cSJoe Wallwork /*@
243cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
244cc2eee55SJoe Wallwork 
245cc2eee55SJoe Wallwork   Input parameters:
246cc2eee55SJoe Wallwork + dm     - The DM
247cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
248cc2eee55SJoe Wallwork 
249cc2eee55SJoe Wallwork   Level: beginner
250cc2eee55SJoe Wallwork 
251cb7bfe8cSJoe Wallwork   Notes:
252cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
253cb7bfe8cSJoe Wallwork 
25476f360caSJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf()
255cb7bfe8cSJoe Wallwork @*/
256cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
257cc2eee55SJoe Wallwork {
258cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
259cc2eee55SJoe Wallwork 
260cc2eee55SJoe Wallwork   PetscFunctionBegin;
261*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
262cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
263cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
264cc2eee55SJoe Wallwork }
265cc2eee55SJoe Wallwork 
266cb7bfe8cSJoe Wallwork /*@
267cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
268cc2eee55SJoe Wallwork 
269cc2eee55SJoe Wallwork   Input parameters:
270cc2eee55SJoe Wallwork . dm     - The DM
271cc2eee55SJoe Wallwork 
272cc2eee55SJoe Wallwork   Output parameters:
273cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
274cc2eee55SJoe Wallwork 
275cc2eee55SJoe Wallwork   Level: beginner
276cc2eee55SJoe Wallwork 
277cb7bfe8cSJoe Wallwork   Notes:
278cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
279cb7bfe8cSJoe Wallwork 
28076f360caSJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf()
281cb7bfe8cSJoe Wallwork @*/
282cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
283cc2eee55SJoe Wallwork {
284cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
285cc2eee55SJoe Wallwork 
286cc2eee55SJoe Wallwork   PetscFunctionBegin;
287*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
288cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
289cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
290cc2eee55SJoe Wallwork }
291cc2eee55SJoe Wallwork 
292cb7bfe8cSJoe Wallwork /*@
293cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
294cc2eee55SJoe Wallwork 
295cc2eee55SJoe Wallwork   Input parameters:
296cc2eee55SJoe Wallwork + dm     - The DM
297cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
298cc2eee55SJoe Wallwork 
299cc2eee55SJoe Wallwork   Level: beginner
300cc2eee55SJoe Wallwork 
301cb7bfe8cSJoe Wallwork   Notes:
302cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
303cb7bfe8cSJoe Wallwork 
30476f360caSJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoSurf()
305cb7bfe8cSJoe Wallwork @*/
306cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
307cc2eee55SJoe Wallwork {
308cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
309cc2eee55SJoe Wallwork 
310cc2eee55SJoe Wallwork   PetscFunctionBegin;
311*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
312cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
313cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
314cc2eee55SJoe Wallwork }
315cc2eee55SJoe Wallwork 
316cb7bfe8cSJoe Wallwork /*@
317cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
318cc2eee55SJoe Wallwork 
319cc2eee55SJoe Wallwork   Input parameters:
320cc2eee55SJoe Wallwork . dm     - The DM
321cc2eee55SJoe Wallwork 
322cc2eee55SJoe Wallwork   Output parameters:
323cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
324cc2eee55SJoe Wallwork 
325cc2eee55SJoe Wallwork   Level: beginner
326cc2eee55SJoe Wallwork 
327cb7bfe8cSJoe Wallwork   Notes:
328cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
329cb7bfe8cSJoe Wallwork 
33076f360caSJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoSurf()
331cb7bfe8cSJoe Wallwork @*/
332cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
333cc2eee55SJoe Wallwork {
334cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
335cc2eee55SJoe Wallwork 
336cc2eee55SJoe Wallwork   PetscFunctionBegin;
337*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
338cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
339cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
340cc2eee55SJoe Wallwork }
341cc2eee55SJoe Wallwork 
342cb7bfe8cSJoe Wallwork /*@
34376f360caSJoe Wallwork   DMPlexMetricSetNoSurf - Should surface modification be turned off?
34476f360caSJoe Wallwork 
34576f360caSJoe Wallwork   Input parameters:
34676f360caSJoe Wallwork + dm     - The DM
34776f360caSJoe Wallwork - noSurf - Should surface modification be turned off?
34876f360caSJoe Wallwork 
34976f360caSJoe Wallwork   Level: beginner
35076f360caSJoe Wallwork 
35176f360caSJoe Wallwork   Notes:
35276f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
35376f360caSJoe Wallwork 
35476f360caSJoe Wallwork .seealso: DMPlexMetricNoSurf(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping()
35576f360caSJoe Wallwork @*/
35676f360caSJoe Wallwork PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
35776f360caSJoe Wallwork {
35876f360caSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
35976f360caSJoe Wallwork 
36076f360caSJoe Wallwork   PetscFunctionBegin;
361*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
36276f360caSJoe Wallwork   plex->metricCtx->noSurf = noSurf;
36376f360caSJoe Wallwork   PetscFunctionReturn(0);
36476f360caSJoe Wallwork }
36576f360caSJoe Wallwork 
36676f360caSJoe Wallwork /*@
36776f360caSJoe Wallwork   DMPlexMetricNoSurf - Is surface modification turned off?
36876f360caSJoe Wallwork 
36976f360caSJoe Wallwork   Input parameters:
37076f360caSJoe Wallwork . dm     - The DM
37176f360caSJoe Wallwork 
37276f360caSJoe Wallwork   Output parameters:
37376f360caSJoe Wallwork . noSurf - Is surface modification turned off?
37476f360caSJoe Wallwork 
37576f360caSJoe Wallwork   Level: beginner
37676f360caSJoe Wallwork 
37776f360caSJoe Wallwork   Notes:
37876f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
37976f360caSJoe Wallwork 
38076f360caSJoe Wallwork .seealso: DMPlexMetricSetNoSurf(), DMPlexMetricNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping()
38176f360caSJoe Wallwork @*/
38276f360caSJoe Wallwork PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
38376f360caSJoe Wallwork {
38476f360caSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
38576f360caSJoe Wallwork 
38676f360caSJoe Wallwork   PetscFunctionBegin;
387*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
38876f360caSJoe Wallwork   *noSurf = plex->metricCtx->noSurf;
38976f360caSJoe Wallwork   PetscFunctionReturn(0);
39076f360caSJoe Wallwork }
39176f360caSJoe Wallwork 
39276f360caSJoe Wallwork /*@
39331b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
39431b70465SJoe Wallwork 
39531b70465SJoe Wallwork   Input parameters:
39631b70465SJoe Wallwork + dm    - The DM
39731b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
39831b70465SJoe Wallwork 
39931b70465SJoe Wallwork   Level: beginner
40031b70465SJoe Wallwork 
40131b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
402cb7bfe8cSJoe Wallwork @*/
40331b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
40431b70465SJoe Wallwork {
40531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
40631b70465SJoe Wallwork 
40731b70465SJoe Wallwork   PetscFunctionBegin;
408*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
409*bc00d7afSJoe Wallwork   PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
41031b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
41131b70465SJoe Wallwork   PetscFunctionReturn(0);
41231b70465SJoe Wallwork }
41331b70465SJoe Wallwork 
414cb7bfe8cSJoe Wallwork /*@
41531b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
41631b70465SJoe Wallwork 
41731b70465SJoe Wallwork   Input parameters:
41831b70465SJoe Wallwork . dm    - The DM
41931b70465SJoe Wallwork 
420cc2eee55SJoe Wallwork   Output parameters:
42131b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
42231b70465SJoe Wallwork 
42331b70465SJoe Wallwork   Level: beginner
42431b70465SJoe Wallwork 
42531b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
426cb7bfe8cSJoe Wallwork @*/
42731b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
42831b70465SJoe Wallwork {
42931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
43031b70465SJoe Wallwork 
43131b70465SJoe Wallwork   PetscFunctionBegin;
432*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
43331b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
43431b70465SJoe Wallwork   PetscFunctionReturn(0);
43531b70465SJoe Wallwork }
43631b70465SJoe Wallwork 
437cb7bfe8cSJoe Wallwork /*@
43831b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
43931b70465SJoe Wallwork 
44031b70465SJoe Wallwork   Input parameters:
44131b70465SJoe Wallwork + dm    - The DM
44231b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44331b70465SJoe Wallwork 
44431b70465SJoe Wallwork   Level: beginner
44531b70465SJoe Wallwork 
44631b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
447cb7bfe8cSJoe Wallwork @*/
44831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
44931b70465SJoe Wallwork {
45031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
45131b70465SJoe Wallwork 
45231b70465SJoe Wallwork   PetscFunctionBegin;
453*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
454*bc00d7afSJoe Wallwork   PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
45531b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
45631b70465SJoe Wallwork   PetscFunctionReturn(0);
45731b70465SJoe Wallwork }
45831b70465SJoe Wallwork 
459cb7bfe8cSJoe Wallwork /*@
46031b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
46131b70465SJoe Wallwork 
46231b70465SJoe Wallwork   Input parameters:
46331b70465SJoe Wallwork . dm    - The DM
46431b70465SJoe Wallwork 
465cc2eee55SJoe Wallwork   Output parameters:
46631b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
46731b70465SJoe Wallwork 
46831b70465SJoe Wallwork   Level: beginner
46931b70465SJoe Wallwork 
47031b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
471cb7bfe8cSJoe Wallwork @*/
47231b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
47331b70465SJoe Wallwork {
47431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
47531b70465SJoe Wallwork 
47631b70465SJoe Wallwork   PetscFunctionBegin;
477*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
47831b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
47931b70465SJoe Wallwork   PetscFunctionReturn(0);
48031b70465SJoe Wallwork }
48131b70465SJoe Wallwork 
482cb7bfe8cSJoe Wallwork /*@
48331b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
48431b70465SJoe Wallwork 
48531b70465SJoe Wallwork   Input parameters:
48631b70465SJoe Wallwork + dm    - The DM
48731b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
48831b70465SJoe Wallwork 
48931b70465SJoe Wallwork   Level: beginner
49031b70465SJoe Wallwork 
49131b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
49231b70465SJoe Wallwork 
49331b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
494cb7bfe8cSJoe Wallwork @*/
49531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
49631b70465SJoe Wallwork {
49731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
49831b70465SJoe Wallwork 
49931b70465SJoe Wallwork   PetscFunctionBegin;
500*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
501*bc00d7afSJoe Wallwork   PetscCheck(a_max >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Anisotropy must be in [1, inf)");
50231b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
50331b70465SJoe Wallwork   PetscFunctionReturn(0);
50431b70465SJoe Wallwork }
50531b70465SJoe Wallwork 
506cb7bfe8cSJoe Wallwork /*@
50731b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
50831b70465SJoe Wallwork 
50931b70465SJoe Wallwork   Input parameters:
51031b70465SJoe Wallwork . dm    - The DM
51131b70465SJoe Wallwork 
512cc2eee55SJoe Wallwork   Output parameters:
51331b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
51431b70465SJoe Wallwork 
51531b70465SJoe Wallwork   Level: beginner
51631b70465SJoe Wallwork 
51731b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
518cb7bfe8cSJoe Wallwork @*/
51931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
52031b70465SJoe Wallwork {
52131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
52231b70465SJoe Wallwork 
52331b70465SJoe Wallwork   PetscFunctionBegin;
524*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
52531b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
52631b70465SJoe Wallwork   PetscFunctionReturn(0);
52731b70465SJoe Wallwork }
52831b70465SJoe Wallwork 
529cb7bfe8cSJoe Wallwork /*@
53031b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
53131b70465SJoe Wallwork 
53231b70465SJoe Wallwork   Input parameters:
53331b70465SJoe Wallwork + dm               - The DM
53431b70465SJoe Wallwork - targetComplexity - The target metric complexity
53531b70465SJoe Wallwork 
53631b70465SJoe Wallwork   Level: beginner
53731b70465SJoe Wallwork 
53831b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
539cb7bfe8cSJoe Wallwork @*/
54031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
54131b70465SJoe Wallwork {
54231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
54331b70465SJoe Wallwork 
54431b70465SJoe Wallwork   PetscFunctionBegin;
545*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
546*bc00d7afSJoe Wallwork   PetscCheck(targetComplexity > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric complexity must be in (0, inf)");
54731b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
54831b70465SJoe Wallwork   PetscFunctionReturn(0);
54931b70465SJoe Wallwork }
55031b70465SJoe Wallwork 
551cb7bfe8cSJoe Wallwork /*@
55231b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
55331b70465SJoe Wallwork 
55431b70465SJoe Wallwork   Input parameters:
55531b70465SJoe Wallwork . dm               - The DM
55631b70465SJoe Wallwork 
557cc2eee55SJoe Wallwork   Output parameters:
55831b70465SJoe Wallwork . targetComplexity - The target metric complexity
55931b70465SJoe Wallwork 
56031b70465SJoe Wallwork   Level: beginner
56131b70465SJoe Wallwork 
56231b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
563cb7bfe8cSJoe Wallwork @*/
56431b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
56531b70465SJoe Wallwork {
56631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
56731b70465SJoe Wallwork 
56831b70465SJoe Wallwork   PetscFunctionBegin;
569*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
57031b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
57131b70465SJoe Wallwork   PetscFunctionReturn(0);
57231b70465SJoe Wallwork }
57331b70465SJoe Wallwork 
574cb7bfe8cSJoe Wallwork /*@
57531b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
57631b70465SJoe Wallwork 
57731b70465SJoe Wallwork   Input parameters:
57831b70465SJoe Wallwork + dm - The DM
57931b70465SJoe Wallwork - p  - The normalization order
58031b70465SJoe Wallwork 
58131b70465SJoe Wallwork   Level: beginner
58231b70465SJoe Wallwork 
58331b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
584cb7bfe8cSJoe Wallwork @*/
58531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
58631b70465SJoe Wallwork {
58731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
58831b70465SJoe Wallwork 
58931b70465SJoe Wallwork   PetscFunctionBegin;
590*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
591*bc00d7afSJoe Wallwork   PetscCheck(p >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Normalization order must be in [1, inf)");
59231b70465SJoe Wallwork   plex->metricCtx->p = p;
59331b70465SJoe Wallwork   PetscFunctionReturn(0);
59431b70465SJoe Wallwork }
59531b70465SJoe Wallwork 
596cb7bfe8cSJoe Wallwork /*@
59731b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
59831b70465SJoe Wallwork 
59931b70465SJoe Wallwork   Input parameters:
60031b70465SJoe Wallwork . dm - The DM
60131b70465SJoe Wallwork 
602cc2eee55SJoe Wallwork   Output parameters:
60331b70465SJoe Wallwork . p - The normalization order
60431b70465SJoe Wallwork 
60531b70465SJoe Wallwork   Level: beginner
60631b70465SJoe Wallwork 
60731b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
608cb7bfe8cSJoe Wallwork @*/
60931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
61031b70465SJoe Wallwork {
61131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
61231b70465SJoe Wallwork 
61331b70465SJoe Wallwork   PetscFunctionBegin;
614*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
61531b70465SJoe Wallwork   *p = plex->metricCtx->p;
61631b70465SJoe Wallwork   PetscFunctionReturn(0);
61731b70465SJoe Wallwork }
61831b70465SJoe Wallwork 
619cb7bfe8cSJoe Wallwork /*@
620cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
621cc2eee55SJoe Wallwork 
622cc2eee55SJoe Wallwork   Input parameters:
623cc2eee55SJoe Wallwork + dm   - The DM
624cc2eee55SJoe Wallwork - beta - The metric gradation factor
625cc2eee55SJoe Wallwork 
626cc2eee55SJoe Wallwork   Level: beginner
627cc2eee55SJoe Wallwork 
628cc2eee55SJoe Wallwork   Notes:
629cc2eee55SJoe Wallwork 
630cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
631cc2eee55SJoe Wallwork 
632cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
633cc2eee55SJoe Wallwork 
634cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
635cb7bfe8cSJoe Wallwork 
636ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber()
637cb7bfe8cSJoe Wallwork @*/
638cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
639cc2eee55SJoe Wallwork {
640cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
641cc2eee55SJoe Wallwork 
642cc2eee55SJoe Wallwork   PetscFunctionBegin;
643*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
644*bc00d7afSJoe Wallwork   PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)");
645cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
646cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
647cc2eee55SJoe Wallwork }
648cc2eee55SJoe Wallwork 
649cb7bfe8cSJoe Wallwork /*@
650cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
651cc2eee55SJoe Wallwork 
652cc2eee55SJoe Wallwork   Input parameters:
653cc2eee55SJoe Wallwork . dm   - The DM
654cc2eee55SJoe Wallwork 
655cc2eee55SJoe Wallwork   Output parameters:
656cc2eee55SJoe Wallwork . beta - The metric gradation factor
657cc2eee55SJoe Wallwork 
658cc2eee55SJoe Wallwork   Level: beginner
659cc2eee55SJoe Wallwork 
660cb7bfe8cSJoe Wallwork   Notes:
661cb7bfe8cSJoe Wallwork 
662cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
663cb7bfe8cSJoe Wallwork 
664cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
665cb7bfe8cSJoe Wallwork 
666cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
667cc2eee55SJoe Wallwork 
668ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber()
669cb7bfe8cSJoe Wallwork @*/
670cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
671cc2eee55SJoe Wallwork {
672cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
673cc2eee55SJoe Wallwork 
674cc2eee55SJoe Wallwork   PetscFunctionBegin;
675*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
676cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
677cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
678cc2eee55SJoe Wallwork }
679cc2eee55SJoe Wallwork 
680cb7bfe8cSJoe Wallwork /*@
681ae8b063eSJoe Wallwork   DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
682ae8b063eSJoe Wallwork 
683ae8b063eSJoe Wallwork   Input parameters:
684ae8b063eSJoe Wallwork + dm    - The DM
685ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number
686ae8b063eSJoe Wallwork 
687ae8b063eSJoe Wallwork   Level: beginner
688ae8b063eSJoe Wallwork 
689ae8b063eSJoe Wallwork   Notes:
690ae8b063eSJoe Wallwork 
691ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
692ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
693ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
694ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
695ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
696ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
697ae8b063eSJoe Wallwork 
698ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
699ae8b063eSJoe Wallwork 
700ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber()
701ae8b063eSJoe Wallwork @*/
702ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
703ae8b063eSJoe Wallwork {
704ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
705ae8b063eSJoe Wallwork 
706ae8b063eSJoe Wallwork   PetscFunctionBegin;
707*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
708*bc00d7afSJoe Wallwork   PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)");
709ae8b063eSJoe Wallwork   plex->metricCtx->hausdorffNumber = hausd;
710ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
711ae8b063eSJoe Wallwork }
712ae8b063eSJoe Wallwork 
713ae8b063eSJoe Wallwork /*@
714ae8b063eSJoe Wallwork   DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
715ae8b063eSJoe Wallwork 
716ae8b063eSJoe Wallwork   Input parameters:
717ae8b063eSJoe Wallwork . dm    - The DM
718ae8b063eSJoe Wallwork 
719ae8b063eSJoe Wallwork   Output parameters:
720ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number
721ae8b063eSJoe Wallwork 
722ae8b063eSJoe Wallwork   Level: beginner
723ae8b063eSJoe Wallwork 
724ae8b063eSJoe Wallwork   Notes:
725ae8b063eSJoe Wallwork 
726ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
727ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
728ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
729ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
730ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
731ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
732ae8b063eSJoe Wallwork 
733ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
734ae8b063eSJoe Wallwork 
735ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber()
736ae8b063eSJoe Wallwork @*/
737ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
738ae8b063eSJoe Wallwork {
739ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
740ae8b063eSJoe Wallwork 
741ae8b063eSJoe Wallwork   PetscFunctionBegin;
742*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
743ae8b063eSJoe Wallwork   *hausd = plex->metricCtx->hausdorffNumber;
744ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
745ae8b063eSJoe Wallwork }
746ae8b063eSJoe Wallwork 
747ae8b063eSJoe Wallwork /*@
748cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
749cc2eee55SJoe Wallwork 
750cc2eee55SJoe Wallwork   Input parameters:
751cc2eee55SJoe Wallwork + dm        - The DM
752cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
753cc2eee55SJoe Wallwork 
754cb7bfe8cSJoe Wallwork   Level: beginner
755cb7bfe8cSJoe Wallwork 
756cb7bfe8cSJoe Wallwork   Notes:
757cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
758cb7bfe8cSJoe Wallwork 
759cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations()
760cb7bfe8cSJoe Wallwork @*/
761cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
762cc2eee55SJoe Wallwork {
763cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
764cc2eee55SJoe Wallwork 
765cc2eee55SJoe Wallwork   PetscFunctionBegin;
766*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
767cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
768cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
769cc2eee55SJoe Wallwork }
770cc2eee55SJoe Wallwork 
771cb7bfe8cSJoe Wallwork /*@
772cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
773cc2eee55SJoe Wallwork 
774cc2eee55SJoe Wallwork   Input parameters:
775cc2eee55SJoe Wallwork . dm        - The DM
776cc2eee55SJoe Wallwork 
777cc2eee55SJoe Wallwork   Output parameters:
778cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
779cc2eee55SJoe Wallwork 
780cb7bfe8cSJoe Wallwork   Level: beginner
781cb7bfe8cSJoe Wallwork 
782cb7bfe8cSJoe Wallwork   Notes:
783cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
784cb7bfe8cSJoe Wallwork 
785cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
786cb7bfe8cSJoe Wallwork @*/
787cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
788cc2eee55SJoe Wallwork {
789cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
790cc2eee55SJoe Wallwork 
791cc2eee55SJoe Wallwork   PetscFunctionBegin;
792*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
793cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
794cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
795cc2eee55SJoe Wallwork }
796cc2eee55SJoe Wallwork 
797cb7bfe8cSJoe Wallwork /*@
798cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
799cc2eee55SJoe Wallwork 
800cc2eee55SJoe Wallwork   Input parameters:
801cc2eee55SJoe Wallwork + dm      - The DM
802cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
803cc2eee55SJoe Wallwork 
804cb7bfe8cSJoe Wallwork   Level: beginner
805cb7bfe8cSJoe Wallwork 
806cb7bfe8cSJoe Wallwork   Notes:
807cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
808cc2eee55SJoe Wallwork 
809cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
810cb7bfe8cSJoe Wallwork @*/
811cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
812cc2eee55SJoe Wallwork {
813cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
814cc2eee55SJoe Wallwork 
815cc2eee55SJoe Wallwork   PetscFunctionBegin;
816*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
817cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
818cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
819cc2eee55SJoe Wallwork }
820cc2eee55SJoe Wallwork 
821cb7bfe8cSJoe Wallwork /*@
822cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
823cc2eee55SJoe Wallwork 
824cc2eee55SJoe Wallwork   Input parameters:
825cc2eee55SJoe Wallwork . dm      - The DM
826cc2eee55SJoe Wallwork 
827cc2eee55SJoe Wallwork   Output parameters:
828cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
829cc2eee55SJoe Wallwork 
830cb7bfe8cSJoe Wallwork   Level: beginner
831cb7bfe8cSJoe Wallwork 
832cb7bfe8cSJoe Wallwork   Notes:
833cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
834cc2eee55SJoe Wallwork 
835cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity()
836cb7bfe8cSJoe Wallwork @*/
837cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
838cc2eee55SJoe Wallwork {
839cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
840cc2eee55SJoe Wallwork 
841cc2eee55SJoe Wallwork   PetscFunctionBegin;
842*bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
843cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
844cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
845cc2eee55SJoe Wallwork }
846cc2eee55SJoe Wallwork 
84720282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
84820282da2SJoe Wallwork {
84920282da2SJoe Wallwork   MPI_Comm       comm;
85020282da2SJoe Wallwork   PetscFE        fe;
85120282da2SJoe Wallwork   PetscInt       dim;
85220282da2SJoe Wallwork 
85320282da2SJoe Wallwork   PetscFunctionBegin;
85420282da2SJoe Wallwork 
85520282da2SJoe Wallwork   /* Extract metadata from dm */
8569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
8579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
85820282da2SJoe Wallwork 
85920282da2SJoe Wallwork   /* Create a P1 field of the requested size */
8609566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
8619566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
8629566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
8639566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
8649566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm, metric));
86520282da2SJoe Wallwork 
86620282da2SJoe Wallwork   PetscFunctionReturn(0);
86720282da2SJoe Wallwork }
86820282da2SJoe Wallwork 
869cb7bfe8cSJoe Wallwork /*@
87020282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
87120282da2SJoe Wallwork 
87220282da2SJoe Wallwork   Input parameters:
87320282da2SJoe Wallwork + dm     - The DM
87420282da2SJoe Wallwork - f      - The field number to use
87520282da2SJoe Wallwork 
87620282da2SJoe Wallwork   Output parameter:
87720282da2SJoe Wallwork . metric - The metric
87820282da2SJoe Wallwork 
87920282da2SJoe Wallwork   Level: beginner
88020282da2SJoe Wallwork 
88131b70465SJoe Wallwork   Notes:
88231b70465SJoe Wallwork 
88331b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
88431b70465SJoe Wallwork 
88531b70465SJoe Wallwork   Command line options for Riemannian metrics:
88631b70465SJoe Wallwork 
887cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
88893520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
889cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
890cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
891cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
892cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
893cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
89467b8a455SSatish Balay - -dm_plex_metric_target_complexity         - Target metric complexity
895cb7bfe8cSJoe Wallwork 
896cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
897cb7bfe8cSJoe Wallwork 
89867b8a455SSatish Balay . -dm_adaptor <pragmatic/mmg/parmmg>        - specify dm adaptor to use
899cb7bfe8cSJoe Wallwork 
900cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
901cb7bfe8cSJoe Wallwork 
902cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
903cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
904cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
905cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
906cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
907cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
90820282da2SJoe Wallwork 
90920282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
910cb7bfe8cSJoe Wallwork @*/
91120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
91220282da2SJoe Wallwork {
91331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
91493520041SJoe Wallwork   PetscBool      isotropic, uniform;
91520282da2SJoe Wallwork   PetscInt       coordDim, Nd;
91620282da2SJoe Wallwork 
91720282da2SJoe Wallwork   PetscFunctionBegin;
91831b70465SJoe Wallwork   if (!plex->metricCtx) {
9199566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
9209566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
92131b70465SJoe Wallwork   }
9229566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
92320282da2SJoe Wallwork   Nd = coordDim*coordDim;
9249566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
9259566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
92693520041SJoe Wallwork   if (uniform) {
92793520041SJoe Wallwork     MPI_Comm comm;
92893520041SJoe Wallwork 
9299566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
930*bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
9319566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, metric));
9329566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
9339566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(*metric));
93493520041SJoe Wallwork   } else if (isotropic) {
9359566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
93693520041SJoe Wallwork   } else {
9379566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
93893520041SJoe Wallwork   }
93920282da2SJoe Wallwork   PetscFunctionReturn(0);
94020282da2SJoe Wallwork }
94120282da2SJoe Wallwork 
942cb7bfe8cSJoe Wallwork /*@
94320282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
94420282da2SJoe Wallwork 
94520282da2SJoe Wallwork   Input parameters:
94620282da2SJoe Wallwork + dm     - The DM
94720282da2SJoe Wallwork . f      - The field number to use
94820282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
94920282da2SJoe Wallwork 
95020282da2SJoe Wallwork   Output parameter:
95120282da2SJoe Wallwork . metric - The uniform metric
95220282da2SJoe Wallwork 
95320282da2SJoe Wallwork   Level: beginner
95420282da2SJoe Wallwork 
95593520041SJoe Wallwork   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
95620282da2SJoe Wallwork 
95720282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
958cb7bfe8cSJoe Wallwork @*/
95920282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
96020282da2SJoe Wallwork {
96120282da2SJoe Wallwork   PetscFunctionBegin;
9629566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
9639566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9642c71b3e2SJacob Faibussowitsch   PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
965*bc00d7afSJoe Wallwork   PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
9669566063dSJacob Faibussowitsch   PetscCall(VecSet(*metric, alpha));
9679566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*metric));
9689566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*metric));
96920282da2SJoe Wallwork   PetscFunctionReturn(0);
97020282da2SJoe Wallwork }
97120282da2SJoe Wallwork 
97293520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
97393520041SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
97493520041SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97593520041SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
97693520041SJoe Wallwork {
97793520041SJoe Wallwork   f0[0] = u[0];
97893520041SJoe Wallwork }
97993520041SJoe Wallwork 
980cb7bfe8cSJoe Wallwork /*@
98120282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
98220282da2SJoe Wallwork 
98320282da2SJoe Wallwork   Input parameters:
98420282da2SJoe Wallwork + dm        - The DM
98520282da2SJoe Wallwork . f         - The field number to use
98620282da2SJoe Wallwork - indicator - The error indicator
98720282da2SJoe Wallwork 
98820282da2SJoe Wallwork   Output parameter:
98920282da2SJoe Wallwork . metric    - The isotropic metric
99020282da2SJoe Wallwork 
99120282da2SJoe Wallwork   Level: beginner
99220282da2SJoe Wallwork 
99320282da2SJoe Wallwork   Notes:
99420282da2SJoe Wallwork 
99520282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
99620282da2SJoe Wallwork 
99793520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
99820282da2SJoe Wallwork 
99920282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
1000cb7bfe8cSJoe Wallwork @*/
100120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
100220282da2SJoe Wallwork {
100393520041SJoe Wallwork   PetscInt       m, n;
100420282da2SJoe Wallwork 
100520282da2SJoe Wallwork   PetscFunctionBegin;
10069566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
10079566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
10089566063dSJacob Faibussowitsch   PetscCall(VecGetSize(indicator, &m));
10099566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metric, &n));
10109566063dSJacob Faibussowitsch   if (m == n) PetscCall(VecCopy(indicator, *metric));
101193520041SJoe Wallwork   else {
101293520041SJoe Wallwork     DM     dmIndi;
101393520041SJoe Wallwork     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
101493520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101593520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101693520041SJoe Wallwork                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
101793520041SJoe Wallwork 
10189566063dSJacob Faibussowitsch     PetscCall(VecGetDM(indicator, &dmIndi));
101993520041SJoe Wallwork     funcs[0] = identity;
10209566063dSJacob Faibussowitsch     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
102120282da2SJoe Wallwork   }
102220282da2SJoe Wallwork   PetscFunctionReturn(0);
102320282da2SJoe Wallwork }
102420282da2SJoe Wallwork 
102582490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
102682490253SJoe Wallwork {
102782490253SJoe Wallwork   PetscInt i, j;
102882490253SJoe Wallwork 
102982490253SJoe Wallwork   PetscFunctionBegin;
103082490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
103182490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
103282490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
103382490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
103482490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
103563a3b9bcSJacob Faibussowitsch       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i*dim+j]));
103663a3b9bcSJacob Faibussowitsch       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i*dim+j]));
103782490253SJoe Wallwork     }
103882490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
103982490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
104082490253SJoe Wallwork   }
104182490253SJoe Wallwork   PetscFunctionReturn(0);
104282490253SJoe Wallwork }
104382490253SJoe Wallwork 
10443f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
104520282da2SJoe Wallwork {
104620282da2SJoe Wallwork   PetscInt       i, j, k;
104720282da2SJoe 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);
104820282da2SJoe Wallwork   PetscScalar   *Mpos;
104920282da2SJoe Wallwork 
105020282da2SJoe Wallwork   PetscFunctionBegin;
10519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim*dim, &Mpos, dim, &eigs));
105220282da2SJoe Wallwork 
105320282da2SJoe Wallwork   /* Symmetrize */
105420282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
105520282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
105620282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
105720282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
105820282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
105920282da2SJoe Wallwork     }
106020282da2SJoe Wallwork   }
106120282da2SJoe Wallwork 
106220282da2SJoe Wallwork   /* Compute eigendecomposition */
106393520041SJoe Wallwork   if (dim == 1) {
106493520041SJoe Wallwork 
106593520041SJoe Wallwork     /* Isotropic case */
106693520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
106793520041SJoe Wallwork     Mpos[0] = 1.0;
106893520041SJoe Wallwork   } else {
106993520041SJoe Wallwork 
107093520041SJoe Wallwork     /* Anisotropic case */
107120282da2SJoe Wallwork     PetscScalar  *work;
107220282da2SJoe Wallwork     PetscBLASInt lwork;
107320282da2SJoe Wallwork 
107420282da2SJoe Wallwork     lwork = 5*dim;
10759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*dim, &work));
107620282da2SJoe Wallwork     {
107720282da2SJoe Wallwork       PetscBLASInt lierr;
107820282da2SJoe Wallwork       PetscBLASInt nb;
107920282da2SJoe Wallwork 
10809566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
10819566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
108220282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
108320282da2SJoe Wallwork       {
108420282da2SJoe Wallwork         PetscReal *rwork;
10859566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
108620282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
10879566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
108820282da2SJoe Wallwork       }
108920282da2SJoe Wallwork #else
109020282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
109120282da2SJoe Wallwork #endif
109282490253SJoe Wallwork       if (lierr) {
109382490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
109482490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
109582490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
109682490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
109782490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
109882490253SJoe Wallwork           }
109982490253SJoe Wallwork         }
110082490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
110198921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
110282490253SJoe Wallwork       }
11039566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
110420282da2SJoe Wallwork     }
11059566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
110620282da2SJoe Wallwork   }
110720282da2SJoe Wallwork 
110820282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
110920282da2SJoe Wallwork   max_eig = 0.0;
111020282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111120282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
111220282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
111320282da2SJoe Wallwork   }
111420282da2SJoe Wallwork 
11153f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11163f07679eSJoe Wallwork   *detMp = 1.0;
111720282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111820282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
11193f07679eSJoe Wallwork     *detMp *= eigs[i];
112020282da2SJoe Wallwork   }
112120282da2SJoe Wallwork 
112220282da2SJoe Wallwork   /* Reconstruct Hessian */
112320282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
112420282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
112520282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
112620282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
112720282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
112820282da2SJoe Wallwork       }
112920282da2SJoe Wallwork     }
113020282da2SJoe Wallwork   }
11319566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Mpos, eigs));
113220282da2SJoe Wallwork 
113320282da2SJoe Wallwork   PetscFunctionReturn(0);
113420282da2SJoe Wallwork }
113520282da2SJoe Wallwork 
1136cb7bfe8cSJoe Wallwork /*@
113720282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
113820282da2SJoe Wallwork 
113920282da2SJoe Wallwork   Input parameters:
114020282da2SJoe Wallwork + dm                 - The DM
11413f07679eSJoe Wallwork . metricIn           - The metric
114299abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11433f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
114420282da2SJoe Wallwork 
114520282da2SJoe Wallwork   Output parameter:
11463f07679eSJoe Wallwork + metricOut          - The metric
11473f07679eSJoe Wallwork - determinant        - Its determinant
114820282da2SJoe Wallwork 
114920282da2SJoe Wallwork   Level: beginner
115020282da2SJoe Wallwork 
1151cb7bfe8cSJoe Wallwork   Notes:
1152cb7bfe8cSJoe Wallwork 
1153cb7bfe8cSJoe Wallwork   Relevant command line options:
1154cb7bfe8cSJoe Wallwork 
115593520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
115693520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
115793520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1158cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1159cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1160cb7bfe8cSJoe Wallwork 
116120282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
1162cb7bfe8cSJoe Wallwork @*/
11633f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant)
116420282da2SJoe Wallwork {
11653f07679eSJoe Wallwork   DM             dmDet;
116693520041SJoe Wallwork   PetscBool      isotropic, uniform;
116720282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
11683f07679eSJoe Wallwork   PetscScalar   *met, *det;
116920282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
117020282da2SJoe Wallwork 
117120282da2SJoe Wallwork   PetscFunctionBegin;
11729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0));
117320282da2SJoe Wallwork 
117420282da2SJoe Wallwork   /* Extract metadata from dm */
11759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
117620282da2SJoe Wallwork   if (restrictSizes) {
11779566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
11789566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
117999abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
118099abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1181*bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
118299abec2bSJoe Wallwork   }
118399abec2bSJoe Wallwork   if (restrictAnisotropy) {
11849566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
118599abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
118620282da2SJoe Wallwork   }
118720282da2SJoe Wallwork 
118893520041SJoe Wallwork   /* Setup output metric */
11899566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, 0, metricOut));
11909566063dSJacob Faibussowitsch   PetscCall(VecCopy(metricIn, *metricOut));
119193520041SJoe Wallwork 
119293520041SJoe Wallwork   /* Enforce SPD and extract determinant */
11939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*metricOut, &met));
11949566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
11959566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
119693520041SJoe Wallwork   if (uniform) {
1197d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
119893520041SJoe Wallwork 
119993520041SJoe Wallwork     /* Uniform case */
12009566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(metricIn, determinant));
12019566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*determinant, &det));
12029566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
12039566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*determinant, &det));
120493520041SJoe Wallwork   } else {
120593520041SJoe Wallwork 
120693520041SJoe Wallwork     /* Spatially varying case */
120793520041SJoe Wallwork     PetscInt nrow;
120893520041SJoe Wallwork 
120993520041SJoe Wallwork     if (isotropic) nrow = 1;
121093520041SJoe Wallwork     else nrow = dim;
12119566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmDet));
12129566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant));
12139566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12149566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*determinant, &det));
121520282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12163f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
12179566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
12189566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
12199566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
122020282da2SJoe Wallwork     }
12219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*determinant, &det));
122293520041SJoe Wallwork   }
12239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*metricOut, &met));
1224fe902aceSJoe Wallwork 
12259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0));
122620282da2SJoe Wallwork   PetscFunctionReturn(0);
122720282da2SJoe Wallwork }
122820282da2SJoe Wallwork 
122920282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
123020282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
123120282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
123220282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
123320282da2SJoe Wallwork {
123420282da2SJoe Wallwork   const PetscScalar p = constants[0];
123520282da2SJoe Wallwork 
1236985ad86eSJose E. Roman   f0[0] = PetscPowScalar(u[0], p/(2.0*p + dim));
123720282da2SJoe Wallwork }
123820282da2SJoe Wallwork 
1239cb7bfe8cSJoe Wallwork /*@
124020282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
124120282da2SJoe Wallwork 
124220282da2SJoe Wallwork   Input parameters:
124320282da2SJoe Wallwork + dm                 - The DM
124420282da2SJoe Wallwork . metricIn           - The unnormalized metric
124516522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
124616522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
124720282da2SJoe Wallwork 
124820282da2SJoe Wallwork   Output parameter:
124920282da2SJoe Wallwork . metricOut          - The normalized metric
125020282da2SJoe Wallwork 
125120282da2SJoe Wallwork   Level: beginner
125220282da2SJoe Wallwork 
1253cb7bfe8cSJoe Wallwork   Notes:
1254cb7bfe8cSJoe Wallwork 
1255cb7bfe8cSJoe Wallwork   Relevant command line options:
1256cb7bfe8cSJoe Wallwork 
125793520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
125893520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
125993520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1260cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1261cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1262cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1263cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1264cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1265cb7bfe8cSJoe Wallwork 
126620282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
1267cb7bfe8cSJoe Wallwork @*/
126816522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut)
126920282da2SJoe Wallwork {
12703f07679eSJoe Wallwork   DM               dmDet;
127120282da2SJoe Wallwork   MPI_Comm         comm;
127293520041SJoe Wallwork   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
127320282da2SJoe Wallwork   PetscDS          ds;
127420282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
12753f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
127693520041SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
12773f07679eSJoe Wallwork   Vec              determinant;
127820282da2SJoe Wallwork 
127920282da2SJoe Wallwork   PetscFunctionBegin;
12809566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0));
128120282da2SJoe Wallwork 
128220282da2SJoe Wallwork   /* Extract metadata from dm */
12839566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
12849566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
12859566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12869566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
128793520041SJoe Wallwork   if (isotropic) Nd = 1;
128893520041SJoe Wallwork   else Nd = dim*dim;
128920282da2SJoe Wallwork 
129020282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
12919566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
12929566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant));
129320282da2SJoe Wallwork 
129420282da2SJoe Wallwork   /* Compute global normalization factor */
12959566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
12969566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
129716522936SJoe Wallwork   constants[0] = p;
129893520041SJoe Wallwork   if (uniform) {
1299*bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
130093520041SJoe Wallwork     DM  dmTmp;
130193520041SJoe Wallwork     Vec tmp;
130293520041SJoe Wallwork 
13039566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmTmp));
13049566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
13059566063dSJacob Faibussowitsch     PetscCall(VecGetArray(determinant, &det));
13069566063dSJacob Faibussowitsch     PetscCall(VecSet(tmp, det[0]));
13079566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(determinant, &det));
13089566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmTmp, &ds));
13099566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13109566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13119566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
13129566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tmp));
13139566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmTmp));
131493520041SJoe Wallwork   } else {
13159566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
13169566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmDet, &ds));
13179566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13189566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13199566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
132093520041SJoe Wallwork   }
13213f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
1322*bc00d7afSJoe Wallwork   PetscCheck(realIntegral > 1.0e-30, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Global metric normalization factor must be in (0, inf). Is the input metric positive-definite?");
13233f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
132420282da2SJoe Wallwork 
132520282da2SJoe Wallwork   /* Apply local scaling */
132620282da2SJoe Wallwork   if (restrictSizes) {
13279566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
13289566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
132916522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
133016522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1331*bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
133216522936SJoe Wallwork   }
133316522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
13349566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
133516522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
133620282da2SJoe Wallwork   }
13379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*metricOut, &met));
13389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(determinant, &det));
133993520041SJoe Wallwork   if (uniform) {
134093520041SJoe Wallwork 
134193520041SJoe Wallwork     /* Uniform case */
134293520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
13439566063dSJacob Faibussowitsch     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
134493520041SJoe Wallwork   } else {
134593520041SJoe Wallwork 
134693520041SJoe Wallwork     /* Spatially varying case */
134793520041SJoe Wallwork     PetscInt nrow;
134893520041SJoe Wallwork 
134993520041SJoe Wallwork     if (isotropic) nrow = 1;
135093520041SJoe Wallwork     else nrow = dim;
13519566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
135220282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13533f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
135420282da2SJoe Wallwork 
13559566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
13569566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
13573f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
135820282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
13599566063dSJacob Faibussowitsch       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
136093520041SJoe Wallwork     }
136120282da2SJoe Wallwork   }
13629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(determinant, &det));
13639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*metricOut, &met));
13649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&determinant));
13659566063dSJacob Faibussowitsch   if (!uniform) PetscCall(DMDestroy(&dmDet));
136620282da2SJoe Wallwork 
13679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0));
136820282da2SJoe Wallwork   PetscFunctionReturn(0);
136920282da2SJoe Wallwork }
137020282da2SJoe Wallwork 
1371cb7bfe8cSJoe Wallwork /*@
137220282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
137320282da2SJoe Wallwork 
1374f1a722f8SMatthew G. Knepley   Input Parameters:
137520282da2SJoe Wallwork + dm         - The DM
137620282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
137720282da2SJoe Wallwork . weights    - Weights for the average
137820282da2SJoe Wallwork - metrics    - The metrics to be averaged
137920282da2SJoe Wallwork 
138020282da2SJoe Wallwork   Output Parameter:
138120282da2SJoe Wallwork . metricAvg  - The averaged metric
138220282da2SJoe Wallwork 
138320282da2SJoe Wallwork   Level: beginner
138420282da2SJoe Wallwork 
138520282da2SJoe Wallwork   Notes:
138620282da2SJoe Wallwork   The weights should sum to unity.
138720282da2SJoe Wallwork 
138820282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
138920282da2SJoe Wallwork 
139020282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
1391cb7bfe8cSJoe Wallwork @*/
139220282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
139320282da2SJoe Wallwork {
139420282da2SJoe Wallwork   PetscBool haveWeights = PETSC_TRUE;
139593520041SJoe Wallwork   PetscInt  i, m, n;
139620282da2SJoe Wallwork   PetscReal sum = 0.0, tol = 1.0e-10;
139720282da2SJoe Wallwork 
139820282da2SJoe Wallwork   PetscFunctionBegin;
13999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0));
140063a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
14019566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, 0, metricAvg));
14029566063dSJacob Faibussowitsch   PetscCall(VecSet(*metricAvg, 0.0));
14039566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metricAvg, &m));
140493520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
14059566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
14065f80ce2aSJacob Faibussowitsch     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
140793520041SJoe Wallwork   }
140820282da2SJoe Wallwork 
140920282da2SJoe Wallwork   /* Default to the unweighted case */
141020282da2SJoe Wallwork   if (!weights) {
14119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numMetrics, &weights));
141220282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
141320282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
141420282da2SJoe Wallwork   }
141520282da2SJoe Wallwork 
141620282da2SJoe Wallwork   /* Check weights sum to unity */
141793520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
14185f80ce2aSJacob Faibussowitsch   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
141920282da2SJoe Wallwork 
142020282da2SJoe Wallwork   /* Compute metric average */
14219566063dSJacob Faibussowitsch   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(*metricAvg, weights[i], metrics[i]));
14229566063dSJacob Faibussowitsch   if (!haveWeights) PetscCall(PetscFree(weights));
1423fe902aceSJoe Wallwork 
14249566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0));
142520282da2SJoe Wallwork   PetscFunctionReturn(0);
142620282da2SJoe Wallwork }
142720282da2SJoe Wallwork 
1428cb7bfe8cSJoe Wallwork /*@
142920282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
143020282da2SJoe Wallwork 
1431f1a722f8SMatthew G. Knepley   Input Parameters:
143220282da2SJoe Wallwork + dm         - The DM
143320282da2SJoe Wallwork . metric1    - The first metric to be averaged
143420282da2SJoe Wallwork - metric2    - The second metric to be averaged
143520282da2SJoe Wallwork 
143620282da2SJoe Wallwork   Output Parameter:
143720282da2SJoe Wallwork . metricAvg  - The averaged metric
143820282da2SJoe Wallwork 
143920282da2SJoe Wallwork   Level: beginner
144020282da2SJoe Wallwork 
144120282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
1442cb7bfe8cSJoe Wallwork @*/
144320282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
144420282da2SJoe Wallwork {
144520282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
144620282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
144720282da2SJoe Wallwork 
144820282da2SJoe Wallwork   PetscFunctionBegin;
14499566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
145020282da2SJoe Wallwork   PetscFunctionReturn(0);
145120282da2SJoe Wallwork }
145220282da2SJoe Wallwork 
1453cb7bfe8cSJoe Wallwork /*@
145420282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
145520282da2SJoe Wallwork 
1456f1a722f8SMatthew G. Knepley   Input Parameters:
145720282da2SJoe Wallwork + dm         - The DM
145820282da2SJoe Wallwork . metric1    - The first metric to be averaged
145920282da2SJoe Wallwork . metric2    - The second metric to be averaged
146020282da2SJoe Wallwork - metric3    - The third metric to be averaged
146120282da2SJoe Wallwork 
146220282da2SJoe Wallwork   Output Parameter:
146320282da2SJoe Wallwork . metricAvg  - The averaged metric
146420282da2SJoe Wallwork 
146520282da2SJoe Wallwork   Level: beginner
146620282da2SJoe Wallwork 
146720282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
1468cb7bfe8cSJoe Wallwork @*/
146920282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
147020282da2SJoe Wallwork {
147120282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
147220282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
147320282da2SJoe Wallwork 
147420282da2SJoe Wallwork   PetscFunctionBegin;
14759566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
147620282da2SJoe Wallwork   PetscFunctionReturn(0);
147720282da2SJoe Wallwork }
147820282da2SJoe Wallwork 
147920282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
148020282da2SJoe Wallwork {
148120282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
148220282da2SJoe Wallwork   PetscReal     *evals, *evals1;
148320282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
148420282da2SJoe Wallwork 
148520282da2SJoe Wallwork   PetscFunctionBegin;
148693520041SJoe Wallwork 
148793520041SJoe Wallwork   /* Isotropic case */
148893520041SJoe Wallwork   if (dim == 1) {
148993520041SJoe Wallwork     M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
149093520041SJoe Wallwork     PetscFunctionReturn(0);
149193520041SJoe Wallwork   }
149293520041SJoe Wallwork 
149393520041SJoe Wallwork   /* Anisotropic case */
14949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1));
149520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
149620282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
149720282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
149820282da2SJoe Wallwork     }
149920282da2SJoe Wallwork   }
150020282da2SJoe Wallwork   {
150120282da2SJoe Wallwork     PetscScalar *work;
150220282da2SJoe Wallwork     PetscBLASInt lwork;
150320282da2SJoe Wallwork 
150420282da2SJoe Wallwork     lwork = 5*dim;
15059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*dim, &work));
150620282da2SJoe Wallwork     {
150720282da2SJoe Wallwork       PetscBLASInt lierr, nb;
150820282da2SJoe Wallwork       PetscReal    sqrtk;
150920282da2SJoe Wallwork 
151020282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
15119566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
15129566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
151320282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
151420282da2SJoe Wallwork       {
151520282da2SJoe Wallwork         PetscReal *rwork;
15169566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
151720282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
15189566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
151920282da2SJoe Wallwork       }
152020282da2SJoe Wallwork #else
152120282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
152220282da2SJoe Wallwork #endif
152382490253SJoe Wallwork       if (lierr) {
152482490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
152598921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
152682490253SJoe Wallwork       }
15279566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
152820282da2SJoe Wallwork 
152920282da2SJoe Wallwork       /* Compute square root and reciprocal */
153020282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
153120282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
153220282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
153320282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
153420282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
153520282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
153620282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
153720282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
153820282da2SJoe Wallwork           }
153920282da2SJoe Wallwork         }
154020282da2SJoe Wallwork       }
154120282da2SJoe Wallwork 
154220282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
154320282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
154420282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
154520282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
154620282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
154720282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
154820282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
154920282da2SJoe Wallwork             }
155020282da2SJoe Wallwork           }
155120282da2SJoe Wallwork         }
155220282da2SJoe Wallwork       }
155320282da2SJoe Wallwork 
155420282da2SJoe Wallwork       /* Compute eigendecomposition */
15559566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
155620282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
155720282da2SJoe Wallwork       {
155820282da2SJoe Wallwork         PetscReal *rwork;
15599566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
156020282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15619566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
156220282da2SJoe Wallwork       }
156320282da2SJoe Wallwork #else
156420282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
156520282da2SJoe Wallwork #endif
156682490253SJoe Wallwork       if (lierr) {
156782490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
156882490253SJoe Wallwork           for (j = 0; j < dim; ++j) {
156982490253SJoe Wallwork             evecs[i*dim+j] = 0.0;
157082490253SJoe Wallwork             for (k = 0; k < dim; ++k) {
157182490253SJoe Wallwork               for (l = 0; l < dim; ++l) {
157282490253SJoe Wallwork                 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
157382490253SJoe Wallwork               }
157482490253SJoe Wallwork             }
157582490253SJoe Wallwork           }
157682490253SJoe Wallwork         }
157782490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
157898921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
157982490253SJoe Wallwork       }
15809566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
158120282da2SJoe Wallwork 
158220282da2SJoe Wallwork       /* Modify eigenvalues */
158320282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
158420282da2SJoe Wallwork 
158520282da2SJoe Wallwork       /* Map back to get the intersection */
158620282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
158720282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
158820282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
158920282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
159020282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
159120282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
159220282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
159320282da2SJoe Wallwork               }
159420282da2SJoe Wallwork             }
159520282da2SJoe Wallwork           }
159620282da2SJoe Wallwork         }
159720282da2SJoe Wallwork       }
159820282da2SJoe Wallwork     }
15999566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
160020282da2SJoe Wallwork   }
16019566063dSJacob Faibussowitsch   PetscCall(PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1));
160220282da2SJoe Wallwork   PetscFunctionReturn(0);
160320282da2SJoe Wallwork }
160420282da2SJoe Wallwork 
1605cb7bfe8cSJoe Wallwork /*@
160620282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
160720282da2SJoe Wallwork 
1608f1a722f8SMatthew G. Knepley   Input Parameters:
160920282da2SJoe Wallwork + dm         - The DM
161020282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
161120282da2SJoe Wallwork - metrics    - The metrics to be intersected
161220282da2SJoe Wallwork 
161320282da2SJoe Wallwork   Output Parameter:
161420282da2SJoe Wallwork . metricInt  - The intersected metric
161520282da2SJoe Wallwork 
161620282da2SJoe Wallwork   Level: beginner
161720282da2SJoe Wallwork 
161820282da2SJoe Wallwork   Notes:
161920282da2SJoe Wallwork 
162020282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
162120282da2SJoe Wallwork 
162220282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
162320282da2SJoe Wallwork 
162420282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
1625cb7bfe8cSJoe Wallwork @*/
162620282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
162720282da2SJoe Wallwork {
162893520041SJoe Wallwork   PetscBool      isotropic, uniform;
162993520041SJoe Wallwork   PetscInt       v, i, m, n;
163093520041SJoe Wallwork   PetscScalar   *met, *meti;
163120282da2SJoe Wallwork 
163220282da2SJoe Wallwork   PetscFunctionBegin;
16339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0));
163463a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
163520282da2SJoe Wallwork 
163620282da2SJoe Wallwork   /* Copy over the first metric */
16379566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, 0, metricInt));
16389566063dSJacob Faibussowitsch   PetscCall(VecCopy(metrics[0], *metricInt));
163993520041SJoe Wallwork   if (numMetrics == 1) PetscFunctionReturn(0);
16409566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metricInt, &m));
164193520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
16429566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
164308401ef6SPierre Jolivet     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
164493520041SJoe Wallwork   }
164520282da2SJoe Wallwork 
164620282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
16479566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
16489566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
164993520041SJoe Wallwork   if (uniform) {
165093520041SJoe Wallwork 
165193520041SJoe Wallwork     /* Uniform case */
16529566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*metricInt, &met));
165393520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16549566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
16559566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
16569566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
165793520041SJoe Wallwork     }
16589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*metricInt, &met));
165993520041SJoe Wallwork   } else {
166093520041SJoe Wallwork 
166193520041SJoe Wallwork     /* Spatially varying case */
166293520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
166393520041SJoe Wallwork     PetscScalar *M, *Mi;
166493520041SJoe Wallwork 
16659566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
166693520041SJoe Wallwork     if (isotropic) nrow = 1;
166793520041SJoe Wallwork     else nrow = dim;
16689566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
16699566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*metricInt, &met));
167020282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16719566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
167220282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
16739566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
16749566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
16759566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
167620282da2SJoe Wallwork       }
16779566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
167820282da2SJoe Wallwork     }
16799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*metricInt, &met));
168020282da2SJoe Wallwork   }
1681fe902aceSJoe Wallwork 
16829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0));
168320282da2SJoe Wallwork   PetscFunctionReturn(0);
168420282da2SJoe Wallwork }
168520282da2SJoe Wallwork 
1686cb7bfe8cSJoe Wallwork /*@
168720282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
168820282da2SJoe Wallwork 
1689f1a722f8SMatthew G. Knepley   Input Parameters:
169020282da2SJoe Wallwork + dm        - The DM
169120282da2SJoe Wallwork . metric1   - The first metric to be intersected
169220282da2SJoe Wallwork - metric2   - The second metric to be intersected
169320282da2SJoe Wallwork 
169420282da2SJoe Wallwork   Output Parameter:
169520282da2SJoe Wallwork . metricInt - The intersected metric
169620282da2SJoe Wallwork 
169720282da2SJoe Wallwork   Level: beginner
169820282da2SJoe Wallwork 
169920282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
1700cb7bfe8cSJoe Wallwork @*/
170120282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
170220282da2SJoe Wallwork {
170320282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
170420282da2SJoe Wallwork 
170520282da2SJoe Wallwork   PetscFunctionBegin;
17069566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
170720282da2SJoe Wallwork   PetscFunctionReturn(0);
170820282da2SJoe Wallwork }
170920282da2SJoe Wallwork 
1710cb7bfe8cSJoe Wallwork /*@
171120282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
171220282da2SJoe Wallwork 
1713f1a722f8SMatthew G. Knepley   Input Parameters:
171420282da2SJoe Wallwork + dm        - The DM
171520282da2SJoe Wallwork . metric1   - The first metric to be intersected
171620282da2SJoe Wallwork . metric2   - The second metric to be intersected
171720282da2SJoe Wallwork - metric3   - The third metric to be intersected
171820282da2SJoe Wallwork 
171920282da2SJoe Wallwork   Output Parameter:
172020282da2SJoe Wallwork . metricInt - The intersected metric
172120282da2SJoe Wallwork 
172220282da2SJoe Wallwork   Level: beginner
172320282da2SJoe Wallwork 
172420282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
1725cb7bfe8cSJoe Wallwork @*/
172620282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
172720282da2SJoe Wallwork {
172820282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
172920282da2SJoe Wallwork 
173020282da2SJoe Wallwork   PetscFunctionBegin;
17319566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
173220282da2SJoe Wallwork   PetscFunctionReturn(0);
173320282da2SJoe Wallwork }
1734