xref: /petsc/src/dm/impls/plex/plexmetric.c (revision d0609ced746bc51b019815ca91d747429db24893)
120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
220282da2SJoe Wallwork #include <petscblaslapack.h>
320282da2SJoe Wallwork 
4fe902aceSJoe Wallwork PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection;
5fe902aceSJoe Wallwork 
631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
731b70465SJoe Wallwork {
831b70465SJoe Wallwork   MPI_Comm       comm;
993520041SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
1076f360caSJoe Wallwork   PetscBool      noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE;
11cc2eee55SJoe Wallwork   PetscInt       verbosity = -1, numIter = 3;
12ae8b063eSJoe 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;
1331b70465SJoe Wallwork 
1431b70465SJoe Wallwork   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
16*d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL));
189566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, isotropic));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL));
209566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, uniform));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL));
229566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL));
249566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoInsertion(dm, noInsert));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL));
269566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoSwapping(dm, noSwap));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL));
289566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoMovement(dm, noMove));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL));
309566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoSurf(dm, noSurf));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0));
329566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNumIterations(dm, numIter));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10));
349566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetVerbosity(dm, verbosity));
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL));
369566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMinimumMagnitude(dm, h_min));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL));
389566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMaximumMagnitude(dm, h_max));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL));
409566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMaximumAnisotropy(dm, a_max));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL));
429566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNormalizationOrder(dm, p));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL));
449566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetTargetComplexity(dm, target));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL));
469566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetGradationFactor(dm, beta));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL));
489566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetHausdorffNumber(dm, hausd));
49*d0609cedSBarry Smith   PetscOptionsEnd();
5031b70465SJoe Wallwork   PetscFunctionReturn(0);
5131b70465SJoe Wallwork }
5231b70465SJoe Wallwork 
53cb7bfe8cSJoe Wallwork /*@
5431b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
5531b70465SJoe Wallwork 
5631b70465SJoe Wallwork   Input parameters:
5731b70465SJoe Wallwork + dm        - The DM
5831b70465SJoe Wallwork - isotropic - Is the metric isotropic?
5931b70465SJoe Wallwork 
6031b70465SJoe Wallwork   Level: beginner
6131b70465SJoe Wallwork 
6293520041SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst()
63cb7bfe8cSJoe Wallwork @*/
6431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
6531b70465SJoe Wallwork {
6631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6731b70465SJoe Wallwork 
6831b70465SJoe Wallwork   PetscFunctionBegin;
6931b70465SJoe Wallwork   if (!plex->metricCtx) {
709566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
719566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
7231b70465SJoe Wallwork   }
7331b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
7431b70465SJoe Wallwork   PetscFunctionReturn(0);
7531b70465SJoe Wallwork }
7631b70465SJoe Wallwork 
77cb7bfe8cSJoe Wallwork /*@
7893520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7931b70465SJoe Wallwork 
8031b70465SJoe Wallwork   Input parameters:
8131b70465SJoe Wallwork . dm        - The DM
8231b70465SJoe Wallwork 
8331b70465SJoe Wallwork   Output parameters:
8431b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8531b70465SJoe Wallwork 
8631b70465SJoe Wallwork   Level: beginner
8731b70465SJoe Wallwork 
8893520041SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst()
89cb7bfe8cSJoe Wallwork @*/
9031b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
9131b70465SJoe Wallwork {
9231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
9331b70465SJoe Wallwork 
9431b70465SJoe Wallwork   PetscFunctionBegin;
9531b70465SJoe Wallwork   if (!plex->metricCtx) {
969566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
979566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
9831b70465SJoe Wallwork   }
9931b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
10031b70465SJoe Wallwork   PetscFunctionReturn(0);
10131b70465SJoe Wallwork }
10231b70465SJoe Wallwork 
103cb7bfe8cSJoe Wallwork /*@
10493520041SJoe Wallwork   DMPlexMetricSetUniform - Record whether a metric is uniform
10593520041SJoe Wallwork 
10693520041SJoe Wallwork   Input parameters:
10793520041SJoe Wallwork + dm      - The DM
10893520041SJoe Wallwork - uniform - Is the metric uniform?
10993520041SJoe Wallwork 
11093520041SJoe Wallwork   Level: beginner
11193520041SJoe Wallwork 
11293520041SJoe Wallwork   Notes:
11393520041SJoe Wallwork 
11493520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11593520041SJoe Wallwork 
11693520041SJoe Wallwork .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
11793520041SJoe Wallwork @*/
11893520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
11993520041SJoe Wallwork {
12093520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
12193520041SJoe Wallwork 
12293520041SJoe Wallwork   PetscFunctionBegin;
12393520041SJoe Wallwork   if (!plex->metricCtx) {
1249566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
1259566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
12693520041SJoe Wallwork   }
12793520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
12893520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
12993520041SJoe Wallwork   PetscFunctionReturn(0);
13093520041SJoe Wallwork }
13193520041SJoe Wallwork 
13293520041SJoe Wallwork /*@
13393520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
13493520041SJoe Wallwork 
13593520041SJoe Wallwork   Input parameters:
13693520041SJoe Wallwork . dm      - The DM
13793520041SJoe Wallwork 
13893520041SJoe Wallwork   Output parameters:
13993520041SJoe Wallwork . uniform - Is the metric uniform?
14093520041SJoe Wallwork 
14193520041SJoe Wallwork   Level: beginner
14293520041SJoe Wallwork 
14393520041SJoe Wallwork .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
14493520041SJoe Wallwork @*/
14593520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
14693520041SJoe Wallwork {
14793520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
14893520041SJoe Wallwork 
14993520041SJoe Wallwork   PetscFunctionBegin;
15093520041SJoe Wallwork   if (!plex->metricCtx) {
1519566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
1529566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
15393520041SJoe Wallwork   }
15493520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
15593520041SJoe Wallwork   PetscFunctionReturn(0);
15693520041SJoe Wallwork }
15793520041SJoe Wallwork 
15893520041SJoe Wallwork /*@
15931b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
16031b70465SJoe Wallwork 
16131b70465SJoe Wallwork   Input parameters:
16231b70465SJoe Wallwork + dm                      - The DM
16331b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
16431b70465SJoe Wallwork 
16531b70465SJoe Wallwork   Level: beginner
16631b70465SJoe Wallwork 
16731b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
168cb7bfe8cSJoe Wallwork @*/
16931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
17031b70465SJoe Wallwork {
17131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
17231b70465SJoe Wallwork 
17331b70465SJoe Wallwork   PetscFunctionBegin;
17431b70465SJoe Wallwork   if (!plex->metricCtx) {
1759566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
1769566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
17731b70465SJoe Wallwork   }
17831b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
17931b70465SJoe Wallwork   PetscFunctionReturn(0);
18031b70465SJoe Wallwork }
18131b70465SJoe Wallwork 
182cb7bfe8cSJoe Wallwork /*@
18331b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
18431b70465SJoe Wallwork 
18531b70465SJoe Wallwork   Input parameters:
18631b70465SJoe Wallwork . dm                      - The DM
18731b70465SJoe Wallwork 
18831b70465SJoe Wallwork   Output parameters:
18931b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
19031b70465SJoe Wallwork 
19131b70465SJoe Wallwork   Level: beginner
19231b70465SJoe Wallwork 
19331b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
194cb7bfe8cSJoe Wallwork @*/
19531b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
19631b70465SJoe Wallwork {
19731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
19831b70465SJoe Wallwork 
19931b70465SJoe Wallwork   PetscFunctionBegin;
20031b70465SJoe Wallwork   if (!plex->metricCtx) {
2019566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
2029566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
20331b70465SJoe Wallwork   }
20431b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
20531b70465SJoe Wallwork   PetscFunctionReturn(0);
20631b70465SJoe Wallwork }
20731b70465SJoe Wallwork 
208cb7bfe8cSJoe Wallwork /*@
209cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
210cc2eee55SJoe Wallwork 
211cc2eee55SJoe Wallwork   Input parameters:
212cc2eee55SJoe Wallwork + dm       - The DM
213cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
214cc2eee55SJoe Wallwork 
215cc2eee55SJoe Wallwork   Level: beginner
216cc2eee55SJoe Wallwork 
217cb7bfe8cSJoe Wallwork   Notes:
218cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
219cb7bfe8cSJoe Wallwork 
22076f360caSJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf()
221cb7bfe8cSJoe Wallwork @*/
222cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
223cc2eee55SJoe Wallwork {
224cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
225cc2eee55SJoe Wallwork 
226cc2eee55SJoe Wallwork   PetscFunctionBegin;
227cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
2289566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
2299566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
230cc2eee55SJoe Wallwork   }
231cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
232cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
233cc2eee55SJoe Wallwork }
234cc2eee55SJoe Wallwork 
235cb7bfe8cSJoe Wallwork /*@
236cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
237cc2eee55SJoe Wallwork 
238cc2eee55SJoe Wallwork   Input parameters:
239cc2eee55SJoe Wallwork . dm       - The DM
240cc2eee55SJoe Wallwork 
241cc2eee55SJoe Wallwork   Output parameters:
242cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
243cc2eee55SJoe Wallwork 
244cc2eee55SJoe Wallwork   Level: beginner
245cc2eee55SJoe Wallwork 
246cb7bfe8cSJoe Wallwork   Notes:
247cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
248cb7bfe8cSJoe Wallwork 
24976f360caSJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf()
250cb7bfe8cSJoe Wallwork @*/
251cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
252cc2eee55SJoe Wallwork {
253cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
254cc2eee55SJoe Wallwork 
255cc2eee55SJoe Wallwork   PetscFunctionBegin;
256cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
2579566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
2589566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
259cc2eee55SJoe Wallwork   }
260cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
261cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
262cc2eee55SJoe Wallwork }
263cc2eee55SJoe Wallwork 
264cb7bfe8cSJoe Wallwork /*@
265cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
266cc2eee55SJoe Wallwork 
267cc2eee55SJoe Wallwork   Input parameters:
268cc2eee55SJoe Wallwork + dm     - The DM
269cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
270cc2eee55SJoe Wallwork 
271cc2eee55SJoe Wallwork   Level: beginner
272cc2eee55SJoe Wallwork 
273cb7bfe8cSJoe Wallwork   Notes:
274cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
275cb7bfe8cSJoe Wallwork 
27676f360caSJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf()
277cb7bfe8cSJoe Wallwork @*/
278cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
279cc2eee55SJoe Wallwork {
280cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
281cc2eee55SJoe Wallwork 
282cc2eee55SJoe Wallwork   PetscFunctionBegin;
283cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
2849566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
2859566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
286cc2eee55SJoe Wallwork   }
287cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
288cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
289cc2eee55SJoe Wallwork }
290cc2eee55SJoe Wallwork 
291cb7bfe8cSJoe Wallwork /*@
292cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
293cc2eee55SJoe Wallwork 
294cc2eee55SJoe Wallwork   Input parameters:
295cc2eee55SJoe Wallwork . dm     - The DM
296cc2eee55SJoe Wallwork 
297cc2eee55SJoe Wallwork   Output parameters:
298cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
299cc2eee55SJoe Wallwork 
300cc2eee55SJoe Wallwork   Level: beginner
301cc2eee55SJoe Wallwork 
302cb7bfe8cSJoe Wallwork   Notes:
303cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
304cb7bfe8cSJoe Wallwork 
30576f360caSJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf()
306cb7bfe8cSJoe Wallwork @*/
307cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
308cc2eee55SJoe Wallwork {
309cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
310cc2eee55SJoe Wallwork 
311cc2eee55SJoe Wallwork   PetscFunctionBegin;
312cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
3139566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
3149566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
315cc2eee55SJoe Wallwork   }
316cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
317cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
318cc2eee55SJoe Wallwork }
319cc2eee55SJoe Wallwork 
320cb7bfe8cSJoe Wallwork /*@
321cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
322cc2eee55SJoe Wallwork 
323cc2eee55SJoe Wallwork   Input parameters:
324cc2eee55SJoe Wallwork + dm     - The DM
325cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
326cc2eee55SJoe Wallwork 
327cc2eee55SJoe Wallwork   Level: beginner
328cc2eee55SJoe Wallwork 
329cb7bfe8cSJoe Wallwork   Notes:
330cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
331cb7bfe8cSJoe Wallwork 
33276f360caSJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoSurf()
333cb7bfe8cSJoe Wallwork @*/
334cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
335cc2eee55SJoe Wallwork {
336cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
337cc2eee55SJoe Wallwork 
338cc2eee55SJoe Wallwork   PetscFunctionBegin;
339cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
3409566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
3419566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
342cc2eee55SJoe Wallwork   }
343cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
344cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
345cc2eee55SJoe Wallwork }
346cc2eee55SJoe Wallwork 
347cb7bfe8cSJoe Wallwork /*@
348cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
349cc2eee55SJoe Wallwork 
350cc2eee55SJoe Wallwork   Input parameters:
351cc2eee55SJoe Wallwork . dm     - The DM
352cc2eee55SJoe Wallwork 
353cc2eee55SJoe Wallwork   Output parameters:
354cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
355cc2eee55SJoe Wallwork 
356cc2eee55SJoe Wallwork   Level: beginner
357cc2eee55SJoe Wallwork 
358cb7bfe8cSJoe Wallwork   Notes:
359cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
360cb7bfe8cSJoe Wallwork 
36176f360caSJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoSurf()
362cb7bfe8cSJoe Wallwork @*/
363cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
364cc2eee55SJoe Wallwork {
365cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
366cc2eee55SJoe Wallwork 
367cc2eee55SJoe Wallwork   PetscFunctionBegin;
368cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
3699566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
3709566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
371cc2eee55SJoe Wallwork   }
372cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
373cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
374cc2eee55SJoe Wallwork }
375cc2eee55SJoe Wallwork 
376cb7bfe8cSJoe Wallwork /*@
37776f360caSJoe Wallwork   DMPlexMetricSetNoSurf - Should surface modification be turned off?
37876f360caSJoe Wallwork 
37976f360caSJoe Wallwork   Input parameters:
38076f360caSJoe Wallwork + dm     - The DM
38176f360caSJoe Wallwork - noSurf - Should surface modification be turned off?
38276f360caSJoe Wallwork 
38376f360caSJoe Wallwork   Level: beginner
38476f360caSJoe Wallwork 
38576f360caSJoe Wallwork   Notes:
38676f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
38776f360caSJoe Wallwork 
38876f360caSJoe Wallwork .seealso: DMPlexMetricNoSurf(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping()
38976f360caSJoe Wallwork @*/
39076f360caSJoe Wallwork PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
39176f360caSJoe Wallwork {
39276f360caSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
39376f360caSJoe Wallwork 
39476f360caSJoe Wallwork   PetscFunctionBegin;
39576f360caSJoe Wallwork   if (!plex->metricCtx) {
3969566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
3979566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
39876f360caSJoe Wallwork   }
39976f360caSJoe Wallwork   plex->metricCtx->noSurf = noSurf;
40076f360caSJoe Wallwork   PetscFunctionReturn(0);
40176f360caSJoe Wallwork }
40276f360caSJoe Wallwork 
40376f360caSJoe Wallwork /*@
40476f360caSJoe Wallwork   DMPlexMetricNoSurf - Is surface modification turned off?
40576f360caSJoe Wallwork 
40676f360caSJoe Wallwork   Input parameters:
40776f360caSJoe Wallwork . dm     - The DM
40876f360caSJoe Wallwork 
40976f360caSJoe Wallwork   Output parameters:
41076f360caSJoe Wallwork . noSurf - Is surface modification turned off?
41176f360caSJoe Wallwork 
41276f360caSJoe Wallwork   Level: beginner
41376f360caSJoe Wallwork 
41476f360caSJoe Wallwork   Notes:
41576f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
41676f360caSJoe Wallwork 
41776f360caSJoe Wallwork .seealso: DMPlexMetricSetNoSurf(), DMPlexMetricNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping()
41876f360caSJoe Wallwork @*/
41976f360caSJoe Wallwork PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
42076f360caSJoe Wallwork {
42176f360caSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
42276f360caSJoe Wallwork 
42376f360caSJoe Wallwork   PetscFunctionBegin;
42476f360caSJoe Wallwork   if (!plex->metricCtx) {
4259566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
4269566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
42776f360caSJoe Wallwork   }
42876f360caSJoe Wallwork   *noSurf = plex->metricCtx->noSurf;
42976f360caSJoe Wallwork   PetscFunctionReturn(0);
43076f360caSJoe Wallwork }
43176f360caSJoe Wallwork 
43276f360caSJoe Wallwork /*@
43331b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
43431b70465SJoe Wallwork 
43531b70465SJoe Wallwork   Input parameters:
43631b70465SJoe Wallwork + dm    - The DM
43731b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
43831b70465SJoe Wallwork 
43931b70465SJoe Wallwork   Level: beginner
44031b70465SJoe Wallwork 
44131b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
442cb7bfe8cSJoe Wallwork @*/
44331b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
44431b70465SJoe Wallwork {
44531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
44631b70465SJoe Wallwork 
44731b70465SJoe Wallwork   PetscFunctionBegin;
44831b70465SJoe Wallwork   if (!plex->metricCtx) {
4499566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
4509566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
45131b70465SJoe Wallwork   }
45208401ef6SPierre Jolivet   PetscCheck(h_min > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min);
45331b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
45431b70465SJoe Wallwork   PetscFunctionReturn(0);
45531b70465SJoe Wallwork }
45631b70465SJoe Wallwork 
457cb7bfe8cSJoe Wallwork /*@
45831b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
45931b70465SJoe Wallwork 
46031b70465SJoe Wallwork   Input parameters:
46131b70465SJoe Wallwork . dm    - The DM
46231b70465SJoe Wallwork 
463cc2eee55SJoe Wallwork   Output parameters:
46431b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
46531b70465SJoe Wallwork 
46631b70465SJoe Wallwork   Level: beginner
46731b70465SJoe Wallwork 
46831b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
469cb7bfe8cSJoe Wallwork @*/
47031b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
47131b70465SJoe Wallwork {
47231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
47331b70465SJoe Wallwork 
47431b70465SJoe Wallwork   PetscFunctionBegin;
47531b70465SJoe Wallwork   if (!plex->metricCtx) {
4769566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
4779566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
47831b70465SJoe Wallwork   }
47931b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
48031b70465SJoe Wallwork   PetscFunctionReturn(0);
48131b70465SJoe Wallwork }
48231b70465SJoe Wallwork 
483cb7bfe8cSJoe Wallwork /*@
48431b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
48531b70465SJoe Wallwork 
48631b70465SJoe Wallwork   Input parameters:
48731b70465SJoe Wallwork + dm    - The DM
48831b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
48931b70465SJoe Wallwork 
49031b70465SJoe Wallwork   Level: beginner
49131b70465SJoe Wallwork 
49231b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
493cb7bfe8cSJoe Wallwork @*/
49431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
49531b70465SJoe Wallwork {
49631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
49731b70465SJoe Wallwork 
49831b70465SJoe Wallwork   PetscFunctionBegin;
49931b70465SJoe Wallwork   if (!plex->metricCtx) {
5009566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
5019566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
50231b70465SJoe Wallwork   }
50308401ef6SPierre Jolivet   PetscCheck(h_max > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max);
50431b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
50531b70465SJoe Wallwork   PetscFunctionReturn(0);
50631b70465SJoe Wallwork }
50731b70465SJoe Wallwork 
508cb7bfe8cSJoe Wallwork /*@
50931b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
51031b70465SJoe Wallwork 
51131b70465SJoe Wallwork   Input parameters:
51231b70465SJoe Wallwork . dm    - The DM
51331b70465SJoe Wallwork 
514cc2eee55SJoe Wallwork   Output parameters:
51531b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
51631b70465SJoe Wallwork 
51731b70465SJoe Wallwork   Level: beginner
51831b70465SJoe Wallwork 
51931b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
520cb7bfe8cSJoe Wallwork @*/
52131b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
52231b70465SJoe Wallwork {
52331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
52431b70465SJoe Wallwork 
52531b70465SJoe Wallwork   PetscFunctionBegin;
52631b70465SJoe Wallwork   if (!plex->metricCtx) {
5279566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
5289566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
52931b70465SJoe Wallwork   }
53031b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
53131b70465SJoe Wallwork   PetscFunctionReturn(0);
53231b70465SJoe Wallwork }
53331b70465SJoe Wallwork 
534cb7bfe8cSJoe Wallwork /*@
53531b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
53631b70465SJoe Wallwork 
53731b70465SJoe Wallwork   Input parameters:
53831b70465SJoe Wallwork + dm    - The DM
53931b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
54031b70465SJoe Wallwork 
54131b70465SJoe Wallwork   Level: beginner
54231b70465SJoe Wallwork 
54331b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
54431b70465SJoe Wallwork 
54531b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
546cb7bfe8cSJoe Wallwork @*/
54731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
54831b70465SJoe Wallwork {
54931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
55031b70465SJoe Wallwork 
55131b70465SJoe Wallwork   PetscFunctionBegin;
55231b70465SJoe Wallwork   if (!plex->metricCtx) {
5539566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
5549566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
55531b70465SJoe Wallwork   }
55608401ef6SPierre Jolivet   PetscCheck(a_max >= 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max);
55731b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
55831b70465SJoe Wallwork   PetscFunctionReturn(0);
55931b70465SJoe Wallwork }
56031b70465SJoe Wallwork 
561cb7bfe8cSJoe Wallwork /*@
56231b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
56331b70465SJoe Wallwork 
56431b70465SJoe Wallwork   Input parameters:
56531b70465SJoe Wallwork . dm    - The DM
56631b70465SJoe Wallwork 
567cc2eee55SJoe Wallwork   Output parameters:
56831b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
56931b70465SJoe Wallwork 
57031b70465SJoe Wallwork   Level: beginner
57131b70465SJoe Wallwork 
57231b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
573cb7bfe8cSJoe Wallwork @*/
57431b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
57531b70465SJoe Wallwork {
57631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
57731b70465SJoe Wallwork 
57831b70465SJoe Wallwork   PetscFunctionBegin;
57931b70465SJoe Wallwork   if (!plex->metricCtx) {
5809566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
5819566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
58231b70465SJoe Wallwork   }
58331b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
58431b70465SJoe Wallwork   PetscFunctionReturn(0);
58531b70465SJoe Wallwork }
58631b70465SJoe Wallwork 
587cb7bfe8cSJoe Wallwork /*@
58831b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
58931b70465SJoe Wallwork 
59031b70465SJoe Wallwork   Input parameters:
59131b70465SJoe Wallwork + dm               - The DM
59231b70465SJoe Wallwork - targetComplexity - The target metric complexity
59331b70465SJoe Wallwork 
59431b70465SJoe Wallwork   Level: beginner
59531b70465SJoe Wallwork 
59631b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
597cb7bfe8cSJoe Wallwork @*/
59831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
59931b70465SJoe Wallwork {
60031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
60131b70465SJoe Wallwork 
60231b70465SJoe Wallwork   PetscFunctionBegin;
60331b70465SJoe Wallwork   if (!plex->metricCtx) {
6049566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
6059566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
60631b70465SJoe Wallwork   }
60708401ef6SPierre Jolivet   PetscCheck(targetComplexity > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive");
60831b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
60931b70465SJoe Wallwork   PetscFunctionReturn(0);
61031b70465SJoe Wallwork }
61131b70465SJoe Wallwork 
612cb7bfe8cSJoe Wallwork /*@
61331b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
61431b70465SJoe Wallwork 
61531b70465SJoe Wallwork   Input parameters:
61631b70465SJoe Wallwork . dm               - The DM
61731b70465SJoe Wallwork 
618cc2eee55SJoe Wallwork   Output parameters:
61931b70465SJoe Wallwork . targetComplexity - The target metric complexity
62031b70465SJoe Wallwork 
62131b70465SJoe Wallwork   Level: beginner
62231b70465SJoe Wallwork 
62331b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
624cb7bfe8cSJoe Wallwork @*/
62531b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
62631b70465SJoe Wallwork {
62731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
62831b70465SJoe Wallwork 
62931b70465SJoe Wallwork   PetscFunctionBegin;
63031b70465SJoe Wallwork   if (!plex->metricCtx) {
6319566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
6329566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
63331b70465SJoe Wallwork   }
63431b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
63531b70465SJoe Wallwork   PetscFunctionReturn(0);
63631b70465SJoe Wallwork }
63731b70465SJoe Wallwork 
638cb7bfe8cSJoe Wallwork /*@
63931b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
64031b70465SJoe Wallwork 
64131b70465SJoe Wallwork   Input parameters:
64231b70465SJoe Wallwork + dm - The DM
64331b70465SJoe Wallwork - p  - The normalization order
64431b70465SJoe Wallwork 
64531b70465SJoe Wallwork   Level: beginner
64631b70465SJoe Wallwork 
64731b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
648cb7bfe8cSJoe Wallwork @*/
64931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
65031b70465SJoe Wallwork {
65131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
65231b70465SJoe Wallwork 
65331b70465SJoe Wallwork   PetscFunctionBegin;
65431b70465SJoe Wallwork   if (!plex->metricCtx) {
6559566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
6569566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
65731b70465SJoe Wallwork   }
65808401ef6SPierre Jolivet   PetscCheck(p >= 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater");
65931b70465SJoe Wallwork   plex->metricCtx->p = p;
66031b70465SJoe Wallwork   PetscFunctionReturn(0);
66131b70465SJoe Wallwork }
66231b70465SJoe Wallwork 
663cb7bfe8cSJoe Wallwork /*@
66431b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
66531b70465SJoe Wallwork 
66631b70465SJoe Wallwork   Input parameters:
66731b70465SJoe Wallwork . dm - The DM
66831b70465SJoe Wallwork 
669cc2eee55SJoe Wallwork   Output parameters:
67031b70465SJoe Wallwork . p - The normalization order
67131b70465SJoe Wallwork 
67231b70465SJoe Wallwork   Level: beginner
67331b70465SJoe Wallwork 
67431b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
675cb7bfe8cSJoe Wallwork @*/
67631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
67731b70465SJoe Wallwork {
67831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
67931b70465SJoe Wallwork 
68031b70465SJoe Wallwork   PetscFunctionBegin;
68131b70465SJoe Wallwork   if (!plex->metricCtx) {
6829566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
6839566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
68431b70465SJoe Wallwork   }
68531b70465SJoe Wallwork   *p = plex->metricCtx->p;
68631b70465SJoe Wallwork   PetscFunctionReturn(0);
68731b70465SJoe Wallwork }
68831b70465SJoe Wallwork 
689cb7bfe8cSJoe Wallwork /*@
690cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
691cc2eee55SJoe Wallwork 
692cc2eee55SJoe Wallwork   Input parameters:
693cc2eee55SJoe Wallwork + dm   - The DM
694cc2eee55SJoe Wallwork - beta - The metric gradation factor
695cc2eee55SJoe Wallwork 
696cc2eee55SJoe Wallwork   Level: beginner
697cc2eee55SJoe Wallwork 
698cc2eee55SJoe Wallwork   Notes:
699cc2eee55SJoe Wallwork 
700cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
701cc2eee55SJoe Wallwork 
702cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
703cc2eee55SJoe Wallwork 
704cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
705cb7bfe8cSJoe Wallwork 
706ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber()
707cb7bfe8cSJoe Wallwork @*/
708cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
709cc2eee55SJoe Wallwork {
710cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
711cc2eee55SJoe Wallwork 
712cc2eee55SJoe Wallwork   PetscFunctionBegin;
713cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
7149566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
7159566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
716cc2eee55SJoe Wallwork   }
717cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
718cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
719cc2eee55SJoe Wallwork }
720cc2eee55SJoe Wallwork 
721cb7bfe8cSJoe Wallwork /*@
722cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
723cc2eee55SJoe Wallwork 
724cc2eee55SJoe Wallwork   Input parameters:
725cc2eee55SJoe Wallwork . dm   - The DM
726cc2eee55SJoe Wallwork 
727cc2eee55SJoe Wallwork   Output parameters:
728cc2eee55SJoe Wallwork . beta - The metric gradation factor
729cc2eee55SJoe Wallwork 
730cc2eee55SJoe Wallwork   Level: beginner
731cc2eee55SJoe Wallwork 
732cb7bfe8cSJoe Wallwork   Notes:
733cb7bfe8cSJoe Wallwork 
734cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
735cb7bfe8cSJoe Wallwork 
736cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
737cb7bfe8cSJoe Wallwork 
738cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
739cc2eee55SJoe Wallwork 
740ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber()
741cb7bfe8cSJoe Wallwork @*/
742cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
743cc2eee55SJoe Wallwork {
744cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
745cc2eee55SJoe Wallwork 
746cc2eee55SJoe Wallwork   PetscFunctionBegin;
747cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
7489566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
7499566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
750cc2eee55SJoe Wallwork   }
751cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
752cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
753cc2eee55SJoe Wallwork }
754cc2eee55SJoe Wallwork 
755cb7bfe8cSJoe Wallwork /*@
756ae8b063eSJoe Wallwork   DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
757ae8b063eSJoe Wallwork 
758ae8b063eSJoe Wallwork   Input parameters:
759ae8b063eSJoe Wallwork + dm    - The DM
760ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number
761ae8b063eSJoe Wallwork 
762ae8b063eSJoe Wallwork   Level: beginner
763ae8b063eSJoe Wallwork 
764ae8b063eSJoe Wallwork   Notes:
765ae8b063eSJoe Wallwork 
766ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
767ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
768ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
769ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
770ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
771ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
772ae8b063eSJoe Wallwork 
773ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
774ae8b063eSJoe Wallwork 
775ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber()
776ae8b063eSJoe Wallwork @*/
777ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
778ae8b063eSJoe Wallwork {
779ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
780ae8b063eSJoe Wallwork 
781ae8b063eSJoe Wallwork   PetscFunctionBegin;
782ae8b063eSJoe Wallwork   if (!plex->metricCtx) {
7839566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
7849566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
785ae8b063eSJoe Wallwork   }
786ae8b063eSJoe Wallwork   plex->metricCtx->hausdorffNumber = hausd;
787ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
788ae8b063eSJoe Wallwork }
789ae8b063eSJoe Wallwork 
790ae8b063eSJoe Wallwork /*@
791ae8b063eSJoe Wallwork   DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
792ae8b063eSJoe Wallwork 
793ae8b063eSJoe Wallwork   Input parameters:
794ae8b063eSJoe Wallwork . dm    - The DM
795ae8b063eSJoe Wallwork 
796ae8b063eSJoe Wallwork   Output parameters:
797ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number
798ae8b063eSJoe Wallwork 
799ae8b063eSJoe Wallwork   Level: beginner
800ae8b063eSJoe Wallwork 
801ae8b063eSJoe Wallwork   Notes:
802ae8b063eSJoe Wallwork 
803ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
804ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
805ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
806ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
807ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
808ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
809ae8b063eSJoe Wallwork 
810ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
811ae8b063eSJoe Wallwork 
812ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber()
813ae8b063eSJoe Wallwork @*/
814ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
815ae8b063eSJoe Wallwork {
816ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
817ae8b063eSJoe Wallwork 
818ae8b063eSJoe Wallwork   PetscFunctionBegin;
819ae8b063eSJoe Wallwork   if (!plex->metricCtx) {
8209566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
8219566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
822ae8b063eSJoe Wallwork   }
823ae8b063eSJoe Wallwork   *hausd = plex->metricCtx->hausdorffNumber;
824ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
825ae8b063eSJoe Wallwork }
826ae8b063eSJoe Wallwork 
827ae8b063eSJoe Wallwork /*@
828cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
829cc2eee55SJoe Wallwork 
830cc2eee55SJoe Wallwork   Input parameters:
831cc2eee55SJoe Wallwork + dm        - The DM
832cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
833cc2eee55SJoe Wallwork 
834cb7bfe8cSJoe Wallwork   Level: beginner
835cb7bfe8cSJoe Wallwork 
836cb7bfe8cSJoe Wallwork   Notes:
837cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
838cb7bfe8cSJoe Wallwork 
839cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations()
840cb7bfe8cSJoe Wallwork @*/
841cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
842cc2eee55SJoe Wallwork {
843cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
844cc2eee55SJoe Wallwork 
845cc2eee55SJoe Wallwork   PetscFunctionBegin;
846cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
8479566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
8489566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
849cc2eee55SJoe Wallwork   }
850cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
851cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
852cc2eee55SJoe Wallwork }
853cc2eee55SJoe Wallwork 
854cb7bfe8cSJoe Wallwork /*@
855cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
856cc2eee55SJoe Wallwork 
857cc2eee55SJoe Wallwork   Input parameters:
858cc2eee55SJoe Wallwork . dm        - The DM
859cc2eee55SJoe Wallwork 
860cc2eee55SJoe Wallwork   Output parameters:
861cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
862cc2eee55SJoe Wallwork 
863cb7bfe8cSJoe Wallwork   Level: beginner
864cb7bfe8cSJoe Wallwork 
865cb7bfe8cSJoe Wallwork   Notes:
866cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
867cb7bfe8cSJoe Wallwork 
868cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
869cb7bfe8cSJoe Wallwork @*/
870cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
871cc2eee55SJoe Wallwork {
872cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
873cc2eee55SJoe Wallwork 
874cc2eee55SJoe Wallwork   PetscFunctionBegin;
875cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
8769566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
8779566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
878cc2eee55SJoe Wallwork   }
879cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
880cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
881cc2eee55SJoe Wallwork }
882cc2eee55SJoe Wallwork 
883cb7bfe8cSJoe Wallwork /*@
884cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
885cc2eee55SJoe Wallwork 
886cc2eee55SJoe Wallwork   Input parameters:
887cc2eee55SJoe Wallwork + dm      - The DM
888cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
889cc2eee55SJoe Wallwork 
890cb7bfe8cSJoe Wallwork   Level: beginner
891cb7bfe8cSJoe Wallwork 
892cb7bfe8cSJoe Wallwork   Notes:
893cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
894cc2eee55SJoe Wallwork 
895cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
896cb7bfe8cSJoe Wallwork @*/
897cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
898cc2eee55SJoe Wallwork {
899cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
900cc2eee55SJoe Wallwork 
901cc2eee55SJoe Wallwork   PetscFunctionBegin;
902cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
9039566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
9049566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
905cc2eee55SJoe Wallwork   }
906cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
907cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
908cc2eee55SJoe Wallwork }
909cc2eee55SJoe Wallwork 
910cb7bfe8cSJoe Wallwork /*@
911cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
912cc2eee55SJoe Wallwork 
913cc2eee55SJoe Wallwork   Input parameters:
914cc2eee55SJoe Wallwork . dm      - The DM
915cc2eee55SJoe Wallwork 
916cc2eee55SJoe Wallwork   Output parameters:
917cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
918cc2eee55SJoe Wallwork 
919cb7bfe8cSJoe Wallwork   Level: beginner
920cb7bfe8cSJoe Wallwork 
921cb7bfe8cSJoe Wallwork   Notes:
922cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
923cc2eee55SJoe Wallwork 
924cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity()
925cb7bfe8cSJoe Wallwork @*/
926cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
927cc2eee55SJoe Wallwork {
928cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
929cc2eee55SJoe Wallwork 
930cc2eee55SJoe Wallwork   PetscFunctionBegin;
931cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
9329566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
9339566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
934cc2eee55SJoe Wallwork   }
935cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
936cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
937cc2eee55SJoe Wallwork }
938cc2eee55SJoe Wallwork 
93920282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
94020282da2SJoe Wallwork {
94120282da2SJoe Wallwork   MPI_Comm       comm;
94220282da2SJoe Wallwork   PetscFE        fe;
94320282da2SJoe Wallwork   PetscInt       dim;
94420282da2SJoe Wallwork 
94520282da2SJoe Wallwork   PetscFunctionBegin;
94620282da2SJoe Wallwork 
94720282da2SJoe Wallwork   /* Extract metadata from dm */
9489566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
9499566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
95020282da2SJoe Wallwork 
95120282da2SJoe Wallwork   /* Create a P1 field of the requested size */
9529566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
9539566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
9549566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
9559566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
9569566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm, metric));
95720282da2SJoe Wallwork 
95820282da2SJoe Wallwork   PetscFunctionReturn(0);
95920282da2SJoe Wallwork }
96020282da2SJoe Wallwork 
961cb7bfe8cSJoe Wallwork /*@
96220282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
96320282da2SJoe Wallwork 
96420282da2SJoe Wallwork   Input parameters:
96520282da2SJoe Wallwork + dm     - The DM
96620282da2SJoe Wallwork - f      - The field number to use
96720282da2SJoe Wallwork 
96820282da2SJoe Wallwork   Output parameter:
96920282da2SJoe Wallwork . metric - The metric
97020282da2SJoe Wallwork 
97120282da2SJoe Wallwork   Level: beginner
97220282da2SJoe Wallwork 
97331b70465SJoe Wallwork   Notes:
97431b70465SJoe Wallwork 
97531b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
97631b70465SJoe Wallwork 
97731b70465SJoe Wallwork   Command line options for Riemannian metrics:
97831b70465SJoe Wallwork 
979cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
98093520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
981cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
982cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
983cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
984cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
985cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
98667b8a455SSatish Balay - -dm_plex_metric_target_complexity         - Target metric complexity
987cb7bfe8cSJoe Wallwork 
988cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
989cb7bfe8cSJoe Wallwork 
99067b8a455SSatish Balay . -dm_adaptor <pragmatic/mmg/parmmg>        - specify dm adaptor to use
991cb7bfe8cSJoe Wallwork 
992cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
993cb7bfe8cSJoe Wallwork 
994cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
995cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
996cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
997cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
998cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
999cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
100020282da2SJoe Wallwork 
100120282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
1002cb7bfe8cSJoe Wallwork @*/
100320282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
100420282da2SJoe Wallwork {
100531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
100693520041SJoe Wallwork   PetscBool      isotropic, uniform;
100720282da2SJoe Wallwork   PetscInt       coordDim, Nd;
100820282da2SJoe Wallwork 
100920282da2SJoe Wallwork   PetscFunctionBegin;
101031b70465SJoe Wallwork   if (!plex->metricCtx) {
10119566063dSJacob Faibussowitsch     PetscCall(PetscNew(&plex->metricCtx));
10129566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricSetFromOptions(dm));
101331b70465SJoe Wallwork   }
10149566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
101520282da2SJoe Wallwork   Nd = coordDim*coordDim;
10169566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
10179566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
101893520041SJoe Wallwork   if (uniform) {
101993520041SJoe Wallwork     MPI_Comm comm;
102093520041SJoe Wallwork 
10219566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
102228b400f6SJacob Faibussowitsch     PetscCheck(isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
10239566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, metric));
10249566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
10259566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(*metric));
102693520041SJoe Wallwork   } else if (isotropic) {
10279566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
102893520041SJoe Wallwork   } else {
10299566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
103093520041SJoe Wallwork   }
103120282da2SJoe Wallwork   PetscFunctionReturn(0);
103220282da2SJoe Wallwork }
103320282da2SJoe Wallwork 
1034cb7bfe8cSJoe Wallwork /*@
103520282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
103620282da2SJoe Wallwork 
103720282da2SJoe Wallwork   Input parameters:
103820282da2SJoe Wallwork + dm     - The DM
103920282da2SJoe Wallwork . f      - The field number to use
104020282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
104120282da2SJoe Wallwork 
104220282da2SJoe Wallwork   Output parameter:
104320282da2SJoe Wallwork . metric - The uniform metric
104420282da2SJoe Wallwork 
104520282da2SJoe Wallwork   Level: beginner
104620282da2SJoe Wallwork 
104793520041SJoe Wallwork   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
104820282da2SJoe Wallwork 
104920282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
1050cb7bfe8cSJoe Wallwork @*/
105120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
105220282da2SJoe Wallwork {
105320282da2SJoe Wallwork   PetscFunctionBegin;
10549566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
10559566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
10562c71b3e2SJacob Faibussowitsch   PetscCheck(alpha,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
10572c71b3e2SJacob Faibussowitsch   PetscCheck(alpha >= 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
10589566063dSJacob Faibussowitsch   PetscCall(VecSet(*metric, alpha));
10599566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*metric));
10609566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*metric));
106120282da2SJoe Wallwork   PetscFunctionReturn(0);
106220282da2SJoe Wallwork }
106320282da2SJoe Wallwork 
106493520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
106593520041SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
106693520041SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106793520041SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
106893520041SJoe Wallwork {
106993520041SJoe Wallwork   f0[0] = u[0];
107093520041SJoe Wallwork }
107193520041SJoe Wallwork 
1072cb7bfe8cSJoe Wallwork /*@
107320282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
107420282da2SJoe Wallwork 
107520282da2SJoe Wallwork   Input parameters:
107620282da2SJoe Wallwork + dm        - The DM
107720282da2SJoe Wallwork . f         - The field number to use
107820282da2SJoe Wallwork - indicator - The error indicator
107920282da2SJoe Wallwork 
108020282da2SJoe Wallwork   Output parameter:
108120282da2SJoe Wallwork . metric    - The isotropic metric
108220282da2SJoe Wallwork 
108320282da2SJoe Wallwork   Level: beginner
108420282da2SJoe Wallwork 
108520282da2SJoe Wallwork   Notes:
108620282da2SJoe Wallwork 
108720282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
108820282da2SJoe Wallwork 
108993520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
109020282da2SJoe Wallwork 
109120282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
1092cb7bfe8cSJoe Wallwork @*/
109320282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
109420282da2SJoe Wallwork {
109593520041SJoe Wallwork   PetscInt       m, n;
109620282da2SJoe Wallwork 
109720282da2SJoe Wallwork   PetscFunctionBegin;
10989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
10999566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
11009566063dSJacob Faibussowitsch   PetscCall(VecGetSize(indicator, &m));
11019566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metric, &n));
11029566063dSJacob Faibussowitsch   if (m == n) PetscCall(VecCopy(indicator, *metric));
110393520041SJoe Wallwork   else {
110493520041SJoe Wallwork     DM     dmIndi;
110593520041SJoe Wallwork     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
110693520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
110793520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
110893520041SJoe Wallwork                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
110993520041SJoe Wallwork 
11109566063dSJacob Faibussowitsch     PetscCall(VecGetDM(indicator, &dmIndi));
111193520041SJoe Wallwork     funcs[0] = identity;
11129566063dSJacob Faibussowitsch     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
111320282da2SJoe Wallwork   }
111420282da2SJoe Wallwork   PetscFunctionReturn(0);
111520282da2SJoe Wallwork }
111620282da2SJoe Wallwork 
111782490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
111882490253SJoe Wallwork {
111982490253SJoe Wallwork   PetscInt i, j;
112082490253SJoe Wallwork 
112182490253SJoe Wallwork   PetscFunctionBegin;
112282490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
112382490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
112482490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
112582490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
112682490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
112782490253SJoe Wallwork       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]);
112882490253SJoe Wallwork       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]);
112982490253SJoe Wallwork     }
113082490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
113182490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
113282490253SJoe Wallwork   }
113382490253SJoe Wallwork   PetscFunctionReturn(0);
113482490253SJoe Wallwork }
113582490253SJoe Wallwork 
11363f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
113720282da2SJoe Wallwork {
113820282da2SJoe Wallwork   PetscInt       i, j, k;
113920282da2SJoe 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);
114020282da2SJoe Wallwork   PetscScalar   *Mpos;
114120282da2SJoe Wallwork 
114220282da2SJoe Wallwork   PetscFunctionBegin;
11439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim*dim, &Mpos, dim, &eigs));
114420282da2SJoe Wallwork 
114520282da2SJoe Wallwork   /* Symmetrize */
114620282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
114720282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
114820282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
114920282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
115020282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
115120282da2SJoe Wallwork     }
115220282da2SJoe Wallwork   }
115320282da2SJoe Wallwork 
115420282da2SJoe Wallwork   /* Compute eigendecomposition */
115593520041SJoe Wallwork   if (dim == 1) {
115693520041SJoe Wallwork 
115793520041SJoe Wallwork     /* Isotropic case */
115893520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
115993520041SJoe Wallwork     Mpos[0] = 1.0;
116093520041SJoe Wallwork   } else {
116193520041SJoe Wallwork 
116293520041SJoe Wallwork     /* Anisotropic case */
116320282da2SJoe Wallwork     PetscScalar  *work;
116420282da2SJoe Wallwork     PetscBLASInt lwork;
116520282da2SJoe Wallwork 
116620282da2SJoe Wallwork     lwork = 5*dim;
11679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*dim, &work));
116820282da2SJoe Wallwork     {
116920282da2SJoe Wallwork       PetscBLASInt lierr;
117020282da2SJoe Wallwork       PetscBLASInt nb;
117120282da2SJoe Wallwork 
11729566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
11739566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
117420282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
117520282da2SJoe Wallwork       {
117620282da2SJoe Wallwork         PetscReal *rwork;
11779566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
117820282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
11799566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
118020282da2SJoe Wallwork       }
118120282da2SJoe Wallwork #else
118220282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
118320282da2SJoe Wallwork #endif
118482490253SJoe Wallwork       if (lierr) {
118582490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
118682490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
118782490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
118882490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
118982490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
119082490253SJoe Wallwork           }
119182490253SJoe Wallwork         }
119282490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
119398921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
119482490253SJoe Wallwork       }
11959566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
119620282da2SJoe Wallwork     }
11979566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
119820282da2SJoe Wallwork   }
119920282da2SJoe Wallwork 
120020282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
120120282da2SJoe Wallwork   max_eig = 0.0;
120220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
120320282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
120420282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
120520282da2SJoe Wallwork   }
120620282da2SJoe Wallwork 
12073f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
12083f07679eSJoe Wallwork   *detMp = 1.0;
120920282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
121020282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
12113f07679eSJoe Wallwork     *detMp *= eigs[i];
121220282da2SJoe Wallwork   }
121320282da2SJoe Wallwork 
121420282da2SJoe Wallwork   /* Reconstruct Hessian */
121520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
121620282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
121720282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
121820282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
121920282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
122020282da2SJoe Wallwork       }
122120282da2SJoe Wallwork     }
122220282da2SJoe Wallwork   }
12239566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Mpos, eigs));
122420282da2SJoe Wallwork 
122520282da2SJoe Wallwork   PetscFunctionReturn(0);
122620282da2SJoe Wallwork }
122720282da2SJoe Wallwork 
1228cb7bfe8cSJoe Wallwork /*@
122920282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
123020282da2SJoe Wallwork 
123120282da2SJoe Wallwork   Input parameters:
123220282da2SJoe Wallwork + dm                 - The DM
12333f07679eSJoe Wallwork . metricIn           - The metric
123499abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
12353f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
123620282da2SJoe Wallwork 
123720282da2SJoe Wallwork   Output parameter:
12383f07679eSJoe Wallwork + metricOut          - The metric
12393f07679eSJoe Wallwork - determinant        - Its determinant
124020282da2SJoe Wallwork 
124120282da2SJoe Wallwork   Level: beginner
124220282da2SJoe Wallwork 
1243cb7bfe8cSJoe Wallwork   Notes:
1244cb7bfe8cSJoe Wallwork 
1245cb7bfe8cSJoe Wallwork   Relevant command line options:
1246cb7bfe8cSJoe Wallwork 
124793520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
124893520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
124993520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1250cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1251cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1252cb7bfe8cSJoe Wallwork 
125320282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
1254cb7bfe8cSJoe Wallwork @*/
12553f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant)
125620282da2SJoe Wallwork {
12573f07679eSJoe Wallwork   DM             dmDet;
125893520041SJoe Wallwork   PetscBool      isotropic, uniform;
125920282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
12603f07679eSJoe Wallwork   PetscScalar   *met, *det;
126120282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
126220282da2SJoe Wallwork 
126320282da2SJoe Wallwork   PetscFunctionBegin;
12649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0));
126520282da2SJoe Wallwork 
126620282da2SJoe Wallwork   /* Extract metadata from dm */
12679566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
126820282da2SJoe Wallwork   if (restrictSizes) {
12699566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
12709566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
127199abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
127299abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
127354c59aa7SJacob Faibussowitsch     PetscCheck(h_min < h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
127499abec2bSJoe Wallwork   }
127599abec2bSJoe Wallwork   if (restrictAnisotropy) {
12769566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
127799abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
127820282da2SJoe Wallwork   }
127920282da2SJoe Wallwork 
128093520041SJoe Wallwork   /* Setup output metric */
12819566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, 0, metricOut));
12829566063dSJacob Faibussowitsch   PetscCall(VecCopy(metricIn, *metricOut));
128393520041SJoe Wallwork 
128493520041SJoe Wallwork   /* Enforce SPD and extract determinant */
12859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*metricOut, &met));
12869566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12879566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
128893520041SJoe Wallwork   if (uniform) {
1289d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
129093520041SJoe Wallwork 
129193520041SJoe Wallwork     /* Uniform case */
12929566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(metricIn, determinant));
12939566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*determinant, &det));
12949566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
12959566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*determinant, &det));
129693520041SJoe Wallwork   } else {
129793520041SJoe Wallwork 
129893520041SJoe Wallwork     /* Spatially varying case */
129993520041SJoe Wallwork     PetscInt nrow;
130093520041SJoe Wallwork 
130193520041SJoe Wallwork     if (isotropic) nrow = 1;
130293520041SJoe Wallwork     else nrow = dim;
13039566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmDet));
13049566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant));
13059566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
13069566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*determinant, &det));
130720282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13083f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
13099566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
13109566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
13119566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
131220282da2SJoe Wallwork     }
13139566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*determinant, &det));
131493520041SJoe Wallwork   }
13159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*metricOut, &met));
1316fe902aceSJoe Wallwork 
13179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0));
131820282da2SJoe Wallwork   PetscFunctionReturn(0);
131920282da2SJoe Wallwork }
132020282da2SJoe Wallwork 
132120282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
132220282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
132320282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
132420282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
132520282da2SJoe Wallwork {
132620282da2SJoe Wallwork   const PetscScalar p = constants[0];
132720282da2SJoe Wallwork 
1328985ad86eSJose E. Roman   f0[0] = PetscPowScalar(u[0], p/(2.0*p + dim));
132920282da2SJoe Wallwork }
133020282da2SJoe Wallwork 
1331cb7bfe8cSJoe Wallwork /*@
133220282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
133320282da2SJoe Wallwork 
133420282da2SJoe Wallwork   Input parameters:
133520282da2SJoe Wallwork + dm                 - The DM
133620282da2SJoe Wallwork . metricIn           - The unnormalized metric
133716522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
133816522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
133920282da2SJoe Wallwork 
134020282da2SJoe Wallwork   Output parameter:
134120282da2SJoe Wallwork . metricOut          - The normalized metric
134220282da2SJoe Wallwork 
134320282da2SJoe Wallwork   Level: beginner
134420282da2SJoe Wallwork 
1345cb7bfe8cSJoe Wallwork   Notes:
1346cb7bfe8cSJoe Wallwork 
1347cb7bfe8cSJoe Wallwork   Relevant command line options:
1348cb7bfe8cSJoe Wallwork 
134993520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
135093520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
135193520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1352cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1353cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1354cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1355cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1356cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1357cb7bfe8cSJoe Wallwork 
135820282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
1359cb7bfe8cSJoe Wallwork @*/
136016522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut)
136120282da2SJoe Wallwork {
13623f07679eSJoe Wallwork   DM               dmDet;
136320282da2SJoe Wallwork   MPI_Comm         comm;
136493520041SJoe Wallwork   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
136520282da2SJoe Wallwork   PetscDS          ds;
136620282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
13673f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
136893520041SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
13693f07679eSJoe Wallwork   Vec              determinant;
137020282da2SJoe Wallwork 
137120282da2SJoe Wallwork   PetscFunctionBegin;
13729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0));
137320282da2SJoe Wallwork 
137420282da2SJoe Wallwork   /* Extract metadata from dm */
13759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
13769566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13779566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
13789566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
137993520041SJoe Wallwork   if (isotropic) Nd = 1;
138093520041SJoe Wallwork   else Nd = dim*dim;
138120282da2SJoe Wallwork 
138220282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
13839566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
13849566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant));
138520282da2SJoe Wallwork 
138620282da2SJoe Wallwork   /* Compute global normalization factor */
13879566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
13889566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
138916522936SJoe Wallwork   constants[0] = p;
139093520041SJoe Wallwork   if (uniform) {
139128b400f6SJacob Faibussowitsch     PetscCheck(isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported");
139293520041SJoe Wallwork     DM  dmTmp;
139393520041SJoe Wallwork     Vec tmp;
139493520041SJoe Wallwork 
13959566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmTmp));
13969566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
13979566063dSJacob Faibussowitsch     PetscCall(VecGetArray(determinant, &det));
13989566063dSJacob Faibussowitsch     PetscCall(VecSet(tmp, det[0]));
13999566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(determinant, &det));
14009566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmTmp, &ds));
14019566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
14029566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
14039566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
14049566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tmp));
14059566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmTmp));
140693520041SJoe Wallwork   } else {
14079566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
14089566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmDet, &ds));
14099566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
14109566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
14119566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
141293520041SJoe Wallwork   }
14133f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
141408401ef6SPierre Jolivet   PetscCheck(realIntegral >= 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global metric normalization factor should be strictly positive, not %.4e Is the input metric positive-definite?", realIntegral);
14153f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
141620282da2SJoe Wallwork 
141720282da2SJoe Wallwork   /* Apply local scaling */
141820282da2SJoe Wallwork   if (restrictSizes) {
14199566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
14209566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
142116522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
142216522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
142308401ef6SPierre Jolivet     PetscCheck(h_min < h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
142416522936SJoe Wallwork   }
142516522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
14269566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
142716522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
142820282da2SJoe Wallwork   }
14299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*metricOut, &met));
14309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(determinant, &det));
143193520041SJoe Wallwork   if (uniform) {
143293520041SJoe Wallwork 
143393520041SJoe Wallwork     /* Uniform case */
143493520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
14359566063dSJacob Faibussowitsch     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
143693520041SJoe Wallwork   } else {
143793520041SJoe Wallwork 
143893520041SJoe Wallwork     /* Spatially varying case */
143993520041SJoe Wallwork     PetscInt nrow;
144093520041SJoe Wallwork 
144193520041SJoe Wallwork     if (isotropic) nrow = 1;
144293520041SJoe Wallwork     else nrow = dim;
14439566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
144420282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
14453f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
144620282da2SJoe Wallwork 
14479566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
14489566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
14493f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
145020282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
14519566063dSJacob Faibussowitsch       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
145293520041SJoe Wallwork     }
145320282da2SJoe Wallwork   }
14549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(determinant, &det));
14559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*metricOut, &met));
14569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&determinant));
14579566063dSJacob Faibussowitsch   if (!uniform) PetscCall(DMDestroy(&dmDet));
145820282da2SJoe Wallwork 
14599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0));
146020282da2SJoe Wallwork   PetscFunctionReturn(0);
146120282da2SJoe Wallwork }
146220282da2SJoe Wallwork 
1463cb7bfe8cSJoe Wallwork /*@
146420282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
146520282da2SJoe Wallwork 
1466f1a722f8SMatthew G. Knepley   Input Parameters:
146720282da2SJoe Wallwork + dm         - The DM
146820282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
146920282da2SJoe Wallwork . weights    - Weights for the average
147020282da2SJoe Wallwork - metrics    - The metrics to be averaged
147120282da2SJoe Wallwork 
147220282da2SJoe Wallwork   Output Parameter:
147320282da2SJoe Wallwork . metricAvg  - The averaged metric
147420282da2SJoe Wallwork 
147520282da2SJoe Wallwork   Level: beginner
147620282da2SJoe Wallwork 
147720282da2SJoe Wallwork   Notes:
147820282da2SJoe Wallwork   The weights should sum to unity.
147920282da2SJoe Wallwork 
148020282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
148120282da2SJoe Wallwork 
148220282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
1483cb7bfe8cSJoe Wallwork @*/
148420282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
148520282da2SJoe Wallwork {
148620282da2SJoe Wallwork   PetscBool haveWeights = PETSC_TRUE;
148793520041SJoe Wallwork   PetscInt  i, m, n;
148820282da2SJoe Wallwork   PetscReal sum = 0.0, tol = 1.0e-10;
148920282da2SJoe Wallwork 
149020282da2SJoe Wallwork   PetscFunctionBegin;
14919566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0));
14922c71b3e2SJacob Faibussowitsch   PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics);
14939566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, 0, metricAvg));
14949566063dSJacob Faibussowitsch   PetscCall(VecSet(*metricAvg, 0.0));
14959566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metricAvg, &m));
149693520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
14979566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
14985f80ce2aSJacob Faibussowitsch     PetscCheck(m == n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
149993520041SJoe Wallwork   }
150020282da2SJoe Wallwork 
150120282da2SJoe Wallwork   /* Default to the unweighted case */
150220282da2SJoe Wallwork   if (!weights) {
15039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numMetrics, &weights));
150420282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
150520282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
150620282da2SJoe Wallwork   }
150720282da2SJoe Wallwork 
150820282da2SJoe Wallwork   /* Check weights sum to unity */
150993520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
15105f80ce2aSJacob Faibussowitsch   PetscCheck(PetscAbsReal(sum - 1) <= tol,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
151120282da2SJoe Wallwork 
151220282da2SJoe Wallwork   /* Compute metric average */
15139566063dSJacob Faibussowitsch   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(*metricAvg, weights[i], metrics[i]));
15149566063dSJacob Faibussowitsch   if (!haveWeights) PetscCall(PetscFree(weights));
1515fe902aceSJoe Wallwork 
15169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0));
151720282da2SJoe Wallwork   PetscFunctionReturn(0);
151820282da2SJoe Wallwork }
151920282da2SJoe Wallwork 
1520cb7bfe8cSJoe Wallwork /*@
152120282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
152220282da2SJoe Wallwork 
1523f1a722f8SMatthew G. Knepley   Input Parameters:
152420282da2SJoe Wallwork + dm         - The DM
152520282da2SJoe Wallwork . metric1    - The first metric to be averaged
152620282da2SJoe Wallwork - metric2    - The second metric to be averaged
152720282da2SJoe Wallwork 
152820282da2SJoe Wallwork   Output Parameter:
152920282da2SJoe Wallwork . metricAvg  - The averaged metric
153020282da2SJoe Wallwork 
153120282da2SJoe Wallwork   Level: beginner
153220282da2SJoe Wallwork 
153320282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
1534cb7bfe8cSJoe Wallwork @*/
153520282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
153620282da2SJoe Wallwork {
153720282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
153820282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
153920282da2SJoe Wallwork 
154020282da2SJoe Wallwork   PetscFunctionBegin;
15419566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
154220282da2SJoe Wallwork   PetscFunctionReturn(0);
154320282da2SJoe Wallwork }
154420282da2SJoe Wallwork 
1545cb7bfe8cSJoe Wallwork /*@
154620282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
154720282da2SJoe Wallwork 
1548f1a722f8SMatthew G. Knepley   Input Parameters:
154920282da2SJoe Wallwork + dm         - The DM
155020282da2SJoe Wallwork . metric1    - The first metric to be averaged
155120282da2SJoe Wallwork . metric2    - The second metric to be averaged
155220282da2SJoe Wallwork - metric3    - The third metric to be averaged
155320282da2SJoe Wallwork 
155420282da2SJoe Wallwork   Output Parameter:
155520282da2SJoe Wallwork . metricAvg  - The averaged metric
155620282da2SJoe Wallwork 
155720282da2SJoe Wallwork   Level: beginner
155820282da2SJoe Wallwork 
155920282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
1560cb7bfe8cSJoe Wallwork @*/
156120282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
156220282da2SJoe Wallwork {
156320282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
156420282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
156520282da2SJoe Wallwork 
156620282da2SJoe Wallwork   PetscFunctionBegin;
15679566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
156820282da2SJoe Wallwork   PetscFunctionReturn(0);
156920282da2SJoe Wallwork }
157020282da2SJoe Wallwork 
157120282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
157220282da2SJoe Wallwork {
157320282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
157420282da2SJoe Wallwork   PetscReal     *evals, *evals1;
157520282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
157620282da2SJoe Wallwork 
157720282da2SJoe Wallwork   PetscFunctionBegin;
157893520041SJoe Wallwork 
157993520041SJoe Wallwork   /* Isotropic case */
158093520041SJoe Wallwork   if (dim == 1) {
158193520041SJoe Wallwork     M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
158293520041SJoe Wallwork     PetscFunctionReturn(0);
158393520041SJoe Wallwork   }
158493520041SJoe Wallwork 
158593520041SJoe Wallwork   /* Anisotropic case */
15869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1));
158720282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
158820282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
158920282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
159020282da2SJoe Wallwork     }
159120282da2SJoe Wallwork   }
159220282da2SJoe Wallwork   {
159320282da2SJoe Wallwork     PetscScalar *work;
159420282da2SJoe Wallwork     PetscBLASInt lwork;
159520282da2SJoe Wallwork 
159620282da2SJoe Wallwork     lwork = 5*dim;
15979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*dim, &work));
159820282da2SJoe Wallwork     {
159920282da2SJoe Wallwork       PetscBLASInt lierr, nb;
160020282da2SJoe Wallwork       PetscReal    sqrtk;
160120282da2SJoe Wallwork 
160220282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
16039566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
16049566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
160520282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
160620282da2SJoe Wallwork       {
160720282da2SJoe Wallwork         PetscReal *rwork;
16089566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
160920282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
16109566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
161120282da2SJoe Wallwork       }
161220282da2SJoe Wallwork #else
161320282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
161420282da2SJoe Wallwork #endif
161582490253SJoe Wallwork       if (lierr) {
161682490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
161798921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
161882490253SJoe Wallwork       }
16199566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
162020282da2SJoe Wallwork 
162120282da2SJoe Wallwork       /* Compute square root and reciprocal */
162220282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
162320282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
162420282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
162520282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
162620282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
162720282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
162820282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
162920282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
163020282da2SJoe Wallwork           }
163120282da2SJoe Wallwork         }
163220282da2SJoe Wallwork       }
163320282da2SJoe Wallwork 
163420282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
163520282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
163620282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
163720282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
163820282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
163920282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
164020282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
164120282da2SJoe Wallwork             }
164220282da2SJoe Wallwork           }
164320282da2SJoe Wallwork         }
164420282da2SJoe Wallwork       }
164520282da2SJoe Wallwork 
164620282da2SJoe Wallwork       /* Compute eigendecomposition */
16479566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
164820282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
164920282da2SJoe Wallwork       {
165020282da2SJoe Wallwork         PetscReal *rwork;
16519566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
165220282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
16539566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
165420282da2SJoe Wallwork       }
165520282da2SJoe Wallwork #else
165620282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
165720282da2SJoe Wallwork #endif
165882490253SJoe Wallwork       if (lierr) {
165982490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
166082490253SJoe Wallwork           for (j = 0; j < dim; ++j) {
166182490253SJoe Wallwork             evecs[i*dim+j] = 0.0;
166282490253SJoe Wallwork             for (k = 0; k < dim; ++k) {
166382490253SJoe Wallwork               for (l = 0; l < dim; ++l) {
166482490253SJoe Wallwork                 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
166582490253SJoe Wallwork               }
166682490253SJoe Wallwork             }
166782490253SJoe Wallwork           }
166882490253SJoe Wallwork         }
166982490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
167098921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
167182490253SJoe Wallwork       }
16729566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
167320282da2SJoe Wallwork 
167420282da2SJoe Wallwork       /* Modify eigenvalues */
167520282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
167620282da2SJoe Wallwork 
167720282da2SJoe Wallwork       /* Map back to get the intersection */
167820282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
167920282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
168020282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
168120282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
168220282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
168320282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
168420282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
168520282da2SJoe Wallwork               }
168620282da2SJoe Wallwork             }
168720282da2SJoe Wallwork           }
168820282da2SJoe Wallwork         }
168920282da2SJoe Wallwork       }
169020282da2SJoe Wallwork     }
16919566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
169220282da2SJoe Wallwork   }
16939566063dSJacob Faibussowitsch   PetscCall(PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1));
169420282da2SJoe Wallwork   PetscFunctionReturn(0);
169520282da2SJoe Wallwork }
169620282da2SJoe Wallwork 
1697cb7bfe8cSJoe Wallwork /*@
169820282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
169920282da2SJoe Wallwork 
1700f1a722f8SMatthew G. Knepley   Input Parameters:
170120282da2SJoe Wallwork + dm         - The DM
170220282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
170320282da2SJoe Wallwork - metrics    - The metrics to be intersected
170420282da2SJoe Wallwork 
170520282da2SJoe Wallwork   Output Parameter:
170620282da2SJoe Wallwork . metricInt  - The intersected metric
170720282da2SJoe Wallwork 
170820282da2SJoe Wallwork   Level: beginner
170920282da2SJoe Wallwork 
171020282da2SJoe Wallwork   Notes:
171120282da2SJoe Wallwork 
171220282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
171320282da2SJoe Wallwork 
171420282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
171520282da2SJoe Wallwork 
171620282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
1717cb7bfe8cSJoe Wallwork @*/
171820282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
171920282da2SJoe Wallwork {
172093520041SJoe Wallwork   PetscBool      isotropic, uniform;
172193520041SJoe Wallwork   PetscInt       v, i, m, n;
172293520041SJoe Wallwork   PetscScalar   *met, *meti;
172320282da2SJoe Wallwork 
172420282da2SJoe Wallwork   PetscFunctionBegin;
17259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0));
17262c71b3e2SJacob Faibussowitsch   PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics);
172720282da2SJoe Wallwork 
172820282da2SJoe Wallwork   /* Copy over the first metric */
17299566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, 0, metricInt));
17309566063dSJacob Faibussowitsch   PetscCall(VecCopy(metrics[0], *metricInt));
173193520041SJoe Wallwork   if (numMetrics == 1) PetscFunctionReturn(0);
17329566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metricInt, &m));
173393520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
17349566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
173508401ef6SPierre Jolivet     PetscCheck(m == n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
173693520041SJoe Wallwork   }
173720282da2SJoe Wallwork 
173820282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
17399566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
17409566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
174193520041SJoe Wallwork   if (uniform) {
174293520041SJoe Wallwork 
174393520041SJoe Wallwork     /* Uniform case */
17449566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*metricInt, &met));
174593520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
17469566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
17479566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
17489566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
174993520041SJoe Wallwork     }
17509566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*metricInt, &met));
175193520041SJoe Wallwork   } else {
175293520041SJoe Wallwork 
175393520041SJoe Wallwork     /* Spatially varying case */
175493520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
175593520041SJoe Wallwork     PetscScalar *M, *Mi;
175693520041SJoe Wallwork 
17579566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
175893520041SJoe Wallwork     if (isotropic) nrow = 1;
175993520041SJoe Wallwork     else nrow = dim;
17609566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
17619566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*metricInt, &met));
176220282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
17639566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
176420282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
17659566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
17669566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
17679566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
176820282da2SJoe Wallwork       }
17699566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
177020282da2SJoe Wallwork     }
17719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*metricInt, &met));
177220282da2SJoe Wallwork   }
1773fe902aceSJoe Wallwork 
17749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0));
177520282da2SJoe Wallwork   PetscFunctionReturn(0);
177620282da2SJoe Wallwork }
177720282da2SJoe Wallwork 
1778cb7bfe8cSJoe Wallwork /*@
177920282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
178020282da2SJoe Wallwork 
1781f1a722f8SMatthew G. Knepley   Input Parameters:
178220282da2SJoe Wallwork + dm        - The DM
178320282da2SJoe Wallwork . metric1   - The first metric to be intersected
178420282da2SJoe Wallwork - metric2   - The second metric to be intersected
178520282da2SJoe Wallwork 
178620282da2SJoe Wallwork   Output Parameter:
178720282da2SJoe Wallwork . metricInt - The intersected metric
178820282da2SJoe Wallwork 
178920282da2SJoe Wallwork   Level: beginner
179020282da2SJoe Wallwork 
179120282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
1792cb7bfe8cSJoe Wallwork @*/
179320282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
179420282da2SJoe Wallwork {
179520282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
179620282da2SJoe Wallwork 
179720282da2SJoe Wallwork   PetscFunctionBegin;
17989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
179920282da2SJoe Wallwork   PetscFunctionReturn(0);
180020282da2SJoe Wallwork }
180120282da2SJoe Wallwork 
1802cb7bfe8cSJoe Wallwork /*@
180320282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
180420282da2SJoe Wallwork 
1805f1a722f8SMatthew G. Knepley   Input Parameters:
180620282da2SJoe Wallwork + dm        - The DM
180720282da2SJoe Wallwork . metric1   - The first metric to be intersected
180820282da2SJoe Wallwork . metric2   - The second metric to be intersected
180920282da2SJoe Wallwork - metric3   - The third metric to be intersected
181020282da2SJoe Wallwork 
181120282da2SJoe Wallwork   Output Parameter:
181220282da2SJoe Wallwork . metricInt - The intersected metric
181320282da2SJoe Wallwork 
181420282da2SJoe Wallwork   Level: beginner
181520282da2SJoe Wallwork 
181620282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
1817cb7bfe8cSJoe Wallwork @*/
181820282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
181920282da2SJoe Wallwork {
182020282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
182120282da2SJoe Wallwork 
182220282da2SJoe Wallwork   PetscFunctionBegin;
18239566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
182420282da2SJoe Wallwork   PetscFunctionReturn(0);
182520282da2SJoe Wallwork }
1826