xref: /petsc/src/dm/impls/plex/plexmetric.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
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 {
8bc00d7afSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
931b70465SJoe Wallwork   MPI_Comm       comm;
1093520041SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
1176f360caSJoe Wallwork   PetscBool      noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE;
12cc2eee55SJoe Wallwork   PetscInt       verbosity = -1, numIter = 3;
13ae8b063eSJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0, beta = 1.3, hausd = 0.01;
1431b70465SJoe Wallwork 
1531b70465SJoe Wallwork   PetscFunctionBegin;
16bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(PetscNew(&plex->metricCtx));
179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
18d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL));
209566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, isotropic));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL));
229566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, uniform));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL));
249566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL));
269566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoInsertion(dm, noInsert));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL));
289566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoSwapping(dm, noSwap));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL));
309566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoMovement(dm, noMove));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL));
329566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNoSurf(dm, noSurf));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0));
349566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNumIterations(dm, numIter));
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10));
369566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetVerbosity(dm, verbosity));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL));
389566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMinimumMagnitude(dm, h_min));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL));
409566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMaximumMagnitude(dm, h_max));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL));
429566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetMaximumAnisotropy(dm, a_max));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL));
449566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetNormalizationOrder(dm, p));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL));
469566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetTargetComplexity(dm, target));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL));
489566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetGradationFactor(dm, beta));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL));
509566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetHausdorffNumber(dm, hausd));
51d0609cedSBarry Smith   PetscOptionsEnd();
5231b70465SJoe Wallwork   PetscFunctionReturn(0);
5331b70465SJoe Wallwork }
5431b70465SJoe Wallwork 
55cb7bfe8cSJoe Wallwork /*@
5631b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
5731b70465SJoe Wallwork 
5831b70465SJoe Wallwork   Input parameters:
5931b70465SJoe Wallwork + dm        - The DM
6031b70465SJoe Wallwork - isotropic - Is the metric isotropic?
6131b70465SJoe Wallwork 
6231b70465SJoe Wallwork   Level: beginner
6331b70465SJoe Wallwork 
64db781477SPatrick Sanan .seealso: `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetUniform()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
65cb7bfe8cSJoe Wallwork @*/
6631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
6731b70465SJoe Wallwork {
6831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6931b70465SJoe Wallwork 
7031b70465SJoe Wallwork   PetscFunctionBegin;
71bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
7231b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
7331b70465SJoe Wallwork   PetscFunctionReturn(0);
7431b70465SJoe Wallwork }
7531b70465SJoe Wallwork 
76cb7bfe8cSJoe Wallwork /*@
7793520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7831b70465SJoe Wallwork 
7931b70465SJoe Wallwork   Input parameters:
8031b70465SJoe Wallwork . dm        - The DM
8131b70465SJoe Wallwork 
8231b70465SJoe Wallwork   Output parameters:
8331b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8431b70465SJoe Wallwork 
8531b70465SJoe Wallwork   Level: beginner
8631b70465SJoe Wallwork 
87db781477SPatrick Sanan .seealso: `DMPlexMetricSetIsotropic()`, `DMPlexMetricIsUniform()`, `DMPlexMetricRestrictAnisotropyFirst()`
88cb7bfe8cSJoe Wallwork @*/
8931b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
9031b70465SJoe Wallwork {
9131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
9231b70465SJoe Wallwork 
9331b70465SJoe Wallwork   PetscFunctionBegin;
94bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
9531b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
9631b70465SJoe Wallwork   PetscFunctionReturn(0);
9731b70465SJoe Wallwork }
9831b70465SJoe Wallwork 
99cb7bfe8cSJoe Wallwork /*@
10093520041SJoe Wallwork   DMPlexMetricSetUniform - Record whether a metric is uniform
10193520041SJoe Wallwork 
10293520041SJoe Wallwork   Input parameters:
10393520041SJoe Wallwork + dm      - The DM
10493520041SJoe Wallwork - uniform - Is the metric uniform?
10593520041SJoe Wallwork 
10693520041SJoe Wallwork   Level: beginner
10793520041SJoe Wallwork 
10893520041SJoe Wallwork   Notes:
10993520041SJoe Wallwork 
11093520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11193520041SJoe Wallwork 
112db781477SPatrick Sanan .seealso: `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
11393520041SJoe Wallwork @*/
11493520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
11593520041SJoe Wallwork {
11693520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
11793520041SJoe Wallwork 
11893520041SJoe Wallwork   PetscFunctionBegin;
119bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
12093520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
12193520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
12293520041SJoe Wallwork   PetscFunctionReturn(0);
12393520041SJoe Wallwork }
12493520041SJoe Wallwork 
12593520041SJoe Wallwork /*@
12693520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
12793520041SJoe Wallwork 
12893520041SJoe Wallwork   Input parameters:
12993520041SJoe Wallwork . dm      - The DM
13093520041SJoe Wallwork 
13193520041SJoe Wallwork   Output parameters:
13293520041SJoe Wallwork . uniform - Is the metric uniform?
13393520041SJoe Wallwork 
13493520041SJoe Wallwork   Level: beginner
13593520041SJoe Wallwork 
136db781477SPatrick Sanan .seealso: `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
13793520041SJoe Wallwork @*/
13893520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
13993520041SJoe Wallwork {
14093520041SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
14193520041SJoe Wallwork 
14293520041SJoe Wallwork   PetscFunctionBegin;
143bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
14493520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
14593520041SJoe Wallwork   PetscFunctionReturn(0);
14693520041SJoe Wallwork }
14793520041SJoe Wallwork 
14893520041SJoe Wallwork /*@
14931b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
15031b70465SJoe Wallwork 
15131b70465SJoe Wallwork   Input parameters:
15231b70465SJoe Wallwork + dm                      - The DM
15331b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
15431b70465SJoe Wallwork 
15531b70465SJoe Wallwork   Level: beginner
15631b70465SJoe Wallwork 
157db781477SPatrick Sanan .seealso: `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
158cb7bfe8cSJoe Wallwork @*/
15931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
16031b70465SJoe Wallwork {
16131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
16231b70465SJoe Wallwork 
16331b70465SJoe Wallwork   PetscFunctionBegin;
164bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
16531b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
16631b70465SJoe Wallwork   PetscFunctionReturn(0);
16731b70465SJoe Wallwork }
16831b70465SJoe Wallwork 
169cb7bfe8cSJoe Wallwork /*@
17031b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
17131b70465SJoe Wallwork 
17231b70465SJoe Wallwork   Input parameters:
17331b70465SJoe Wallwork . dm                      - The DM
17431b70465SJoe Wallwork 
17531b70465SJoe Wallwork   Output parameters:
17631b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
17731b70465SJoe Wallwork 
17831b70465SJoe Wallwork   Level: beginner
17931b70465SJoe Wallwork 
180db781477SPatrick Sanan .seealso: `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
181cb7bfe8cSJoe Wallwork @*/
18231b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
18331b70465SJoe Wallwork {
18431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
18531b70465SJoe Wallwork 
18631b70465SJoe Wallwork   PetscFunctionBegin;
187bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
18831b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
18931b70465SJoe Wallwork   PetscFunctionReturn(0);
19031b70465SJoe Wallwork }
19131b70465SJoe Wallwork 
192cb7bfe8cSJoe Wallwork /*@
193cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
194cc2eee55SJoe Wallwork 
195cc2eee55SJoe Wallwork   Input parameters:
196cc2eee55SJoe Wallwork + dm       - The DM
197cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
198cc2eee55SJoe Wallwork 
199cc2eee55SJoe Wallwork   Level: beginner
200cc2eee55SJoe Wallwork 
201cb7bfe8cSJoe Wallwork   Notes:
202cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
203cb7bfe8cSJoe Wallwork 
204db781477SPatrick Sanan .seealso: `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
205cb7bfe8cSJoe Wallwork @*/
206cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
207cc2eee55SJoe Wallwork {
208cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
209cc2eee55SJoe Wallwork 
210cc2eee55SJoe Wallwork   PetscFunctionBegin;
211bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
212cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
213cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
214cc2eee55SJoe Wallwork }
215cc2eee55SJoe Wallwork 
216cb7bfe8cSJoe Wallwork /*@
217cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
218cc2eee55SJoe Wallwork 
219cc2eee55SJoe Wallwork   Input parameters:
220cc2eee55SJoe Wallwork . dm       - The DM
221cc2eee55SJoe Wallwork 
222cc2eee55SJoe Wallwork   Output parameters:
223cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
224cc2eee55SJoe Wallwork 
225cc2eee55SJoe Wallwork   Level: beginner
226cc2eee55SJoe Wallwork 
227cb7bfe8cSJoe Wallwork   Notes:
228cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
229cb7bfe8cSJoe Wallwork 
230db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
231cb7bfe8cSJoe Wallwork @*/
232cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
233cc2eee55SJoe Wallwork {
234cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
235cc2eee55SJoe Wallwork 
236cc2eee55SJoe Wallwork   PetscFunctionBegin;
237bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
238cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
239cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
240cc2eee55SJoe Wallwork }
241cc2eee55SJoe Wallwork 
242cb7bfe8cSJoe Wallwork /*@
243cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
244cc2eee55SJoe Wallwork 
245cc2eee55SJoe Wallwork   Input parameters:
246cc2eee55SJoe Wallwork + dm     - The DM
247cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
248cc2eee55SJoe Wallwork 
249cc2eee55SJoe Wallwork   Level: beginner
250cc2eee55SJoe Wallwork 
251cb7bfe8cSJoe Wallwork   Notes:
252cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
253cb7bfe8cSJoe Wallwork 
254db781477SPatrick Sanan .seealso: `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
255cb7bfe8cSJoe Wallwork @*/
256cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
257cc2eee55SJoe Wallwork {
258cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
259cc2eee55SJoe Wallwork 
260cc2eee55SJoe Wallwork   PetscFunctionBegin;
261bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
262cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
263cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
264cc2eee55SJoe Wallwork }
265cc2eee55SJoe Wallwork 
266cb7bfe8cSJoe Wallwork /*@
267cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
268cc2eee55SJoe Wallwork 
269cc2eee55SJoe Wallwork   Input parameters:
270cc2eee55SJoe Wallwork . dm     - The DM
271cc2eee55SJoe Wallwork 
272cc2eee55SJoe Wallwork   Output parameters:
273cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
274cc2eee55SJoe Wallwork 
275cc2eee55SJoe Wallwork   Level: beginner
276cc2eee55SJoe Wallwork 
277cb7bfe8cSJoe Wallwork   Notes:
278cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
279cb7bfe8cSJoe Wallwork 
280db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
281cb7bfe8cSJoe Wallwork @*/
282cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
283cc2eee55SJoe Wallwork {
284cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
285cc2eee55SJoe Wallwork 
286cc2eee55SJoe Wallwork   PetscFunctionBegin;
287bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
288cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
289cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
290cc2eee55SJoe Wallwork }
291cc2eee55SJoe Wallwork 
292cb7bfe8cSJoe Wallwork /*@
293cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
294cc2eee55SJoe Wallwork 
295cc2eee55SJoe Wallwork   Input parameters:
296cc2eee55SJoe Wallwork + dm     - The DM
297cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
298cc2eee55SJoe Wallwork 
299cc2eee55SJoe Wallwork   Level: beginner
300cc2eee55SJoe Wallwork 
301cb7bfe8cSJoe Wallwork   Notes:
302cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
303cb7bfe8cSJoe Wallwork 
304db781477SPatrick Sanan .seealso: `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()`
305cb7bfe8cSJoe Wallwork @*/
306cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
307cc2eee55SJoe Wallwork {
308cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
309cc2eee55SJoe Wallwork 
310cc2eee55SJoe Wallwork   PetscFunctionBegin;
311bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
312cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
313cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
314cc2eee55SJoe Wallwork }
315cc2eee55SJoe Wallwork 
316cb7bfe8cSJoe Wallwork /*@
317cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
318cc2eee55SJoe Wallwork 
319cc2eee55SJoe Wallwork   Input parameters:
320cc2eee55SJoe Wallwork . dm     - The DM
321cc2eee55SJoe Wallwork 
322cc2eee55SJoe Wallwork   Output parameters:
323cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
324cc2eee55SJoe Wallwork 
325cc2eee55SJoe Wallwork   Level: beginner
326cc2eee55SJoe Wallwork 
327cb7bfe8cSJoe Wallwork   Notes:
328cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
329cb7bfe8cSJoe Wallwork 
330db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()`
331cb7bfe8cSJoe Wallwork @*/
332cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
333cc2eee55SJoe Wallwork {
334cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
335cc2eee55SJoe Wallwork 
336cc2eee55SJoe Wallwork   PetscFunctionBegin;
337bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
338cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
339cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
340cc2eee55SJoe Wallwork }
341cc2eee55SJoe Wallwork 
342cb7bfe8cSJoe Wallwork /*@
34376f360caSJoe Wallwork   DMPlexMetricSetNoSurf - Should surface modification be turned off?
34476f360caSJoe Wallwork 
34576f360caSJoe Wallwork   Input parameters:
34676f360caSJoe Wallwork + dm     - The DM
34776f360caSJoe Wallwork - noSurf - Should surface modification be turned off?
34876f360caSJoe Wallwork 
34976f360caSJoe Wallwork   Level: beginner
35076f360caSJoe Wallwork 
35176f360caSJoe Wallwork   Notes:
35276f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
35376f360caSJoe Wallwork 
354db781477SPatrick Sanan .seealso: `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`
35576f360caSJoe Wallwork @*/
35676f360caSJoe Wallwork PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
35776f360caSJoe Wallwork {
35876f360caSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
35976f360caSJoe Wallwork 
36076f360caSJoe Wallwork   PetscFunctionBegin;
361bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
36276f360caSJoe Wallwork   plex->metricCtx->noSurf = noSurf;
36376f360caSJoe Wallwork   PetscFunctionReturn(0);
36476f360caSJoe Wallwork }
36576f360caSJoe Wallwork 
36676f360caSJoe Wallwork /*@
36776f360caSJoe Wallwork   DMPlexMetricNoSurf - Is surface modification turned off?
36876f360caSJoe Wallwork 
36976f360caSJoe Wallwork   Input parameters:
37076f360caSJoe Wallwork . dm     - The DM
37176f360caSJoe Wallwork 
37276f360caSJoe Wallwork   Output parameters:
37376f360caSJoe Wallwork . noSurf - Is surface modification turned off?
37476f360caSJoe Wallwork 
37576f360caSJoe Wallwork   Level: beginner
37676f360caSJoe Wallwork 
37776f360caSJoe Wallwork   Notes:
37876f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
37976f360caSJoe Wallwork 
380db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`
38176f360caSJoe Wallwork @*/
38276f360caSJoe Wallwork PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
38376f360caSJoe Wallwork {
38476f360caSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
38576f360caSJoe Wallwork 
38676f360caSJoe Wallwork   PetscFunctionBegin;
387bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
38876f360caSJoe Wallwork   *noSurf = plex->metricCtx->noSurf;
38976f360caSJoe Wallwork   PetscFunctionReturn(0);
39076f360caSJoe Wallwork }
39176f360caSJoe Wallwork 
39276f360caSJoe Wallwork /*@
39331b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
39431b70465SJoe Wallwork 
39531b70465SJoe Wallwork   Input parameters:
39631b70465SJoe Wallwork + dm    - The DM
39731b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
39831b70465SJoe Wallwork 
39931b70465SJoe Wallwork   Level: beginner
40031b70465SJoe Wallwork 
401db781477SPatrick Sanan .seealso: `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()`
402cb7bfe8cSJoe Wallwork @*/
40331b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
40431b70465SJoe Wallwork {
40531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
40631b70465SJoe Wallwork 
40731b70465SJoe Wallwork   PetscFunctionBegin;
408bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
409bc00d7afSJoe Wallwork   PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
41031b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
41131b70465SJoe Wallwork   PetscFunctionReturn(0);
41231b70465SJoe Wallwork }
41331b70465SJoe Wallwork 
414cb7bfe8cSJoe Wallwork /*@
41531b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
41631b70465SJoe Wallwork 
41731b70465SJoe Wallwork   Input parameters:
41831b70465SJoe Wallwork . dm    - The DM
41931b70465SJoe Wallwork 
420cc2eee55SJoe Wallwork   Output parameters:
42131b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
42231b70465SJoe Wallwork 
42331b70465SJoe Wallwork   Level: beginner
42431b70465SJoe Wallwork 
425db781477SPatrick Sanan .seealso: `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()`
426cb7bfe8cSJoe Wallwork @*/
42731b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
42831b70465SJoe Wallwork {
42931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
43031b70465SJoe Wallwork 
43131b70465SJoe Wallwork   PetscFunctionBegin;
432bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
43331b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
43431b70465SJoe Wallwork   PetscFunctionReturn(0);
43531b70465SJoe Wallwork }
43631b70465SJoe Wallwork 
437cb7bfe8cSJoe Wallwork /*@
43831b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
43931b70465SJoe Wallwork 
44031b70465SJoe Wallwork   Input parameters:
44131b70465SJoe Wallwork + dm    - The DM
44231b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44331b70465SJoe Wallwork 
44431b70465SJoe Wallwork   Level: beginner
44531b70465SJoe Wallwork 
446db781477SPatrick Sanan .seealso: `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()`
447cb7bfe8cSJoe Wallwork @*/
44831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
44931b70465SJoe Wallwork {
45031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
45131b70465SJoe Wallwork 
45231b70465SJoe Wallwork   PetscFunctionBegin;
453bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
454bc00d7afSJoe Wallwork   PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
45531b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
45631b70465SJoe Wallwork   PetscFunctionReturn(0);
45731b70465SJoe Wallwork }
45831b70465SJoe Wallwork 
459cb7bfe8cSJoe Wallwork /*@
46031b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
46131b70465SJoe Wallwork 
46231b70465SJoe Wallwork   Input parameters:
46331b70465SJoe Wallwork . dm    - The DM
46431b70465SJoe Wallwork 
465cc2eee55SJoe Wallwork   Output parameters:
46631b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
46731b70465SJoe Wallwork 
46831b70465SJoe Wallwork   Level: beginner
46931b70465SJoe Wallwork 
470db781477SPatrick Sanan .seealso: `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()`
471cb7bfe8cSJoe Wallwork @*/
47231b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
47331b70465SJoe Wallwork {
47431b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
47531b70465SJoe Wallwork 
47631b70465SJoe Wallwork   PetscFunctionBegin;
477bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
47831b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
47931b70465SJoe Wallwork   PetscFunctionReturn(0);
48031b70465SJoe Wallwork }
48131b70465SJoe Wallwork 
482cb7bfe8cSJoe Wallwork /*@
48331b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
48431b70465SJoe Wallwork 
48531b70465SJoe Wallwork   Input parameters:
48631b70465SJoe Wallwork + dm    - The DM
48731b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
48831b70465SJoe Wallwork 
48931b70465SJoe Wallwork   Level: beginner
49031b70465SJoe Wallwork 
49131b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
49231b70465SJoe Wallwork 
493db781477SPatrick Sanan .seealso: `DMPlexMetricGetMaximumAnisotropy()`, `DMPlexMetricSetMaximumMagnitude()`
494cb7bfe8cSJoe Wallwork @*/
49531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
49631b70465SJoe Wallwork {
49731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
49831b70465SJoe Wallwork 
49931b70465SJoe Wallwork   PetscFunctionBegin;
500bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
501bc00d7afSJoe Wallwork   PetscCheck(a_max >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Anisotropy must be in [1, inf)");
50231b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
50331b70465SJoe Wallwork   PetscFunctionReturn(0);
50431b70465SJoe Wallwork }
50531b70465SJoe Wallwork 
506cb7bfe8cSJoe Wallwork /*@
50731b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
50831b70465SJoe Wallwork 
50931b70465SJoe Wallwork   Input parameters:
51031b70465SJoe Wallwork . dm    - The DM
51131b70465SJoe Wallwork 
512cc2eee55SJoe Wallwork   Output parameters:
51331b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
51431b70465SJoe Wallwork 
51531b70465SJoe Wallwork   Level: beginner
51631b70465SJoe Wallwork 
517db781477SPatrick Sanan .seealso: `DMPlexMetricSetMaximumAnisotropy()`, `DMPlexMetricGetMaximumMagnitude()`
518cb7bfe8cSJoe Wallwork @*/
51931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
52031b70465SJoe Wallwork {
52131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
52231b70465SJoe Wallwork 
52331b70465SJoe Wallwork   PetscFunctionBegin;
524bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
52531b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
52631b70465SJoe Wallwork   PetscFunctionReturn(0);
52731b70465SJoe Wallwork }
52831b70465SJoe Wallwork 
529cb7bfe8cSJoe Wallwork /*@
53031b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
53131b70465SJoe Wallwork 
53231b70465SJoe Wallwork   Input parameters:
53331b70465SJoe Wallwork + dm               - The DM
53431b70465SJoe Wallwork - targetComplexity - The target metric complexity
53531b70465SJoe Wallwork 
53631b70465SJoe Wallwork   Level: beginner
53731b70465SJoe Wallwork 
538db781477SPatrick Sanan .seealso: `DMPlexMetricGetTargetComplexity()`, `DMPlexMetricSetNormalizationOrder()`
539cb7bfe8cSJoe Wallwork @*/
54031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
54131b70465SJoe Wallwork {
54231b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
54331b70465SJoe Wallwork 
54431b70465SJoe Wallwork   PetscFunctionBegin;
545bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
546bc00d7afSJoe Wallwork   PetscCheck(targetComplexity > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric complexity must be in (0, inf)");
54731b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
54831b70465SJoe Wallwork   PetscFunctionReturn(0);
54931b70465SJoe Wallwork }
55031b70465SJoe Wallwork 
551cb7bfe8cSJoe Wallwork /*@
55231b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
55331b70465SJoe Wallwork 
55431b70465SJoe Wallwork   Input parameters:
55531b70465SJoe Wallwork . dm               - The DM
55631b70465SJoe Wallwork 
557cc2eee55SJoe Wallwork   Output parameters:
55831b70465SJoe Wallwork . targetComplexity - The target metric complexity
55931b70465SJoe Wallwork 
56031b70465SJoe Wallwork   Level: beginner
56131b70465SJoe Wallwork 
562db781477SPatrick Sanan .seealso: `DMPlexMetricSetTargetComplexity()`, `DMPlexMetricGetNormalizationOrder()`
563cb7bfe8cSJoe Wallwork @*/
56431b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
56531b70465SJoe Wallwork {
56631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
56731b70465SJoe Wallwork 
56831b70465SJoe Wallwork   PetscFunctionBegin;
569bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
57031b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
57131b70465SJoe Wallwork   PetscFunctionReturn(0);
57231b70465SJoe Wallwork }
57331b70465SJoe Wallwork 
574cb7bfe8cSJoe Wallwork /*@
57531b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
57631b70465SJoe Wallwork 
57731b70465SJoe Wallwork   Input parameters:
57831b70465SJoe Wallwork + dm - The DM
57931b70465SJoe Wallwork - p  - The normalization order
58031b70465SJoe Wallwork 
58131b70465SJoe Wallwork   Level: beginner
58231b70465SJoe Wallwork 
583db781477SPatrick Sanan .seealso: `DMPlexMetricGetNormalizationOrder()`, `DMPlexMetricSetTargetComplexity()`
584cb7bfe8cSJoe Wallwork @*/
58531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
58631b70465SJoe Wallwork {
58731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
58831b70465SJoe Wallwork 
58931b70465SJoe Wallwork   PetscFunctionBegin;
590bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
591bc00d7afSJoe Wallwork   PetscCheck(p >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Normalization order must be in [1, inf)");
59231b70465SJoe Wallwork   plex->metricCtx->p = p;
59331b70465SJoe Wallwork   PetscFunctionReturn(0);
59431b70465SJoe Wallwork }
59531b70465SJoe Wallwork 
596cb7bfe8cSJoe Wallwork /*@
59731b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
59831b70465SJoe Wallwork 
59931b70465SJoe Wallwork   Input parameters:
60031b70465SJoe Wallwork . dm - The DM
60131b70465SJoe Wallwork 
602cc2eee55SJoe Wallwork   Output parameters:
60331b70465SJoe Wallwork . p - The normalization order
60431b70465SJoe Wallwork 
60531b70465SJoe Wallwork   Level: beginner
60631b70465SJoe Wallwork 
607db781477SPatrick Sanan .seealso: `DMPlexMetricSetNormalizationOrder()`, `DMPlexMetricGetTargetComplexity()`
608cb7bfe8cSJoe Wallwork @*/
60931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
61031b70465SJoe Wallwork {
61131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
61231b70465SJoe Wallwork 
61331b70465SJoe Wallwork   PetscFunctionBegin;
614bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
61531b70465SJoe Wallwork   *p = plex->metricCtx->p;
61631b70465SJoe Wallwork   PetscFunctionReturn(0);
61731b70465SJoe Wallwork }
61831b70465SJoe Wallwork 
619cb7bfe8cSJoe Wallwork /*@
620cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
621cc2eee55SJoe Wallwork 
622cc2eee55SJoe Wallwork   Input parameters:
623cc2eee55SJoe Wallwork + dm   - The DM
624cc2eee55SJoe Wallwork - beta - The metric gradation factor
625cc2eee55SJoe Wallwork 
626cc2eee55SJoe Wallwork   Level: beginner
627cc2eee55SJoe Wallwork 
628cc2eee55SJoe Wallwork   Notes:
629cc2eee55SJoe Wallwork 
630cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
631cc2eee55SJoe Wallwork 
632cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
633cc2eee55SJoe Wallwork 
634cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
635cb7bfe8cSJoe Wallwork 
636db781477SPatrick Sanan .seealso: `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
637cb7bfe8cSJoe Wallwork @*/
638cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
639cc2eee55SJoe Wallwork {
640cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
641cc2eee55SJoe Wallwork 
642cc2eee55SJoe Wallwork   PetscFunctionBegin;
643bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
644bc00d7afSJoe Wallwork   PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)");
645cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
646cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
647cc2eee55SJoe Wallwork }
648cc2eee55SJoe Wallwork 
649cb7bfe8cSJoe Wallwork /*@
650cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
651cc2eee55SJoe Wallwork 
652cc2eee55SJoe Wallwork   Input parameters:
653cc2eee55SJoe Wallwork . dm   - The DM
654cc2eee55SJoe Wallwork 
655cc2eee55SJoe Wallwork   Output parameters:
656cc2eee55SJoe Wallwork . beta - The metric gradation factor
657cc2eee55SJoe Wallwork 
658cc2eee55SJoe Wallwork   Level: beginner
659cc2eee55SJoe Wallwork 
660cb7bfe8cSJoe Wallwork   Notes:
661cb7bfe8cSJoe Wallwork 
662cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
663cb7bfe8cSJoe Wallwork 
664cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
665cb7bfe8cSJoe Wallwork 
666cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
667cc2eee55SJoe Wallwork 
668db781477SPatrick Sanan .seealso: `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
669cb7bfe8cSJoe Wallwork @*/
670cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
671cc2eee55SJoe Wallwork {
672cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
673cc2eee55SJoe Wallwork 
674cc2eee55SJoe Wallwork   PetscFunctionBegin;
675bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
676cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
677cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
678cc2eee55SJoe Wallwork }
679cc2eee55SJoe Wallwork 
680cb7bfe8cSJoe Wallwork /*@
681ae8b063eSJoe Wallwork   DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
682ae8b063eSJoe Wallwork 
683ae8b063eSJoe Wallwork   Input parameters:
684ae8b063eSJoe Wallwork + dm    - The DM
685ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number
686ae8b063eSJoe Wallwork 
687ae8b063eSJoe Wallwork   Level: beginner
688ae8b063eSJoe Wallwork 
689ae8b063eSJoe Wallwork   Notes:
690ae8b063eSJoe Wallwork 
691ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
692ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
693ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
694ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
695ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
696ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
697ae8b063eSJoe Wallwork 
698ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
699ae8b063eSJoe Wallwork 
700db781477SPatrick Sanan .seealso: `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
701ae8b063eSJoe Wallwork @*/
702ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
703ae8b063eSJoe Wallwork {
704ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
705ae8b063eSJoe Wallwork 
706ae8b063eSJoe Wallwork   PetscFunctionBegin;
707bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
708bc00d7afSJoe Wallwork   PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)");
709ae8b063eSJoe Wallwork   plex->metricCtx->hausdorffNumber = hausd;
710ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
711ae8b063eSJoe Wallwork }
712ae8b063eSJoe Wallwork 
713ae8b063eSJoe Wallwork /*@
714ae8b063eSJoe Wallwork   DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
715ae8b063eSJoe Wallwork 
716ae8b063eSJoe Wallwork   Input parameters:
717ae8b063eSJoe Wallwork . dm    - The DM
718ae8b063eSJoe Wallwork 
719ae8b063eSJoe Wallwork   Output parameters:
720ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number
721ae8b063eSJoe Wallwork 
722ae8b063eSJoe Wallwork   Level: beginner
723ae8b063eSJoe Wallwork 
724ae8b063eSJoe Wallwork   Notes:
725ae8b063eSJoe Wallwork 
726ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
727ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
728ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
729ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
730ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
731ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
732ae8b063eSJoe Wallwork 
733ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
734ae8b063eSJoe Wallwork 
735db781477SPatrick Sanan .seealso: `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
736ae8b063eSJoe Wallwork @*/
737ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
738ae8b063eSJoe Wallwork {
739ae8b063eSJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
740ae8b063eSJoe Wallwork 
741ae8b063eSJoe Wallwork   PetscFunctionBegin;
742bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
743ae8b063eSJoe Wallwork   *hausd = plex->metricCtx->hausdorffNumber;
744ae8b063eSJoe Wallwork   PetscFunctionReturn(0);
745ae8b063eSJoe Wallwork }
746ae8b063eSJoe Wallwork 
747ae8b063eSJoe Wallwork /*@
748cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
749cc2eee55SJoe Wallwork 
750cc2eee55SJoe Wallwork   Input parameters:
751cc2eee55SJoe Wallwork + dm        - The DM
752cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
753cc2eee55SJoe Wallwork 
754cb7bfe8cSJoe Wallwork   Level: beginner
755cb7bfe8cSJoe Wallwork 
756cb7bfe8cSJoe Wallwork   Notes:
757cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
758cb7bfe8cSJoe Wallwork 
759db781477SPatrick Sanan .seealso: `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()`
760cb7bfe8cSJoe Wallwork @*/
761cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
762cc2eee55SJoe Wallwork {
763cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
764cc2eee55SJoe Wallwork 
765cc2eee55SJoe Wallwork   PetscFunctionBegin;
766bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
767cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
768cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
769cc2eee55SJoe Wallwork }
770cc2eee55SJoe Wallwork 
771cb7bfe8cSJoe Wallwork /*@
772cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
773cc2eee55SJoe Wallwork 
774cc2eee55SJoe Wallwork   Input parameters:
775cc2eee55SJoe Wallwork . dm        - The DM
776cc2eee55SJoe Wallwork 
777cc2eee55SJoe Wallwork   Output parameters:
778cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
779cc2eee55SJoe Wallwork 
780cb7bfe8cSJoe Wallwork   Level: beginner
781cb7bfe8cSJoe Wallwork 
782cb7bfe8cSJoe Wallwork   Notes:
783cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
784cb7bfe8cSJoe Wallwork 
785db781477SPatrick Sanan .seealso: `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
786cb7bfe8cSJoe Wallwork @*/
787cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
788cc2eee55SJoe Wallwork {
789cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
790cc2eee55SJoe Wallwork 
791cc2eee55SJoe Wallwork   PetscFunctionBegin;
792bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
793cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
794cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
795cc2eee55SJoe Wallwork }
796cc2eee55SJoe Wallwork 
797cb7bfe8cSJoe Wallwork /*@
798cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
799cc2eee55SJoe Wallwork 
800cc2eee55SJoe Wallwork   Input parameters:
801cc2eee55SJoe Wallwork + dm      - The DM
802cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
803cc2eee55SJoe Wallwork 
804cb7bfe8cSJoe Wallwork   Level: beginner
805cb7bfe8cSJoe Wallwork 
806cb7bfe8cSJoe Wallwork   Notes:
807cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
808cc2eee55SJoe Wallwork 
809db781477SPatrick Sanan .seealso: `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
810cb7bfe8cSJoe Wallwork @*/
811cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
812cc2eee55SJoe Wallwork {
813cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
814cc2eee55SJoe Wallwork 
815cc2eee55SJoe Wallwork   PetscFunctionBegin;
816bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
817cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
818cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
819cc2eee55SJoe Wallwork }
820cc2eee55SJoe Wallwork 
821cb7bfe8cSJoe Wallwork /*@
822cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
823cc2eee55SJoe Wallwork 
824cc2eee55SJoe Wallwork   Input parameters:
825cc2eee55SJoe Wallwork . dm      - The DM
826cc2eee55SJoe Wallwork 
827cc2eee55SJoe Wallwork   Output parameters:
828cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
829cc2eee55SJoe Wallwork 
830cb7bfe8cSJoe Wallwork   Level: beginner
831cb7bfe8cSJoe Wallwork 
832cb7bfe8cSJoe Wallwork   Notes:
833cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
834cc2eee55SJoe Wallwork 
835db781477SPatrick Sanan .seealso: `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()`
836cb7bfe8cSJoe Wallwork @*/
837cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
838cc2eee55SJoe Wallwork {
839cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
840cc2eee55SJoe Wallwork 
841cc2eee55SJoe Wallwork   PetscFunctionBegin;
842bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
843cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
844cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
845cc2eee55SJoe Wallwork }
846cc2eee55SJoe Wallwork 
84720282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
84820282da2SJoe Wallwork {
84920282da2SJoe Wallwork   MPI_Comm       comm;
85020282da2SJoe Wallwork   PetscFE        fe;
85120282da2SJoe Wallwork   PetscInt       dim;
85220282da2SJoe Wallwork 
85320282da2SJoe Wallwork   PetscFunctionBegin;
85420282da2SJoe Wallwork 
85520282da2SJoe Wallwork   /* Extract metadata from dm */
8569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
8579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
85820282da2SJoe Wallwork 
85920282da2SJoe Wallwork   /* Create a P1 field of the requested size */
8609566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
8619566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
8629566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
8639566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
8649566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm, metric));
86520282da2SJoe Wallwork 
86620282da2SJoe Wallwork   PetscFunctionReturn(0);
86720282da2SJoe Wallwork }
86820282da2SJoe Wallwork 
869cb7bfe8cSJoe Wallwork /*@
87020282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
87120282da2SJoe Wallwork 
87220282da2SJoe Wallwork   Input parameters:
87320282da2SJoe Wallwork + dm     - The DM
87420282da2SJoe Wallwork - f      - The field number to use
87520282da2SJoe Wallwork 
87620282da2SJoe Wallwork   Output parameter:
87720282da2SJoe Wallwork . metric - The metric
87820282da2SJoe Wallwork 
87920282da2SJoe Wallwork   Level: beginner
88020282da2SJoe Wallwork 
88131b70465SJoe Wallwork   Notes:
88231b70465SJoe Wallwork 
88331b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
88431b70465SJoe Wallwork 
88531b70465SJoe Wallwork   Command line options for Riemannian metrics:
88631b70465SJoe Wallwork 
887cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
88893520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
889cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
890cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
891cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
892cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
893cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
89467b8a455SSatish Balay - -dm_plex_metric_target_complexity         - Target metric complexity
895cb7bfe8cSJoe Wallwork 
896cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
897cb7bfe8cSJoe Wallwork 
89867b8a455SSatish Balay . -dm_adaptor <pragmatic/mmg/parmmg>        - specify dm adaptor to use
899cb7bfe8cSJoe Wallwork 
900cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
901cb7bfe8cSJoe Wallwork 
902cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
903cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
904cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
905cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
906cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
907cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
90820282da2SJoe Wallwork 
909db781477SPatrick Sanan .seealso: `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`
910cb7bfe8cSJoe Wallwork @*/
91120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
91220282da2SJoe Wallwork {
91393520041SJoe Wallwork   PetscBool isotropic, uniform;
91420282da2SJoe Wallwork   PetscInt  coordDim, Nd;
91520282da2SJoe Wallwork 
91620282da2SJoe Wallwork   PetscFunctionBegin;
9179566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
91820282da2SJoe Wallwork   Nd = coordDim*coordDim;
9199566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
9209566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
92193520041SJoe Wallwork   if (uniform) {
92293520041SJoe Wallwork     MPI_Comm comm;
92393520041SJoe Wallwork 
9249566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
925bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
9269566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, metric));
9279566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
9289566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(*metric));
929f6db075bSJoe Wallwork   } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
930f6db075bSJoe Wallwork   else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
93120282da2SJoe Wallwork   PetscFunctionReturn(0);
93220282da2SJoe Wallwork }
93320282da2SJoe Wallwork 
934cb7bfe8cSJoe Wallwork /*@
93520282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
93620282da2SJoe Wallwork 
93720282da2SJoe Wallwork   Input parameters:
93820282da2SJoe Wallwork + dm     - The DM
93920282da2SJoe Wallwork . f      - The field number to use
94020282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
94120282da2SJoe Wallwork 
94220282da2SJoe Wallwork   Output parameter:
94320282da2SJoe Wallwork . metric - The uniform metric
94420282da2SJoe Wallwork 
94520282da2SJoe Wallwork   Level: beginner
94620282da2SJoe Wallwork 
94793520041SJoe Wallwork   Note: In this case, the metric is constant in space and so is only specified for a single vertex.
94820282da2SJoe Wallwork 
949db781477SPatrick Sanan .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
950cb7bfe8cSJoe Wallwork @*/
95120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
95220282da2SJoe Wallwork {
95320282da2SJoe Wallwork   PetscFunctionBegin;
9549566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
9559566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9562c71b3e2SJacob Faibussowitsch   PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
957bc00d7afSJoe Wallwork   PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
9589566063dSJacob Faibussowitsch   PetscCall(VecSet(*metric, alpha));
9599566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*metric));
9609566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*metric));
96120282da2SJoe Wallwork   PetscFunctionReturn(0);
96220282da2SJoe Wallwork }
96320282da2SJoe Wallwork 
96493520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
96593520041SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96693520041SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96793520041SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
96893520041SJoe Wallwork {
96993520041SJoe Wallwork   f0[0] = u[0];
97093520041SJoe Wallwork }
97193520041SJoe Wallwork 
972cb7bfe8cSJoe Wallwork /*@
97320282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
97420282da2SJoe Wallwork 
97520282da2SJoe Wallwork   Input parameters:
97620282da2SJoe Wallwork + dm        - The DM
97720282da2SJoe Wallwork . f         - The field number to use
97820282da2SJoe Wallwork - indicator - The error indicator
97920282da2SJoe Wallwork 
98020282da2SJoe Wallwork   Output parameter:
98120282da2SJoe Wallwork . metric    - The isotropic metric
98220282da2SJoe Wallwork 
98320282da2SJoe Wallwork   Level: beginner
98420282da2SJoe Wallwork 
98520282da2SJoe Wallwork   Notes:
98620282da2SJoe Wallwork 
98720282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
98820282da2SJoe Wallwork 
98993520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
99020282da2SJoe Wallwork 
991db781477SPatrick Sanan .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
992cb7bfe8cSJoe Wallwork @*/
99320282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
99420282da2SJoe Wallwork {
99593520041SJoe Wallwork   PetscInt       m, n;
99620282da2SJoe Wallwork 
99720282da2SJoe Wallwork   PetscFunctionBegin;
9989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
9999566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
10009566063dSJacob Faibussowitsch   PetscCall(VecGetSize(indicator, &m));
10019566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metric, &n));
10029566063dSJacob Faibussowitsch   if (m == n) PetscCall(VecCopy(indicator, *metric));
100393520041SJoe Wallwork   else {
100493520041SJoe Wallwork     DM     dmIndi;
100593520041SJoe Wallwork     void (*funcs[1])(PetscInt, PetscInt, PetscInt,
100693520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
100793520041SJoe Wallwork                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
100893520041SJoe Wallwork                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
100993520041SJoe Wallwork 
10109566063dSJacob Faibussowitsch     PetscCall(VecGetDM(indicator, &dmIndi));
101193520041SJoe Wallwork     funcs[0] = identity;
10129566063dSJacob Faibussowitsch     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
101320282da2SJoe Wallwork   }
101420282da2SJoe Wallwork   PetscFunctionReturn(0);
101520282da2SJoe Wallwork }
101620282da2SJoe Wallwork 
1017f6db075bSJoe Wallwork /*@
1018f6db075bSJoe Wallwork   DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1019f6db075bSJoe Wallwork 
1020f6db075bSJoe Wallwork   Input parameters:
1021f6db075bSJoe Wallwork + dm          - The DM of the metric field
1022f6db075bSJoe Wallwork - f           - The field number to use
1023f6db075bSJoe Wallwork 
1024f6db075bSJoe Wallwork   Output parameter:
1025f6db075bSJoe Wallwork + determinant - The determinant field
1026f6db075bSJoe Wallwork - dmDet       - The corresponding DM
1027f6db075bSJoe Wallwork 
1028f6db075bSJoe Wallwork   Level: beginner
1029f6db075bSJoe Wallwork 
1030f6db075bSJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic(), DMPlexMetricCreate()
1031f6db075bSJoe Wallwork @*/
1032f6db075bSJoe Wallwork PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1033f6db075bSJoe Wallwork {
1034f6db075bSJoe Wallwork   PetscBool uniform;
1035f6db075bSJoe Wallwork 
1036f6db075bSJoe Wallwork   PetscFunctionBegin;
1037f6db075bSJoe Wallwork   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1038f6db075bSJoe Wallwork   PetscCall(DMClone(dm, dmDet));
1039f6db075bSJoe Wallwork   if (uniform) {
1040f6db075bSJoe Wallwork     MPI_Comm comm;
1041f6db075bSJoe Wallwork 
1042f6db075bSJoe Wallwork     PetscCall(PetscObjectGetComm((PetscObject) *dmDet, &comm));
1043f6db075bSJoe Wallwork     PetscCall(VecCreate(comm, determinant));
1044f6db075bSJoe Wallwork     PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1045f6db075bSJoe Wallwork     PetscCall(VecSetFromOptions(*determinant));
1046f6db075bSJoe Wallwork   } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
1047f6db075bSJoe Wallwork   PetscFunctionReturn(0);
1048f6db075bSJoe Wallwork }
1049f6db075bSJoe Wallwork 
105082490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
105182490253SJoe Wallwork {
105282490253SJoe Wallwork   PetscInt i, j;
105382490253SJoe Wallwork 
105482490253SJoe Wallwork   PetscFunctionBegin;
105582490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
105682490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
105782490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
105882490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
105982490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
106063a3b9bcSJacob Faibussowitsch       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i*dim+j]));
106163a3b9bcSJacob Faibussowitsch       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i*dim+j]));
106282490253SJoe Wallwork     }
106382490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
106482490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
106582490253SJoe Wallwork   }
106682490253SJoe Wallwork   PetscFunctionReturn(0);
106782490253SJoe Wallwork }
106882490253SJoe Wallwork 
10693f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
107020282da2SJoe Wallwork {
107120282da2SJoe Wallwork   PetscInt       i, j, k;
107220282da2SJoe 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);
107320282da2SJoe Wallwork   PetscScalar   *Mpos;
107420282da2SJoe Wallwork 
107520282da2SJoe Wallwork   PetscFunctionBegin;
10769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim*dim, &Mpos, dim, &eigs));
107720282da2SJoe Wallwork 
107820282da2SJoe Wallwork   /* Symmetrize */
107920282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
108020282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
108120282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
108220282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
108320282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
108420282da2SJoe Wallwork     }
108520282da2SJoe Wallwork   }
108620282da2SJoe Wallwork 
108720282da2SJoe Wallwork   /* Compute eigendecomposition */
108893520041SJoe Wallwork   if (dim == 1) {
108993520041SJoe Wallwork 
109093520041SJoe Wallwork     /* Isotropic case */
109193520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
109293520041SJoe Wallwork     Mpos[0] = 1.0;
109393520041SJoe Wallwork   } else {
109493520041SJoe Wallwork 
109593520041SJoe Wallwork     /* Anisotropic case */
109620282da2SJoe Wallwork     PetscScalar  *work;
109720282da2SJoe Wallwork     PetscBLASInt lwork;
109820282da2SJoe Wallwork 
109920282da2SJoe Wallwork     lwork = 5*dim;
11009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*dim, &work));
110120282da2SJoe Wallwork     {
110220282da2SJoe Wallwork       PetscBLASInt lierr;
110320282da2SJoe Wallwork       PetscBLASInt nb;
110420282da2SJoe Wallwork 
11059566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
11069566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
110720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
110820282da2SJoe Wallwork       {
110920282da2SJoe Wallwork         PetscReal *rwork;
11109566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
1111*792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
11129566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
111320282da2SJoe Wallwork       }
111420282da2SJoe Wallwork #else
1115*792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
111620282da2SJoe Wallwork #endif
111782490253SJoe Wallwork       if (lierr) {
111882490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
111982490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
112082490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
112182490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
112282490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
112382490253SJoe Wallwork           }
112482490253SJoe Wallwork         }
112582490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
112698921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
112782490253SJoe Wallwork       }
11289566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
112920282da2SJoe Wallwork     }
11309566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
113120282da2SJoe Wallwork   }
113220282da2SJoe Wallwork 
113320282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
113420282da2SJoe Wallwork   max_eig = 0.0;
113520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
113620282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
113720282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
113820282da2SJoe Wallwork   }
113920282da2SJoe Wallwork 
11403f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11413f07679eSJoe Wallwork   *detMp = 1.0;
114220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
114320282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
11443f07679eSJoe Wallwork     *detMp *= eigs[i];
114520282da2SJoe Wallwork   }
114620282da2SJoe Wallwork 
114720282da2SJoe Wallwork   /* Reconstruct Hessian */
114820282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
114920282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
115020282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
115120282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
115220282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
115320282da2SJoe Wallwork       }
115420282da2SJoe Wallwork     }
115520282da2SJoe Wallwork   }
11569566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Mpos, eigs));
115720282da2SJoe Wallwork 
115820282da2SJoe Wallwork   PetscFunctionReturn(0);
115920282da2SJoe Wallwork }
116020282da2SJoe Wallwork 
1161cb7bfe8cSJoe Wallwork /*@
116220282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
116320282da2SJoe Wallwork 
116420282da2SJoe Wallwork   Input parameters:
116520282da2SJoe Wallwork + dm                 - The DM
11663f07679eSJoe Wallwork . metricIn           - The metric
116799abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11683f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
116920282da2SJoe Wallwork 
117020282da2SJoe Wallwork   Output parameter:
11713f07679eSJoe Wallwork + metricOut          - The metric
11723f07679eSJoe Wallwork - determinant        - Its determinant
117320282da2SJoe Wallwork 
117420282da2SJoe Wallwork   Level: beginner
117520282da2SJoe Wallwork 
1176cb7bfe8cSJoe Wallwork   Notes:
1177cb7bfe8cSJoe Wallwork 
1178cb7bfe8cSJoe Wallwork   Relevant command line options:
1179cb7bfe8cSJoe Wallwork 
118093520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
118193520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
118293520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1183cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1184cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1185cb7bfe8cSJoe Wallwork 
1186db781477SPatrick Sanan .seealso: `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1187cb7bfe8cSJoe Wallwork @*/
11885241a91fSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
118920282da2SJoe Wallwork {
11903f07679eSJoe Wallwork   DM             dmDet;
119193520041SJoe Wallwork   PetscBool      isotropic, uniform;
119220282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
11933f07679eSJoe Wallwork   PetscScalar   *met, *det;
119420282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
119520282da2SJoe Wallwork 
119620282da2SJoe Wallwork   PetscFunctionBegin;
11979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0));
119820282da2SJoe Wallwork 
119920282da2SJoe Wallwork   /* Extract metadata from dm */
12009566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
120120282da2SJoe Wallwork   if (restrictSizes) {
12029566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
12039566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
120499abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
120599abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1206bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
120799abec2bSJoe Wallwork   }
120899abec2bSJoe Wallwork   if (restrictAnisotropy) {
12099566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
121099abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
121120282da2SJoe Wallwork   }
121220282da2SJoe Wallwork 
121393520041SJoe Wallwork   /* Setup output metric */
12145241a91fSJoe Wallwork   PetscCall(VecCopy(metricIn, metricOut));
121593520041SJoe Wallwork 
121693520041SJoe Wallwork   /* Enforce SPD and extract determinant */
12175241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
12189566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12199566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
122093520041SJoe Wallwork   if (uniform) {
1221d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
122293520041SJoe Wallwork 
122393520041SJoe Wallwork     /* Uniform case */
12245241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
12259566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
12265241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
122793520041SJoe Wallwork   } else {
122893520041SJoe Wallwork 
122993520041SJoe Wallwork     /* Spatially varying case */
123093520041SJoe Wallwork     PetscInt nrow;
123193520041SJoe Wallwork 
123293520041SJoe Wallwork     if (isotropic) nrow = 1;
123393520041SJoe Wallwork     else nrow = dim;
12345241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
12359566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12365241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
123720282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12383f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
12399566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
12409566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
12419566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
124220282da2SJoe Wallwork     }
12435241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
124493520041SJoe Wallwork   }
12455241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
1246fe902aceSJoe Wallwork 
12479566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0));
124820282da2SJoe Wallwork   PetscFunctionReturn(0);
124920282da2SJoe Wallwork }
125020282da2SJoe Wallwork 
125120282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
125220282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
125320282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
125420282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
125520282da2SJoe Wallwork {
125620282da2SJoe Wallwork   const PetscScalar p = constants[0];
125720282da2SJoe Wallwork 
1258985ad86eSJose E. Roman   f0[0] = PetscPowScalar(u[0], p/(2.0*p + dim));
125920282da2SJoe Wallwork }
126020282da2SJoe Wallwork 
1261cb7bfe8cSJoe Wallwork /*@
126220282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
126320282da2SJoe Wallwork 
126420282da2SJoe Wallwork   Input parameters:
126520282da2SJoe Wallwork + dm                 - The DM
126620282da2SJoe Wallwork . metricIn           - The unnormalized metric
126716522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
126816522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
126920282da2SJoe Wallwork 
127020282da2SJoe Wallwork   Output parameter:
127120282da2SJoe Wallwork . metricOut          - The normalized metric
127220282da2SJoe Wallwork 
127320282da2SJoe Wallwork   Level: beginner
127420282da2SJoe Wallwork 
1275cb7bfe8cSJoe Wallwork   Notes:
1276cb7bfe8cSJoe Wallwork 
1277cb7bfe8cSJoe Wallwork   Relevant command line options:
1278cb7bfe8cSJoe Wallwork 
127993520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
128093520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
128193520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1282cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1283cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1284cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1285cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1286cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1287cb7bfe8cSJoe Wallwork 
1288db781477SPatrick Sanan .seealso: `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1289cb7bfe8cSJoe Wallwork @*/
12905241a91fSJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
129120282da2SJoe Wallwork {
12923f07679eSJoe Wallwork   DM               dmDet;
129320282da2SJoe Wallwork   MPI_Comm         comm;
129493520041SJoe Wallwork   PetscBool        restrictAnisotropyFirst, isotropic, uniform;
129520282da2SJoe Wallwork   PetscDS          ds;
129620282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
12973f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
129893520041SJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
129920282da2SJoe Wallwork 
130020282da2SJoe Wallwork   PetscFunctionBegin;
13019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0));
130220282da2SJoe Wallwork 
130320282da2SJoe Wallwork   /* Extract metadata from dm */
13049566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
13059566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13069566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
13079566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
130893520041SJoe Wallwork   if (isotropic) Nd = 1;
130993520041SJoe Wallwork   else Nd = dim*dim;
131020282da2SJoe Wallwork 
131120282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
13129566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
13135241a91fSJoe Wallwork   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
131420282da2SJoe Wallwork 
131520282da2SJoe Wallwork   /* Compute global normalization factor */
13169566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
13179566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
131816522936SJoe Wallwork   constants[0] = p;
131993520041SJoe Wallwork   if (uniform) {
1320bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
132193520041SJoe Wallwork     DM  dmTmp;
132293520041SJoe Wallwork     Vec tmp;
132393520041SJoe Wallwork 
13249566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmTmp));
13259566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
13269566063dSJacob Faibussowitsch     PetscCall(VecGetArray(determinant, &det));
13279566063dSJacob Faibussowitsch     PetscCall(VecSet(tmp, det[0]));
13289566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(determinant, &det));
13299566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmTmp, &ds));
13309566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13319566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13329566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
13339566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tmp));
13349566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmTmp));
133593520041SJoe Wallwork   } else {
13369566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
13379566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmDet, &ds));
13389566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13399566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13409566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
134193520041SJoe Wallwork   }
13423f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
1343bc00d7afSJoe Wallwork   PetscCheck(realIntegral > 1.0e-30, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Global metric normalization factor must be in (0, inf). Is the input metric positive-definite?");
13443f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
134520282da2SJoe Wallwork 
134620282da2SJoe Wallwork   /* Apply local scaling */
134720282da2SJoe Wallwork   if (restrictSizes) {
13489566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
13499566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
135016522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
135116522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1352bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
135316522936SJoe Wallwork   }
135416522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
13559566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
135616522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
135720282da2SJoe Wallwork   }
13585241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
13599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(determinant, &det));
136093520041SJoe Wallwork   if (uniform) {
136193520041SJoe Wallwork 
136293520041SJoe Wallwork     /* Uniform case */
136393520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim));
13649566063dSJacob Faibussowitsch     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
136593520041SJoe Wallwork   } else {
136693520041SJoe Wallwork 
136793520041SJoe Wallwork     /* Spatially varying case */
136893520041SJoe Wallwork     PetscInt nrow;
136993520041SJoe Wallwork 
137093520041SJoe Wallwork     if (isotropic) nrow = 1;
137193520041SJoe Wallwork     else nrow = dim;
13729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
13735241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
137420282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13753f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
137620282da2SJoe Wallwork 
13779566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
13789566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
13793f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
138020282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
13819566063dSJacob Faibussowitsch       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
138293520041SJoe Wallwork     }
138320282da2SJoe Wallwork   }
13849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(determinant, &det));
13855241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
138620282da2SJoe Wallwork 
13879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0));
138820282da2SJoe Wallwork   PetscFunctionReturn(0);
138920282da2SJoe Wallwork }
139020282da2SJoe Wallwork 
1391cb7bfe8cSJoe Wallwork /*@
139220282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
139320282da2SJoe Wallwork 
1394f1a722f8SMatthew G. Knepley   Input Parameters:
139520282da2SJoe Wallwork + dm         - The DM
139620282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
139720282da2SJoe Wallwork . weights    - Weights for the average
139820282da2SJoe Wallwork - metrics    - The metrics to be averaged
139920282da2SJoe Wallwork 
140020282da2SJoe Wallwork   Output Parameter:
140120282da2SJoe Wallwork . metricAvg  - The averaged metric
140220282da2SJoe Wallwork 
140320282da2SJoe Wallwork   Level: beginner
140420282da2SJoe Wallwork 
140520282da2SJoe Wallwork   Notes:
140620282da2SJoe Wallwork   The weights should sum to unity.
140720282da2SJoe Wallwork 
140820282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
140920282da2SJoe Wallwork 
1410db781477SPatrick Sanan .seealso: `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1411cb7bfe8cSJoe Wallwork @*/
14125241a91fSJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
141320282da2SJoe Wallwork {
141420282da2SJoe Wallwork   PetscBool haveWeights = PETSC_TRUE;
141593520041SJoe Wallwork   PetscInt  i, m, n;
141620282da2SJoe Wallwork   PetscReal sum = 0.0, tol = 1.0e-10;
141720282da2SJoe Wallwork 
141820282da2SJoe Wallwork   PetscFunctionBegin;
14199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0));
142063a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
14215241a91fSJoe Wallwork   PetscCall(VecSet(metricAvg, 0.0));
14225241a91fSJoe Wallwork   PetscCall(VecGetSize(metricAvg, &m));
142393520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
14249566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
14255f80ce2aSJacob Faibussowitsch     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
142693520041SJoe Wallwork   }
142720282da2SJoe Wallwork 
142820282da2SJoe Wallwork   /* Default to the unweighted case */
142920282da2SJoe Wallwork   if (!weights) {
14309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numMetrics, &weights));
143120282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
143220282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
143320282da2SJoe Wallwork   }
143420282da2SJoe Wallwork 
143520282da2SJoe Wallwork   /* Check weights sum to unity */
143693520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
14375f80ce2aSJacob Faibussowitsch   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
143820282da2SJoe Wallwork 
143920282da2SJoe Wallwork   /* Compute metric average */
14405241a91fSJoe Wallwork   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
14419566063dSJacob Faibussowitsch   if (!haveWeights) PetscCall(PetscFree(weights));
1442fe902aceSJoe Wallwork 
14439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0));
144420282da2SJoe Wallwork   PetscFunctionReturn(0);
144520282da2SJoe Wallwork }
144620282da2SJoe Wallwork 
1447cb7bfe8cSJoe Wallwork /*@
144820282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
144920282da2SJoe Wallwork 
1450f1a722f8SMatthew G. Knepley   Input Parameters:
145120282da2SJoe Wallwork + dm         - The DM
145220282da2SJoe Wallwork . metric1    - The first metric to be averaged
145320282da2SJoe Wallwork - metric2    - The second metric to be averaged
145420282da2SJoe Wallwork 
145520282da2SJoe Wallwork   Output Parameter:
145620282da2SJoe Wallwork . metricAvg  - The averaged metric
145720282da2SJoe Wallwork 
145820282da2SJoe Wallwork   Level: beginner
145920282da2SJoe Wallwork 
1460db781477SPatrick Sanan .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1461cb7bfe8cSJoe Wallwork @*/
14625241a91fSJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
146320282da2SJoe Wallwork {
146420282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
146520282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
146620282da2SJoe Wallwork 
146720282da2SJoe Wallwork   PetscFunctionBegin;
14689566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
146920282da2SJoe Wallwork   PetscFunctionReturn(0);
147020282da2SJoe Wallwork }
147120282da2SJoe Wallwork 
1472cb7bfe8cSJoe Wallwork /*@
147320282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
147420282da2SJoe Wallwork 
1475f1a722f8SMatthew G. Knepley   Input Parameters:
147620282da2SJoe Wallwork + dm         - The DM
147720282da2SJoe Wallwork . metric1    - The first metric to be averaged
147820282da2SJoe Wallwork . metric2    - The second metric to be averaged
147920282da2SJoe Wallwork - metric3    - The third metric to be averaged
148020282da2SJoe Wallwork 
148120282da2SJoe Wallwork   Output Parameter:
148220282da2SJoe Wallwork . metricAvg  - The averaged metric
148320282da2SJoe Wallwork 
148420282da2SJoe Wallwork   Level: beginner
148520282da2SJoe Wallwork 
1486db781477SPatrick Sanan .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1487cb7bfe8cSJoe Wallwork @*/
14885241a91fSJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
148920282da2SJoe Wallwork {
149020282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
149120282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
149220282da2SJoe Wallwork 
149320282da2SJoe Wallwork   PetscFunctionBegin;
14949566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
149520282da2SJoe Wallwork   PetscFunctionReturn(0);
149620282da2SJoe Wallwork }
149720282da2SJoe Wallwork 
149820282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
149920282da2SJoe Wallwork {
150020282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
1501e2606525SJoe Wallwork   PetscReal     *evals;
150220282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
150320282da2SJoe Wallwork 
150420282da2SJoe Wallwork   PetscFunctionBegin;
150593520041SJoe Wallwork 
150693520041SJoe Wallwork   /* Isotropic case */
150793520041SJoe Wallwork   if (dim == 1) {
15086f71502aSJoe Wallwork     M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
150993520041SJoe Wallwork     PetscFunctionReturn(0);
151093520041SJoe Wallwork   }
151193520041SJoe Wallwork 
151293520041SJoe Wallwork   /* Anisotropic case */
1513e2606525SJoe Wallwork   PetscCall(PetscMalloc4(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals));
151420282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
151520282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
151620282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
151720282da2SJoe Wallwork     }
151820282da2SJoe Wallwork   }
151920282da2SJoe Wallwork   {
152020282da2SJoe Wallwork     PetscScalar *work;
152120282da2SJoe Wallwork     PetscBLASInt lwork;
152220282da2SJoe Wallwork 
152320282da2SJoe Wallwork     lwork = 5*dim;
15249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*dim, &work));
152520282da2SJoe Wallwork     {
152620282da2SJoe Wallwork       PetscBLASInt lierr, nb;
1527e2606525SJoe Wallwork       PetscReal    sqrtj;
152820282da2SJoe Wallwork 
152920282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
15309566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
15319566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
153220282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
153320282da2SJoe Wallwork       {
153420282da2SJoe Wallwork         PetscReal *rwork;
15359566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
1536*792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15379566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
153820282da2SJoe Wallwork       }
153920282da2SJoe Wallwork #else
1540*792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
154120282da2SJoe Wallwork #endif
154282490253SJoe Wallwork       if (lierr) {
154382490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
154498921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
154582490253SJoe Wallwork       }
15469566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
154720282da2SJoe Wallwork 
1548e2606525SJoe Wallwork       /* Compute square root and the reciprocal thereof */
154920282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
155020282da2SJoe Wallwork         for (k = 0; k < dim; ++k) {
1551e2606525SJoe Wallwork           sqrtM1[i*dim+k] = 0.0;
1552e2606525SJoe Wallwork           isqrtM1[i*dim+k] = 0.0;
1553e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1554e2606525SJoe Wallwork             sqrtj = PetscSqrtReal(evals[j]);
1555e2606525SJoe Wallwork             sqrtM1[i*dim+k] += evecs[j*dim+i] * sqrtj * evecs[j*dim+k];
1556e2606525SJoe Wallwork             isqrtM1[i*dim+k] += evecs[j*dim+i] * (1.0/sqrtj) * evecs[j*dim+k];
155720282da2SJoe Wallwork           }
155820282da2SJoe Wallwork         }
155920282da2SJoe Wallwork       }
156020282da2SJoe Wallwork 
1561e2606525SJoe Wallwork       /* Map M2 into the space spanned by the eigenvectors of M1 */
156220282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
156320282da2SJoe Wallwork         for (l = 0; l < dim; ++l) {
1564e2606525SJoe Wallwork           evecs[i*dim+l] = 0.0;
1565e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1566e2606525SJoe Wallwork             for (k = 0; k < dim; ++k) {
1567e2606525SJoe Wallwork               evecs[i*dim+l] += isqrtM1[j*dim+i] * M2[j*dim+k] * isqrtM1[k*dim+l];
156820282da2SJoe Wallwork             }
156920282da2SJoe Wallwork           }
157020282da2SJoe Wallwork         }
157120282da2SJoe Wallwork       }
157220282da2SJoe Wallwork 
157320282da2SJoe Wallwork       /* Compute eigendecomposition */
15749566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
157520282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
157620282da2SJoe Wallwork       {
157720282da2SJoe Wallwork         PetscReal *rwork;
15789566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3*dim, &rwork));
1579*792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15809566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
158120282da2SJoe Wallwork       }
158220282da2SJoe Wallwork #else
1583*792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
158420282da2SJoe Wallwork #endif
158582490253SJoe Wallwork       if (lierr) {
158682490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
158782490253SJoe Wallwork           for (l = 0; l < dim; ++l) {
1588e2606525SJoe Wallwork             evecs[i*dim+l] = 0.0;
1589e2606525SJoe Wallwork             for (j = 0; j < dim; ++j) {
1590e2606525SJoe Wallwork               for (k = 0; k < dim; ++k) {
1591e2606525SJoe Wallwork                 evecs[i*dim+l] += isqrtM1[j*dim+i] * M2[j*dim+k] * isqrtM1[k*dim+l];
159282490253SJoe Wallwork               }
159382490253SJoe Wallwork             }
159482490253SJoe Wallwork           }
159582490253SJoe Wallwork         }
159682490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
159798921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
159882490253SJoe Wallwork       }
15999566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
160020282da2SJoe Wallwork 
160120282da2SJoe Wallwork       /* Modify eigenvalues */
1602e2606525SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
160320282da2SJoe Wallwork 
160420282da2SJoe Wallwork       /* Map back to get the intersection */
160520282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
1606e2606525SJoe Wallwork         for (m = 0; m < dim; ++m) {
1607e2606525SJoe Wallwork           M2[i*dim+m] = 0.0;
160820282da2SJoe Wallwork           for (j = 0; j < dim; ++j) {
160920282da2SJoe Wallwork             for (k = 0; k < dim; ++k) {
161020282da2SJoe Wallwork               for (l = 0; l < dim; ++l) {
1611e2606525SJoe Wallwork                 M2[i*dim+m] += sqrtM1[j*dim+i] * evecs[j*dim+k] * evals[k] * evecs[l*dim+k] * sqrtM1[l*dim+m];
161220282da2SJoe Wallwork               }
161320282da2SJoe Wallwork             }
161420282da2SJoe Wallwork           }
161520282da2SJoe Wallwork         }
161620282da2SJoe Wallwork       }
161720282da2SJoe Wallwork     }
16189566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
161920282da2SJoe Wallwork   }
1620e2606525SJoe Wallwork   PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
162120282da2SJoe Wallwork   PetscFunctionReturn(0);
162220282da2SJoe Wallwork }
162320282da2SJoe Wallwork 
1624cb7bfe8cSJoe Wallwork /*@
162520282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
162620282da2SJoe Wallwork 
1627f1a722f8SMatthew G. Knepley   Input Parameters:
162820282da2SJoe Wallwork + dm         - The DM
162920282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
163020282da2SJoe Wallwork - metrics    - The metrics to be intersected
163120282da2SJoe Wallwork 
163220282da2SJoe Wallwork   Output Parameter:
163320282da2SJoe Wallwork . metricInt  - The intersected metric
163420282da2SJoe Wallwork 
163520282da2SJoe Wallwork   Level: beginner
163620282da2SJoe Wallwork 
163720282da2SJoe Wallwork   Notes:
163820282da2SJoe Wallwork 
1639e2606525SJoe Wallwork   The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
164020282da2SJoe Wallwork 
1641e2606525SJoe Wallwork   The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
164220282da2SJoe Wallwork 
1643db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1644cb7bfe8cSJoe Wallwork @*/
16455241a91fSJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
164620282da2SJoe Wallwork {
164793520041SJoe Wallwork   PetscBool      isotropic, uniform;
164893520041SJoe Wallwork   PetscInt       v, i, m, n;
164993520041SJoe Wallwork   PetscScalar   *met, *meti;
165020282da2SJoe Wallwork 
165120282da2SJoe Wallwork   PetscFunctionBegin;
16529566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0));
165363a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
165420282da2SJoe Wallwork 
165520282da2SJoe Wallwork   /* Copy over the first metric */
16565241a91fSJoe Wallwork   PetscCall(VecCopy(metrics[0], metricInt));
165793520041SJoe Wallwork   if (numMetrics == 1) PetscFunctionReturn(0);
16585241a91fSJoe Wallwork   PetscCall(VecGetSize(metricInt, &m));
165993520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
16609566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
166108401ef6SPierre Jolivet     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
166293520041SJoe Wallwork   }
166320282da2SJoe Wallwork 
166420282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
16659566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
16669566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
166793520041SJoe Wallwork   if (uniform) {
166893520041SJoe Wallwork 
166993520041SJoe Wallwork     /* Uniform case */
16705241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
167193520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16729566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
16739566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
16749566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
167593520041SJoe Wallwork     }
16765241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
167793520041SJoe Wallwork   } else {
167893520041SJoe Wallwork 
167993520041SJoe Wallwork     /* Spatially varying case */
168093520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
168193520041SJoe Wallwork     PetscScalar *M, *Mi;
168293520041SJoe Wallwork 
16839566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
168493520041SJoe Wallwork     if (isotropic) nrow = 1;
168593520041SJoe Wallwork     else nrow = dim;
16869566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
16875241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
168820282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16899566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
169020282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
16919566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
16929566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
16939566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
169420282da2SJoe Wallwork       }
16959566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
169620282da2SJoe Wallwork     }
16975241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
169820282da2SJoe Wallwork   }
1699fe902aceSJoe Wallwork 
17009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0));
170120282da2SJoe Wallwork   PetscFunctionReturn(0);
170220282da2SJoe Wallwork }
170320282da2SJoe Wallwork 
1704cb7bfe8cSJoe Wallwork /*@
170520282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
170620282da2SJoe Wallwork 
1707f1a722f8SMatthew G. Knepley   Input Parameters:
170820282da2SJoe Wallwork + dm        - The DM
170920282da2SJoe Wallwork . metric1   - The first metric to be intersected
171020282da2SJoe Wallwork - metric2   - The second metric to be intersected
171120282da2SJoe Wallwork 
171220282da2SJoe Wallwork   Output Parameter:
171320282da2SJoe Wallwork . metricInt - The intersected metric
171420282da2SJoe Wallwork 
171520282da2SJoe Wallwork   Level: beginner
171620282da2SJoe Wallwork 
1717db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1718cb7bfe8cSJoe Wallwork @*/
17195241a91fSJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
172020282da2SJoe Wallwork {
172120282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
172220282da2SJoe Wallwork 
172320282da2SJoe Wallwork   PetscFunctionBegin;
17249566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
172520282da2SJoe Wallwork   PetscFunctionReturn(0);
172620282da2SJoe Wallwork }
172720282da2SJoe Wallwork 
1728cb7bfe8cSJoe Wallwork /*@
172920282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
173020282da2SJoe Wallwork 
1731f1a722f8SMatthew G. Knepley   Input Parameters:
173220282da2SJoe Wallwork + dm        - The DM
173320282da2SJoe Wallwork . metric1   - The first metric to be intersected
173420282da2SJoe Wallwork . metric2   - The second metric to be intersected
173520282da2SJoe Wallwork - metric3   - The third metric to be intersected
173620282da2SJoe Wallwork 
173720282da2SJoe Wallwork   Output Parameter:
173820282da2SJoe Wallwork . metricInt - The intersected metric
173920282da2SJoe Wallwork 
174020282da2SJoe Wallwork   Level: beginner
174120282da2SJoe Wallwork 
1742db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1743cb7bfe8cSJoe Wallwork @*/
17445241a91fSJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
174520282da2SJoe Wallwork {
174620282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
174720282da2SJoe Wallwork 
174820282da2SJoe Wallwork   PetscFunctionBegin;
17499566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
175020282da2SJoe Wallwork   PetscFunctionReturn(0);
175120282da2SJoe Wallwork }
1752