xref: /petsc/src/dm/impls/plex/plexmetric.c (revision be37439ebbbdb2f81c3420c175a94aa72e59929c)
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 
DMPlexMetricSetFromOptions(DM dm)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 
5860225df5SJacob Faibussowitsch   Input Parameters:
592fe279fdSBarry Smith + dm        - The `DM`
6031b70465SJoe Wallwork - isotropic - Is the metric isotropic?
6131b70465SJoe Wallwork 
6231b70465SJoe Wallwork   Level: beginner
6331b70465SJoe Wallwork 
642fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetUniform()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
65cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetIsotropic(DM dm,PetscBool isotropic)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 
7960225df5SJacob Faibussowitsch   Input Parameters:
802fe279fdSBarry Smith . dm - The `DM`
8131b70465SJoe Wallwork 
8260225df5SJacob Faibussowitsch   Output Parameters:
8331b70465SJoe Wallwork . isotropic - Is the metric isotropic?
8431b70465SJoe Wallwork 
8531b70465SJoe Wallwork   Level: beginner
8631b70465SJoe Wallwork 
872fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricIsUniform()`, `DMPlexMetricRestrictAnisotropyFirst()`
88cb7bfe8cSJoe Wallwork @*/
DMPlexMetricIsIsotropic(DM dm,PetscBool * isotropic)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 
10260225df5SJacob Faibussowitsch   Input Parameters:
1032fe279fdSBarry Smith + dm      - The `DM`
10493520041SJoe Wallwork - uniform - Is the metric uniform?
10593520041SJoe Wallwork 
10693520041SJoe Wallwork   Level: beginner
10793520041SJoe Wallwork 
1082fe279fdSBarry Smith   Note:
10993520041SJoe Wallwork   If the metric is specified as uniform then it is assumed to be isotropic, too.
11093520041SJoe Wallwork 
1112fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
11293520041SJoe Wallwork @*/
DMPlexMetricSetUniform(DM dm,PetscBool uniform)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 
12760225df5SJacob Faibussowitsch   Input Parameters:
1282fe279fdSBarry Smith . dm - The `DM`
12993520041SJoe Wallwork 
13060225df5SJacob Faibussowitsch   Output Parameters:
13193520041SJoe Wallwork . uniform - Is the metric uniform?
13293520041SJoe Wallwork 
13393520041SJoe Wallwork   Level: beginner
13493520041SJoe Wallwork 
1352fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
13693520041SJoe Wallwork @*/
DMPlexMetricIsUniform(DM dm,PetscBool * uniform)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 
15060225df5SJacob Faibussowitsch   Input Parameters:
1512fe279fdSBarry Smith + dm                      - The `DM`
15231b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
15331b70465SJoe Wallwork 
15431b70465SJoe Wallwork   Level: beginner
15531b70465SJoe Wallwork 
1562fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
157cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetRestrictAnisotropyFirst(DM dm,PetscBool restrictAnisotropyFirst)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 
17160225df5SJacob Faibussowitsch   Input Parameters:
1722fe279fdSBarry Smith . dm - The `DM`
17331b70465SJoe Wallwork 
17460225df5SJacob Faibussowitsch   Output Parameters:
17531b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
17631b70465SJoe Wallwork 
17731b70465SJoe Wallwork   Level: beginner
17831b70465SJoe Wallwork 
1792fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
180cb7bfe8cSJoe Wallwork @*/
DMPlexMetricRestrictAnisotropyFirst(DM dm,PetscBool * restrictAnisotropyFirst)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 
19460225df5SJacob Faibussowitsch   Input Parameters:
1952fe279fdSBarry Smith + dm       - The `DM`
196cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
197cc2eee55SJoe Wallwork 
198cc2eee55SJoe Wallwork   Level: beginner
199cc2eee55SJoe Wallwork 
2002fe279fdSBarry Smith   Note:
201cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
202cb7bfe8cSJoe Wallwork 
2032fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
204cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetNoInsertion(DM dm,PetscBool noInsert)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 
21860225df5SJacob Faibussowitsch   Input Parameters:
2192fe279fdSBarry Smith . dm - The `DM`
220cc2eee55SJoe Wallwork 
22160225df5SJacob Faibussowitsch   Output Parameters:
222cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
223cc2eee55SJoe Wallwork 
224cc2eee55SJoe Wallwork   Level: beginner
225cc2eee55SJoe Wallwork 
2262fe279fdSBarry Smith   Note:
227cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
228cb7bfe8cSJoe Wallwork 
2292fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
230cb7bfe8cSJoe Wallwork @*/
DMPlexMetricNoInsertion(DM dm,PetscBool * noInsert)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 
24460225df5SJacob Faibussowitsch   Input Parameters:
2452fe279fdSBarry Smith + dm     - The `DM`
246cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
247cc2eee55SJoe Wallwork 
248cc2eee55SJoe Wallwork   Level: beginner
249cc2eee55SJoe Wallwork 
2502fe279fdSBarry Smith   Note:
251cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
252cb7bfe8cSJoe Wallwork 
2532fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
254cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetNoSwapping(DM dm,PetscBool noSwap)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 
26860225df5SJacob Faibussowitsch   Input Parameters:
2692fe279fdSBarry Smith . dm - The `DM`
270cc2eee55SJoe Wallwork 
27160225df5SJacob Faibussowitsch   Output Parameters:
272cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
273cc2eee55SJoe Wallwork 
274cc2eee55SJoe Wallwork   Level: beginner
275cc2eee55SJoe Wallwork 
2762fe279fdSBarry Smith   Note:
277cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
278cb7bfe8cSJoe Wallwork 
2792fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
280cb7bfe8cSJoe Wallwork @*/
DMPlexMetricNoSwapping(DM dm,PetscBool * noSwap)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 
29460225df5SJacob Faibussowitsch   Input Parameters:
2952fe279fdSBarry Smith + dm     - The `DM`
296cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
297cc2eee55SJoe Wallwork 
298cc2eee55SJoe Wallwork   Level: beginner
299cc2eee55SJoe Wallwork 
3002fe279fdSBarry Smith   Note:
301cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
302cb7bfe8cSJoe Wallwork 
3032fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()`
304cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetNoMovement(DM dm,PetscBool noMove)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 
31860225df5SJacob Faibussowitsch   Input Parameters:
3192fe279fdSBarry Smith . dm - The `DM`
320cc2eee55SJoe Wallwork 
32160225df5SJacob Faibussowitsch   Output Parameters:
322cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
323cc2eee55SJoe Wallwork 
324cc2eee55SJoe Wallwork   Level: beginner
325cc2eee55SJoe Wallwork 
3262fe279fdSBarry Smith   Note:
327cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
328cb7bfe8cSJoe Wallwork 
3292fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()`
330cb7bfe8cSJoe Wallwork @*/
DMPlexMetricNoMovement(DM dm,PetscBool * noMove)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 
34460225df5SJacob Faibussowitsch   Input Parameters:
3452fe279fdSBarry Smith + dm     - The `DM`
34676f360caSJoe Wallwork - noSurf - Should surface modification be turned off?
34776f360caSJoe Wallwork 
34876f360caSJoe Wallwork   Level: beginner
34976f360caSJoe Wallwork 
3502fe279fdSBarry Smith   Note:
35176f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
35276f360caSJoe Wallwork 
3532fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`
35476f360caSJoe Wallwork @*/
DMPlexMetricSetNoSurf(DM dm,PetscBool noSurf)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 
36860225df5SJacob Faibussowitsch   Input Parameters:
3692fe279fdSBarry Smith . dm - The `DM`
37076f360caSJoe Wallwork 
37160225df5SJacob Faibussowitsch   Output Parameters:
37276f360caSJoe Wallwork . noSurf - Is surface modification turned off?
37376f360caSJoe Wallwork 
37476f360caSJoe Wallwork   Level: beginner
37576f360caSJoe Wallwork 
3762fe279fdSBarry Smith   Note:
37776f360caSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
37876f360caSJoe Wallwork 
3792fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`
38076f360caSJoe Wallwork @*/
DMPlexMetricNoSurf(DM dm,PetscBool * noSurf)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 
39460225df5SJacob Faibussowitsch   Input Parameters:
3952fe279fdSBarry Smith + dm    - The `DM`
39631b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
39731b70465SJoe Wallwork 
39831b70465SJoe Wallwork   Level: beginner
39931b70465SJoe Wallwork 
4002fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()`
401cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetMinimumMagnitude(DM dm,PetscReal h_min)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 
41660225df5SJacob Faibussowitsch   Input Parameters:
4172fe279fdSBarry Smith . dm - The `DM`
41831b70465SJoe Wallwork 
41960225df5SJacob Faibussowitsch   Output Parameters:
42031b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
42131b70465SJoe Wallwork 
42231b70465SJoe Wallwork   Level: beginner
42331b70465SJoe Wallwork 
4242fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()`
425cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetMinimumMagnitude(DM dm,PetscReal * h_min)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 
43960225df5SJacob Faibussowitsch   Input Parameters:
4402fe279fdSBarry Smith + dm    - The `DM`
44131b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
44231b70465SJoe Wallwork 
44331b70465SJoe Wallwork   Level: beginner
44431b70465SJoe Wallwork 
4452fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()`
446cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetMaximumMagnitude(DM dm,PetscReal h_max)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 
46160225df5SJacob Faibussowitsch   Input Parameters:
4622fe279fdSBarry Smith . dm - The `DM`
46331b70465SJoe Wallwork 
46460225df5SJacob Faibussowitsch   Output Parameters:
46531b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
46631b70465SJoe Wallwork 
46731b70465SJoe Wallwork   Level: beginner
46831b70465SJoe Wallwork 
4692fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()`
470cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetMaximumMagnitude(DM dm,PetscReal * h_max)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 
48460225df5SJacob Faibussowitsch   Input Parameters:
4852fe279fdSBarry Smith + dm    - The `DM`
48631b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
48731b70465SJoe Wallwork 
48831b70465SJoe Wallwork   Level: beginner
48931b70465SJoe Wallwork 
4902fe279fdSBarry Smith   Note:
4912fe279fdSBarry Smith   If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
49231b70465SJoe Wallwork 
4932fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetMaximumAnisotropy()`, `DMPlexMetricSetMaximumMagnitude()`
494cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetMaximumAnisotropy(DM dm,PetscReal a_max)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 
50960225df5SJacob Faibussowitsch   Input Parameters:
5102fe279fdSBarry Smith . dm - The `DM`
51131b70465SJoe Wallwork 
51260225df5SJacob Faibussowitsch   Output Parameters:
51331b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
51431b70465SJoe Wallwork 
51531b70465SJoe Wallwork   Level: beginner
51631b70465SJoe Wallwork 
5172fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetMaximumAnisotropy()`, `DMPlexMetricGetMaximumMagnitude()`
518cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetMaximumAnisotropy(DM dm,PetscReal * a_max)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 
53260225df5SJacob Faibussowitsch   Input Parameters:
5332fe279fdSBarry Smith + dm               - The `DM`
53431b70465SJoe Wallwork - targetComplexity - The target metric complexity
53531b70465SJoe Wallwork 
53631b70465SJoe Wallwork   Level: beginner
53731b70465SJoe Wallwork 
5382fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetTargetComplexity()`, `DMPlexMetricSetNormalizationOrder()`
539cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetTargetComplexity(DM dm,PetscReal targetComplexity)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 
55460225df5SJacob Faibussowitsch   Input Parameters:
5552fe279fdSBarry Smith . dm - The `DM`
55631b70465SJoe Wallwork 
55760225df5SJacob Faibussowitsch   Output Parameters:
55831b70465SJoe Wallwork . targetComplexity - The target metric complexity
55931b70465SJoe Wallwork 
56031b70465SJoe Wallwork   Level: beginner
56131b70465SJoe Wallwork 
5622fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetTargetComplexity()`, `DMPlexMetricGetNormalizationOrder()`
563cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetTargetComplexity(DM dm,PetscReal * targetComplexity)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 
57760225df5SJacob Faibussowitsch   Input Parameters:
5782fe279fdSBarry Smith + dm - The `DM`
57931b70465SJoe Wallwork - p  - The normalization order
58031b70465SJoe Wallwork 
58131b70465SJoe Wallwork   Level: beginner
58231b70465SJoe Wallwork 
5832fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetNormalizationOrder()`, `DMPlexMetricSetTargetComplexity()`
584cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetNormalizationOrder(DM dm,PetscReal p)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 
59960225df5SJacob Faibussowitsch   Input Parameters:
6002fe279fdSBarry Smith . dm - The `DM`
60131b70465SJoe Wallwork 
60260225df5SJacob Faibussowitsch   Output Parameters:
60331b70465SJoe Wallwork . p - The normalization order
60431b70465SJoe Wallwork 
60531b70465SJoe Wallwork   Level: beginner
60631b70465SJoe Wallwork 
6072fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNormalizationOrder()`, `DMPlexMetricGetTargetComplexity()`
608cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetNormalizationOrder(DM dm,PetscReal * p)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 
62260225df5SJacob Faibussowitsch   Input Parameters:
6232fe279fdSBarry 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 
6352fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
636cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetGradationFactor(DM dm,PetscReal beta)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 
65160225df5SJacob Faibussowitsch   Input Parameters:
6522fe279fdSBarry Smith . dm - The `DM`
653cc2eee55SJoe Wallwork 
65460225df5SJacob Faibussowitsch   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 
6662fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
667cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetGradationFactor(DM dm,PetscReal * beta)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 
68160225df5SJacob Faibussowitsch   Input Parameters:
6822fe279fdSBarry 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 
6972fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
698ae8b063eSJoe Wallwork @*/
DMPlexMetricSetHausdorffNumber(DM dm,PetscReal hausd)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 
71360225df5SJacob Faibussowitsch   Input Parameters:
7142fe279fdSBarry Smith . dm - The `DM`
715ae8b063eSJoe Wallwork 
71660225df5SJacob Faibussowitsch   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 
7312fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
732ae8b063eSJoe Wallwork @*/
DMPlexMetricGetHausdorffNumber(DM dm,PetscReal * hausd)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 
74660225df5SJacob Faibussowitsch   Input Parameters:
7472fe279fdSBarry 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 
7522fe279fdSBarry Smith   Note:
753cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
754cb7bfe8cSJoe Wallwork 
7552fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()`
756cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetVerbosity(DM dm,PetscInt verbosity)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 
77060225df5SJacob Faibussowitsch   Input Parameters:
7712fe279fdSBarry Smith . dm - The `DM`
772cc2eee55SJoe Wallwork 
77360225df5SJacob Faibussowitsch   Output Parameters:
774cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
775cc2eee55SJoe Wallwork 
776cb7bfe8cSJoe Wallwork   Level: beginner
777cb7bfe8cSJoe Wallwork 
7782fe279fdSBarry Smith   Note:
779cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
780cb7bfe8cSJoe Wallwork 
7812fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
782cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetVerbosity(DM dm,PetscInt * verbosity)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 
79660225df5SJacob Faibussowitsch   Input Parameters:
7972fe279fdSBarry Smith + dm      - The `DM`
798cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
799cc2eee55SJoe Wallwork 
800cb7bfe8cSJoe Wallwork   Level: beginner
801cb7bfe8cSJoe Wallwork 
8022fe279fdSBarry Smith   Note:
803cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
804cc2eee55SJoe Wallwork 
8052fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
806cb7bfe8cSJoe Wallwork @*/
DMPlexMetricSetNumIterations(DM dm,PetscInt numIter)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 
82060225df5SJacob Faibussowitsch   Input Parameters:
8212fe279fdSBarry Smith . dm - The `DM`
822cc2eee55SJoe Wallwork 
82360225df5SJacob Faibussowitsch   Output Parameters:
824cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
825cc2eee55SJoe Wallwork 
826cb7bfe8cSJoe Wallwork   Level: beginner
827cb7bfe8cSJoe Wallwork 
8282fe279fdSBarry Smith   Note:
829cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
830cc2eee55SJoe Wallwork 
8312fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()`
832cb7bfe8cSJoe Wallwork @*/
DMPlexMetricGetNumIterations(DM dm,PetscInt * numIter)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 
DMPlexP1FieldCreate_Private(DM dm,PetscInt f,PetscInt size,Vec * metric)84366976f2fSJacob Faibussowitsch static 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   /* Extract metadata from dm */
8519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
8529566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
85320282da2SJoe Wallwork 
85420282da2SJoe Wallwork   /* Create a P1 field of the requested size */
8559566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
8569566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
8579566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
8589566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
8599566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm, metric));
8603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86120282da2SJoe Wallwork }
86220282da2SJoe Wallwork 
863cb7bfe8cSJoe Wallwork /*@
86420282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
86520282da2SJoe Wallwork 
86660225df5SJacob Faibussowitsch   Input Parameters:
8672fe279fdSBarry Smith + dm - The `DM`
86820282da2SJoe Wallwork - f  - The field number to use
86920282da2SJoe Wallwork 
87060225df5SJacob Faibussowitsch   Output Parameter:
87120282da2SJoe Wallwork . metric - The metric
87220282da2SJoe Wallwork 
8732fe279fdSBarry Smith   Options Database Key:
8742fe279fdSBarry Smith . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use
87520282da2SJoe Wallwork 
8762fe279fdSBarry Smith   Options Database Keys for Mmg and ParMmg:
8772fe279fdSBarry Smith + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation
8782fe279fdSBarry Smith . -dm_plex_metric_num_iterations   - Number of parallel mesh adaptation iterations for ParMmg
8792fe279fdSBarry Smith . -dm_plex_metric_no_insert        - Should node insertion/deletion be turned off?
8802fe279fdSBarry Smith . -dm_plex_metric_no_swap          - Should facet swapping be turned off?
8812fe279fdSBarry Smith . -dm_plex_metric_no_move          - Should node movement be turned off?
8822fe279fdSBarry Smith - -dm_plex_metric_verbosity        - Choose a verbosity level from -1 (silent) to 10 (maximum).
88331b70465SJoe Wallwork 
8842fe279fdSBarry Smith   Options Database Keys for Riemannian metrics:
885cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
88693520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
887cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
888cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
889cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
890cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
891cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
89267b8a455SSatish Balay - -dm_plex_metric_target_complexity         - Target metric complexity
893cb7bfe8cSJoe Wallwork 
8942fe279fdSBarry Smith   Level: beginner
895cb7bfe8cSJoe Wallwork 
8962fe279fdSBarry Smith   Note:
8972fe279fdSBarry Smith   It is assumed that the `DM` is comprised of simplices.
898cb7bfe8cSJoe Wallwork 
8992fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`
900cb7bfe8cSJoe Wallwork @*/
DMPlexMetricCreate(DM dm,PetscInt f,Vec * metric)901d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
902d71ae5a4SJacob Faibussowitsch {
90393520041SJoe Wallwork   PetscBool isotropic, uniform;
90420282da2SJoe Wallwork   PetscInt  coordDim, Nd;
90520282da2SJoe Wallwork 
90620282da2SJoe Wallwork   PetscFunctionBegin;
9079566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
90820282da2SJoe Wallwork   Nd = coordDim * coordDim;
9099566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
9109566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
91193520041SJoe Wallwork   if (uniform) {
91293520041SJoe Wallwork     MPI_Comm comm;
91393520041SJoe Wallwork 
9149566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
915bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
9169566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, metric));
9179566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
9189566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(*metric));
919f6db075bSJoe Wallwork   } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
920f6db075bSJoe Wallwork   else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
9213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92220282da2SJoe Wallwork }
92320282da2SJoe Wallwork 
924cb7bfe8cSJoe Wallwork /*@
92520282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
92620282da2SJoe Wallwork 
92760225df5SJacob Faibussowitsch   Input Parameters:
9282fe279fdSBarry Smith + dm    - The `DM`
92920282da2SJoe Wallwork . f     - The field number to use
93020282da2SJoe Wallwork - alpha - Scaling parameter for the diagonal
93120282da2SJoe Wallwork 
93260225df5SJacob Faibussowitsch   Output Parameter:
93320282da2SJoe Wallwork . metric - The uniform metric
93420282da2SJoe Wallwork 
93520282da2SJoe Wallwork   Level: beginner
93620282da2SJoe Wallwork 
9372fe279fdSBarry Smith   Note:
9382fe279fdSBarry Smith   In this case, the metric is constant in space and so is only specified for a single vertex.
93920282da2SJoe Wallwork 
9402fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
941cb7bfe8cSJoe Wallwork @*/
DMPlexMetricCreateUniform(DM dm,PetscInt f,PetscReal alpha,Vec * metric)942d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
943d71ae5a4SJacob Faibussowitsch {
94420282da2SJoe Wallwork   PetscFunctionBegin;
9459566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
9469566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9472c71b3e2SJacob Faibussowitsch   PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
948bc00d7afSJoe Wallwork   PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
9499566063dSJacob Faibussowitsch   PetscCall(VecSet(*metric, alpha));
9509566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*metric));
9519566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*metric));
9523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95320282da2SJoe Wallwork }
95420282da2SJoe Wallwork 
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[])955d71ae5a4SJacob 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[])
956d71ae5a4SJacob Faibussowitsch {
95793520041SJoe Wallwork   f0[0] = u[0];
95893520041SJoe Wallwork }
95993520041SJoe Wallwork 
960cb7bfe8cSJoe Wallwork /*@
96120282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
96220282da2SJoe Wallwork 
96360225df5SJacob Faibussowitsch   Input Parameters:
9642fe279fdSBarry Smith + dm        - The `DM`
96520282da2SJoe Wallwork . f         - The field number to use
96620282da2SJoe Wallwork - indicator - The error indicator
96720282da2SJoe Wallwork 
96860225df5SJacob Faibussowitsch   Output Parameter:
96920282da2SJoe Wallwork . metric - The isotropic metric
97020282da2SJoe Wallwork 
97120282da2SJoe Wallwork   Level: beginner
97220282da2SJoe Wallwork 
97320282da2SJoe Wallwork   Notes:
9742fe279fdSBarry Smith   It is assumed that the `DM` is comprised of simplices.
97520282da2SJoe Wallwork 
97693520041SJoe Wallwork   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
97720282da2SJoe Wallwork 
9782fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
979cb7bfe8cSJoe Wallwork @*/
DMPlexMetricCreateIsotropic(DM dm,PetscInt f,Vec indicator,Vec * metric)980d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
981d71ae5a4SJacob Faibussowitsch {
98293520041SJoe Wallwork   PetscInt m, n;
98320282da2SJoe Wallwork 
98420282da2SJoe Wallwork   PetscFunctionBegin;
9859566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
9869566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dm, f, metric));
9879566063dSJacob Faibussowitsch   PetscCall(VecGetSize(indicator, &m));
9889566063dSJacob Faibussowitsch   PetscCall(VecGetSize(*metric, &n));
9899566063dSJacob Faibussowitsch   if (m == n) PetscCall(VecCopy(indicator, *metric));
99093520041SJoe Wallwork   else {
99193520041SJoe Wallwork     DM dmIndi;
9929371c9d4SSatish 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[]);
99393520041SJoe Wallwork 
9949566063dSJacob Faibussowitsch     PetscCall(VecGetDM(indicator, &dmIndi));
99593520041SJoe Wallwork     funcs[0] = identity;
9969566063dSJacob Faibussowitsch     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
99720282da2SJoe Wallwork   }
9983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99920282da2SJoe Wallwork }
100020282da2SJoe Wallwork 
1001f6db075bSJoe Wallwork /*@
1002f6db075bSJoe Wallwork   DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1003f6db075bSJoe Wallwork 
100460225df5SJacob Faibussowitsch   Input Parameters:
10052fe279fdSBarry Smith + dm - The `DM` of the metric field
1006f6db075bSJoe Wallwork - f  - The field number to use
1007f6db075bSJoe Wallwork 
100860225df5SJacob Faibussowitsch   Output Parameters:
1009f6db075bSJoe Wallwork + determinant - The determinant field
10102fe279fdSBarry Smith - dmDet       - The corresponding `DM`
1011f6db075bSJoe Wallwork 
1012f6db075bSJoe Wallwork   Level: beginner
1013f6db075bSJoe Wallwork 
10142fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`, `DMPlexMetricCreate()`
1015f6db075bSJoe Wallwork @*/
DMPlexMetricDeterminantCreate(DM dm,PetscInt f,Vec * determinant,DM * dmDet)1016d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1017d71ae5a4SJacob Faibussowitsch {
1018f6db075bSJoe Wallwork   PetscBool uniform;
1019f6db075bSJoe Wallwork 
1020f6db075bSJoe Wallwork   PetscFunctionBegin;
1021f6db075bSJoe Wallwork   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1022f6db075bSJoe Wallwork   PetscCall(DMClone(dm, dmDet));
1023f6db075bSJoe Wallwork   if (uniform) {
1024f6db075bSJoe Wallwork     MPI_Comm comm;
1025f6db075bSJoe Wallwork 
1026f6db075bSJoe Wallwork     PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm));
1027f6db075bSJoe Wallwork     PetscCall(VecCreate(comm, determinant));
1028f6db075bSJoe Wallwork     PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1029f6db075bSJoe Wallwork     PetscCall(VecSetFromOptions(*determinant));
1030f6db075bSJoe Wallwork   } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
10313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1032f6db075bSJoe Wallwork }
1033f6db075bSJoe Wallwork 
LAPACKsyevFail(PetscInt dim,PetscScalar Mpos[])1034d71ae5a4SJacob Faibussowitsch static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1035d71ae5a4SJacob Faibussowitsch {
103682490253SJoe Wallwork   PetscInt i, j;
103782490253SJoe Wallwork 
103882490253SJoe Wallwork   PetscFunctionBegin;
10393ba16761SJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"));
104082490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
10413ba16761SJacob Faibussowitsch     if (i == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    [["));
10423ba16761SJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "     ["));
104382490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
10443ba16761SJacob Faibussowitsch       if (j < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])));
10453ba16761SJacob Faibussowitsch       else PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])));
104682490253SJoe Wallwork     }
10473ba16761SJacob Faibussowitsch     if (i < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
10483ba16761SJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "]]\n"));
104982490253SJoe Wallwork   }
10503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105182490253SJoe Wallwork }
105282490253SJoe Wallwork 
DMPlexMetricModify_Private(PetscInt dim,PetscReal h_min,PetscReal h_max,PetscReal a_max,PetscScalar Mp[],PetscScalar * detMp)1053d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1054d71ae5a4SJacob Faibussowitsch {
105520282da2SJoe Wallwork   PetscInt     i, j, k;
105620282da2SJoe 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);
105720282da2SJoe Wallwork   PetscScalar *Mpos;
105820282da2SJoe Wallwork 
105920282da2SJoe Wallwork   PetscFunctionBegin;
10609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs));
106120282da2SJoe Wallwork 
106220282da2SJoe Wallwork   /* Symmetrize */
106320282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
106420282da2SJoe Wallwork     Mpos[i * dim + i] = Mp[i * dim + i];
106520282da2SJoe Wallwork     for (j = i + 1; j < dim; ++j) {
106620282da2SJoe Wallwork       Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
106720282da2SJoe Wallwork       Mpos[j * dim + i] = Mpos[i * dim + j];
106820282da2SJoe Wallwork     }
106920282da2SJoe Wallwork   }
107020282da2SJoe Wallwork 
107120282da2SJoe Wallwork   /* Compute eigendecomposition */
107293520041SJoe Wallwork   if (dim == 1) {
107393520041SJoe Wallwork     /* Isotropic case */
107493520041SJoe Wallwork     eigs[0] = PetscRealPart(Mpos[0]);
107593520041SJoe Wallwork     Mpos[0] = 1.0;
107693520041SJoe Wallwork   } else {
107793520041SJoe Wallwork     /* Anisotropic case */
107820282da2SJoe Wallwork     PetscScalar *work;
1079*835f2295SStefano Zampini     PetscBLASInt lwork;
108020282da2SJoe Wallwork 
1081*835f2295SStefano Zampini     PetscCall(PetscBLASIntCast(5 * dim, &lwork));
10829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * dim, &work));
108320282da2SJoe Wallwork     {
108420282da2SJoe Wallwork       PetscBLASInt lierr;
108520282da2SJoe Wallwork       PetscBLASInt nb;
108620282da2SJoe Wallwork 
10879566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
10889566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
108920282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
109020282da2SJoe Wallwork       {
109120282da2SJoe Wallwork         PetscReal *rwork;
10929566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1093792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr));
10949566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
109520282da2SJoe Wallwork       }
109620282da2SJoe Wallwork #else
1097792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr));
109820282da2SJoe Wallwork #endif
109982490253SJoe Wallwork       if (lierr) {
110082490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
110182490253SJoe Wallwork           Mpos[i * dim + i] = Mp[i * dim + i];
110282490253SJoe Wallwork           for (j = i + 1; j < dim; ++j) {
110382490253SJoe Wallwork             Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
110482490253SJoe Wallwork             Mpos[j * dim + i] = Mpos[i * dim + j];
110582490253SJoe Wallwork           }
110682490253SJoe Wallwork         }
11073ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, Mpos));
1108*835f2295SStefano Zampini         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %" PetscBLASInt_FMT, lierr);
110982490253SJoe Wallwork       }
11109566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
111120282da2SJoe Wallwork     }
11129566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
111320282da2SJoe Wallwork   }
111420282da2SJoe Wallwork 
111520282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
111620282da2SJoe Wallwork   max_eig = 0.0;
111720282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
111820282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
111920282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
112020282da2SJoe Wallwork   }
112120282da2SJoe Wallwork 
11223f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
11233f07679eSJoe Wallwork   *detMp = 1.0;
112420282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
1125ac314011SJoe Wallwork     if (a_max >= 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min);
11263f07679eSJoe Wallwork     *detMp *= eigs[i];
112720282da2SJoe Wallwork   }
112820282da2SJoe Wallwork 
112920282da2SJoe Wallwork   /* Reconstruct Hessian */
113020282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
113120282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
113220282da2SJoe Wallwork       Mp[i * dim + j] = 0.0;
1133ad540459SPierre Jolivet       for (k = 0; k < dim; ++k) Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j];
113420282da2SJoe Wallwork     }
113520282da2SJoe Wallwork   }
11369566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Mpos, eigs));
11373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113820282da2SJoe Wallwork }
113920282da2SJoe Wallwork 
1140cb7bfe8cSJoe Wallwork /*@
114120282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
114220282da2SJoe Wallwork 
114360225df5SJacob Faibussowitsch   Input Parameters:
11442fe279fdSBarry Smith + dm                 - The `DM`
11453f07679eSJoe Wallwork . metricIn           - The metric
114699abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
11473f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
114820282da2SJoe Wallwork 
114960225df5SJacob Faibussowitsch   Output Parameters:
11503f07679eSJoe Wallwork + metricOut   - The metric
11513f07679eSJoe Wallwork - determinant - Its determinant
115220282da2SJoe Wallwork 
11532fe279fdSBarry Smith   Options Database Keys:
115493520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic?
115593520041SJoe Wallwork . -dm_plex_metric_uniform   - Is the metric uniform?
115693520041SJoe Wallwork . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1157cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1158cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1159cb7bfe8cSJoe Wallwork 
11602fe279fdSBarry Smith   Level: beginner
11612fe279fdSBarry Smith 
11622fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1163cb7bfe8cSJoe Wallwork @*/
DMPlexMetricEnforceSPD(DM dm,Vec metricIn,PetscBool restrictSizes,PetscBool restrictAnisotropy,Vec metricOut,Vec determinant)1164d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1165d71ae5a4SJacob Faibussowitsch {
11663f07679eSJoe Wallwork   DM           dmDet;
116793520041SJoe Wallwork   PetscBool    isotropic, uniform;
116820282da2SJoe Wallwork   PetscInt     dim, vStart, vEnd, v;
11693f07679eSJoe Wallwork   PetscScalar *met, *det;
117061451c10SMatthew G. Knepley   PetscReal    h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+15;
117120282da2SJoe Wallwork 
117220282da2SJoe Wallwork   PetscFunctionBegin;
11739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
117420282da2SJoe Wallwork 
117520282da2SJoe Wallwork   /* Extract metadata from dm */
11769566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
117720282da2SJoe Wallwork   if (restrictSizes) {
11789566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
11799566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
118099abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
118199abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1182bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
118399abec2bSJoe Wallwork   }
118499abec2bSJoe Wallwork   if (restrictAnisotropy) {
11859566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
118699abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
118720282da2SJoe Wallwork   }
118820282da2SJoe Wallwork 
118993520041SJoe Wallwork   /* Setup output metric */
11905241a91fSJoe Wallwork   PetscCall(VecCopy(metricIn, metricOut));
119193520041SJoe Wallwork 
119293520041SJoe Wallwork   /* Enforce SPD and extract determinant */
11935241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
11949566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
11959566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
119693520041SJoe Wallwork   if (uniform) {
1197d1a5b989SMatthew G. Knepley     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
119893520041SJoe Wallwork 
119993520041SJoe Wallwork     /* Uniform case */
12005241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
12019566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
12025241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
120393520041SJoe Wallwork   } else {
120493520041SJoe Wallwork     /* Spatially varying case */
120593520041SJoe Wallwork     PetscInt nrow;
120693520041SJoe Wallwork 
120793520041SJoe Wallwork     if (isotropic) nrow = 1;
120893520041SJoe Wallwork     else nrow = dim;
12095241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
12109566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12115241a91fSJoe Wallwork     PetscCall(VecGetArray(determinant, &det));
121220282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
12133f07679eSJoe Wallwork       PetscScalar *vmet, *vdet;
12149566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
12159566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
12169566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
121720282da2SJoe Wallwork     }
12185241a91fSJoe Wallwork     PetscCall(VecRestoreArray(determinant, &det));
121993520041SJoe Wallwork   }
12205241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
1221fe902aceSJoe Wallwork 
12229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
12233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122420282da2SJoe Wallwork }
122520282da2SJoe Wallwork 
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[])1226d71ae5a4SJacob 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[])
1227d71ae5a4SJacob Faibussowitsch {
122820282da2SJoe Wallwork   const PetscScalar p = constants[0];
122920282da2SJoe Wallwork 
1230985ad86eSJose E. Roman   f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim));
123120282da2SJoe Wallwork }
123220282da2SJoe Wallwork 
1233cb7bfe8cSJoe Wallwork /*@
123420282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
123520282da2SJoe Wallwork 
123660225df5SJacob Faibussowitsch   Input Parameters:
12372fe279fdSBarry Smith + dm                 - The `DM`
123820282da2SJoe Wallwork . metricIn           - The unnormalized metric
123916522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
124016522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
124120282da2SJoe Wallwork 
124260225df5SJacob Faibussowitsch   Output Parameters:
12432fe279fdSBarry Smith + metricOut   - The normalized metric
12442fe279fdSBarry Smith - determinant - computed determinant
124520282da2SJoe Wallwork 
12462fe279fdSBarry Smith   Options Database Keys:
124793520041SJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
124893520041SJoe Wallwork . -dm_plex_metric_uniform                   - Is the metric uniform?
124993520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1250cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1251cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1252cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1253cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1254cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1255cb7bfe8cSJoe Wallwork 
12562fe279fdSBarry Smith   Level: beginner
12572fe279fdSBarry Smith 
12582fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1259cb7bfe8cSJoe Wallwork @*/
DMPlexMetricNormalize(DM dm,Vec metricIn,PetscBool restrictSizes,PetscBool restrictAnisotropy,Vec metricOut,Vec determinant)1260d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1261d71ae5a4SJacob Faibussowitsch {
12623f07679eSJoe Wallwork   DM           dmDet;
126320282da2SJoe Wallwork   MPI_Comm     comm;
126493520041SJoe Wallwork   PetscBool    restrictAnisotropyFirst, isotropic, uniform;
126520282da2SJoe Wallwork   PetscDS      ds;
126620282da2SJoe Wallwork   PetscInt     dim, Nd, vStart, vEnd, v, i;
12673f07679eSJoe Wallwork   PetscScalar *met, *det, integral, constants[1];
126893520041SJoe Wallwork   PetscReal    p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
126920282da2SJoe Wallwork 
127020282da2SJoe Wallwork   PetscFunctionBegin;
12719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0));
127220282da2SJoe Wallwork 
127320282da2SJoe Wallwork   /* Extract metadata from dm */
12749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
12759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
12769566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
12779566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
127893520041SJoe Wallwork   if (isotropic) Nd = 1;
127993520041SJoe Wallwork   else Nd = dim * dim;
128020282da2SJoe Wallwork 
128120282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
12829566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
12835241a91fSJoe Wallwork   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
128420282da2SJoe Wallwork 
128520282da2SJoe Wallwork   /* Compute global normalization factor */
12869566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
12879566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
128816522936SJoe Wallwork   constants[0] = p;
128993520041SJoe Wallwork   if (uniform) {
1290bc00d7afSJoe Wallwork     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
129193520041SJoe Wallwork     DM  dmTmp;
129293520041SJoe Wallwork     Vec tmp;
129393520041SJoe Wallwork 
12949566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmTmp));
12959566063dSJacob Faibussowitsch     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
12969566063dSJacob Faibussowitsch     PetscCall(VecGetArray(determinant, &det));
12979566063dSJacob Faibussowitsch     PetscCall(VecSet(tmp, det[0]));
12989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(determinant, &det));
12999566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmTmp, &ds));
13009566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13019566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13029566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
13039566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tmp));
13049566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmTmp));
130593520041SJoe Wallwork   } else {
13069566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
13079566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmDet, &ds));
13089566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
13099566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
13109566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
131193520041SJoe Wallwork   }
13123f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
1313bc00d7afSJoe 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?");
13143f07679eSJoe Wallwork   factGlob = PetscPowReal(target / realIntegral, 2.0 / dim);
131520282da2SJoe Wallwork 
131620282da2SJoe Wallwork   /* Apply local scaling */
131720282da2SJoe Wallwork   if (restrictSizes) {
13189566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
13199566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
132016522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
132116522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
1322bc00d7afSJoe Wallwork     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
132316522936SJoe Wallwork   }
132416522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
13259566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
132616522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
132720282da2SJoe Wallwork   }
13285241a91fSJoe Wallwork   PetscCall(VecGetArray(metricOut, &met));
13299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(determinant, &det));
133093520041SJoe Wallwork   if (uniform) {
133193520041SJoe Wallwork     /* Uniform case */
133293520041SJoe Wallwork     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim));
13339566063dSJacob Faibussowitsch     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
133493520041SJoe Wallwork   } else {
133593520041SJoe Wallwork     /* Spatially varying case */
133693520041SJoe Wallwork     PetscInt nrow;
133793520041SJoe Wallwork 
133893520041SJoe Wallwork     if (isotropic) nrow = 1;
133993520041SJoe Wallwork     else nrow = dim;
13409566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
13415241a91fSJoe Wallwork     PetscCall(VecGetDM(determinant, &dmDet));
134220282da2SJoe Wallwork     for (v = vStart; v < vEnd; ++v) {
13433f07679eSJoe Wallwork       PetscScalar *Mp, *detM;
134420282da2SJoe Wallwork 
13459566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
13469566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
13473f07679eSJoe Wallwork       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim));
134820282da2SJoe Wallwork       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
13499566063dSJacob Faibussowitsch       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
135093520041SJoe Wallwork     }
135120282da2SJoe Wallwork   }
13529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(determinant, &det));
13535241a91fSJoe Wallwork   PetscCall(VecRestoreArray(metricOut, &met));
135420282da2SJoe Wallwork 
13559566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0));
13563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135720282da2SJoe Wallwork }
135820282da2SJoe Wallwork 
1359cb7bfe8cSJoe Wallwork /*@
136020282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
136120282da2SJoe Wallwork 
1362f1a722f8SMatthew G. Knepley   Input Parameters:
13632fe279fdSBarry Smith + dm         - The `DM`
136420282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
136520282da2SJoe Wallwork . weights    - Weights for the average
136620282da2SJoe Wallwork - metrics    - The metrics to be averaged
136720282da2SJoe Wallwork 
136820282da2SJoe Wallwork   Output Parameter:
136920282da2SJoe Wallwork . metricAvg - The averaged metric
137020282da2SJoe Wallwork 
137120282da2SJoe Wallwork   Level: beginner
137220282da2SJoe Wallwork 
137320282da2SJoe Wallwork   Notes:
137420282da2SJoe Wallwork   The weights should sum to unity.
137520282da2SJoe Wallwork 
137620282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
137720282da2SJoe Wallwork 
13782fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1379cb7bfe8cSJoe Wallwork @*/
DMPlexMetricAverage(DM dm,PetscInt numMetrics,PetscReal weights[],Vec metrics[],Vec metricAvg)1380d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1381d71ae5a4SJacob Faibussowitsch {
138220282da2SJoe Wallwork   PetscBool haveWeights = PETSC_TRUE;
138393520041SJoe Wallwork   PetscInt  i, m, n;
138420282da2SJoe Wallwork   PetscReal sum = 0.0, tol = 1.0e-10;
138520282da2SJoe Wallwork 
138620282da2SJoe Wallwork   PetscFunctionBegin;
13879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0));
138863a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
13895241a91fSJoe Wallwork   PetscCall(VecSet(metricAvg, 0.0));
13905241a91fSJoe Wallwork   PetscCall(VecGetSize(metricAvg, &m));
139193520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
13929566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
13935f80ce2aSJacob Faibussowitsch     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
139493520041SJoe Wallwork   }
139520282da2SJoe Wallwork 
139620282da2SJoe Wallwork   /* Default to the unweighted case */
139720282da2SJoe Wallwork   if (!weights) {
13989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numMetrics, &weights));
139920282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
1400ad540459SPierre Jolivet     for (i = 0; i < numMetrics; ++i) weights[i] = 1.0 / numMetrics;
140120282da2SJoe Wallwork   }
140220282da2SJoe Wallwork 
140320282da2SJoe Wallwork   /* Check weights sum to unity */
140493520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) sum += weights[i];
14055f80ce2aSJacob Faibussowitsch   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
140620282da2SJoe Wallwork 
140720282da2SJoe Wallwork   /* Compute metric average */
14085241a91fSJoe Wallwork   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
14099566063dSJacob Faibussowitsch   if (!haveWeights) PetscCall(PetscFree(weights));
1410fe902aceSJoe Wallwork 
14119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0));
14123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141320282da2SJoe Wallwork }
141420282da2SJoe Wallwork 
1415cb7bfe8cSJoe Wallwork /*@
141620282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
141720282da2SJoe Wallwork 
1418f1a722f8SMatthew G. Knepley   Input Parameters:
14192fe279fdSBarry Smith + dm      - The `DM`
142020282da2SJoe Wallwork . metric1 - The first metric to be averaged
142120282da2SJoe Wallwork - metric2 - The second metric to be averaged
142220282da2SJoe Wallwork 
142320282da2SJoe Wallwork   Output Parameter:
142420282da2SJoe Wallwork . metricAvg - The averaged metric
142520282da2SJoe Wallwork 
142620282da2SJoe Wallwork   Level: beginner
142720282da2SJoe Wallwork 
14282fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1429cb7bfe8cSJoe Wallwork @*/
DMPlexMetricAverage2(DM dm,Vec metric1,Vec metric2,Vec metricAvg)1430d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1431d71ae5a4SJacob Faibussowitsch {
143220282da2SJoe Wallwork   PetscReal weights[2] = {0.5, 0.5};
143320282da2SJoe Wallwork   Vec       metrics[2] = {metric1, metric2};
143420282da2SJoe Wallwork 
143520282da2SJoe Wallwork   PetscFunctionBegin;
14369566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
14373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143820282da2SJoe Wallwork }
143920282da2SJoe Wallwork 
1440cb7bfe8cSJoe Wallwork /*@
144120282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
144220282da2SJoe Wallwork 
1443f1a722f8SMatthew G. Knepley   Input Parameters:
14442fe279fdSBarry Smith + dm      - The `DM`
144520282da2SJoe Wallwork . metric1 - The first metric to be averaged
144620282da2SJoe Wallwork . metric2 - The second metric to be averaged
144720282da2SJoe Wallwork - metric3 - The third metric to be averaged
144820282da2SJoe Wallwork 
144920282da2SJoe Wallwork   Output Parameter:
145020282da2SJoe Wallwork . metricAvg - The averaged metric
145120282da2SJoe Wallwork 
145220282da2SJoe Wallwork   Level: beginner
145320282da2SJoe Wallwork 
14542fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1455cb7bfe8cSJoe Wallwork @*/
DMPlexMetricAverage3(DM dm,Vec metric1,Vec metric2,Vec metric3,Vec metricAvg)1456d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1457d71ae5a4SJacob Faibussowitsch {
145820282da2SJoe Wallwork   PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0};
145920282da2SJoe Wallwork   Vec       metrics[3] = {metric1, metric2, metric3};
146020282da2SJoe Wallwork 
146120282da2SJoe Wallwork   PetscFunctionBegin;
14629566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
14633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146420282da2SJoe Wallwork }
146520282da2SJoe Wallwork 
DMPlexMetricIntersection_Private(PetscInt dim,PetscScalar M1[],PetscScalar M2[])1466d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1467d71ae5a4SJacob Faibussowitsch {
146820282da2SJoe Wallwork   PetscInt     i, j, k, l, m;
1469e2606525SJoe Wallwork   PetscReal   *evals;
147020282da2SJoe Wallwork   PetscScalar *evecs, *sqrtM1, *isqrtM1;
147120282da2SJoe Wallwork 
147220282da2SJoe Wallwork   PetscFunctionBegin;
147393520041SJoe Wallwork   /* Isotropic case */
147493520041SJoe Wallwork   if (dim == 1) {
1475*835f2295SStefano Zampini     M2[0] = PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
14763ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
147793520041SJoe Wallwork   }
147893520041SJoe Wallwork 
147993520041SJoe Wallwork   /* Anisotropic case */
1480e2606525SJoe Wallwork   PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals));
148120282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
1482ad540459SPierre Jolivet     for (j = 0; j < dim; ++j) evecs[i * dim + j] = M1[i * dim + j];
148320282da2SJoe Wallwork   }
148420282da2SJoe Wallwork   {
148520282da2SJoe Wallwork     PetscScalar *work;
1486*835f2295SStefano Zampini     PetscBLASInt lwork;
1487*835f2295SStefano Zampini 
1488*835f2295SStefano Zampini     PetscCall(PetscBLASIntCast(5 * dim, &lwork));
14899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * dim, &work));
149020282da2SJoe Wallwork     {
149120282da2SJoe Wallwork       PetscBLASInt lierr, nb;
1492e2606525SJoe Wallwork       PetscReal    sqrtj;
149320282da2SJoe Wallwork 
149420282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
14959566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(dim, &nb));
14969566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
149720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
149820282da2SJoe Wallwork       {
149920282da2SJoe Wallwork         PetscReal *rwork;
15009566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1501792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15029566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
150320282da2SJoe Wallwork       }
150420282da2SJoe Wallwork #else
1505792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
150620282da2SJoe Wallwork #endif
150782490253SJoe Wallwork       if (lierr) {
15083ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, M1));
1509*835f2295SStefano Zampini         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %" PetscBLASInt_FMT, lierr);
151082490253SJoe Wallwork       }
15119566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
151220282da2SJoe Wallwork 
1513e2606525SJoe Wallwork       /* Compute square root and the reciprocal thereof */
151420282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
151520282da2SJoe Wallwork         for (k = 0; k < dim; ++k) {
1516e2606525SJoe Wallwork           sqrtM1[i * dim + k]  = 0.0;
1517e2606525SJoe Wallwork           isqrtM1[i * dim + k] = 0.0;
1518e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1519e2606525SJoe Wallwork             sqrtj = PetscSqrtReal(evals[j]);
1520e2606525SJoe Wallwork             sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k];
1521e2606525SJoe Wallwork             isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k];
152220282da2SJoe Wallwork           }
152320282da2SJoe Wallwork         }
152420282da2SJoe Wallwork       }
152520282da2SJoe Wallwork 
1526e2606525SJoe Wallwork       /* Map M2 into the space spanned by the eigenvectors of M1 */
152720282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
152820282da2SJoe Wallwork         for (l = 0; l < dim; ++l) {
1529e2606525SJoe Wallwork           evecs[i * dim + l] = 0.0;
1530e2606525SJoe Wallwork           for (j = 0; j < dim; ++j) {
1531ad540459SPierre Jolivet             for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
153220282da2SJoe Wallwork           }
153320282da2SJoe Wallwork         }
153420282da2SJoe Wallwork       }
153520282da2SJoe Wallwork 
153620282da2SJoe Wallwork       /* Compute eigendecomposition */
15379566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
153820282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
153920282da2SJoe Wallwork       {
154020282da2SJoe Wallwork         PetscReal *rwork;
15419566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(3 * dim, &rwork));
1542792fecdfSBarry Smith         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
15439566063dSJacob Faibussowitsch         PetscCall(PetscFree(rwork));
154420282da2SJoe Wallwork       }
154520282da2SJoe Wallwork #else
1546792fecdfSBarry Smith       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
154720282da2SJoe Wallwork #endif
154882490253SJoe Wallwork       if (lierr) {
154982490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
155082490253SJoe Wallwork           for (l = 0; l < dim; ++l) {
1551e2606525SJoe Wallwork             evecs[i * dim + l] = 0.0;
1552e2606525SJoe Wallwork             for (j = 0; j < dim; ++j) {
1553ad540459SPierre Jolivet               for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
155482490253SJoe Wallwork             }
155582490253SJoe Wallwork           }
155682490253SJoe Wallwork         }
15573ba16761SJacob Faibussowitsch         PetscCall(LAPACKsyevFail(dim, evecs));
1558*835f2295SStefano Zampini         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %" PetscBLASInt_FMT, lierr);
155982490253SJoe Wallwork       }
15609566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
156120282da2SJoe Wallwork 
156220282da2SJoe Wallwork       /* Modify eigenvalues */
1563e2606525SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
156420282da2SJoe Wallwork 
156520282da2SJoe Wallwork       /* Map back to get the intersection */
156620282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
1567e2606525SJoe Wallwork         for (m = 0; m < dim; ++m) {
1568e2606525SJoe Wallwork           M2[i * dim + m] = 0.0;
156920282da2SJoe Wallwork           for (j = 0; j < dim; ++j) {
157020282da2SJoe Wallwork             for (k = 0; k < dim; ++k) {
1571ad540459SPierre 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];
157220282da2SJoe Wallwork             }
157320282da2SJoe Wallwork           }
157420282da2SJoe Wallwork         }
157520282da2SJoe Wallwork       }
157620282da2SJoe Wallwork     }
15779566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
157820282da2SJoe Wallwork   }
1579e2606525SJoe Wallwork   PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
15803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
158120282da2SJoe Wallwork }
158220282da2SJoe Wallwork 
1583cb7bfe8cSJoe Wallwork /*@
158420282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
158520282da2SJoe Wallwork 
1586f1a722f8SMatthew G. Knepley   Input Parameters:
15872fe279fdSBarry Smith + dm         - The `DM`
158820282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
158920282da2SJoe Wallwork - metrics    - The metrics to be intersected
159020282da2SJoe Wallwork 
159120282da2SJoe Wallwork   Output Parameter:
159220282da2SJoe Wallwork . metricInt - The intersected metric
159320282da2SJoe Wallwork 
159420282da2SJoe Wallwork   Level: beginner
159520282da2SJoe Wallwork 
159620282da2SJoe Wallwork   Notes:
1597e2606525SJoe Wallwork   The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
159820282da2SJoe Wallwork 
1599e2606525SJoe Wallwork   The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
160020282da2SJoe Wallwork 
16012fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1602cb7bfe8cSJoe Wallwork @*/
DMPlexMetricIntersection(DM dm,PetscInt numMetrics,Vec metrics[],Vec metricInt)1603d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1604d71ae5a4SJacob Faibussowitsch {
160593520041SJoe Wallwork   PetscBool    isotropic, uniform;
160693520041SJoe Wallwork   PetscInt     v, i, m, n;
160793520041SJoe Wallwork   PetscScalar *met, *meti;
160820282da2SJoe Wallwork 
160920282da2SJoe Wallwork   PetscFunctionBegin;
16109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0));
161163a3b9bcSJacob Faibussowitsch   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
161220282da2SJoe Wallwork 
161320282da2SJoe Wallwork   /* Copy over the first metric */
16145241a91fSJoe Wallwork   PetscCall(VecCopy(metrics[0], metricInt));
16153ba16761SJacob Faibussowitsch   if (numMetrics == 1) PetscFunctionReturn(PETSC_SUCCESS);
16165241a91fSJoe Wallwork   PetscCall(VecGetSize(metricInt, &m));
161793520041SJoe Wallwork   for (i = 0; i < numMetrics; ++i) {
16189566063dSJacob Faibussowitsch     PetscCall(VecGetSize(metrics[i], &n));
161908401ef6SPierre Jolivet     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
162093520041SJoe Wallwork   }
162120282da2SJoe Wallwork 
162220282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
16239566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
16249566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
162593520041SJoe Wallwork   if (uniform) {
162693520041SJoe Wallwork     /* Uniform case */
16275241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
162893520041SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16299566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
16309566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
16319566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
163293520041SJoe Wallwork     }
16335241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
163493520041SJoe Wallwork   } else {
163593520041SJoe Wallwork     /* Spatially varying case */
163693520041SJoe Wallwork     PetscInt     dim, vStart, vEnd, nrow;
163793520041SJoe Wallwork     PetscScalar *M, *Mi;
163893520041SJoe Wallwork 
16399566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
164093520041SJoe Wallwork     if (isotropic) nrow = 1;
164193520041SJoe Wallwork     else nrow = dim;
16429566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
16435241a91fSJoe Wallwork     PetscCall(VecGetArray(metricInt, &met));
164420282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
16459566063dSJacob Faibussowitsch       PetscCall(VecGetArray(metrics[i], &meti));
164620282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
16479566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
16489566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
16499566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
165020282da2SJoe Wallwork       }
16519566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(metrics[i], &meti));
165220282da2SJoe Wallwork     }
16535241a91fSJoe Wallwork     PetscCall(VecRestoreArray(metricInt, &met));
165420282da2SJoe Wallwork   }
1655fe902aceSJoe Wallwork 
16569566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0));
16573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165820282da2SJoe Wallwork }
165920282da2SJoe Wallwork 
1660cb7bfe8cSJoe Wallwork /*@
166120282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
166220282da2SJoe Wallwork 
1663f1a722f8SMatthew G. Knepley   Input Parameters:
16642fe279fdSBarry Smith + dm      - The `DM`
166520282da2SJoe Wallwork . metric1 - The first metric to be intersected
166620282da2SJoe Wallwork - metric2 - The second metric to be intersected
166720282da2SJoe Wallwork 
166820282da2SJoe Wallwork   Output Parameter:
166920282da2SJoe Wallwork . metricInt - The intersected metric
167020282da2SJoe Wallwork 
167120282da2SJoe Wallwork   Level: beginner
167220282da2SJoe Wallwork 
16732fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1674cb7bfe8cSJoe Wallwork @*/
DMPlexMetricIntersection2(DM dm,Vec metric1,Vec metric2,Vec metricInt)1675d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1676d71ae5a4SJacob Faibussowitsch {
167720282da2SJoe Wallwork   Vec metrics[2] = {metric1, metric2};
167820282da2SJoe Wallwork 
167920282da2SJoe Wallwork   PetscFunctionBegin;
16809566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
16813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168220282da2SJoe Wallwork }
168320282da2SJoe Wallwork 
1684cb7bfe8cSJoe Wallwork /*@
168520282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
168620282da2SJoe Wallwork 
1687f1a722f8SMatthew G. Knepley   Input Parameters:
16882fe279fdSBarry Smith + dm      - The `DM`
168920282da2SJoe Wallwork . metric1 - The first metric to be intersected
169020282da2SJoe Wallwork . metric2 - The second metric to be intersected
169120282da2SJoe Wallwork - metric3 - The third metric to be intersected
169220282da2SJoe Wallwork 
169320282da2SJoe Wallwork   Output Parameter:
169420282da2SJoe Wallwork . metricInt - The intersected metric
169520282da2SJoe Wallwork 
169620282da2SJoe Wallwork   Level: beginner
169720282da2SJoe Wallwork 
16982fe279fdSBarry Smith .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1699cb7bfe8cSJoe Wallwork @*/
DMPlexMetricIntersection3(DM dm,Vec metric1,Vec metric2,Vec metric3,Vec metricInt)1700d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1701d71ae5a4SJacob Faibussowitsch {
170220282da2SJoe Wallwork   Vec metrics[3] = {metric1, metric2, metric3};
170320282da2SJoe Wallwork 
170420282da2SJoe Wallwork   PetscFunctionBegin;
17059566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
17063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170720282da2SJoe Wallwork }
1708