xref: /petsc/src/dm/impls/plex/plexmetric.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
6d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
7d71ae5a4SJacob Faibussowitsch {
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();
52*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
66d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
67d71ae5a4SJacob Faibussowitsch {
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;
73*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
89d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
90d71ae5a4SJacob Faibussowitsch {
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;
96*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
114d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
115d71ae5a4SJacob Faibussowitsch {
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;
122*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
138d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
139d71ae5a4SJacob Faibussowitsch {
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;
145*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
159d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
160d71ae5a4SJacob Faibussowitsch {
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;
166*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
182d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
183d71ae5a4SJacob Faibussowitsch {
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;
189*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
206d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
207d71ae5a4SJacob Faibussowitsch {
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;
213*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
232d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
233d71ae5a4SJacob Faibussowitsch {
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;
239*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
256d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
257d71ae5a4SJacob Faibussowitsch {
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;
263*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
282d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
283d71ae5a4SJacob Faibussowitsch {
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;
289*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
306d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
307d71ae5a4SJacob Faibussowitsch {
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;
313*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
332d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
333d71ae5a4SJacob Faibussowitsch {
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;
339*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
356d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
357d71ae5a4SJacob Faibussowitsch {
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;
363*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
382d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
383d71ae5a4SJacob Faibussowitsch {
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;
389*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
403d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
404d71ae5a4SJacob Faibussowitsch {
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;
411*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
427d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
428d71ae5a4SJacob Faibussowitsch {
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;
434*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
448d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
449d71ae5a4SJacob Faibussowitsch {
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;
456*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
472d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
473d71ae5a4SJacob Faibussowitsch {
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;
479*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
495d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
496d71ae5a4SJacob Faibussowitsch {
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;
503*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
519d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
520d71ae5a4SJacob Faibussowitsch {
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;
526*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
540d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
541d71ae5a4SJacob Faibussowitsch {
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;
548*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
564d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
565d71ae5a4SJacob Faibussowitsch {
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;
571*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
585d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
586d71ae5a4SJacob Faibussowitsch {
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;
593*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
609d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
610d71ae5a4SJacob Faibussowitsch {
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;
616*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
638d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
639d71ae5a4SJacob Faibussowitsch {
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;
646*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
670d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
671d71ae5a4SJacob Faibussowitsch {
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;
677*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
702d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
703d71ae5a4SJacob Faibussowitsch {
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;
710*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
737d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
738d71ae5a4SJacob Faibussowitsch {
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;
744*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
761d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
762d71ae5a4SJacob Faibussowitsch {
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;
768*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
787d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
788d71ae5a4SJacob Faibussowitsch {
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;
794*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
811d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
812d71ae5a4SJacob Faibussowitsch {
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;
818*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
837d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
838d71ae5a4SJacob Faibussowitsch {
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;
844*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
845cc2eee55SJoe Wallwork }
846cc2eee55SJoe Wallwork 
847d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
848d71ae5a4SJacob Faibussowitsch {
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 
866*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
911d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
912d71ae5a4SJacob Faibussowitsch {
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));
931*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 @*/
951d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
952d71ae5a4SJacob Faibussowitsch {
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));
961*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96220282da2SJoe Wallwork }
96320282da2SJoe Wallwork 
964d71ae5a4SJacob Faibussowitsch static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
965d71ae5a4SJacob Faibussowitsch {
96693520041SJoe Wallwork   f0[0] = u[0];
96793520041SJoe Wallwork }
96893520041SJoe Wallwork 
969cb7bfe8cSJoe Wallwork /*@
97020282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
97120282da2SJoe Wallwork 
97220282da2SJoe Wallwork   Input parameters:
97320282da2SJoe Wallwork + dm        - The DM
97420282da2SJoe Wallwork . f         - The field number to use
97520282da2SJoe Wallwork - indicator - The error indicator
97620282da2SJoe Wallwork 
97720282da2SJoe Wallwork   Output parameter:
97820282da2SJoe Wallwork . metric    - The isotropic metric
97920282da2SJoe Wallwork 
98020282da2SJoe Wallwork   Level: beginner
98120282da2SJoe Wallwork 
98220282da2SJoe Wallwork   Notes:
98320282da2SJoe Wallwork 
98420282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
98520282da2SJoe Wallwork 
98693520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
98720282da2SJoe Wallwork 
988db781477SPatrick Sanan .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
989cb7bfe8cSJoe Wallwork @*/
990d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
991d71ae5a4SJacob Faibussowitsch {
99293520041SJoe Wallwork   PetscInt m, n;
99320282da2SJoe Wallwork 
99420282da2SJoe Wallwork   PetscFunctionBegin;
9959566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
9969566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9979566063dSJacob Faibussowitsch   PetscCall(VecGetSize(indicator, &m));
9989566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metric, &n));
9999566063dSJacob Faibussowitsch   if (m == n) PetscCall(VecCopy(indicator, *metric));
100093520041SJoe Wallwork   else {
100193520041SJoe Wallwork     DM dmIndi;
10029371c9d4SSatish Balay     void (*funcs[1])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
100393520041SJoe Wallwork 
10049566063dSJacob Faibussowitsch     PetscCall(VecGetDM(indicator, &dmIndi));
100593520041SJoe Wallwork     funcs[0] = identity;
10069566063dSJacob Faibussowitsch     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
100720282da2SJoe Wallwork   }
1008*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100920282da2SJoe Wallwork }
101020282da2SJoe Wallwork 
1011f6db075bSJoe Wallwork /*@
1012f6db075bSJoe Wallwork   DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1013f6db075bSJoe Wallwork 
1014f6db075bSJoe Wallwork   Input parameters:
1015f6db075bSJoe Wallwork + dm          - The DM of the metric field
1016f6db075bSJoe Wallwork - f           - The field number to use
1017f6db075bSJoe Wallwork 
1018f6db075bSJoe Wallwork   Output parameter:
1019f6db075bSJoe Wallwork + determinant - The determinant field
1020f6db075bSJoe Wallwork - dmDet       - The corresponding DM
1021f6db075bSJoe Wallwork 
1022f6db075bSJoe Wallwork   Level: beginner
1023f6db075bSJoe Wallwork 
1024f6db075bSJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic(), DMPlexMetricCreate()
1025f6db075bSJoe Wallwork @*/
1026d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1027d71ae5a4SJacob Faibussowitsch {
1028f6db075bSJoe Wallwork   PetscBool uniform;
1029f6db075bSJoe Wallwork 
1030f6db075bSJoe Wallwork   PetscFunctionBegin;
1031f6db075bSJoe Wallwork   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1032f6db075bSJoe Wallwork   PetscCall(DMClone(dm, dmDet));
1033f6db075bSJoe Wallwork   if (uniform) {
1034f6db075bSJoe Wallwork     MPI_Comm comm;
1035f6db075bSJoe Wallwork 
1036f6db075bSJoe Wallwork     PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm));
1037f6db075bSJoe Wallwork     PetscCall(VecCreate(comm, determinant));
1038f6db075bSJoe Wallwork     PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1039f6db075bSJoe Wallwork     PetscCall(VecSetFromOptions(*determinant));
1040f6db075bSJoe Wallwork   } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
1041*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1042f6db075bSJoe Wallwork }
1043f6db075bSJoe Wallwork 
1044d71ae5a4SJacob Faibussowitsch static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1045d71ae5a4SJacob Faibussowitsch {
104682490253SJoe Wallwork   PetscInt i, j;
104782490253SJoe Wallwork 
104882490253SJoe Wallwork   PetscFunctionBegin;
1049*3ba16761SJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"));
105082490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
1051*3ba16761SJacob Faibussowitsch     if (i == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    [["));
1052*3ba16761SJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "     ["));
105382490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
1054*3ba16761SJacob Faibussowitsch       if (j < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])));
1055*3ba16761SJacob Faibussowitsch       else PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])));
105682490253SJoe Wallwork     }
1057*3ba16761SJacob Faibussowitsch     if (i < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
1058*3ba16761SJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "]]\n"));
105982490253SJoe Wallwork   }
1060*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106182490253SJoe Wallwork }
106282490253SJoe Wallwork 
1063d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1064d71ae5a4SJacob Faibussowitsch {
106520282da2SJoe Wallwork   PetscInt     i, j, k;
106620282da2SJoe 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);
106720282da2SJoe Wallwork   PetscScalar *Mpos;
106820282da2SJoe Wallwork 
106920282da2SJoe Wallwork   PetscFunctionBegin;
10709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs));
107120282da2SJoe Wallwork 
107220282da2SJoe Wallwork   /* Symmetrize */
107320282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
107420282da2SJoe Wallwork     Mpos[i * dim + i] = Mp[i * dim + i];
107520282da2SJoe Wallwork     for (j = i + 1; j < dim; ++j) {
107620282da2SJoe Wallwork       Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
107720282da2SJoe Wallwork       Mpos[j * dim + i] = Mpos[i * dim + j];
107820282da2SJoe Wallwork     }
107920282da2SJoe Wallwork   }
108020282da2SJoe Wallwork 
108120282da2SJoe Wallwork   /* Compute eigendecomposition */
108293520041SJoe Wallwork   if (dim == 1) {
108393520041SJoe Wallwork     /* Isotropic case */
108493520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
108593520041SJoe Wallwork     Mpos[0] = 1.0;
108693520041SJoe Wallwork   } else {
108793520041SJoe Wallwork     /* Anisotropic case */
108820282da2SJoe Wallwork     PetscScalar *work;
108920282da2SJoe Wallwork     PetscBLASInt lwork;
109020282da2SJoe Wallwork 
109120282da2SJoe Wallwork     lwork = 5 * dim;
10929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * dim, &work));
109320282da2SJoe Wallwork     {
109420282da2SJoe Wallwork       PetscBLASInt lierr;
109520282da2SJoe Wallwork       PetscBLASInt nb;
109620282da2SJoe Wallwork 
10979566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
10989566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
109920282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
110020282da2SJoe Wallwork       {
110120282da2SJoe Wallwork         PetscReal *rwork;
11029566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1103792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr));
11049566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
110520282da2SJoe Wallwork       }
110620282da2SJoe Wallwork #else
1107792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr));
110820282da2SJoe Wallwork #endif
110982490253SJoe Wallwork       if (lierr) {
111082490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
111182490253SJoe Wallwork           Mpos[i * dim + i] = Mp[i * dim + i];
111282490253SJoe Wallwork           for (j = i + 1; j < dim; ++j) {
111382490253SJoe Wallwork             Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
111482490253SJoe Wallwork             Mpos[j * dim + i] = Mpos[i * dim + j];
111582490253SJoe Wallwork           }
111682490253SJoe Wallwork         }
1117*3ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, Mpos));
111898921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
111982490253SJoe Wallwork       }
11209566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
112120282da2SJoe Wallwork     }
11229566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
112320282da2SJoe Wallwork   }
112420282da2SJoe Wallwork 
112520282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
112620282da2SJoe Wallwork   max_eig = 0.0;
112720282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
112820282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
112920282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
113020282da2SJoe Wallwork   }
113120282da2SJoe Wallwork 
11323f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11333f07679eSJoe Wallwork   *detMp = 1.0;
113420282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
113520282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min);
11363f07679eSJoe Wallwork     *detMp *= eigs[i];
113720282da2SJoe Wallwork   }
113820282da2SJoe Wallwork 
113920282da2SJoe Wallwork   /* Reconstruct Hessian */
114020282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
114120282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
114220282da2SJoe Wallwork       Mp[i * dim + j] = 0.0;
1143ad540459SPierre Jolivet       for (k = 0; k < dim; ++k) Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j];
114420282da2SJoe Wallwork     }
114520282da2SJoe Wallwork   }
11469566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Mpos, eigs));
114720282da2SJoe Wallwork 
1148*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114920282da2SJoe Wallwork }
115020282da2SJoe Wallwork 
1151cb7bfe8cSJoe Wallwork /*@
115220282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
115320282da2SJoe Wallwork 
115420282da2SJoe Wallwork   Input parameters:
115520282da2SJoe Wallwork + dm                 - The DM
11563f07679eSJoe Wallwork . metricIn           - The metric
115799abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11583f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
115920282da2SJoe Wallwork 
116020282da2SJoe Wallwork   Output parameter:
11613f07679eSJoe Wallwork + metricOut          - The metric
11623f07679eSJoe Wallwork - determinant        - Its determinant
116320282da2SJoe Wallwork 
116420282da2SJoe Wallwork   Level: beginner
116520282da2SJoe Wallwork 
1166cb7bfe8cSJoe Wallwork   Notes:
1167cb7bfe8cSJoe Wallwork 
1168cb7bfe8cSJoe Wallwork   Relevant command line options:
1169cb7bfe8cSJoe Wallwork 
117093520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
117193520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
117293520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1173cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1174cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1175cb7bfe8cSJoe Wallwork 
1176db781477SPatrick Sanan .seealso: `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1177cb7bfe8cSJoe Wallwork @*/
1178d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1179d71ae5a4SJacob Faibussowitsch {
11803f07679eSJoe Wallwork   DM           dmDet;
118193520041SJoe Wallwork   PetscBool    isotropic, uniform;
118220282da2SJoe Wallwork   PetscInt     dim, vStart, vEnd, v;
11833f07679eSJoe Wallwork   PetscScalar *met, *det;
118420282da2SJoe Wallwork   PetscReal    h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
118520282da2SJoe Wallwork 
118620282da2SJoe Wallwork   PetscFunctionBegin;
11879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
118820282da2SJoe Wallwork 
118920282da2SJoe Wallwork   /* Extract metadata from dm */
11909566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
119120282da2SJoe Wallwork   if (restrictSizes) {
11929566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
11939566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
119499abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
119599abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1196bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
119799abec2bSJoe Wallwork   }
119899abec2bSJoe Wallwork   if (restrictAnisotropy) {
11999566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
120099abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
120120282da2SJoe Wallwork   }
120220282da2SJoe Wallwork 
120393520041SJoe Wallwork   /* Setup output metric */
12045241a91fSJoe Wallwork   PetscCall(VecCopy(metricIn, metricOut));
120593520041SJoe Wallwork 
120693520041SJoe Wallwork   /* Enforce SPD and extract determinant */
12075241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
12089566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12099566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
121093520041SJoe Wallwork   if (uniform) {
1211d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
121293520041SJoe Wallwork 
121393520041SJoe Wallwork     /* Uniform case */
12145241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
12159566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
12165241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
121793520041SJoe Wallwork   } else {
121893520041SJoe Wallwork     /* Spatially varying case */
121993520041SJoe Wallwork     PetscInt nrow;
122093520041SJoe Wallwork 
122193520041SJoe Wallwork     if (isotropic) nrow = 1;
122293520041SJoe Wallwork     else nrow = dim;
12235241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
12249566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12255241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
122620282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12273f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
12289566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
12299566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
12309566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
123120282da2SJoe Wallwork     }
12325241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
123393520041SJoe Wallwork   }
12345241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
1235fe902aceSJoe Wallwork 
12369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1237*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123820282da2SJoe Wallwork }
123920282da2SJoe Wallwork 
1240d71ae5a4SJacob Faibussowitsch static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1241d71ae5a4SJacob Faibussowitsch {
124220282da2SJoe Wallwork   const PetscScalar p = constants[0];
124320282da2SJoe Wallwork 
1244985ad86eSJose E. Roman   f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim));
124520282da2SJoe Wallwork }
124620282da2SJoe Wallwork 
1247cb7bfe8cSJoe Wallwork /*@
124820282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
124920282da2SJoe Wallwork 
125020282da2SJoe Wallwork   Input parameters:
125120282da2SJoe Wallwork + dm                 - The DM
125220282da2SJoe Wallwork . metricIn           - The unnormalized metric
125316522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
125416522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
125520282da2SJoe Wallwork 
125620282da2SJoe Wallwork   Output parameter:
125720282da2SJoe Wallwork . metricOut          - The normalized metric
125820282da2SJoe Wallwork 
125920282da2SJoe Wallwork   Level: beginner
126020282da2SJoe Wallwork 
1261cb7bfe8cSJoe Wallwork   Notes:
1262cb7bfe8cSJoe Wallwork 
1263cb7bfe8cSJoe Wallwork   Relevant command line options:
1264cb7bfe8cSJoe Wallwork 
126593520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
126693520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
126793520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1268cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1269cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1270cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1271cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1272cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1273cb7bfe8cSJoe Wallwork 
1274db781477SPatrick Sanan .seealso: `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1275cb7bfe8cSJoe Wallwork @*/
1276d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1277d71ae5a4SJacob Faibussowitsch {
12783f07679eSJoe Wallwork   DM           dmDet;
127920282da2SJoe Wallwork   MPI_Comm     comm;
128093520041SJoe Wallwork   PetscBool    restrictAnisotropyFirst, isotropic, uniform;
128120282da2SJoe Wallwork   PetscDS      ds;
128220282da2SJoe Wallwork   PetscInt     dim, Nd, vStart, vEnd, v, i;
12833f07679eSJoe Wallwork   PetscScalar *met, *det, integral, constants[1];
128493520041SJoe Wallwork   PetscReal    p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
128520282da2SJoe Wallwork 
128620282da2SJoe Wallwork   PetscFunctionBegin;
12879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0));
128820282da2SJoe Wallwork 
128920282da2SJoe Wallwork   /* Extract metadata from dm */
12909566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
12919566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
12929566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12939566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
129493520041SJoe Wallwork   if (isotropic) Nd = 1;
129593520041SJoe Wallwork   else Nd = dim * dim;
129620282da2SJoe Wallwork 
129720282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
12989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
12995241a91fSJoe Wallwork   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
130020282da2SJoe Wallwork 
130120282da2SJoe Wallwork   /* Compute global normalization factor */
13029566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
13039566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
130416522936SJoe Wallwork   constants[0] = p;
130593520041SJoe Wallwork   if (uniform) {
1306bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
130793520041SJoe Wallwork     DM  dmTmp;
130893520041SJoe Wallwork     Vec tmp;
130993520041SJoe Wallwork 
13109566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmTmp));
13119566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
13129566063dSJacob Faibussowitsch     PetscCall(VecGetArray(determinant, &det));
13139566063dSJacob Faibussowitsch     PetscCall(VecSet(tmp, det[0]));
13149566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(determinant, &det));
13159566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmTmp, &ds));
13169566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13179566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13189566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
13199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tmp));
13209566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmTmp));
132193520041SJoe Wallwork   } else {
13229566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
13239566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmDet, &ds));
13249566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13259566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13269566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
132793520041SJoe Wallwork   }
13283f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
1329bc00d7afSJoe 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?");
13303f07679eSJoe Wallwork   factGlob = PetscPowReal(target / realIntegral, 2.0 / dim);
133120282da2SJoe Wallwork 
133220282da2SJoe Wallwork   /* Apply local scaling */
133320282da2SJoe Wallwork   if (restrictSizes) {
13349566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
13359566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
133616522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
133716522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1338bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
133916522936SJoe Wallwork   }
134016522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
13419566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
134216522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
134320282da2SJoe Wallwork   }
13445241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
13459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(determinant, &det));
134693520041SJoe Wallwork   if (uniform) {
134793520041SJoe Wallwork     /* Uniform case */
134893520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim));
13499566063dSJacob Faibussowitsch     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
135093520041SJoe Wallwork   } else {
135193520041SJoe Wallwork     /* Spatially varying case */
135293520041SJoe Wallwork     PetscInt nrow;
135393520041SJoe Wallwork 
135493520041SJoe Wallwork     if (isotropic) nrow = 1;
135593520041SJoe Wallwork     else nrow = dim;
13569566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
13575241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
135820282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13593f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
136020282da2SJoe Wallwork 
13619566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
13629566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
13633f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim));
136420282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
13659566063dSJacob Faibussowitsch       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
136693520041SJoe Wallwork     }
136720282da2SJoe Wallwork   }
13689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(determinant, &det));
13695241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
137020282da2SJoe Wallwork 
13719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1372*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
137320282da2SJoe Wallwork }
137420282da2SJoe Wallwork 
1375cb7bfe8cSJoe Wallwork /*@
137620282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
137720282da2SJoe Wallwork 
1378f1a722f8SMatthew G. Knepley   Input Parameters:
137920282da2SJoe Wallwork + dm         - The DM
138020282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
138120282da2SJoe Wallwork . weights    - Weights for the average
138220282da2SJoe Wallwork - metrics    - The metrics to be averaged
138320282da2SJoe Wallwork 
138420282da2SJoe Wallwork   Output Parameter:
138520282da2SJoe Wallwork . metricAvg  - The averaged metric
138620282da2SJoe Wallwork 
138720282da2SJoe Wallwork   Level: beginner
138820282da2SJoe Wallwork 
138920282da2SJoe Wallwork   Notes:
139020282da2SJoe Wallwork   The weights should sum to unity.
139120282da2SJoe Wallwork 
139220282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
139320282da2SJoe Wallwork 
1394db781477SPatrick Sanan .seealso: `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1395cb7bfe8cSJoe Wallwork @*/
1396d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1397d71ae5a4SJacob Faibussowitsch {
139820282da2SJoe Wallwork   PetscBool haveWeights = PETSC_TRUE;
139993520041SJoe Wallwork   PetscInt  i, m, n;
140020282da2SJoe Wallwork   PetscReal sum = 0.0, tol = 1.0e-10;
140120282da2SJoe Wallwork 
140220282da2SJoe Wallwork   PetscFunctionBegin;
14039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0));
140463a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
14055241a91fSJoe Wallwork   PetscCall(VecSet(metricAvg, 0.0));
14065241a91fSJoe Wallwork   PetscCall(VecGetSize(metricAvg, &m));
140793520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
14089566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
14095f80ce2aSJacob Faibussowitsch     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
141093520041SJoe Wallwork   }
141120282da2SJoe Wallwork 
141220282da2SJoe Wallwork   /* Default to the unweighted case */
141320282da2SJoe Wallwork   if (!weights) {
14149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numMetrics, &weights));
141520282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
1416ad540459SPierre Jolivet     for (i = 0; i < numMetrics; ++i) weights[i] = 1.0 / numMetrics;
141720282da2SJoe Wallwork   }
141820282da2SJoe Wallwork 
141920282da2SJoe Wallwork   /* Check weights sum to unity */
142093520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
14215f80ce2aSJacob Faibussowitsch   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
142220282da2SJoe Wallwork 
142320282da2SJoe Wallwork   /* Compute metric average */
14245241a91fSJoe Wallwork   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
14259566063dSJacob Faibussowitsch   if (!haveWeights) PetscCall(PetscFree(weights));
1426fe902aceSJoe Wallwork 
14279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0));
1428*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142920282da2SJoe Wallwork }
143020282da2SJoe Wallwork 
1431cb7bfe8cSJoe Wallwork /*@
143220282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
143320282da2SJoe Wallwork 
1434f1a722f8SMatthew G. Knepley   Input Parameters:
143520282da2SJoe Wallwork + dm         - The DM
143620282da2SJoe Wallwork . metric1    - The first metric to be averaged
143720282da2SJoe Wallwork - metric2    - The second metric to be averaged
143820282da2SJoe Wallwork 
143920282da2SJoe Wallwork   Output Parameter:
144020282da2SJoe Wallwork . metricAvg  - The averaged metric
144120282da2SJoe Wallwork 
144220282da2SJoe Wallwork   Level: beginner
144320282da2SJoe Wallwork 
1444db781477SPatrick Sanan .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1445cb7bfe8cSJoe Wallwork @*/
1446d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1447d71ae5a4SJacob Faibussowitsch {
144820282da2SJoe Wallwork   PetscReal weights[2] = {0.5, 0.5};
144920282da2SJoe Wallwork   Vec       metrics[2] = {metric1, metric2};
145020282da2SJoe Wallwork 
145120282da2SJoe Wallwork   PetscFunctionBegin;
14529566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
1453*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145420282da2SJoe Wallwork }
145520282da2SJoe Wallwork 
1456cb7bfe8cSJoe Wallwork /*@
145720282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
145820282da2SJoe Wallwork 
1459f1a722f8SMatthew G. Knepley   Input Parameters:
146020282da2SJoe Wallwork + dm         - The DM
146120282da2SJoe Wallwork . metric1    - The first metric to be averaged
146220282da2SJoe Wallwork . metric2    - The second metric to be averaged
146320282da2SJoe Wallwork - metric3    - The third metric to be averaged
146420282da2SJoe Wallwork 
146520282da2SJoe Wallwork   Output Parameter:
146620282da2SJoe Wallwork . metricAvg  - The averaged metric
146720282da2SJoe Wallwork 
146820282da2SJoe Wallwork   Level: beginner
146920282da2SJoe Wallwork 
1470db781477SPatrick Sanan .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1471cb7bfe8cSJoe Wallwork @*/
1472d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1473d71ae5a4SJacob Faibussowitsch {
147420282da2SJoe Wallwork   PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0};
147520282da2SJoe Wallwork   Vec       metrics[3] = {metric1, metric2, metric3};
147620282da2SJoe Wallwork 
147720282da2SJoe Wallwork   PetscFunctionBegin;
14789566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
1479*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148020282da2SJoe Wallwork }
148120282da2SJoe Wallwork 
1482d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1483d71ae5a4SJacob Faibussowitsch {
148420282da2SJoe Wallwork   PetscInt     i, j, k, l, m;
1485e2606525SJoe Wallwork   PetscReal   *evals;
148620282da2SJoe Wallwork   PetscScalar *evecs, *sqrtM1, *isqrtM1;
148720282da2SJoe Wallwork 
148820282da2SJoe Wallwork   PetscFunctionBegin;
148993520041SJoe Wallwork 
149093520041SJoe Wallwork   /* Isotropic case */
149193520041SJoe Wallwork   if (dim == 1) {
14926f71502aSJoe Wallwork     M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
1493*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
149493520041SJoe Wallwork   }
149593520041SJoe Wallwork 
149693520041SJoe Wallwork   /* Anisotropic case */
1497e2606525SJoe Wallwork   PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals));
149820282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
1499ad540459SPierre Jolivet     for (j = 0; j < dim; ++j) evecs[i * dim + j] = M1[i * dim + j];
150020282da2SJoe Wallwork   }
150120282da2SJoe Wallwork   {
150220282da2SJoe Wallwork     PetscScalar *work;
150320282da2SJoe Wallwork     PetscBLASInt lwork;
150420282da2SJoe Wallwork 
150520282da2SJoe Wallwork     lwork = 5 * dim;
15069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * dim, &work));
150720282da2SJoe Wallwork     {
150820282da2SJoe Wallwork       PetscBLASInt lierr, nb;
1509e2606525SJoe Wallwork       PetscReal    sqrtj;
151020282da2SJoe Wallwork 
151120282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
15129566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
15139566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
151420282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
151520282da2SJoe Wallwork       {
151620282da2SJoe Wallwork         PetscReal *rwork;
15179566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1518792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15199566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
152020282da2SJoe Wallwork       }
152120282da2SJoe Wallwork #else
1522792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
152320282da2SJoe Wallwork #endif
152482490253SJoe Wallwork       if (lierr) {
1525*3ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, M1));
152698921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
152782490253SJoe Wallwork       }
15289566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
152920282da2SJoe Wallwork 
1530e2606525SJoe Wallwork       /* Compute square root and the reciprocal thereof */
153120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
153220282da2SJoe Wallwork         for (k = 0; k < dim; ++k) {
1533e2606525SJoe Wallwork           sqrtM1[i * dim + k]  = 0.0;
1534e2606525SJoe Wallwork           isqrtM1[i * dim + k] = 0.0;
1535e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1536e2606525SJoe Wallwork             sqrtj = PetscSqrtReal(evals[j]);
1537e2606525SJoe Wallwork             sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k];
1538e2606525SJoe Wallwork             isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k];
153920282da2SJoe Wallwork           }
154020282da2SJoe Wallwork         }
154120282da2SJoe Wallwork       }
154220282da2SJoe Wallwork 
1543e2606525SJoe Wallwork       /* Map M2 into the space spanned by the eigenvectors of M1 */
154420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
154520282da2SJoe Wallwork         for (l = 0; l < dim; ++l) {
1546e2606525SJoe Wallwork           evecs[i * dim + l] = 0.0;
1547e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1548ad540459SPierre Jolivet             for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
154920282da2SJoe Wallwork           }
155020282da2SJoe Wallwork         }
155120282da2SJoe Wallwork       }
155220282da2SJoe Wallwork 
155320282da2SJoe Wallwork       /* Compute eigendecomposition */
15549566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
155520282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
155620282da2SJoe Wallwork       {
155720282da2SJoe Wallwork         PetscReal *rwork;
15589566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1559792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15609566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
156120282da2SJoe Wallwork       }
156220282da2SJoe Wallwork #else
1563792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
156420282da2SJoe Wallwork #endif
156582490253SJoe Wallwork       if (lierr) {
156682490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
156782490253SJoe Wallwork           for (l = 0; l < dim; ++l) {
1568e2606525SJoe Wallwork             evecs[i * dim + l] = 0.0;
1569e2606525SJoe Wallwork             for (j = 0; j < dim; ++j) {
1570ad540459SPierre Jolivet               for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
157182490253SJoe Wallwork             }
157282490253SJoe Wallwork           }
157382490253SJoe Wallwork         }
1574*3ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, evecs));
157598921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
157682490253SJoe Wallwork       }
15779566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
157820282da2SJoe Wallwork 
157920282da2SJoe Wallwork       /* Modify eigenvalues */
1580e2606525SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
158120282da2SJoe Wallwork 
158220282da2SJoe Wallwork       /* Map back to get the intersection */
158320282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
1584e2606525SJoe Wallwork         for (m = 0; m < dim; ++m) {
1585e2606525SJoe Wallwork           M2[i * dim + m] = 0.0;
158620282da2SJoe Wallwork           for (j = 0; j < dim; ++j) {
158720282da2SJoe Wallwork             for (k = 0; k < dim; ++k) {
1588ad540459SPierre Jolivet               for (l = 0; l < dim; ++l) M2[i * dim + m] += sqrtM1[j * dim + i] * evecs[j * dim + k] * evals[k] * evecs[l * dim + k] * sqrtM1[l * dim + m];
158920282da2SJoe Wallwork             }
159020282da2SJoe Wallwork           }
159120282da2SJoe Wallwork         }
159220282da2SJoe Wallwork       }
159320282da2SJoe Wallwork     }
15949566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
159520282da2SJoe Wallwork   }
1596e2606525SJoe Wallwork   PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
1597*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159820282da2SJoe Wallwork }
159920282da2SJoe Wallwork 
1600cb7bfe8cSJoe Wallwork /*@
160120282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
160220282da2SJoe Wallwork 
1603f1a722f8SMatthew G. Knepley   Input Parameters:
160420282da2SJoe Wallwork + dm         - The DM
160520282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
160620282da2SJoe Wallwork - metrics    - The metrics to be intersected
160720282da2SJoe Wallwork 
160820282da2SJoe Wallwork   Output Parameter:
160920282da2SJoe Wallwork . metricInt  - The intersected metric
161020282da2SJoe Wallwork 
161120282da2SJoe Wallwork   Level: beginner
161220282da2SJoe Wallwork 
161320282da2SJoe Wallwork   Notes:
161420282da2SJoe Wallwork 
1615e2606525SJoe Wallwork   The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
161620282da2SJoe Wallwork 
1617e2606525SJoe Wallwork   The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
161820282da2SJoe Wallwork 
1619db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1620cb7bfe8cSJoe Wallwork @*/
1621d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1622d71ae5a4SJacob Faibussowitsch {
162393520041SJoe Wallwork   PetscBool    isotropic, uniform;
162493520041SJoe Wallwork   PetscInt     v, i, m, n;
162593520041SJoe Wallwork   PetscScalar *met, *meti;
162620282da2SJoe Wallwork 
162720282da2SJoe Wallwork   PetscFunctionBegin;
16289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0));
162963a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
163020282da2SJoe Wallwork 
163120282da2SJoe Wallwork   /* Copy over the first metric */
16325241a91fSJoe Wallwork   PetscCall(VecCopy(metrics[0], metricInt));
1633*3ba16761SJacob Faibussowitsch   if (numMetrics == 1) PetscFunctionReturn(PETSC_SUCCESS);
16345241a91fSJoe Wallwork   PetscCall(VecGetSize(metricInt, &m));
163593520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
16369566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
163708401ef6SPierre Jolivet     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
163893520041SJoe Wallwork   }
163920282da2SJoe Wallwork 
164020282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
16419566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
16429566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
164393520041SJoe Wallwork   if (uniform) {
164493520041SJoe Wallwork     /* Uniform case */
16455241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
164693520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16479566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
16489566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
16499566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
165093520041SJoe Wallwork     }
16515241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
165293520041SJoe Wallwork   } else {
165393520041SJoe Wallwork     /* Spatially varying case */
165493520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
165593520041SJoe Wallwork     PetscScalar *M, *Mi;
165693520041SJoe Wallwork 
16579566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
165893520041SJoe Wallwork     if (isotropic) nrow = 1;
165993520041SJoe Wallwork     else nrow = dim;
16609566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
16615241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
166220282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16639566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
166420282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
16659566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
16669566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
16679566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
166820282da2SJoe Wallwork       }
16699566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
167020282da2SJoe Wallwork     }
16715241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
167220282da2SJoe Wallwork   }
1673fe902aceSJoe Wallwork 
16749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1675*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167620282da2SJoe Wallwork }
167720282da2SJoe Wallwork 
1678cb7bfe8cSJoe Wallwork /*@
167920282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
168020282da2SJoe Wallwork 
1681f1a722f8SMatthew G. Knepley   Input Parameters:
168220282da2SJoe Wallwork + dm        - The DM
168320282da2SJoe Wallwork . metric1   - The first metric to be intersected
168420282da2SJoe Wallwork - metric2   - The second metric to be intersected
168520282da2SJoe Wallwork 
168620282da2SJoe Wallwork   Output Parameter:
168720282da2SJoe Wallwork . metricInt - The intersected metric
168820282da2SJoe Wallwork 
168920282da2SJoe Wallwork   Level: beginner
169020282da2SJoe Wallwork 
1691db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1692cb7bfe8cSJoe Wallwork @*/
1693d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1694d71ae5a4SJacob Faibussowitsch {
169520282da2SJoe Wallwork   Vec metrics[2] = {metric1, metric2};
169620282da2SJoe Wallwork 
169720282da2SJoe Wallwork   PetscFunctionBegin;
16989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
1699*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170020282da2SJoe Wallwork }
170120282da2SJoe Wallwork 
1702cb7bfe8cSJoe Wallwork /*@
170320282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
170420282da2SJoe Wallwork 
1705f1a722f8SMatthew G. Knepley   Input Parameters:
170620282da2SJoe Wallwork + dm        - The DM
170720282da2SJoe Wallwork . metric1   - The first metric to be intersected
170820282da2SJoe Wallwork . metric2   - The second metric to be intersected
170920282da2SJoe Wallwork - metric3   - The third metric to be intersected
171020282da2SJoe Wallwork 
171120282da2SJoe Wallwork   Output Parameter:
171220282da2SJoe Wallwork . metricInt - The intersected metric
171320282da2SJoe Wallwork 
171420282da2SJoe Wallwork   Level: beginner
171520282da2SJoe Wallwork 
1716db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1717cb7bfe8cSJoe Wallwork @*/
1718d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1719d71ae5a4SJacob Faibussowitsch {
172020282da2SJoe Wallwork   Vec metrics[3] = {metric1, metric2, metric3};
172120282da2SJoe Wallwork 
172220282da2SJoe Wallwork   PetscFunctionBegin;
17239566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
1724*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172520282da2SJoe Wallwork }
1726