xref: /petsc/src/dm/impls/plex/plexmetric.c (revision 2fe279fdf3e687a416e4eadb7d3c7a82d60442c6)
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();
523ba16761SJacob 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:
59*2fe279fdSBarry Smith + dm        - The `DM`
6031b70465SJoe Wallwork - isotropic - Is the metric isotropic?
6131b70465SJoe Wallwork 
6231b70465SJoe Wallwork   Level: beginner
6331b70465SJoe Wallwork 
64*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7431b70465SJoe Wallwork }
7531b70465SJoe Wallwork 
76cb7bfe8cSJoe Wallwork /*@
7793520041SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric isotropic?
7831b70465SJoe Wallwork 
7931b70465SJoe Wallwork   Input parameters:
80*2fe279fdSBarry Smith . dm        - The `DM`
8131b70465SJoe Wallwork 
8231b70465SJoe Wallwork   Output parameters:
8331b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8431b70465SJoe Wallwork 
8531b70465SJoe Wallwork   Level: beginner
8631b70465SJoe Wallwork 
87*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
963ba16761SJacob 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:
103*2fe279fdSBarry Smith + dm      - The `DM`
10493520041SJoe Wallwork - uniform - Is the metric uniform?
10593520041SJoe Wallwork 
10693520041SJoe Wallwork   Level: beginner
10793520041SJoe Wallwork 
108*2fe279fdSBarry Smith   Note:
10993520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11093520041SJoe Wallwork 
111*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
11293520041SJoe Wallwork @*/
113d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
114d71ae5a4SJacob Faibussowitsch {
11593520041SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
11693520041SJoe Wallwork 
11793520041SJoe Wallwork   PetscFunctionBegin;
118bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
11993520041SJoe Wallwork   plex->metricCtx->uniform = uniform;
12093520041SJoe Wallwork   if (uniform) plex->metricCtx->isotropic = uniform;
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12293520041SJoe Wallwork }
12393520041SJoe Wallwork 
12493520041SJoe Wallwork /*@
12593520041SJoe Wallwork   DMPlexMetricIsUniform - Is a metric uniform?
12693520041SJoe Wallwork 
12793520041SJoe Wallwork   Input parameters:
128*2fe279fdSBarry Smith . dm      - The `DM`
12993520041SJoe Wallwork 
13093520041SJoe Wallwork   Output parameters:
13193520041SJoe Wallwork . uniform - Is the metric uniform?
13293520041SJoe Wallwork 
13393520041SJoe Wallwork   Level: beginner
13493520041SJoe Wallwork 
135*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
13693520041SJoe Wallwork @*/
137d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
138d71ae5a4SJacob Faibussowitsch {
13993520041SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
14093520041SJoe Wallwork 
14193520041SJoe Wallwork   PetscFunctionBegin;
142bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
14393520041SJoe Wallwork   *uniform = plex->metricCtx->uniform;
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14593520041SJoe Wallwork }
14693520041SJoe Wallwork 
14793520041SJoe Wallwork /*@
14831b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
14931b70465SJoe Wallwork 
15031b70465SJoe Wallwork   Input parameters:
151*2fe279fdSBarry Smith + dm                      - The `DM`
15231b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
15331b70465SJoe Wallwork 
15431b70465SJoe Wallwork   Level: beginner
15531b70465SJoe Wallwork 
156*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
157cb7bfe8cSJoe Wallwork @*/
158d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
159d71ae5a4SJacob Faibussowitsch {
16031b70465SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
16131b70465SJoe Wallwork 
16231b70465SJoe Wallwork   PetscFunctionBegin;
163bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
16431b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16631b70465SJoe Wallwork }
16731b70465SJoe Wallwork 
168cb7bfe8cSJoe Wallwork /*@
16931b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
17031b70465SJoe Wallwork 
17131b70465SJoe Wallwork   Input parameters:
172*2fe279fdSBarry Smith . dm                      - The `DM`
17331b70465SJoe Wallwork 
17431b70465SJoe Wallwork   Output parameters:
17531b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
17631b70465SJoe Wallwork 
17731b70465SJoe Wallwork   Level: beginner
17831b70465SJoe Wallwork 
179*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
180cb7bfe8cSJoe Wallwork @*/
181d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
182d71ae5a4SJacob Faibussowitsch {
18331b70465SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
18431b70465SJoe Wallwork 
18531b70465SJoe Wallwork   PetscFunctionBegin;
186bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
18731b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18931b70465SJoe Wallwork }
19031b70465SJoe Wallwork 
191cb7bfe8cSJoe Wallwork /*@
192cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
193cc2eee55SJoe Wallwork 
194cc2eee55SJoe Wallwork   Input parameters:
195*2fe279fdSBarry Smith + dm       - The `DM`
196cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
197cc2eee55SJoe Wallwork 
198cc2eee55SJoe Wallwork   Level: beginner
199cc2eee55SJoe Wallwork 
200*2fe279fdSBarry Smith   Note:
201cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
202cb7bfe8cSJoe Wallwork 
203*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
204cb7bfe8cSJoe Wallwork @*/
205d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
206d71ae5a4SJacob Faibussowitsch {
207cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
208cc2eee55SJoe Wallwork 
209cc2eee55SJoe Wallwork   PetscFunctionBegin;
210bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
211cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
213cc2eee55SJoe Wallwork }
214cc2eee55SJoe Wallwork 
215cb7bfe8cSJoe Wallwork /*@
216cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
217cc2eee55SJoe Wallwork 
218cc2eee55SJoe Wallwork   Input parameters:
219*2fe279fdSBarry Smith . dm       - The `DM`
220cc2eee55SJoe Wallwork 
221cc2eee55SJoe Wallwork   Output parameters:
222cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
223cc2eee55SJoe Wallwork 
224cc2eee55SJoe Wallwork   Level: beginner
225cc2eee55SJoe Wallwork 
226*2fe279fdSBarry Smith   Note:
227cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
228cb7bfe8cSJoe Wallwork 
229*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
230cb7bfe8cSJoe Wallwork @*/
231d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
232d71ae5a4SJacob Faibussowitsch {
233cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
234cc2eee55SJoe Wallwork 
235cc2eee55SJoe Wallwork   PetscFunctionBegin;
236bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
237cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
2383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
239cc2eee55SJoe Wallwork }
240cc2eee55SJoe Wallwork 
241cb7bfe8cSJoe Wallwork /*@
242cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
243cc2eee55SJoe Wallwork 
244cc2eee55SJoe Wallwork   Input parameters:
245*2fe279fdSBarry Smith + dm     - The `DM`
246cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
247cc2eee55SJoe Wallwork 
248cc2eee55SJoe Wallwork   Level: beginner
249cc2eee55SJoe Wallwork 
250*2fe279fdSBarry Smith   Note:
251cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
252cb7bfe8cSJoe Wallwork 
253*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
254cb7bfe8cSJoe Wallwork @*/
255d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
256d71ae5a4SJacob Faibussowitsch {
257cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
258cc2eee55SJoe Wallwork 
259cc2eee55SJoe Wallwork   PetscFunctionBegin;
260bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
261cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263cc2eee55SJoe Wallwork }
264cc2eee55SJoe Wallwork 
265cb7bfe8cSJoe Wallwork /*@
266cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
267cc2eee55SJoe Wallwork 
268cc2eee55SJoe Wallwork   Input parameters:
269*2fe279fdSBarry Smith . dm     - The `DM`
270cc2eee55SJoe Wallwork 
271cc2eee55SJoe Wallwork   Output parameters:
272cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
273cc2eee55SJoe Wallwork 
274cc2eee55SJoe Wallwork   Level: beginner
275cc2eee55SJoe Wallwork 
276*2fe279fdSBarry Smith   Note:
277cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
278cb7bfe8cSJoe Wallwork 
279*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
280cb7bfe8cSJoe Wallwork @*/
281d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
282d71ae5a4SJacob Faibussowitsch {
283cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
284cc2eee55SJoe Wallwork 
285cc2eee55SJoe Wallwork   PetscFunctionBegin;
286bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
287cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289cc2eee55SJoe Wallwork }
290cc2eee55SJoe Wallwork 
291cb7bfe8cSJoe Wallwork /*@
292cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
293cc2eee55SJoe Wallwork 
294cc2eee55SJoe Wallwork   Input parameters:
295*2fe279fdSBarry Smith + dm     - The `DM`
296cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
297cc2eee55SJoe Wallwork 
298cc2eee55SJoe Wallwork   Level: beginner
299cc2eee55SJoe Wallwork 
300*2fe279fdSBarry Smith   Note:
301cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
302cb7bfe8cSJoe Wallwork 
303*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()`
304cb7bfe8cSJoe Wallwork @*/
305d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
306d71ae5a4SJacob Faibussowitsch {
307cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
308cc2eee55SJoe Wallwork 
309cc2eee55SJoe Wallwork   PetscFunctionBegin;
310bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
311cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
313cc2eee55SJoe Wallwork }
314cc2eee55SJoe Wallwork 
315cb7bfe8cSJoe Wallwork /*@
316cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
317cc2eee55SJoe Wallwork 
318cc2eee55SJoe Wallwork   Input parameters:
319*2fe279fdSBarry Smith . dm     - The `DM`
320cc2eee55SJoe Wallwork 
321cc2eee55SJoe Wallwork   Output parameters:
322cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
323cc2eee55SJoe Wallwork 
324cc2eee55SJoe Wallwork   Level: beginner
325cc2eee55SJoe Wallwork 
326*2fe279fdSBarry Smith   Note:
327cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
328cb7bfe8cSJoe Wallwork 
329*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()`
330cb7bfe8cSJoe Wallwork @*/
331d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
332d71ae5a4SJacob Faibussowitsch {
333cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
334cc2eee55SJoe Wallwork 
335cc2eee55SJoe Wallwork   PetscFunctionBegin;
336bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
337cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339cc2eee55SJoe Wallwork }
340cc2eee55SJoe Wallwork 
341cb7bfe8cSJoe Wallwork /*@
34276f360caSJoe Wallwork   DMPlexMetricSetNoSurf - Should surface modification be turned off?
34376f360caSJoe Wallwork 
34476f360caSJoe Wallwork   Input parameters:
345*2fe279fdSBarry Smith + dm     - The `DM`
34676f360caSJoe Wallwork - noSurf - Should surface modification be turned off?
34776f360caSJoe Wallwork 
34876f360caSJoe Wallwork   Level: beginner
34976f360caSJoe Wallwork 
350*2fe279fdSBarry Smith   Note:
35176f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
35276f360caSJoe Wallwork 
353*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`
35476f360caSJoe Wallwork @*/
355d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
356d71ae5a4SJacob Faibussowitsch {
35776f360caSJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
35876f360caSJoe Wallwork 
35976f360caSJoe Wallwork   PetscFunctionBegin;
360bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
36176f360caSJoe Wallwork   plex->metricCtx->noSurf = noSurf;
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36376f360caSJoe Wallwork }
36476f360caSJoe Wallwork 
36576f360caSJoe Wallwork /*@
36676f360caSJoe Wallwork   DMPlexMetricNoSurf - Is surface modification turned off?
36776f360caSJoe Wallwork 
36876f360caSJoe Wallwork   Input parameters:
369*2fe279fdSBarry Smith . dm     - The `DM`
37076f360caSJoe Wallwork 
37176f360caSJoe Wallwork   Output parameters:
37276f360caSJoe Wallwork . noSurf - Is surface modification turned off?
37376f360caSJoe Wallwork 
37476f360caSJoe Wallwork   Level: beginner
37576f360caSJoe Wallwork 
376*2fe279fdSBarry Smith   Note:
37776f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
37876f360caSJoe Wallwork 
379*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`
38076f360caSJoe Wallwork @*/
381d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
382d71ae5a4SJacob Faibussowitsch {
38376f360caSJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
38476f360caSJoe Wallwork 
38576f360caSJoe Wallwork   PetscFunctionBegin;
386bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
38776f360caSJoe Wallwork   *noSurf = plex->metricCtx->noSurf;
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38976f360caSJoe Wallwork }
39076f360caSJoe Wallwork 
39176f360caSJoe Wallwork /*@
39231b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
39331b70465SJoe Wallwork 
39431b70465SJoe Wallwork   Input parameters:
395*2fe279fdSBarry Smith + dm    - The `DM`
39631b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
39731b70465SJoe Wallwork 
39831b70465SJoe Wallwork   Level: beginner
39931b70465SJoe Wallwork 
400*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()`
401cb7bfe8cSJoe Wallwork @*/
402d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
403d71ae5a4SJacob Faibussowitsch {
40431b70465SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
40531b70465SJoe Wallwork 
40631b70465SJoe Wallwork   PetscFunctionBegin;
407bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
408bc00d7afSJoe Wallwork   PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
40931b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
4103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41131b70465SJoe Wallwork }
41231b70465SJoe Wallwork 
413cb7bfe8cSJoe Wallwork /*@
41431b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
41531b70465SJoe Wallwork 
41631b70465SJoe Wallwork   Input parameters:
417*2fe279fdSBarry Smith . dm    - The `DM`
41831b70465SJoe Wallwork 
419cc2eee55SJoe Wallwork   Output parameters:
42031b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
42131b70465SJoe Wallwork 
42231b70465SJoe Wallwork   Level: beginner
42331b70465SJoe Wallwork 
424*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()`
425cb7bfe8cSJoe Wallwork @*/
426d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
427d71ae5a4SJacob Faibussowitsch {
42831b70465SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
42931b70465SJoe Wallwork 
43031b70465SJoe Wallwork   PetscFunctionBegin;
431bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
43231b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
4333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43431b70465SJoe Wallwork }
43531b70465SJoe Wallwork 
436cb7bfe8cSJoe Wallwork /*@
43731b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
43831b70465SJoe Wallwork 
43931b70465SJoe Wallwork   Input parameters:
440*2fe279fdSBarry Smith + dm    - The `DM`
44131b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44231b70465SJoe Wallwork 
44331b70465SJoe Wallwork   Level: beginner
44431b70465SJoe Wallwork 
445*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()`
446cb7bfe8cSJoe Wallwork @*/
447d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
448d71ae5a4SJacob Faibussowitsch {
44931b70465SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
45031b70465SJoe Wallwork 
45131b70465SJoe Wallwork   PetscFunctionBegin;
452bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
453bc00d7afSJoe Wallwork   PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
45431b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45631b70465SJoe Wallwork }
45731b70465SJoe Wallwork 
458cb7bfe8cSJoe Wallwork /*@
45931b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
46031b70465SJoe Wallwork 
46131b70465SJoe Wallwork   Input parameters:
462*2fe279fdSBarry Smith . dm    - The `DM`
46331b70465SJoe Wallwork 
464cc2eee55SJoe Wallwork   Output parameters:
46531b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
46631b70465SJoe Wallwork 
46731b70465SJoe Wallwork   Level: beginner
46831b70465SJoe Wallwork 
469*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()`
470cb7bfe8cSJoe Wallwork @*/
471d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
472d71ae5a4SJacob Faibussowitsch {
47331b70465SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
47431b70465SJoe Wallwork 
47531b70465SJoe Wallwork   PetscFunctionBegin;
476bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
47731b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
4783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47931b70465SJoe Wallwork }
48031b70465SJoe Wallwork 
481cb7bfe8cSJoe Wallwork /*@
48231b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
48331b70465SJoe Wallwork 
48431b70465SJoe Wallwork   Input parameters:
485*2fe279fdSBarry Smith + dm    - The `DM`
48631b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
48731b70465SJoe Wallwork 
48831b70465SJoe Wallwork   Level: beginner
48931b70465SJoe Wallwork 
490*2fe279fdSBarry Smith   Note:
491*2fe279fdSBarry Smith   If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
49231b70465SJoe Wallwork 
493*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
5033ba16761SJacob 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:
510*2fe279fdSBarry Smith . 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 
517*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52731b70465SJoe Wallwork }
52831b70465SJoe Wallwork 
529cb7bfe8cSJoe Wallwork /*@
53031b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
53131b70465SJoe Wallwork 
53231b70465SJoe Wallwork   Input parameters:
533*2fe279fdSBarry Smith + dm               - The `DM`
53431b70465SJoe Wallwork - targetComplexity - The target metric complexity
53531b70465SJoe Wallwork 
53631b70465SJoe Wallwork   Level: beginner
53731b70465SJoe Wallwork 
538*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54931b70465SJoe Wallwork }
55031b70465SJoe Wallwork 
551cb7bfe8cSJoe Wallwork /*@
55231b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
55331b70465SJoe Wallwork 
55431b70465SJoe Wallwork   Input parameters:
555*2fe279fdSBarry Smith . dm               - The `DM`
55631b70465SJoe Wallwork 
557cc2eee55SJoe Wallwork   Output parameters:
55831b70465SJoe Wallwork . targetComplexity - The target metric complexity
55931b70465SJoe Wallwork 
56031b70465SJoe Wallwork   Level: beginner
56131b70465SJoe Wallwork 
562*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
5713ba16761SJacob 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:
578*2fe279fdSBarry Smith + dm - The `DM`
57931b70465SJoe Wallwork - p  - The normalization order
58031b70465SJoe Wallwork 
58131b70465SJoe Wallwork   Level: beginner
58231b70465SJoe Wallwork 
583*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
5933ba16761SJacob 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:
600*2fe279fdSBarry Smith . dm - The `DM`
60131b70465SJoe Wallwork 
602cc2eee55SJoe Wallwork   Output parameters:
60331b70465SJoe Wallwork . p - The normalization order
60431b70465SJoe Wallwork 
60531b70465SJoe Wallwork   Level: beginner
60631b70465SJoe Wallwork 
607*2fe279fdSBarry Smith .seealso: `DMPLEX`, `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;
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61731b70465SJoe Wallwork }
61831b70465SJoe Wallwork 
619cb7bfe8cSJoe Wallwork /*@
620cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
621cc2eee55SJoe Wallwork 
622cc2eee55SJoe Wallwork   Input parameters:
623*2fe279fdSBarry Smith + dm   - The `DM`
624cc2eee55SJoe Wallwork - beta - The metric gradation factor
625cc2eee55SJoe Wallwork 
626cc2eee55SJoe Wallwork   Level: beginner
627cc2eee55SJoe Wallwork 
628cc2eee55SJoe Wallwork   Notes:
629cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
630cc2eee55SJoe Wallwork 
631cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
632cc2eee55SJoe Wallwork 
633cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
634cb7bfe8cSJoe Wallwork 
635*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
636cb7bfe8cSJoe Wallwork @*/
637d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
638d71ae5a4SJacob Faibussowitsch {
639cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
640cc2eee55SJoe Wallwork 
641cc2eee55SJoe Wallwork   PetscFunctionBegin;
642bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
643bc00d7afSJoe Wallwork   PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)");
644cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
646cc2eee55SJoe Wallwork }
647cc2eee55SJoe Wallwork 
648cb7bfe8cSJoe Wallwork /*@
649cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
650cc2eee55SJoe Wallwork 
651cc2eee55SJoe Wallwork   Input parameters:
652*2fe279fdSBarry Smith . dm   - The `DM`
653cc2eee55SJoe Wallwork 
654cc2eee55SJoe Wallwork   Output parameters:
655cc2eee55SJoe Wallwork . beta - The metric gradation factor
656cc2eee55SJoe Wallwork 
657cc2eee55SJoe Wallwork   Level: beginner
658cc2eee55SJoe Wallwork 
659cb7bfe8cSJoe Wallwork   Notes:
660cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
661cb7bfe8cSJoe Wallwork 
662cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
663cb7bfe8cSJoe Wallwork 
664cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
665cc2eee55SJoe Wallwork 
666*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
667cb7bfe8cSJoe Wallwork @*/
668d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
669d71ae5a4SJacob Faibussowitsch {
670cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
671cc2eee55SJoe Wallwork 
672cc2eee55SJoe Wallwork   PetscFunctionBegin;
673bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
674cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
676cc2eee55SJoe Wallwork }
677cc2eee55SJoe Wallwork 
678cb7bfe8cSJoe Wallwork /*@
679ae8b063eSJoe Wallwork   DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
680ae8b063eSJoe Wallwork 
681ae8b063eSJoe Wallwork   Input parameters:
682*2fe279fdSBarry Smith + dm    - The `DM`
683ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number
684ae8b063eSJoe Wallwork 
685ae8b063eSJoe Wallwork   Level: beginner
686ae8b063eSJoe Wallwork 
687ae8b063eSJoe Wallwork   Notes:
688ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
689ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
690ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
691ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
692ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
693ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
694ae8b063eSJoe Wallwork 
695ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
696ae8b063eSJoe Wallwork 
697*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
698ae8b063eSJoe Wallwork @*/
699d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
700d71ae5a4SJacob Faibussowitsch {
701ae8b063eSJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
702ae8b063eSJoe Wallwork 
703ae8b063eSJoe Wallwork   PetscFunctionBegin;
704bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
705bc00d7afSJoe Wallwork   PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)");
706ae8b063eSJoe Wallwork   plex->metricCtx->hausdorffNumber = hausd;
7073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
708ae8b063eSJoe Wallwork }
709ae8b063eSJoe Wallwork 
710ae8b063eSJoe Wallwork /*@
711ae8b063eSJoe Wallwork   DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
712ae8b063eSJoe Wallwork 
713ae8b063eSJoe Wallwork   Input parameters:
714*2fe279fdSBarry Smith . dm    - The `DM`
715ae8b063eSJoe Wallwork 
716ae8b063eSJoe Wallwork   Output parameters:
717ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number
718ae8b063eSJoe Wallwork 
719ae8b063eSJoe Wallwork   Level: beginner
720ae8b063eSJoe Wallwork 
721ae8b063eSJoe Wallwork   Notes:
722ae8b063eSJoe Wallwork   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
723ae8b063eSJoe Wallwork   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
724ae8b063eSJoe Wallwork   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
725ae8b063eSJoe Wallwork   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
726ae8b063eSJoe Wallwork   (resp. increase) the Hausdorff parameter. (Taken from
727ae8b063eSJoe Wallwork   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
728ae8b063eSJoe Wallwork 
729ae8b063eSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
730ae8b063eSJoe Wallwork 
731*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
732ae8b063eSJoe Wallwork @*/
733d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
734d71ae5a4SJacob Faibussowitsch {
735ae8b063eSJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
736ae8b063eSJoe Wallwork 
737ae8b063eSJoe Wallwork   PetscFunctionBegin;
738bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
739ae8b063eSJoe Wallwork   *hausd = plex->metricCtx->hausdorffNumber;
7403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
741ae8b063eSJoe Wallwork }
742ae8b063eSJoe Wallwork 
743ae8b063eSJoe Wallwork /*@
744cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
745cc2eee55SJoe Wallwork 
746cc2eee55SJoe Wallwork   Input parameters:
747*2fe279fdSBarry Smith + dm        - The `DM`
748cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
749cc2eee55SJoe Wallwork 
750cb7bfe8cSJoe Wallwork   Level: beginner
751cb7bfe8cSJoe Wallwork 
752*2fe279fdSBarry Smith   Note:
753cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
754cb7bfe8cSJoe Wallwork 
755*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()`
756cb7bfe8cSJoe Wallwork @*/
757d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
758d71ae5a4SJacob Faibussowitsch {
759cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
760cc2eee55SJoe Wallwork 
761cc2eee55SJoe Wallwork   PetscFunctionBegin;
762bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
763cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
7643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
765cc2eee55SJoe Wallwork }
766cc2eee55SJoe Wallwork 
767cb7bfe8cSJoe Wallwork /*@
768cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
769cc2eee55SJoe Wallwork 
770cc2eee55SJoe Wallwork   Input parameters:
771*2fe279fdSBarry Smith . dm        - The `DM`
772cc2eee55SJoe Wallwork 
773cc2eee55SJoe Wallwork   Output parameters:
774cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
775cc2eee55SJoe Wallwork 
776cb7bfe8cSJoe Wallwork   Level: beginner
777cb7bfe8cSJoe Wallwork 
778*2fe279fdSBarry Smith   Note:
779cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
780cb7bfe8cSJoe Wallwork 
781*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
782cb7bfe8cSJoe Wallwork @*/
783d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
784d71ae5a4SJacob Faibussowitsch {
785cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
786cc2eee55SJoe Wallwork 
787cc2eee55SJoe Wallwork   PetscFunctionBegin;
788bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
789cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
7903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
791cc2eee55SJoe Wallwork }
792cc2eee55SJoe Wallwork 
793cb7bfe8cSJoe Wallwork /*@
794cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
795cc2eee55SJoe Wallwork 
796cc2eee55SJoe Wallwork   Input parameters:
797*2fe279fdSBarry Smith + dm      - The `DM`
798cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
799cc2eee55SJoe Wallwork 
800cb7bfe8cSJoe Wallwork   Level: beginner
801cb7bfe8cSJoe Wallwork 
802*2fe279fdSBarry Smith   Note:
803cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
804cc2eee55SJoe Wallwork 
805*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
806cb7bfe8cSJoe Wallwork @*/
807d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
808d71ae5a4SJacob Faibussowitsch {
809cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
810cc2eee55SJoe Wallwork 
811cc2eee55SJoe Wallwork   PetscFunctionBegin;
812bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
813cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
8143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
815cc2eee55SJoe Wallwork }
816cc2eee55SJoe Wallwork 
817cb7bfe8cSJoe Wallwork /*@
818cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
819cc2eee55SJoe Wallwork 
820cc2eee55SJoe Wallwork   Input parameters:
821*2fe279fdSBarry Smith . dm      - The `DM`
822cc2eee55SJoe Wallwork 
823cc2eee55SJoe Wallwork   Output parameters:
824cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
825cc2eee55SJoe Wallwork 
826cb7bfe8cSJoe Wallwork   Level: beginner
827cb7bfe8cSJoe Wallwork 
828*2fe279fdSBarry Smith   Note:
829cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
830cc2eee55SJoe Wallwork 
831*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()`
832cb7bfe8cSJoe Wallwork @*/
833d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
834d71ae5a4SJacob Faibussowitsch {
835cc2eee55SJoe Wallwork   DM_Plex *plex = (DM_Plex *)dm->data;
836cc2eee55SJoe Wallwork 
837cc2eee55SJoe Wallwork   PetscFunctionBegin;
838bc00d7afSJoe Wallwork   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
839cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
8403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
841cc2eee55SJoe Wallwork }
842cc2eee55SJoe Wallwork 
843d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
844d71ae5a4SJacob Faibussowitsch {
84520282da2SJoe Wallwork   MPI_Comm comm;
84620282da2SJoe Wallwork   PetscFE  fe;
84720282da2SJoe Wallwork   PetscInt dim;
84820282da2SJoe Wallwork 
84920282da2SJoe Wallwork   PetscFunctionBegin;
85020282da2SJoe Wallwork 
85120282da2SJoe Wallwork   /* Extract metadata from dm */
8529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
8539566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
85420282da2SJoe Wallwork 
85520282da2SJoe Wallwork   /* Create a P1 field of the requested size */
8569566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
8579566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
8589566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
8599566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
8609566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm, metric));
86120282da2SJoe Wallwork 
8623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86320282da2SJoe Wallwork }
86420282da2SJoe Wallwork 
865cb7bfe8cSJoe Wallwork /*@
86620282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
86720282da2SJoe Wallwork 
86820282da2SJoe Wallwork   Input parameters:
869*2fe279fdSBarry Smith + dm     - The `DM`
87020282da2SJoe Wallwork - f      - The field number to use
87120282da2SJoe Wallwork 
87220282da2SJoe Wallwork   Output parameter:
87320282da2SJoe Wallwork . metric - The metric
87420282da2SJoe Wallwork 
875*2fe279fdSBarry Smith   Options Database Key:
876*2fe279fdSBarry Smith . -dm_adaptor <pragmatic/mmg/parmmg>        - specify dm adaptor to use
87720282da2SJoe Wallwork 
878*2fe279fdSBarry Smith   Options Database Keys for Mmg and ParMmg:
879*2fe279fdSBarry Smith + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
880*2fe279fdSBarry Smith . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
881*2fe279fdSBarry Smith . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
882*2fe279fdSBarry Smith . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
883*2fe279fdSBarry Smith . -dm_plex_metric_no_move                   - Should node movement be turned off?
884*2fe279fdSBarry Smith - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
88531b70465SJoe Wallwork 
886*2fe279fdSBarry Smith   Options Database Keys for Riemannian metrics:
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 
896*2fe279fdSBarry Smith   Level: beginner
897cb7bfe8cSJoe Wallwork 
898*2fe279fdSBarry Smith   Note:
899*2fe279fdSBarry Smith   It is assumed that the `DM` is comprised of simplices.
900cb7bfe8cSJoe Wallwork 
901*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`
902cb7bfe8cSJoe Wallwork @*/
903d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
904d71ae5a4SJacob Faibussowitsch {
90593520041SJoe Wallwork   PetscBool isotropic, uniform;
90620282da2SJoe Wallwork   PetscInt  coordDim, Nd;
90720282da2SJoe Wallwork 
90820282da2SJoe Wallwork   PetscFunctionBegin;
9099566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
91020282da2SJoe Wallwork   Nd = coordDim * coordDim;
9119566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
9129566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
91393520041SJoe Wallwork   if (uniform) {
91493520041SJoe Wallwork     MPI_Comm comm;
91593520041SJoe Wallwork 
9169566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
917bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
9189566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, metric));
9199566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
9209566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(*metric));
921f6db075bSJoe Wallwork   } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
922f6db075bSJoe Wallwork   else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
9233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92420282da2SJoe Wallwork }
92520282da2SJoe Wallwork 
926cb7bfe8cSJoe Wallwork /*@
92720282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
92820282da2SJoe Wallwork 
92920282da2SJoe Wallwork   Input parameters:
930*2fe279fdSBarry Smith + dm     - The `DM`
93120282da2SJoe Wallwork . f      - The field number to use
93220282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
93320282da2SJoe Wallwork 
93420282da2SJoe Wallwork   Output parameter:
93520282da2SJoe Wallwork . metric - The uniform metric
93620282da2SJoe Wallwork 
93720282da2SJoe Wallwork   Level: beginner
93820282da2SJoe Wallwork 
939*2fe279fdSBarry Smith   Note:
940*2fe279fdSBarry Smith   In this case, the metric is constant in space and so is only specified for a single vertex.
94120282da2SJoe Wallwork 
942*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
943cb7bfe8cSJoe Wallwork @*/
944d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
945d71ae5a4SJacob Faibussowitsch {
94620282da2SJoe Wallwork   PetscFunctionBegin;
9479566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
9489566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9492c71b3e2SJacob Faibussowitsch   PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
950bc00d7afSJoe Wallwork   PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
9519566063dSJacob Faibussowitsch   PetscCall(VecSet(*metric, alpha));
9529566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*metric));
9539566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*metric));
9543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95520282da2SJoe Wallwork }
95620282da2SJoe Wallwork 
957d71ae5a4SJacob 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[])
958d71ae5a4SJacob Faibussowitsch {
95993520041SJoe Wallwork   f0[0] = u[0];
96093520041SJoe Wallwork }
96193520041SJoe Wallwork 
962cb7bfe8cSJoe Wallwork /*@
96320282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
96420282da2SJoe Wallwork 
96520282da2SJoe Wallwork   Input parameters:
966*2fe279fdSBarry Smith + dm        - The `DM`
96720282da2SJoe Wallwork . f         - The field number to use
96820282da2SJoe Wallwork - indicator - The error indicator
96920282da2SJoe Wallwork 
97020282da2SJoe Wallwork   Output parameter:
97120282da2SJoe Wallwork . metric    - The isotropic metric
97220282da2SJoe Wallwork 
97320282da2SJoe Wallwork   Level: beginner
97420282da2SJoe Wallwork 
97520282da2SJoe Wallwork   Notes:
976*2fe279fdSBarry Smith   It is assumed that the `DM` is comprised of simplices.
97720282da2SJoe Wallwork 
97893520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
97920282da2SJoe Wallwork 
980*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
981cb7bfe8cSJoe Wallwork @*/
982d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
983d71ae5a4SJacob Faibussowitsch {
98493520041SJoe Wallwork   PetscInt m, n;
98520282da2SJoe Wallwork 
98620282da2SJoe Wallwork   PetscFunctionBegin;
9879566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
9889566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9899566063dSJacob Faibussowitsch   PetscCall(VecGetSize(indicator, &m));
9909566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metric, &n));
9919566063dSJacob Faibussowitsch   if (m == n) PetscCall(VecCopy(indicator, *metric));
99293520041SJoe Wallwork   else {
99393520041SJoe Wallwork     DM dmIndi;
9949371c9d4SSatish 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[]);
99593520041SJoe Wallwork 
9969566063dSJacob Faibussowitsch     PetscCall(VecGetDM(indicator, &dmIndi));
99793520041SJoe Wallwork     funcs[0] = identity;
9989566063dSJacob Faibussowitsch     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
99920282da2SJoe Wallwork   }
10003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100120282da2SJoe Wallwork }
100220282da2SJoe Wallwork 
1003f6db075bSJoe Wallwork /*@
1004f6db075bSJoe Wallwork   DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1005f6db075bSJoe Wallwork 
1006f6db075bSJoe Wallwork   Input parameters:
1007*2fe279fdSBarry Smith + dm          - The `DM` of the metric field
1008f6db075bSJoe Wallwork - f           - The field number to use
1009f6db075bSJoe Wallwork 
1010f6db075bSJoe Wallwork   Output parameter:
1011f6db075bSJoe Wallwork + determinant - The determinant field
1012*2fe279fdSBarry Smith - dmDet       - The corresponding `DM`
1013f6db075bSJoe Wallwork 
1014f6db075bSJoe Wallwork   Level: beginner
1015f6db075bSJoe Wallwork 
1016*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`, `DMPlexMetricCreate()`
1017f6db075bSJoe Wallwork @*/
1018d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1019d71ae5a4SJacob Faibussowitsch {
1020f6db075bSJoe Wallwork   PetscBool uniform;
1021f6db075bSJoe Wallwork 
1022f6db075bSJoe Wallwork   PetscFunctionBegin;
1023f6db075bSJoe Wallwork   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1024f6db075bSJoe Wallwork   PetscCall(DMClone(dm, dmDet));
1025f6db075bSJoe Wallwork   if (uniform) {
1026f6db075bSJoe Wallwork     MPI_Comm comm;
1027f6db075bSJoe Wallwork 
1028f6db075bSJoe Wallwork     PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm));
1029f6db075bSJoe Wallwork     PetscCall(VecCreate(comm, determinant));
1030f6db075bSJoe Wallwork     PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1031f6db075bSJoe Wallwork     PetscCall(VecSetFromOptions(*determinant));
1032f6db075bSJoe Wallwork   } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
10333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1034f6db075bSJoe Wallwork }
1035f6db075bSJoe Wallwork 
1036d71ae5a4SJacob Faibussowitsch static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1037d71ae5a4SJacob Faibussowitsch {
103882490253SJoe Wallwork   PetscInt i, j;
103982490253SJoe Wallwork 
104082490253SJoe Wallwork   PetscFunctionBegin;
10413ba16761SJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"));
104282490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
10433ba16761SJacob Faibussowitsch     if (i == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    [["));
10443ba16761SJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "     ["));
104582490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
10463ba16761SJacob Faibussowitsch       if (j < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])));
10473ba16761SJacob Faibussowitsch       else PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])));
104882490253SJoe Wallwork     }
10493ba16761SJacob Faibussowitsch     if (i < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
10503ba16761SJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "]]\n"));
105182490253SJoe Wallwork   }
10523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105382490253SJoe Wallwork }
105482490253SJoe Wallwork 
1055d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1056d71ae5a4SJacob Faibussowitsch {
105720282da2SJoe Wallwork   PetscInt     i, j, k;
105820282da2SJoe 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);
105920282da2SJoe Wallwork   PetscScalar *Mpos;
106020282da2SJoe Wallwork 
106120282da2SJoe Wallwork   PetscFunctionBegin;
10629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs));
106320282da2SJoe Wallwork 
106420282da2SJoe Wallwork   /* Symmetrize */
106520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
106620282da2SJoe Wallwork     Mpos[i * dim + i] = Mp[i * dim + i];
106720282da2SJoe Wallwork     for (j = i + 1; j < dim; ++j) {
106820282da2SJoe Wallwork       Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
106920282da2SJoe Wallwork       Mpos[j * dim + i] = Mpos[i * dim + j];
107020282da2SJoe Wallwork     }
107120282da2SJoe Wallwork   }
107220282da2SJoe Wallwork 
107320282da2SJoe Wallwork   /* Compute eigendecomposition */
107493520041SJoe Wallwork   if (dim == 1) {
107593520041SJoe Wallwork     /* Isotropic case */
107693520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
107793520041SJoe Wallwork     Mpos[0] = 1.0;
107893520041SJoe Wallwork   } else {
107993520041SJoe Wallwork     /* Anisotropic case */
108020282da2SJoe Wallwork     PetscScalar *work;
108120282da2SJoe Wallwork     PetscBLASInt lwork;
108220282da2SJoe Wallwork 
108320282da2SJoe Wallwork     lwork = 5 * dim;
10849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * dim, &work));
108520282da2SJoe Wallwork     {
108620282da2SJoe Wallwork       PetscBLASInt lierr;
108720282da2SJoe Wallwork       PetscBLASInt nb;
108820282da2SJoe Wallwork 
10899566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
10909566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
109120282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
109220282da2SJoe Wallwork       {
109320282da2SJoe Wallwork         PetscReal *rwork;
10949566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1095792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr));
10969566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
109720282da2SJoe Wallwork       }
109820282da2SJoe Wallwork #else
1099792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr));
110020282da2SJoe Wallwork #endif
110182490253SJoe Wallwork       if (lierr) {
110282490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
110382490253SJoe Wallwork           Mpos[i * dim + i] = Mp[i * dim + i];
110482490253SJoe Wallwork           for (j = i + 1; j < dim; ++j) {
110582490253SJoe Wallwork             Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
110682490253SJoe Wallwork             Mpos[j * dim + i] = Mpos[i * dim + j];
110782490253SJoe Wallwork           }
110882490253SJoe Wallwork         }
11093ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, Mpos));
111098921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
111182490253SJoe Wallwork       }
11129566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
111320282da2SJoe Wallwork     }
11149566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
111520282da2SJoe Wallwork   }
111620282da2SJoe Wallwork 
111720282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
111820282da2SJoe Wallwork   max_eig = 0.0;
111920282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
112020282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
112120282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
112220282da2SJoe Wallwork   }
112320282da2SJoe Wallwork 
11243f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11253f07679eSJoe Wallwork   *detMp = 1.0;
112620282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
112720282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min);
11283f07679eSJoe Wallwork     *detMp *= eigs[i];
112920282da2SJoe Wallwork   }
113020282da2SJoe Wallwork 
113120282da2SJoe Wallwork   /* Reconstruct Hessian */
113220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
113320282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
113420282da2SJoe Wallwork       Mp[i * dim + j] = 0.0;
1135ad540459SPierre Jolivet       for (k = 0; k < dim; ++k) Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j];
113620282da2SJoe Wallwork     }
113720282da2SJoe Wallwork   }
11389566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Mpos, eigs));
113920282da2SJoe Wallwork 
11403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114120282da2SJoe Wallwork }
114220282da2SJoe Wallwork 
1143cb7bfe8cSJoe Wallwork /*@
114420282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
114520282da2SJoe Wallwork 
114620282da2SJoe Wallwork   Input parameters:
1147*2fe279fdSBarry Smith + dm                 - The `DM`
11483f07679eSJoe Wallwork . metricIn           - The metric
114999abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11503f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
115120282da2SJoe Wallwork 
115220282da2SJoe Wallwork   Output parameter:
11533f07679eSJoe Wallwork + metricOut          - The metric
11543f07679eSJoe Wallwork - determinant        - Its determinant
115520282da2SJoe Wallwork 
1156*2fe279fdSBarry Smith   Options Database Keys:
115793520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
115893520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
115993520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1160cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1161cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1162cb7bfe8cSJoe Wallwork 
1163*2fe279fdSBarry Smith   Level: beginner
1164*2fe279fdSBarry Smith 
1165*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1166cb7bfe8cSJoe Wallwork @*/
1167d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1168d71ae5a4SJacob Faibussowitsch {
11693f07679eSJoe Wallwork   DM           dmDet;
117093520041SJoe Wallwork   PetscBool    isotropic, uniform;
117120282da2SJoe Wallwork   PetscInt     dim, vStart, vEnd, v;
11723f07679eSJoe Wallwork   PetscScalar *met, *det;
117320282da2SJoe Wallwork   PetscReal    h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
117420282da2SJoe Wallwork 
117520282da2SJoe Wallwork   PetscFunctionBegin;
11769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
117720282da2SJoe Wallwork 
117820282da2SJoe Wallwork   /* Extract metadata from dm */
11799566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
118020282da2SJoe Wallwork   if (restrictSizes) {
11819566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
11829566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
118399abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
118499abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1185bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
118699abec2bSJoe Wallwork   }
118799abec2bSJoe Wallwork   if (restrictAnisotropy) {
11889566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
118999abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
119020282da2SJoe Wallwork   }
119120282da2SJoe Wallwork 
119293520041SJoe Wallwork   /* Setup output metric */
11935241a91fSJoe Wallwork   PetscCall(VecCopy(metricIn, metricOut));
119493520041SJoe Wallwork 
119593520041SJoe Wallwork   /* Enforce SPD and extract determinant */
11965241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
11979566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
11989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
119993520041SJoe Wallwork   if (uniform) {
1200d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
120193520041SJoe Wallwork 
120293520041SJoe Wallwork     /* Uniform case */
12035241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
12049566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
12055241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
120693520041SJoe Wallwork   } else {
120793520041SJoe Wallwork     /* Spatially varying case */
120893520041SJoe Wallwork     PetscInt nrow;
120993520041SJoe Wallwork 
121093520041SJoe Wallwork     if (isotropic) nrow = 1;
121193520041SJoe Wallwork     else nrow = dim;
12125241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
12139566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12145241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
121520282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12163f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
12179566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
12189566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
12199566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
122020282da2SJoe Wallwork     }
12215241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
122293520041SJoe Wallwork   }
12235241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
1224fe902aceSJoe Wallwork 
12259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
12263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122720282da2SJoe Wallwork }
122820282da2SJoe Wallwork 
1229d71ae5a4SJacob 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[])
1230d71ae5a4SJacob Faibussowitsch {
123120282da2SJoe Wallwork   const PetscScalar p = constants[0];
123220282da2SJoe Wallwork 
1233985ad86eSJose E. Roman   f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim));
123420282da2SJoe Wallwork }
123520282da2SJoe Wallwork 
1236cb7bfe8cSJoe Wallwork /*@
123720282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
123820282da2SJoe Wallwork 
123920282da2SJoe Wallwork   Input parameters:
1240*2fe279fdSBarry Smith + dm                 - The `DM`
124120282da2SJoe Wallwork . metricIn           - The unnormalized metric
124216522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
124316522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
124420282da2SJoe Wallwork 
1245*2fe279fdSBarry Smith   Output parameters:
1246*2fe279fdSBarry Smith + metricOut          - The normalized metric
1247*2fe279fdSBarry Smith - determinant - computed determinant
124820282da2SJoe Wallwork 
1249*2fe279fdSBarry Smith   Options Database Keys:
125093520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
125193520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
125293520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1253cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1254cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1255cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1256cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1257cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1258cb7bfe8cSJoe Wallwork 
1259*2fe279fdSBarry Smith   Level: beginner
1260*2fe279fdSBarry Smith 
1261*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1262cb7bfe8cSJoe Wallwork @*/
1263d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1264d71ae5a4SJacob Faibussowitsch {
12653f07679eSJoe Wallwork   DM           dmDet;
126620282da2SJoe Wallwork   MPI_Comm     comm;
126793520041SJoe Wallwork   PetscBool    restrictAnisotropyFirst, isotropic, uniform;
126820282da2SJoe Wallwork   PetscDS      ds;
126920282da2SJoe Wallwork   PetscInt     dim, Nd, vStart, vEnd, v, i;
12703f07679eSJoe Wallwork   PetscScalar *met, *det, integral, constants[1];
127193520041SJoe Wallwork   PetscReal    p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
127220282da2SJoe Wallwork 
127320282da2SJoe Wallwork   PetscFunctionBegin;
12749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0));
127520282da2SJoe Wallwork 
127620282da2SJoe Wallwork   /* Extract metadata from dm */
12779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
12789566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
12799566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12809566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
128193520041SJoe Wallwork   if (isotropic) Nd = 1;
128293520041SJoe Wallwork   else Nd = dim * dim;
128320282da2SJoe Wallwork 
128420282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
12859566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
12865241a91fSJoe Wallwork   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
128720282da2SJoe Wallwork 
128820282da2SJoe Wallwork   /* Compute global normalization factor */
12899566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
12909566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
129116522936SJoe Wallwork   constants[0] = p;
129293520041SJoe Wallwork   if (uniform) {
1293bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
129493520041SJoe Wallwork     DM  dmTmp;
129593520041SJoe Wallwork     Vec tmp;
129693520041SJoe Wallwork 
12979566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmTmp));
12989566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
12999566063dSJacob Faibussowitsch     PetscCall(VecGetArray(determinant, &det));
13009566063dSJacob Faibussowitsch     PetscCall(VecSet(tmp, det[0]));
13019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(determinant, &det));
13029566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmTmp, &ds));
13039566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13049566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13059566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
13069566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tmp));
13079566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmTmp));
130893520041SJoe Wallwork   } else {
13099566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
13109566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmDet, &ds));
13119566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13129566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13139566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
131493520041SJoe Wallwork   }
13153f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
1316bc00d7afSJoe 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?");
13173f07679eSJoe Wallwork   factGlob = PetscPowReal(target / realIntegral, 2.0 / dim);
131820282da2SJoe Wallwork 
131920282da2SJoe Wallwork   /* Apply local scaling */
132020282da2SJoe Wallwork   if (restrictSizes) {
13219566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
13229566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
132316522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
132416522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1325bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
132616522936SJoe Wallwork   }
132716522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
13289566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
132916522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
133020282da2SJoe Wallwork   }
13315241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
13329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(determinant, &det));
133393520041SJoe Wallwork   if (uniform) {
133493520041SJoe Wallwork     /* Uniform case */
133593520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim));
13369566063dSJacob Faibussowitsch     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
133793520041SJoe Wallwork   } else {
133893520041SJoe Wallwork     /* Spatially varying case */
133993520041SJoe Wallwork     PetscInt nrow;
134093520041SJoe Wallwork 
134193520041SJoe Wallwork     if (isotropic) nrow = 1;
134293520041SJoe Wallwork     else nrow = dim;
13439566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
13445241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
134520282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13463f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
134720282da2SJoe Wallwork 
13489566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
13499566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
13503f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim));
135120282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
13529566063dSJacob Faibussowitsch       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
135393520041SJoe Wallwork     }
135420282da2SJoe Wallwork   }
13559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(determinant, &det));
13565241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
135720282da2SJoe Wallwork 
13589566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0));
13593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136020282da2SJoe Wallwork }
136120282da2SJoe Wallwork 
1362cb7bfe8cSJoe Wallwork /*@
136320282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
136420282da2SJoe Wallwork 
1365f1a722f8SMatthew G. Knepley   Input Parameters:
1366*2fe279fdSBarry Smith + dm         - The `DM`
136720282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
136820282da2SJoe Wallwork . weights    - Weights for the average
136920282da2SJoe Wallwork - metrics    - The metrics to be averaged
137020282da2SJoe Wallwork 
137120282da2SJoe Wallwork   Output Parameter:
137220282da2SJoe Wallwork . metricAvg  - The averaged metric
137320282da2SJoe Wallwork 
137420282da2SJoe Wallwork   Level: beginner
137520282da2SJoe Wallwork 
137620282da2SJoe Wallwork   Notes:
137720282da2SJoe Wallwork   The weights should sum to unity.
137820282da2SJoe Wallwork 
137920282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
138020282da2SJoe Wallwork 
1381*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1382cb7bfe8cSJoe Wallwork @*/
1383d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1384d71ae5a4SJacob Faibussowitsch {
138520282da2SJoe Wallwork   PetscBool haveWeights = PETSC_TRUE;
138693520041SJoe Wallwork   PetscInt  i, m, n;
138720282da2SJoe Wallwork   PetscReal sum = 0.0, tol = 1.0e-10;
138820282da2SJoe Wallwork 
138920282da2SJoe Wallwork   PetscFunctionBegin;
13909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0));
139163a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
13925241a91fSJoe Wallwork   PetscCall(VecSet(metricAvg, 0.0));
13935241a91fSJoe Wallwork   PetscCall(VecGetSize(metricAvg, &m));
139493520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
13959566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
13965f80ce2aSJacob Faibussowitsch     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
139793520041SJoe Wallwork   }
139820282da2SJoe Wallwork 
139920282da2SJoe Wallwork   /* Default to the unweighted case */
140020282da2SJoe Wallwork   if (!weights) {
14019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numMetrics, &weights));
140220282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
1403ad540459SPierre Jolivet     for (i = 0; i < numMetrics; ++i) weights[i] = 1.0 / numMetrics;
140420282da2SJoe Wallwork   }
140520282da2SJoe Wallwork 
140620282da2SJoe Wallwork   /* Check weights sum to unity */
140793520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
14085f80ce2aSJacob Faibussowitsch   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
140920282da2SJoe Wallwork 
141020282da2SJoe Wallwork   /* Compute metric average */
14115241a91fSJoe Wallwork   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
14129566063dSJacob Faibussowitsch   if (!haveWeights) PetscCall(PetscFree(weights));
1413fe902aceSJoe Wallwork 
14149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0));
14153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141620282da2SJoe Wallwork }
141720282da2SJoe Wallwork 
1418cb7bfe8cSJoe Wallwork /*@
141920282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
142020282da2SJoe Wallwork 
1421f1a722f8SMatthew G. Knepley   Input Parameters:
1422*2fe279fdSBarry Smith + dm         - The `DM`
142320282da2SJoe Wallwork . metric1    - The first metric to be averaged
142420282da2SJoe Wallwork - metric2    - The second metric to be averaged
142520282da2SJoe Wallwork 
142620282da2SJoe Wallwork   Output Parameter:
142720282da2SJoe Wallwork . metricAvg  - The averaged metric
142820282da2SJoe Wallwork 
142920282da2SJoe Wallwork   Level: beginner
143020282da2SJoe Wallwork 
1431*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1432cb7bfe8cSJoe Wallwork @*/
1433d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1434d71ae5a4SJacob Faibussowitsch {
143520282da2SJoe Wallwork   PetscReal weights[2] = {0.5, 0.5};
143620282da2SJoe Wallwork   Vec       metrics[2] = {metric1, metric2};
143720282da2SJoe Wallwork 
143820282da2SJoe Wallwork   PetscFunctionBegin;
14399566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
14403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144120282da2SJoe Wallwork }
144220282da2SJoe Wallwork 
1443cb7bfe8cSJoe Wallwork /*@
144420282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
144520282da2SJoe Wallwork 
1446f1a722f8SMatthew G. Knepley   Input Parameters:
1447*2fe279fdSBarry Smith + dm         - The `DM`
144820282da2SJoe Wallwork . metric1    - The first metric to be averaged
144920282da2SJoe Wallwork . metric2    - The second metric to be averaged
145020282da2SJoe Wallwork - metric3    - The third metric to be averaged
145120282da2SJoe Wallwork 
145220282da2SJoe Wallwork   Output Parameter:
145320282da2SJoe Wallwork . metricAvg  - The averaged metric
145420282da2SJoe Wallwork 
145520282da2SJoe Wallwork   Level: beginner
145620282da2SJoe Wallwork 
1457*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1458cb7bfe8cSJoe Wallwork @*/
1459d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1460d71ae5a4SJacob Faibussowitsch {
146120282da2SJoe Wallwork   PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0};
146220282da2SJoe Wallwork   Vec       metrics[3] = {metric1, metric2, metric3};
146320282da2SJoe Wallwork 
146420282da2SJoe Wallwork   PetscFunctionBegin;
14659566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
14663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146720282da2SJoe Wallwork }
146820282da2SJoe Wallwork 
1469d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1470d71ae5a4SJacob Faibussowitsch {
147120282da2SJoe Wallwork   PetscInt     i, j, k, l, m;
1472e2606525SJoe Wallwork   PetscReal   *evals;
147320282da2SJoe Wallwork   PetscScalar *evecs, *sqrtM1, *isqrtM1;
147420282da2SJoe Wallwork 
147520282da2SJoe Wallwork   PetscFunctionBegin;
147693520041SJoe Wallwork 
147793520041SJoe Wallwork   /* Isotropic case */
147893520041SJoe Wallwork   if (dim == 1) {
14796f71502aSJoe Wallwork     M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
14803ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
148193520041SJoe Wallwork   }
148293520041SJoe Wallwork 
148393520041SJoe Wallwork   /* Anisotropic case */
1484e2606525SJoe Wallwork   PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals));
148520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
1486ad540459SPierre Jolivet     for (j = 0; j < dim; ++j) evecs[i * dim + j] = M1[i * dim + j];
148720282da2SJoe Wallwork   }
148820282da2SJoe Wallwork   {
148920282da2SJoe Wallwork     PetscScalar *work;
149020282da2SJoe Wallwork     PetscBLASInt lwork;
149120282da2SJoe Wallwork 
149220282da2SJoe Wallwork     lwork = 5 * dim;
14939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * dim, &work));
149420282da2SJoe Wallwork     {
149520282da2SJoe Wallwork       PetscBLASInt lierr, nb;
1496e2606525SJoe Wallwork       PetscReal    sqrtj;
149720282da2SJoe Wallwork 
149820282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
14999566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
15009566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
150120282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
150220282da2SJoe Wallwork       {
150320282da2SJoe Wallwork         PetscReal *rwork;
15049566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1505792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15069566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
150720282da2SJoe Wallwork       }
150820282da2SJoe Wallwork #else
1509792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
151020282da2SJoe Wallwork #endif
151182490253SJoe Wallwork       if (lierr) {
15123ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, M1));
151398921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
151482490253SJoe Wallwork       }
15159566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
151620282da2SJoe Wallwork 
1517e2606525SJoe Wallwork       /* Compute square root and the reciprocal thereof */
151820282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
151920282da2SJoe Wallwork         for (k = 0; k < dim; ++k) {
1520e2606525SJoe Wallwork           sqrtM1[i * dim + k]  = 0.0;
1521e2606525SJoe Wallwork           isqrtM1[i * dim + k] = 0.0;
1522e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1523e2606525SJoe Wallwork             sqrtj = PetscSqrtReal(evals[j]);
1524e2606525SJoe Wallwork             sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k];
1525e2606525SJoe Wallwork             isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k];
152620282da2SJoe Wallwork           }
152720282da2SJoe Wallwork         }
152820282da2SJoe Wallwork       }
152920282da2SJoe Wallwork 
1530e2606525SJoe Wallwork       /* Map M2 into the space spanned by the eigenvectors of M1 */
153120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
153220282da2SJoe Wallwork         for (l = 0; l < dim; ++l) {
1533e2606525SJoe Wallwork           evecs[i * dim + l] = 0.0;
1534e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1535ad540459SPierre Jolivet             for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
153620282da2SJoe Wallwork           }
153720282da2SJoe Wallwork         }
153820282da2SJoe Wallwork       }
153920282da2SJoe Wallwork 
154020282da2SJoe Wallwork       /* Compute eigendecomposition */
15419566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
154220282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
154320282da2SJoe Wallwork       {
154420282da2SJoe Wallwork         PetscReal *rwork;
15459566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1546792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15479566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
154820282da2SJoe Wallwork       }
154920282da2SJoe Wallwork #else
1550792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
155120282da2SJoe Wallwork #endif
155282490253SJoe Wallwork       if (lierr) {
155382490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
155482490253SJoe Wallwork           for (l = 0; l < dim; ++l) {
1555e2606525SJoe Wallwork             evecs[i * dim + l] = 0.0;
1556e2606525SJoe Wallwork             for (j = 0; j < dim; ++j) {
1557ad540459SPierre Jolivet               for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
155882490253SJoe Wallwork             }
155982490253SJoe Wallwork           }
156082490253SJoe Wallwork         }
15613ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, evecs));
156298921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
156382490253SJoe Wallwork       }
15649566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
156520282da2SJoe Wallwork 
156620282da2SJoe Wallwork       /* Modify eigenvalues */
1567e2606525SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
156820282da2SJoe Wallwork 
156920282da2SJoe Wallwork       /* Map back to get the intersection */
157020282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
1571e2606525SJoe Wallwork         for (m = 0; m < dim; ++m) {
1572e2606525SJoe Wallwork           M2[i * dim + m] = 0.0;
157320282da2SJoe Wallwork           for (j = 0; j < dim; ++j) {
157420282da2SJoe Wallwork             for (k = 0; k < dim; ++k) {
1575ad540459SPierre 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];
157620282da2SJoe Wallwork             }
157720282da2SJoe Wallwork           }
157820282da2SJoe Wallwork         }
157920282da2SJoe Wallwork       }
158020282da2SJoe Wallwork     }
15819566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
158220282da2SJoe Wallwork   }
1583e2606525SJoe Wallwork   PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
15843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
158520282da2SJoe Wallwork }
158620282da2SJoe Wallwork 
1587cb7bfe8cSJoe Wallwork /*@
158820282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
158920282da2SJoe Wallwork 
1590f1a722f8SMatthew G. Knepley   Input Parameters:
1591*2fe279fdSBarry Smith + dm         - The `DM`
159220282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
159320282da2SJoe Wallwork - metrics    - The metrics to be intersected
159420282da2SJoe Wallwork 
159520282da2SJoe Wallwork   Output Parameter:
159620282da2SJoe Wallwork . metricInt  - The intersected metric
159720282da2SJoe Wallwork 
159820282da2SJoe Wallwork   Level: beginner
159920282da2SJoe Wallwork 
160020282da2SJoe Wallwork   Notes:
1601e2606525SJoe Wallwork   The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
160220282da2SJoe Wallwork 
1603e2606525SJoe Wallwork   The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
160420282da2SJoe Wallwork 
1605*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1606cb7bfe8cSJoe Wallwork @*/
1607d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1608d71ae5a4SJacob Faibussowitsch {
160993520041SJoe Wallwork   PetscBool    isotropic, uniform;
161093520041SJoe Wallwork   PetscInt     v, i, m, n;
161193520041SJoe Wallwork   PetscScalar *met, *meti;
161220282da2SJoe Wallwork 
161320282da2SJoe Wallwork   PetscFunctionBegin;
16149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0));
161563a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
161620282da2SJoe Wallwork 
161720282da2SJoe Wallwork   /* Copy over the first metric */
16185241a91fSJoe Wallwork   PetscCall(VecCopy(metrics[0], metricInt));
16193ba16761SJacob Faibussowitsch   if (numMetrics == 1) PetscFunctionReturn(PETSC_SUCCESS);
16205241a91fSJoe Wallwork   PetscCall(VecGetSize(metricInt, &m));
162193520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
16229566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
162308401ef6SPierre Jolivet     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
162493520041SJoe Wallwork   }
162520282da2SJoe Wallwork 
162620282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
16279566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
16289566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
162993520041SJoe Wallwork   if (uniform) {
163093520041SJoe Wallwork     /* Uniform case */
16315241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
163293520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16339566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
16349566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
16359566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
163693520041SJoe Wallwork     }
16375241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
163893520041SJoe Wallwork   } else {
163993520041SJoe Wallwork     /* Spatially varying case */
164093520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
164193520041SJoe Wallwork     PetscScalar *M, *Mi;
164293520041SJoe Wallwork 
16439566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
164493520041SJoe Wallwork     if (isotropic) nrow = 1;
164593520041SJoe Wallwork     else nrow = dim;
16469566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
16475241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
164820282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16499566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
165020282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
16519566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
16529566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
16539566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
165420282da2SJoe Wallwork       }
16559566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
165620282da2SJoe Wallwork     }
16575241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
165820282da2SJoe Wallwork   }
1659fe902aceSJoe Wallwork 
16609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0));
16613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166220282da2SJoe Wallwork }
166320282da2SJoe Wallwork 
1664cb7bfe8cSJoe Wallwork /*@
166520282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
166620282da2SJoe Wallwork 
1667f1a722f8SMatthew G. Knepley   Input Parameters:
1668*2fe279fdSBarry Smith + dm        - The `DM`
166920282da2SJoe Wallwork . metric1   - The first metric to be intersected
167020282da2SJoe Wallwork - metric2   - The second metric to be intersected
167120282da2SJoe Wallwork 
167220282da2SJoe Wallwork   Output Parameter:
167320282da2SJoe Wallwork . metricInt - The intersected metric
167420282da2SJoe Wallwork 
167520282da2SJoe Wallwork   Level: beginner
167620282da2SJoe Wallwork 
1677*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1678cb7bfe8cSJoe Wallwork @*/
1679d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1680d71ae5a4SJacob Faibussowitsch {
168120282da2SJoe Wallwork   Vec metrics[2] = {metric1, metric2};
168220282da2SJoe Wallwork 
168320282da2SJoe Wallwork   PetscFunctionBegin;
16849566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
16853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168620282da2SJoe Wallwork }
168720282da2SJoe Wallwork 
1688cb7bfe8cSJoe Wallwork /*@
168920282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
169020282da2SJoe Wallwork 
1691f1a722f8SMatthew G. Knepley   Input Parameters:
1692*2fe279fdSBarry Smith + dm        - The `DM`
169320282da2SJoe Wallwork . metric1   - The first metric to be intersected
169420282da2SJoe Wallwork . metric2   - The second metric to be intersected
169520282da2SJoe Wallwork - metric3   - The third metric to be intersected
169620282da2SJoe Wallwork 
169720282da2SJoe Wallwork   Output Parameter:
169820282da2SJoe Wallwork . metricInt - The intersected metric
169920282da2SJoe Wallwork 
170020282da2SJoe Wallwork   Level: beginner
170120282da2SJoe Wallwork 
1702*2fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1703cb7bfe8cSJoe Wallwork @*/
1704d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1705d71ae5a4SJacob Faibussowitsch {
170620282da2SJoe Wallwork   Vec metrics[3] = {metric1, metric2, metric3};
170720282da2SJoe Wallwork 
170820282da2SJoe Wallwork   PetscFunctionBegin;
17099566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
17103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171120282da2SJoe Wallwork }
1712