xref: /petsc/src/dm/impls/plex/plexmetric.c (revision cb7bfe8c620965d1e127b83ac87f7429f0b9029f)
120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
220282da2SJoe Wallwork #include <petscblaslapack.h>
320282da2SJoe Wallwork 
431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
531b70465SJoe Wallwork {
631b70465SJoe Wallwork   MPI_Comm       comm;
731b70465SJoe Wallwork   PetscBool      isotropic = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
8cc2eee55SJoe Wallwork   PetscBool      noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE;
931b70465SJoe Wallwork   PetscErrorCode ierr;
10cc2eee55SJoe Wallwork   PetscInt       verbosity = -1, numIter = 3;
11cc2eee55SJoe 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;
1231b70465SJoe Wallwork 
1331b70465SJoe Wallwork   PetscFunctionBegin;
1431b70465SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1531b70465SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr);
1631b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr);
17cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetIsotropic(dm, isotropic);CHKERRQ(ierr);
1831b70465SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr);
19cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst);CHKERRQ(ierr);
20cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL);CHKERRQ(ierr);
21cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoInsertion(dm, noInsert);CHKERRQ(ierr);
22cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL);CHKERRQ(ierr);
23cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoSwapping(dm, noSwap);CHKERRQ(ierr);
24cc2eee55SJoe Wallwork   ierr = PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL);CHKERRQ(ierr);
25cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNoMovement(dm, noMove);CHKERRQ(ierr);
26cc2eee55SJoe Wallwork   ierr = PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0);CHKERRQ(ierr);
27cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetNumIterations(dm, numIter);CHKERRQ(ierr);
28cc2eee55SJoe Wallwork   ierr = PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10);CHKERRQ(ierr);
29cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetVerbosity(dm, verbosity);CHKERRQ(ierr);
3031b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr);
3131b70465SJoe Wallwork   ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr);
3231b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr);
3331b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr);
3431b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr);
3531b70465SJoe Wallwork   ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr);
3631b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr);
3731b70465SJoe Wallwork   ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr);
3831b70465SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr);
3931b70465SJoe Wallwork   ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr);
40cc2eee55SJoe Wallwork   ierr = PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL);CHKERRQ(ierr);
41cc2eee55SJoe Wallwork   ierr = DMPlexMetricSetGradationFactor(dm, beta);CHKERRQ(ierr);
4231b70465SJoe Wallwork   ierr = PetscOptionsEnd();CHKERRQ(ierr);
4331b70465SJoe Wallwork   PetscFunctionReturn(0);
4431b70465SJoe Wallwork }
4531b70465SJoe Wallwork 
46*cb7bfe8cSJoe Wallwork /*@
4731b70465SJoe Wallwork   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
4831b70465SJoe Wallwork 
4931b70465SJoe Wallwork   Input parameters:
5031b70465SJoe Wallwork + dm        - The DM
5131b70465SJoe Wallwork - isotropic - Is the metric isotropic?
5231b70465SJoe Wallwork 
5331b70465SJoe Wallwork   Level: beginner
5431b70465SJoe Wallwork 
5531b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
56*cb7bfe8cSJoe Wallwork @*/
5731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
5831b70465SJoe Wallwork {
5931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
6031b70465SJoe Wallwork   PetscErrorCode ierr;
6131b70465SJoe Wallwork 
6231b70465SJoe Wallwork   PetscFunctionBegin;
6331b70465SJoe Wallwork   if (!plex->metricCtx) {
6431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
6531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
6631b70465SJoe Wallwork   }
6731b70465SJoe Wallwork   plex->metricCtx->isotropic = isotropic;
6831b70465SJoe Wallwork   PetscFunctionReturn(0);
6931b70465SJoe Wallwork }
7031b70465SJoe Wallwork 
71*cb7bfe8cSJoe Wallwork /*@
7231b70465SJoe Wallwork   DMPlexMetricIsIsotropic - Is a metric is isotropic?
7331b70465SJoe Wallwork 
7431b70465SJoe Wallwork   Input parameters:
7531b70465SJoe Wallwork . dm        - The DM
7631b70465SJoe Wallwork 
7731b70465SJoe Wallwork   Output parameters:
7831b70465SJoe Wallwork . isotropic - Is the metric isotropic?
7931b70465SJoe Wallwork 
8031b70465SJoe Wallwork   Level: beginner
8131b70465SJoe Wallwork 
8231b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
83*cb7bfe8cSJoe Wallwork @*/
8431b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
8531b70465SJoe Wallwork {
8631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
8731b70465SJoe Wallwork   PetscErrorCode ierr;
8831b70465SJoe Wallwork 
8931b70465SJoe Wallwork   PetscFunctionBegin;
9031b70465SJoe Wallwork   if (!plex->metricCtx) {
9131b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
9231b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
9331b70465SJoe Wallwork   }
9431b70465SJoe Wallwork   *isotropic = plex->metricCtx->isotropic;
9531b70465SJoe Wallwork   PetscFunctionReturn(0);
9631b70465SJoe Wallwork }
9731b70465SJoe Wallwork 
98*cb7bfe8cSJoe Wallwork /*@
9931b70465SJoe Wallwork   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
10031b70465SJoe Wallwork 
10131b70465SJoe Wallwork   Input parameters:
10231b70465SJoe Wallwork + dm                      - The DM
10331b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first?
10431b70465SJoe Wallwork 
10531b70465SJoe Wallwork   Level: beginner
10631b70465SJoe Wallwork 
10731b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst()
108*cb7bfe8cSJoe Wallwork @*/
10931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
11031b70465SJoe Wallwork {
11131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
11231b70465SJoe Wallwork   PetscErrorCode ierr;
11331b70465SJoe Wallwork 
11431b70465SJoe Wallwork   PetscFunctionBegin;
11531b70465SJoe Wallwork   if (!plex->metricCtx) {
11631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
11731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
11831b70465SJoe Wallwork   }
11931b70465SJoe Wallwork   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
12031b70465SJoe Wallwork   PetscFunctionReturn(0);
12131b70465SJoe Wallwork }
12231b70465SJoe Wallwork 
123*cb7bfe8cSJoe Wallwork /*@
12431b70465SJoe Wallwork   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
12531b70465SJoe Wallwork 
12631b70465SJoe Wallwork   Input parameters:
12731b70465SJoe Wallwork . dm                      - The DM
12831b70465SJoe Wallwork 
12931b70465SJoe Wallwork   Output parameters:
13031b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first?
13131b70465SJoe Wallwork 
13231b70465SJoe Wallwork   Level: beginner
13331b70465SJoe Wallwork 
13431b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst()
135*cb7bfe8cSJoe Wallwork @*/
13631b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
13731b70465SJoe Wallwork {
13831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
13931b70465SJoe Wallwork   PetscErrorCode ierr;
14031b70465SJoe Wallwork 
14131b70465SJoe Wallwork   PetscFunctionBegin;
14231b70465SJoe Wallwork   if (!plex->metricCtx) {
14331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
14431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
14531b70465SJoe Wallwork   }
14631b70465SJoe Wallwork   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
14731b70465SJoe Wallwork   PetscFunctionReturn(0);
14831b70465SJoe Wallwork }
14931b70465SJoe Wallwork 
150*cb7bfe8cSJoe Wallwork /*@
151cc2eee55SJoe Wallwork   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
152cc2eee55SJoe Wallwork 
153cc2eee55SJoe Wallwork   Input parameters:
154cc2eee55SJoe Wallwork + dm       - The DM
155cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off?
156cc2eee55SJoe Wallwork 
157cc2eee55SJoe Wallwork   Level: beginner
158cc2eee55SJoe Wallwork 
159*cb7bfe8cSJoe Wallwork   Notes:
160*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
161*cb7bfe8cSJoe Wallwork 
162cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement()
163*cb7bfe8cSJoe Wallwork @*/
164cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
165cc2eee55SJoe Wallwork {
166cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
167cc2eee55SJoe Wallwork   PetscErrorCode ierr;
168cc2eee55SJoe Wallwork 
169cc2eee55SJoe Wallwork   PetscFunctionBegin;
170cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
171cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
172cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
173cc2eee55SJoe Wallwork   }
174cc2eee55SJoe Wallwork   plex->metricCtx->noInsert = noInsert;
175cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
176cc2eee55SJoe Wallwork }
177cc2eee55SJoe Wallwork 
178*cb7bfe8cSJoe Wallwork /*@
179cc2eee55SJoe Wallwork   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
180cc2eee55SJoe Wallwork 
181cc2eee55SJoe Wallwork   Input parameters:
182cc2eee55SJoe Wallwork . dm       - The DM
183cc2eee55SJoe Wallwork 
184cc2eee55SJoe Wallwork   Output parameters:
185cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off?
186cc2eee55SJoe Wallwork 
187cc2eee55SJoe Wallwork   Level: beginner
188cc2eee55SJoe Wallwork 
189*cb7bfe8cSJoe Wallwork   Notes:
190*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
191*cb7bfe8cSJoe Wallwork 
192cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement()
193*cb7bfe8cSJoe Wallwork @*/
194cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
195cc2eee55SJoe Wallwork {
196cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
197cc2eee55SJoe Wallwork   PetscErrorCode ierr;
198cc2eee55SJoe Wallwork 
199cc2eee55SJoe Wallwork   PetscFunctionBegin;
200cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
201cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
202cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
203cc2eee55SJoe Wallwork   }
204cc2eee55SJoe Wallwork   *noInsert = plex->metricCtx->noInsert;
205cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
206cc2eee55SJoe Wallwork }
207cc2eee55SJoe Wallwork 
208*cb7bfe8cSJoe Wallwork /*@
209cc2eee55SJoe Wallwork   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
210cc2eee55SJoe Wallwork 
211cc2eee55SJoe Wallwork   Input parameters:
212cc2eee55SJoe Wallwork + dm     - The DM
213cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off?
214cc2eee55SJoe Wallwork 
215cc2eee55SJoe Wallwork   Level: beginner
216cc2eee55SJoe Wallwork 
217*cb7bfe8cSJoe Wallwork   Notes:
218*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
219*cb7bfe8cSJoe Wallwork 
220cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement()
221*cb7bfe8cSJoe Wallwork @*/
222cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
223cc2eee55SJoe Wallwork {
224cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
225cc2eee55SJoe Wallwork   PetscErrorCode ierr;
226cc2eee55SJoe Wallwork 
227cc2eee55SJoe Wallwork   PetscFunctionBegin;
228cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
229cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
230cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
231cc2eee55SJoe Wallwork   }
232cc2eee55SJoe Wallwork   plex->metricCtx->noSwap = noSwap;
233cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
234cc2eee55SJoe Wallwork }
235cc2eee55SJoe Wallwork 
236*cb7bfe8cSJoe Wallwork /*@
237cc2eee55SJoe Wallwork   DMPlexMetricNoSwapping - Is facet swapping turned off?
238cc2eee55SJoe Wallwork 
239cc2eee55SJoe Wallwork   Input parameters:
240cc2eee55SJoe Wallwork . dm     - The DM
241cc2eee55SJoe Wallwork 
242cc2eee55SJoe Wallwork   Output parameters:
243cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off?
244cc2eee55SJoe Wallwork 
245cc2eee55SJoe Wallwork   Level: beginner
246cc2eee55SJoe Wallwork 
247*cb7bfe8cSJoe Wallwork   Notes:
248*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
249*cb7bfe8cSJoe Wallwork 
250cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement()
251*cb7bfe8cSJoe Wallwork @*/
252cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
253cc2eee55SJoe Wallwork {
254cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
255cc2eee55SJoe Wallwork   PetscErrorCode ierr;
256cc2eee55SJoe Wallwork 
257cc2eee55SJoe Wallwork   PetscFunctionBegin;
258cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
259cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
260cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
261cc2eee55SJoe Wallwork   }
262cc2eee55SJoe Wallwork   *noSwap = plex->metricCtx->noSwap;
263cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
264cc2eee55SJoe Wallwork }
265cc2eee55SJoe Wallwork 
266*cb7bfe8cSJoe Wallwork /*@
267cc2eee55SJoe Wallwork   DMPlexMetricSetNoMovement - Should node movement be turned off?
268cc2eee55SJoe Wallwork 
269cc2eee55SJoe Wallwork   Input parameters:
270cc2eee55SJoe Wallwork + dm     - The DM
271cc2eee55SJoe Wallwork - noMove - Should node movement be turned off?
272cc2eee55SJoe Wallwork 
273cc2eee55SJoe Wallwork   Level: beginner
274cc2eee55SJoe Wallwork 
275*cb7bfe8cSJoe Wallwork   Notes:
276*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
277*cb7bfe8cSJoe Wallwork 
278cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping()
279*cb7bfe8cSJoe Wallwork @*/
280cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
281cc2eee55SJoe Wallwork {
282cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
283cc2eee55SJoe Wallwork   PetscErrorCode ierr;
284cc2eee55SJoe Wallwork 
285cc2eee55SJoe Wallwork   PetscFunctionBegin;
286cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
287cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
288cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
289cc2eee55SJoe Wallwork   }
290cc2eee55SJoe Wallwork   plex->metricCtx->noMove = noMove;
291cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
292cc2eee55SJoe Wallwork }
293cc2eee55SJoe Wallwork 
294*cb7bfe8cSJoe Wallwork /*@
295cc2eee55SJoe Wallwork   DMPlexMetricNoMovement - Is node movement turned off?
296cc2eee55SJoe Wallwork 
297cc2eee55SJoe Wallwork   Input parameters:
298cc2eee55SJoe Wallwork . dm     - The DM
299cc2eee55SJoe Wallwork 
300cc2eee55SJoe Wallwork   Output parameters:
301cc2eee55SJoe Wallwork . noMove - Is node movement turned off?
302cc2eee55SJoe Wallwork 
303cc2eee55SJoe Wallwork   Level: beginner
304cc2eee55SJoe Wallwork 
305*cb7bfe8cSJoe Wallwork   Notes:
306*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
307*cb7bfe8cSJoe Wallwork 
308cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping()
309*cb7bfe8cSJoe Wallwork @*/
310cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
311cc2eee55SJoe Wallwork {
312cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
313cc2eee55SJoe Wallwork   PetscErrorCode ierr;
314cc2eee55SJoe Wallwork 
315cc2eee55SJoe Wallwork   PetscFunctionBegin;
316cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
317cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
318cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
319cc2eee55SJoe Wallwork   }
320cc2eee55SJoe Wallwork   *noMove = plex->metricCtx->noMove;
321cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
322cc2eee55SJoe Wallwork }
323cc2eee55SJoe Wallwork 
324*cb7bfe8cSJoe Wallwork /*@
32531b70465SJoe Wallwork   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
32631b70465SJoe Wallwork 
32731b70465SJoe Wallwork   Input parameters:
32831b70465SJoe Wallwork + dm    - The DM
32931b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude
33031b70465SJoe Wallwork 
33131b70465SJoe Wallwork   Level: beginner
33231b70465SJoe Wallwork 
33331b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude()
334*cb7bfe8cSJoe Wallwork @*/
33531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
33631b70465SJoe Wallwork {
33731b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
33831b70465SJoe Wallwork   PetscErrorCode ierr;
33931b70465SJoe Wallwork 
34031b70465SJoe Wallwork   PetscFunctionBegin;
34131b70465SJoe Wallwork   if (!plex->metricCtx) {
34231b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
34331b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
34431b70465SJoe Wallwork   }
34531b70465SJoe Wallwork   if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min);
34631b70465SJoe Wallwork   plex->metricCtx->h_min = h_min;
34731b70465SJoe Wallwork   PetscFunctionReturn(0);
34831b70465SJoe Wallwork }
34931b70465SJoe Wallwork 
350*cb7bfe8cSJoe Wallwork /*@
35131b70465SJoe Wallwork   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
35231b70465SJoe Wallwork 
35331b70465SJoe Wallwork   Input parameters:
35431b70465SJoe Wallwork . dm    - The DM
35531b70465SJoe Wallwork 
356cc2eee55SJoe Wallwork   Output parameters:
35731b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude
35831b70465SJoe Wallwork 
35931b70465SJoe Wallwork   Level: beginner
36031b70465SJoe Wallwork 
36131b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude()
362*cb7bfe8cSJoe Wallwork @*/
36331b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
36431b70465SJoe Wallwork {
36531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
36631b70465SJoe Wallwork   PetscErrorCode ierr;
36731b70465SJoe Wallwork 
36831b70465SJoe Wallwork   PetscFunctionBegin;
36931b70465SJoe Wallwork   if (!plex->metricCtx) {
37031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
37131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
37231b70465SJoe Wallwork   }
37331b70465SJoe Wallwork   *h_min = plex->metricCtx->h_min;
37431b70465SJoe Wallwork   PetscFunctionReturn(0);
37531b70465SJoe Wallwork }
37631b70465SJoe Wallwork 
377*cb7bfe8cSJoe Wallwork /*@
37831b70465SJoe Wallwork   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
37931b70465SJoe Wallwork 
38031b70465SJoe Wallwork   Input parameters:
38131b70465SJoe Wallwork + dm    - The DM
38231b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude
38331b70465SJoe Wallwork 
38431b70465SJoe Wallwork   Level: beginner
38531b70465SJoe Wallwork 
38631b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude()
387*cb7bfe8cSJoe Wallwork @*/
38831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
38931b70465SJoe Wallwork {
39031b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
39131b70465SJoe Wallwork   PetscErrorCode ierr;
39231b70465SJoe Wallwork 
39331b70465SJoe Wallwork   PetscFunctionBegin;
39431b70465SJoe Wallwork   if (!plex->metricCtx) {
39531b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
39631b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
39731b70465SJoe Wallwork   }
39831b70465SJoe Wallwork   if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max);
39931b70465SJoe Wallwork   plex->metricCtx->h_max = h_max;
40031b70465SJoe Wallwork   PetscFunctionReturn(0);
40131b70465SJoe Wallwork }
40231b70465SJoe Wallwork 
403*cb7bfe8cSJoe Wallwork /*@
40431b70465SJoe Wallwork   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
40531b70465SJoe Wallwork 
40631b70465SJoe Wallwork   Input parameters:
40731b70465SJoe Wallwork . dm    - The DM
40831b70465SJoe Wallwork 
409cc2eee55SJoe Wallwork   Output parameters:
41031b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude
41131b70465SJoe Wallwork 
41231b70465SJoe Wallwork   Level: beginner
41331b70465SJoe Wallwork 
41431b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude()
415*cb7bfe8cSJoe Wallwork @*/
41631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
41731b70465SJoe Wallwork {
41831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
41931b70465SJoe Wallwork   PetscErrorCode ierr;
42031b70465SJoe Wallwork 
42131b70465SJoe Wallwork   PetscFunctionBegin;
42231b70465SJoe Wallwork   if (!plex->metricCtx) {
42331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
42431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
42531b70465SJoe Wallwork   }
42631b70465SJoe Wallwork   *h_max = plex->metricCtx->h_max;
42731b70465SJoe Wallwork   PetscFunctionReturn(0);
42831b70465SJoe Wallwork }
42931b70465SJoe Wallwork 
430*cb7bfe8cSJoe Wallwork /*@
43131b70465SJoe Wallwork   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
43231b70465SJoe Wallwork 
43331b70465SJoe Wallwork   Input parameters:
43431b70465SJoe Wallwork + dm    - The DM
43531b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy
43631b70465SJoe Wallwork 
43731b70465SJoe Wallwork   Level: beginner
43831b70465SJoe Wallwork 
43931b70465SJoe Wallwork   Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
44031b70465SJoe Wallwork 
44131b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude()
442*cb7bfe8cSJoe Wallwork @*/
44331b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
44431b70465SJoe Wallwork {
44531b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
44631b70465SJoe Wallwork   PetscErrorCode ierr;
44731b70465SJoe Wallwork 
44831b70465SJoe Wallwork   PetscFunctionBegin;
44931b70465SJoe Wallwork   if (!plex->metricCtx) {
45031b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
45131b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
45231b70465SJoe Wallwork   }
45331b70465SJoe Wallwork   if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max);
45431b70465SJoe Wallwork   plex->metricCtx->a_max = a_max;
45531b70465SJoe Wallwork   PetscFunctionReturn(0);
45631b70465SJoe Wallwork }
45731b70465SJoe Wallwork 
458*cb7bfe8cSJoe Wallwork /*@
45931b70465SJoe Wallwork   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
46031b70465SJoe Wallwork 
46131b70465SJoe Wallwork   Input parameters:
46231b70465SJoe Wallwork . dm    - The DM
46331b70465SJoe Wallwork 
464cc2eee55SJoe Wallwork   Output parameters:
46531b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy
46631b70465SJoe Wallwork 
46731b70465SJoe Wallwork   Level: beginner
46831b70465SJoe Wallwork 
46931b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude()
470*cb7bfe8cSJoe Wallwork @*/
47131b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
47231b70465SJoe Wallwork {
47331b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
47431b70465SJoe Wallwork   PetscErrorCode ierr;
47531b70465SJoe Wallwork 
47631b70465SJoe Wallwork   PetscFunctionBegin;
47731b70465SJoe Wallwork   if (!plex->metricCtx) {
47831b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
47931b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
48031b70465SJoe Wallwork   }
48131b70465SJoe Wallwork   *a_max = plex->metricCtx->a_max;
48231b70465SJoe Wallwork   PetscFunctionReturn(0);
48331b70465SJoe Wallwork }
48431b70465SJoe Wallwork 
485*cb7bfe8cSJoe Wallwork /*@
48631b70465SJoe Wallwork   DMPlexMetricSetTargetComplexity - Set the target metric complexity
48731b70465SJoe Wallwork 
48831b70465SJoe Wallwork   Input parameters:
48931b70465SJoe Wallwork + dm               - The DM
49031b70465SJoe Wallwork - targetComplexity - The target metric complexity
49131b70465SJoe Wallwork 
49231b70465SJoe Wallwork   Level: beginner
49331b70465SJoe Wallwork 
49431b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder()
495*cb7bfe8cSJoe Wallwork @*/
49631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
49731b70465SJoe Wallwork {
49831b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
49931b70465SJoe Wallwork   PetscErrorCode ierr;
50031b70465SJoe Wallwork 
50131b70465SJoe Wallwork   PetscFunctionBegin;
50231b70465SJoe Wallwork   if (!plex->metricCtx) {
50331b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
50431b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
50531b70465SJoe Wallwork   }
50631b70465SJoe Wallwork   if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive");
50731b70465SJoe Wallwork   plex->metricCtx->targetComplexity = targetComplexity;
50831b70465SJoe Wallwork   PetscFunctionReturn(0);
50931b70465SJoe Wallwork }
51031b70465SJoe Wallwork 
511*cb7bfe8cSJoe Wallwork /*@
51231b70465SJoe Wallwork   DMPlexMetricGetTargetComplexity - Get the target metric complexity
51331b70465SJoe Wallwork 
51431b70465SJoe Wallwork   Input parameters:
51531b70465SJoe Wallwork . dm               - The DM
51631b70465SJoe Wallwork 
517cc2eee55SJoe Wallwork   Output parameters:
51831b70465SJoe Wallwork . targetComplexity - The target metric complexity
51931b70465SJoe Wallwork 
52031b70465SJoe Wallwork   Level: beginner
52131b70465SJoe Wallwork 
52231b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder()
523*cb7bfe8cSJoe Wallwork @*/
52431b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
52531b70465SJoe Wallwork {
52631b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
52731b70465SJoe Wallwork   PetscErrorCode ierr;
52831b70465SJoe Wallwork 
52931b70465SJoe Wallwork   PetscFunctionBegin;
53031b70465SJoe Wallwork   if (!plex->metricCtx) {
53131b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
53231b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
53331b70465SJoe Wallwork   }
53431b70465SJoe Wallwork   *targetComplexity = plex->metricCtx->targetComplexity;
53531b70465SJoe Wallwork   PetscFunctionReturn(0);
53631b70465SJoe Wallwork }
53731b70465SJoe Wallwork 
538*cb7bfe8cSJoe Wallwork /*@
53931b70465SJoe Wallwork   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
54031b70465SJoe Wallwork 
54131b70465SJoe Wallwork   Input parameters:
54231b70465SJoe Wallwork + dm - The DM
54331b70465SJoe Wallwork - p  - The normalization order
54431b70465SJoe Wallwork 
54531b70465SJoe Wallwork   Level: beginner
54631b70465SJoe Wallwork 
54731b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity()
548*cb7bfe8cSJoe Wallwork @*/
54931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
55031b70465SJoe Wallwork {
55131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
55231b70465SJoe Wallwork   PetscErrorCode ierr;
55331b70465SJoe Wallwork 
55431b70465SJoe Wallwork   PetscFunctionBegin;
55531b70465SJoe Wallwork   if (!plex->metricCtx) {
55631b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
55731b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
55831b70465SJoe Wallwork   }
55931b70465SJoe Wallwork   if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater");
56031b70465SJoe Wallwork   plex->metricCtx->p = p;
56131b70465SJoe Wallwork   PetscFunctionReturn(0);
56231b70465SJoe Wallwork }
56331b70465SJoe Wallwork 
564*cb7bfe8cSJoe Wallwork /*@
56531b70465SJoe Wallwork   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
56631b70465SJoe Wallwork 
56731b70465SJoe Wallwork   Input parameters:
56831b70465SJoe Wallwork . dm - The DM
56931b70465SJoe Wallwork 
570cc2eee55SJoe Wallwork   Output parameters:
57131b70465SJoe Wallwork . p - The normalization order
57231b70465SJoe Wallwork 
57331b70465SJoe Wallwork   Level: beginner
57431b70465SJoe Wallwork 
57531b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity()
576*cb7bfe8cSJoe Wallwork @*/
57731b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
57831b70465SJoe Wallwork {
57931b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
58031b70465SJoe Wallwork   PetscErrorCode ierr;
58131b70465SJoe Wallwork 
58231b70465SJoe Wallwork   PetscFunctionBegin;
58331b70465SJoe Wallwork   if (!plex->metricCtx) {
58431b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
58531b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
58631b70465SJoe Wallwork   }
58731b70465SJoe Wallwork   *p = plex->metricCtx->p;
58831b70465SJoe Wallwork   PetscFunctionReturn(0);
58931b70465SJoe Wallwork }
59031b70465SJoe Wallwork 
591*cb7bfe8cSJoe Wallwork /*@
592cc2eee55SJoe Wallwork   DMPlexMetricSetGradationFactor - Set the metric gradation factor
593cc2eee55SJoe Wallwork 
594cc2eee55SJoe Wallwork   Input parameters:
595cc2eee55SJoe Wallwork + dm   - The DM
596cc2eee55SJoe Wallwork - beta - The metric gradation factor
597cc2eee55SJoe Wallwork 
598cc2eee55SJoe Wallwork   Level: beginner
599cc2eee55SJoe Wallwork 
600cc2eee55SJoe Wallwork   Notes:
601cc2eee55SJoe Wallwork 
602cc2eee55SJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
603cc2eee55SJoe Wallwork 
604cc2eee55SJoe Wallwork   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
605cc2eee55SJoe Wallwork 
606*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
607*cb7bfe8cSJoe Wallwork 
608cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetGradationFactor()
609*cb7bfe8cSJoe Wallwork @*/
610cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
611cc2eee55SJoe Wallwork {
612cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
613cc2eee55SJoe Wallwork   PetscErrorCode ierr;
614cc2eee55SJoe Wallwork 
615cc2eee55SJoe Wallwork   PetscFunctionBegin;
616cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
617cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
618cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
619cc2eee55SJoe Wallwork   }
620cc2eee55SJoe Wallwork   plex->metricCtx->gradationFactor = beta;
621cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
622cc2eee55SJoe Wallwork }
623cc2eee55SJoe Wallwork 
624*cb7bfe8cSJoe Wallwork /*@
625cc2eee55SJoe Wallwork   DMPlexMetricGetGradationFactor - Get the metric gradation factor
626cc2eee55SJoe Wallwork 
627cc2eee55SJoe Wallwork   Input parameters:
628cc2eee55SJoe Wallwork . dm   - The DM
629cc2eee55SJoe Wallwork 
630cc2eee55SJoe Wallwork   Output parameters:
631cc2eee55SJoe Wallwork . beta - The metric gradation factor
632cc2eee55SJoe Wallwork 
633cc2eee55SJoe Wallwork   Level: beginner
634cc2eee55SJoe Wallwork 
635*cb7bfe8cSJoe Wallwork   Notes:
636*cb7bfe8cSJoe Wallwork 
637*cb7bfe8cSJoe Wallwork   The gradation factor is the maximum tolerated length ratio between adjacent edges.
638*cb7bfe8cSJoe Wallwork 
639*cb7bfe8cSJoe Wallwork   The value -1 implies that gradation is turned off.
640*cb7bfe8cSJoe Wallwork 
641*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
642cc2eee55SJoe Wallwork 
643cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetGradationFactor()
644*cb7bfe8cSJoe Wallwork @*/
645cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
646cc2eee55SJoe Wallwork {
647cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
648cc2eee55SJoe Wallwork   PetscErrorCode ierr;
649cc2eee55SJoe Wallwork 
650cc2eee55SJoe Wallwork   PetscFunctionBegin;
651cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
652cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
653cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
654cc2eee55SJoe Wallwork   }
655cc2eee55SJoe Wallwork   *beta = plex->metricCtx->gradationFactor;
656cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
657cc2eee55SJoe Wallwork }
658cc2eee55SJoe Wallwork 
659*cb7bfe8cSJoe Wallwork /*@
660cc2eee55SJoe Wallwork   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
661cc2eee55SJoe Wallwork 
662cc2eee55SJoe Wallwork   Input parameters:
663cc2eee55SJoe Wallwork + dm        - The DM
664cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum
665cc2eee55SJoe Wallwork 
666*cb7bfe8cSJoe Wallwork   Level: beginner
667*cb7bfe8cSJoe Wallwork 
668*cb7bfe8cSJoe Wallwork   Notes:
669*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
670*cb7bfe8cSJoe Wallwork 
671cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations()
672*cb7bfe8cSJoe Wallwork @*/
673cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
674cc2eee55SJoe Wallwork {
675cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
676cc2eee55SJoe Wallwork   PetscErrorCode ierr;
677cc2eee55SJoe Wallwork 
678cc2eee55SJoe Wallwork   PetscFunctionBegin;
679cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
680cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
681cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
682cc2eee55SJoe Wallwork   }
683cc2eee55SJoe Wallwork   plex->metricCtx->verbosity = verbosity;
684cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
685cc2eee55SJoe Wallwork }
686cc2eee55SJoe Wallwork 
687*cb7bfe8cSJoe Wallwork /*@
688cc2eee55SJoe Wallwork   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
689cc2eee55SJoe Wallwork 
690cc2eee55SJoe Wallwork   Input parameters:
691cc2eee55SJoe Wallwork . dm        - The DM
692cc2eee55SJoe Wallwork 
693cc2eee55SJoe Wallwork   Output parameters:
694cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum
695cc2eee55SJoe Wallwork 
696*cb7bfe8cSJoe Wallwork   Level: beginner
697*cb7bfe8cSJoe Wallwork 
698*cb7bfe8cSJoe Wallwork   Notes:
699*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic).
700*cb7bfe8cSJoe Wallwork 
701cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
702*cb7bfe8cSJoe Wallwork @*/
703cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
704cc2eee55SJoe Wallwork {
705cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
706cc2eee55SJoe Wallwork   PetscErrorCode ierr;
707cc2eee55SJoe Wallwork 
708cc2eee55SJoe Wallwork   PetscFunctionBegin;
709cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
710cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
711cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
712cc2eee55SJoe Wallwork   }
713cc2eee55SJoe Wallwork   *verbosity = plex->metricCtx->verbosity;
714cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
715cc2eee55SJoe Wallwork }
716cc2eee55SJoe Wallwork 
717*cb7bfe8cSJoe Wallwork /*@
718cc2eee55SJoe Wallwork   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
719cc2eee55SJoe Wallwork 
720cc2eee55SJoe Wallwork   Input parameters:
721cc2eee55SJoe Wallwork + dm      - The DM
722cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations
723cc2eee55SJoe Wallwork 
724*cb7bfe8cSJoe Wallwork   Level: beginner
725*cb7bfe8cSJoe Wallwork 
726*cb7bfe8cSJoe Wallwork   Notes:
727*cb7bfe8cSJoe Wallwork   This is only used by ParMmg (not Pragmatic or Mmg).
728cc2eee55SJoe Wallwork 
729cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations()
730*cb7bfe8cSJoe Wallwork @*/
731cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
732cc2eee55SJoe Wallwork {
733cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
734cc2eee55SJoe Wallwork   PetscErrorCode ierr;
735cc2eee55SJoe Wallwork 
736cc2eee55SJoe Wallwork   PetscFunctionBegin;
737cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
738cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
739cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
740cc2eee55SJoe Wallwork   }
741cc2eee55SJoe Wallwork   plex->metricCtx->numIter = numIter;
742cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
743cc2eee55SJoe Wallwork }
744cc2eee55SJoe Wallwork 
745*cb7bfe8cSJoe Wallwork /*@
746cc2eee55SJoe Wallwork   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
747cc2eee55SJoe Wallwork 
748cc2eee55SJoe Wallwork   Input parameters:
749cc2eee55SJoe Wallwork . dm      - The DM
750cc2eee55SJoe Wallwork 
751cc2eee55SJoe Wallwork   Output parameters:
752cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations
753cc2eee55SJoe Wallwork 
754*cb7bfe8cSJoe Wallwork   Level: beginner
755*cb7bfe8cSJoe Wallwork 
756*cb7bfe8cSJoe Wallwork   Notes:
757*cb7bfe8cSJoe Wallwork   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
758cc2eee55SJoe Wallwork 
759cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity()
760*cb7bfe8cSJoe Wallwork @*/
761cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
762cc2eee55SJoe Wallwork {
763cc2eee55SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
764cc2eee55SJoe Wallwork   PetscErrorCode ierr;
765cc2eee55SJoe Wallwork 
766cc2eee55SJoe Wallwork   PetscFunctionBegin;
767cc2eee55SJoe Wallwork   if (!plex->metricCtx) {
768cc2eee55SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
769cc2eee55SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
770cc2eee55SJoe Wallwork   }
771cc2eee55SJoe Wallwork   *numIter = plex->metricCtx->numIter;
772cc2eee55SJoe Wallwork   PetscFunctionReturn(0);
773cc2eee55SJoe Wallwork }
774cc2eee55SJoe Wallwork 
77520282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
77620282da2SJoe Wallwork {
77720282da2SJoe Wallwork   MPI_Comm       comm;
77820282da2SJoe Wallwork   PetscErrorCode ierr;
77920282da2SJoe Wallwork   PetscFE        fe;
78020282da2SJoe Wallwork   PetscInt       dim;
78120282da2SJoe Wallwork 
78220282da2SJoe Wallwork   PetscFunctionBegin;
78320282da2SJoe Wallwork 
78420282da2SJoe Wallwork   /* Extract metadata from dm */
78520282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
78620282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
78720282da2SJoe Wallwork 
78820282da2SJoe Wallwork   /* Create a P1 field of the requested size */
78920282da2SJoe Wallwork   ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
79020282da2SJoe Wallwork   ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr);
79120282da2SJoe Wallwork   ierr = DMCreateDS(dm);CHKERRQ(ierr);
79220282da2SJoe Wallwork   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
79320282da2SJoe Wallwork   ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr);
79420282da2SJoe Wallwork 
79520282da2SJoe Wallwork   PetscFunctionReturn(0);
79620282da2SJoe Wallwork }
79720282da2SJoe Wallwork 
798*cb7bfe8cSJoe Wallwork /*@
79920282da2SJoe Wallwork   DMPlexMetricCreate - Create a Riemannian metric field
80020282da2SJoe Wallwork 
80120282da2SJoe Wallwork   Input parameters:
80220282da2SJoe Wallwork + dm     - The DM
80320282da2SJoe Wallwork - f      - The field number to use
80420282da2SJoe Wallwork 
80520282da2SJoe Wallwork   Output parameter:
80620282da2SJoe Wallwork . metric - The metric
80720282da2SJoe Wallwork 
80820282da2SJoe Wallwork   Level: beginner
80920282da2SJoe Wallwork 
81031b70465SJoe Wallwork   Notes:
81131b70465SJoe Wallwork 
81231b70465SJoe Wallwork   It is assumed that the DM is comprised of simplices.
81331b70465SJoe Wallwork 
81431b70465SJoe Wallwork   Command line options for Riemannian metrics:
81531b70465SJoe Wallwork 
816*cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic                 - Is the metric isotropic?
817*cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
818*cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
819*cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
820*cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
821*cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
822*cb7bfe8cSJoe Wallwork . -dm_plex_metric_target_complexity         - Target metric complexity
823*cb7bfe8cSJoe Wallwork 
824*cb7bfe8cSJoe Wallwork   Switching between remeshers can be achieved using
825*cb7bfe8cSJoe Wallwork 
826*cb7bfe8cSJoe Wallwork . -dm_adaptor <pragmatic/mmg/parmmg>
827*cb7bfe8cSJoe Wallwork 
828*cb7bfe8cSJoe Wallwork   Further options that are only relevant to Mmg and ParMmg:
829*cb7bfe8cSJoe Wallwork 
830*cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor          - Maximum ratio by which edge lengths may grow during gradation
831*cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations            - Number of parallel mesh adaptation iterations for ParMmg
832*cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert                 - Should node insertion/deletion be turned off?
833*cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap                   - Should facet swapping be turned off?
834*cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move                   - Should node movement be turned off?
835*cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity                 - Choose a verbosity level from -1 (silent) to 10 (maximum).
83620282da2SJoe Wallwork 
83720282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
838*cb7bfe8cSJoe Wallwork @*/
83920282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
84020282da2SJoe Wallwork {
84131b70465SJoe Wallwork   DM_Plex       *plex = (DM_Plex *) dm->data;
84220282da2SJoe Wallwork   PetscErrorCode ierr;
84320282da2SJoe Wallwork   PetscInt       coordDim, Nd;
84420282da2SJoe Wallwork 
84520282da2SJoe Wallwork   PetscFunctionBegin;
84631b70465SJoe Wallwork   if (!plex->metricCtx) {
84731b70465SJoe Wallwork     ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr);
84831b70465SJoe Wallwork     ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
84931b70465SJoe Wallwork   }
85020282da2SJoe Wallwork   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
85120282da2SJoe Wallwork   Nd = coordDim*coordDim;
85220282da2SJoe Wallwork   ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr);
85320282da2SJoe Wallwork   PetscFunctionReturn(0);
85420282da2SJoe Wallwork }
85520282da2SJoe Wallwork 
85620282da2SJoe Wallwork typedef struct {
85720282da2SJoe Wallwork   PetscReal scaling;  /* Scaling for uniform metric diagonal */
85820282da2SJoe Wallwork } DMPlexMetricUniformCtx;
85920282da2SJoe Wallwork 
86020282da2SJoe Wallwork static PetscErrorCode diagonal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
86120282da2SJoe Wallwork {
86220282da2SJoe Wallwork   DMPlexMetricUniformCtx *user = (DMPlexMetricUniformCtx*)ctx;
86320282da2SJoe Wallwork   PetscInt                i, j;
86420282da2SJoe Wallwork 
86520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
86620282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
86720282da2SJoe Wallwork       if (i == j) u[i+dim*j] = user->scaling;
86820282da2SJoe Wallwork       else u[i+dim*j] = 0.0;
86920282da2SJoe Wallwork     }
87020282da2SJoe Wallwork   }
87120282da2SJoe Wallwork   return 0;
87220282da2SJoe Wallwork }
87320282da2SJoe Wallwork 
874*cb7bfe8cSJoe Wallwork /*@
87520282da2SJoe Wallwork   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
87620282da2SJoe Wallwork 
87720282da2SJoe Wallwork   Input parameters:
87820282da2SJoe Wallwork + dm     - The DM
87920282da2SJoe Wallwork . f      - The field number to use
88020282da2SJoe Wallwork - alpha  - Scaling parameter for the diagonal
88120282da2SJoe Wallwork 
88220282da2SJoe Wallwork   Output parameter:
88320282da2SJoe Wallwork . metric - The uniform metric
88420282da2SJoe Wallwork 
88520282da2SJoe Wallwork   Level: beginner
88620282da2SJoe Wallwork 
88720282da2SJoe Wallwork   Note: It is assumed that the DM is comprised of simplices.
88820282da2SJoe Wallwork 
88920282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
890*cb7bfe8cSJoe Wallwork @*/
89120282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
89220282da2SJoe Wallwork {
89320282da2SJoe Wallwork   DMPlexMetricUniformCtx user;
89420282da2SJoe Wallwork   PetscErrorCode       (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
89520282da2SJoe Wallwork   PetscErrorCode         ierr;
89620282da2SJoe Wallwork   void                  *ctxs[1];
89720282da2SJoe Wallwork 
89820282da2SJoe Wallwork   PetscFunctionBegin;
89920282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
90020282da2SJoe Wallwork   if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
90120282da2SJoe Wallwork   if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
90220282da2SJoe Wallwork   else user.scaling = alpha;
90320282da2SJoe Wallwork   funcs[0] = diagonal;
90420282da2SJoe Wallwork   ctxs[0] = &user;
90520282da2SJoe Wallwork   ierr = DMProjectFunctionLocal(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, *metric);CHKERRQ(ierr);
90620282da2SJoe Wallwork   PetscFunctionReturn(0);
90720282da2SJoe Wallwork }
90820282da2SJoe Wallwork 
909*cb7bfe8cSJoe Wallwork /*@
91020282da2SJoe Wallwork   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
91120282da2SJoe Wallwork 
91220282da2SJoe Wallwork   Input parameters:
91320282da2SJoe Wallwork + dm        - The DM
91420282da2SJoe Wallwork . f         - The field number to use
91520282da2SJoe Wallwork - indicator - The error indicator
91620282da2SJoe Wallwork 
91720282da2SJoe Wallwork   Output parameter:
91820282da2SJoe Wallwork . metric    - The isotropic metric
91920282da2SJoe Wallwork 
92020282da2SJoe Wallwork   Level: beginner
92120282da2SJoe Wallwork 
92220282da2SJoe Wallwork   Notes:
92320282da2SJoe Wallwork 
92420282da2SJoe Wallwork   It is assumed that the DM is comprised of simplices.
92520282da2SJoe Wallwork 
92620282da2SJoe Wallwork   The indicator needs to be a scalar field defined at *vertices*.
92720282da2SJoe Wallwork 
92820282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
929*cb7bfe8cSJoe Wallwork @*/
93020282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
93120282da2SJoe Wallwork {
93220282da2SJoe Wallwork   DM                 dmIndi;
93320282da2SJoe Wallwork   PetscErrorCode     ierr;
93420282da2SJoe Wallwork   PetscInt           dim, d, vStart, vEnd, v;
93520282da2SJoe Wallwork   const PetscScalar *indi;
93620282da2SJoe Wallwork   PetscScalar       *met;
93720282da2SJoe Wallwork 
93820282da2SJoe Wallwork   PetscFunctionBegin;
93920282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
94020282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
94120282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
94220282da2SJoe Wallwork   ierr = VecGetArrayRead(indicator, &indi);CHKERRQ(ierr);
94320282da2SJoe Wallwork   ierr = VecGetArrayWrite(*metric, &met);CHKERRQ(ierr);
94420282da2SJoe Wallwork   ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr);
94520282da2SJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
94620282da2SJoe Wallwork     PetscScalar *vindi, *vmet;
94720282da2SJoe Wallwork     ierr = DMPlexPointLocalRead(dmIndi, v, indi, &vindi);CHKERRQ(ierr);
94820282da2SJoe Wallwork     ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
94920282da2SJoe Wallwork     for (d = 0; d < dim; ++d) vmet[d*(dim+1)] = vindi[0];
95020282da2SJoe Wallwork   }
95120282da2SJoe Wallwork   ierr = VecRestoreArrayWrite(*metric, &met);CHKERRQ(ierr);
95220282da2SJoe Wallwork   ierr = VecRestoreArrayRead(indicator, &indi);CHKERRQ(ierr);
95320282da2SJoe Wallwork   PetscFunctionReturn(0);
95420282da2SJoe Wallwork }
95520282da2SJoe Wallwork 
95682490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
95782490253SJoe Wallwork {
95882490253SJoe Wallwork   PetscInt i, j;
95982490253SJoe Wallwork 
96082490253SJoe Wallwork   PetscFunctionBegin;
96182490253SJoe Wallwork   PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n");
96282490253SJoe Wallwork   for (i = 0; i < dim; ++i) {
96382490253SJoe Wallwork     if (i == 0) PetscPrintf(PETSC_COMM_SELF, "    [[");
96482490253SJoe Wallwork     else        PetscPrintf(PETSC_COMM_SELF, "     [");
96582490253SJoe Wallwork     for (j = 0; j < dim; ++j) {
96682490253SJoe Wallwork       if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]);
96782490253SJoe Wallwork       else           PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]);
96882490253SJoe Wallwork     }
96982490253SJoe Wallwork     if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n");
97082490253SJoe Wallwork     else           PetscPrintf(PETSC_COMM_SELF, "]]\n");
97182490253SJoe Wallwork   }
97282490253SJoe Wallwork   PetscFunctionReturn(0);
97382490253SJoe Wallwork }
97482490253SJoe Wallwork 
9753f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
97620282da2SJoe Wallwork {
97720282da2SJoe Wallwork   PetscErrorCode ierr;
97820282da2SJoe Wallwork   PetscInt       i, j, k;
97920282da2SJoe 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);
98020282da2SJoe Wallwork   PetscScalar   *Mpos;
98120282da2SJoe Wallwork 
98220282da2SJoe Wallwork   PetscFunctionBegin;
98320282da2SJoe Wallwork   ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr);
98420282da2SJoe Wallwork 
98520282da2SJoe Wallwork   /* Symmetrize */
98620282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
98720282da2SJoe Wallwork     Mpos[i*dim+i] = Mp[i*dim+i];
98820282da2SJoe Wallwork     for (j = i+1; j < dim; ++j) {
98920282da2SJoe Wallwork       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
99020282da2SJoe Wallwork       Mpos[j*dim+i] = Mpos[i*dim+j];
99120282da2SJoe Wallwork     }
99220282da2SJoe Wallwork   }
99320282da2SJoe Wallwork 
99420282da2SJoe Wallwork   /* Compute eigendecomposition */
99520282da2SJoe Wallwork   {
99620282da2SJoe Wallwork     PetscScalar  *work;
99720282da2SJoe Wallwork     PetscBLASInt lwork;
99820282da2SJoe Wallwork 
99920282da2SJoe Wallwork     lwork = 5*dim;
100020282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
100120282da2SJoe Wallwork     {
100220282da2SJoe Wallwork       PetscBLASInt lierr;
100320282da2SJoe Wallwork       PetscBLASInt nb;
100420282da2SJoe Wallwork 
100520282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
100620282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
100720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
100820282da2SJoe Wallwork       {
100920282da2SJoe Wallwork         PetscReal *rwork;
101020282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
101120282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
101220282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
101320282da2SJoe Wallwork       }
101420282da2SJoe Wallwork #else
101520282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
101620282da2SJoe Wallwork #endif
101782490253SJoe Wallwork       if (lierr) {
101882490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
101982490253SJoe Wallwork           Mpos[i*dim+i] = Mp[i*dim+i];
102082490253SJoe Wallwork           for (j = i+1; j < dim; ++j) {
102182490253SJoe Wallwork             Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
102282490253SJoe Wallwork             Mpos[j*dim+i] = Mpos[i*dim+j];
102382490253SJoe Wallwork           }
102482490253SJoe Wallwork         }
102582490253SJoe Wallwork         LAPACKsyevFail(dim, Mpos);
102682490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
102782490253SJoe Wallwork       }
102820282da2SJoe Wallwork       ierr = PetscFPTrapPop();CHKERRQ(ierr);
102920282da2SJoe Wallwork     }
103020282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
103120282da2SJoe Wallwork   }
103220282da2SJoe Wallwork 
103320282da2SJoe Wallwork   /* Reflect to positive orthant and enforce maximum and minimum size */
103420282da2SJoe Wallwork   max_eig = 0.0;
103520282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
103620282da2SJoe Wallwork     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
103720282da2SJoe Wallwork     max_eig = PetscMax(eigs[i], max_eig);
103820282da2SJoe Wallwork   }
103920282da2SJoe Wallwork 
10403f07679eSJoe Wallwork   /* Enforce maximum anisotropy and compute determinant */
10413f07679eSJoe Wallwork   *detMp = 1.0;
104220282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
104320282da2SJoe Wallwork     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
10443f07679eSJoe Wallwork     *detMp *= eigs[i];
104520282da2SJoe Wallwork   }
104620282da2SJoe Wallwork 
104720282da2SJoe Wallwork   /* Reconstruct Hessian */
104820282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
104920282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
105020282da2SJoe Wallwork       Mp[i*dim+j] = 0.0;
105120282da2SJoe Wallwork       for (k = 0; k < dim; ++k) {
105220282da2SJoe Wallwork         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
105320282da2SJoe Wallwork       }
105420282da2SJoe Wallwork     }
105520282da2SJoe Wallwork   }
105620282da2SJoe Wallwork   ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr);
105720282da2SJoe Wallwork 
105820282da2SJoe Wallwork   PetscFunctionReturn(0);
105920282da2SJoe Wallwork }
106020282da2SJoe Wallwork 
1061*cb7bfe8cSJoe Wallwork /*@
106220282da2SJoe Wallwork   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
106320282da2SJoe Wallwork 
106420282da2SJoe Wallwork   Input parameters:
106520282da2SJoe Wallwork + dm                 - The DM
10663f07679eSJoe Wallwork . metricIn           - The metric
106799abec2bSJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
10683f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced?
106920282da2SJoe Wallwork 
107020282da2SJoe Wallwork   Output parameter:
10713f07679eSJoe Wallwork + metricOut          - The metric
10723f07679eSJoe Wallwork - determinant        - Its determinant
107320282da2SJoe Wallwork 
107420282da2SJoe Wallwork   Level: beginner
107520282da2SJoe Wallwork 
1076*cb7bfe8cSJoe Wallwork   Notes:
1077*cb7bfe8cSJoe Wallwork 
1078*cb7bfe8cSJoe Wallwork   Relevant command line options:
1079*cb7bfe8cSJoe Wallwork 
1080*cb7bfe8cSJoe Wallwork + -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1081*cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1082*cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1083*cb7bfe8cSJoe Wallwork 
108420282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
1085*cb7bfe8cSJoe Wallwork @*/
10863f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant)
108720282da2SJoe Wallwork {
10883f07679eSJoe Wallwork   DM             dmDet;
108920282da2SJoe Wallwork   PetscErrorCode ierr;
109020282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v;
10913f07679eSJoe Wallwork   PetscScalar   *met, *det;
109220282da2SJoe Wallwork   PetscReal      h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
109320282da2SJoe Wallwork 
109420282da2SJoe Wallwork   PetscFunctionBegin;
109520282da2SJoe Wallwork 
109620282da2SJoe Wallwork   /* Extract metadata from dm */
109720282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
109820282da2SJoe Wallwork   if (restrictSizes) {
109999abec2bSJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
110099abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
110199abec2bSJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
110299abec2bSJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
110399abec2bSJoe Wallwork     if (h_min >= h_max) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
110499abec2bSJoe Wallwork   }
110599abec2bSJoe Wallwork   if (restrictAnisotropy) {
110699abec2bSJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
110799abec2bSJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
110820282da2SJoe Wallwork   }
110920282da2SJoe Wallwork 
11103f07679eSJoe Wallwork   /* Setup output metric and determinant */
11113f07679eSJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr);
11123f07679eSJoe Wallwork   ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr);
11133f07679eSJoe Wallwork   ierr = DMClone(dm, &dmDet);CHKERRQ(ierr);
11143f07679eSJoe Wallwork   ierr = DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant);CHKERRQ(ierr);
11153f07679eSJoe Wallwork 
111620282da2SJoe Wallwork   /* Enforce SPD */
111720282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
11183f07679eSJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
11193f07679eSJoe Wallwork   ierr = VecGetArray(*determinant, &det);CHKERRQ(ierr);
112020282da2SJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
11213f07679eSJoe Wallwork     PetscScalar *vmet, *vdet;
112220282da2SJoe Wallwork     ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
11233f07679eSJoe Wallwork     ierr = DMPlexPointLocalRef(dmDet, v, det, &vdet);CHKERRQ(ierr);
11243f07679eSJoe Wallwork     ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, vmet, vdet);CHKERRQ(ierr);
112520282da2SJoe Wallwork   }
11263f07679eSJoe Wallwork   ierr = VecRestoreArray(*determinant, &det);CHKERRQ(ierr);
11273f07679eSJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
112820282da2SJoe Wallwork 
112920282da2SJoe Wallwork   PetscFunctionReturn(0);
113020282da2SJoe Wallwork }
113120282da2SJoe Wallwork 
113220282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
113320282da2SJoe Wallwork                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
113420282da2SJoe Wallwork                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
113520282da2SJoe Wallwork                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
113620282da2SJoe Wallwork {
113720282da2SJoe Wallwork   const PetscScalar p = constants[0];
113820282da2SJoe Wallwork 
11393f07679eSJoe Wallwork   f0[0] = PetscPowReal(u[0], p/(2.0*p + dim));
114020282da2SJoe Wallwork }
114120282da2SJoe Wallwork 
1142*cb7bfe8cSJoe Wallwork /*@
114320282da2SJoe Wallwork   DMPlexMetricNormalize - Apply L-p normalization to a metric
114420282da2SJoe Wallwork 
114520282da2SJoe Wallwork   Input parameters:
114620282da2SJoe Wallwork + dm                 - The DM
114720282da2SJoe Wallwork . metricIn           - The unnormalized metric
114816522936SJoe Wallwork . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
114916522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced?
115020282da2SJoe Wallwork 
115120282da2SJoe Wallwork   Output parameter:
115220282da2SJoe Wallwork . metricOut          - The normalized metric
115320282da2SJoe Wallwork 
115420282da2SJoe Wallwork   Level: beginner
115520282da2SJoe Wallwork 
1156*cb7bfe8cSJoe Wallwork   Notes:
1157*cb7bfe8cSJoe Wallwork 
1158*cb7bfe8cSJoe Wallwork   Relevant command line options:
1159*cb7bfe8cSJoe Wallwork 
1160*cb7bfe8cSJoe Wallwork + -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1161*cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1162*cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1163*cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1164*cb7bfe8cSJoe Wallwork . -dm_plex_metric_p                         - L-p normalization order
1165*cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity         - Target metric complexity
1166*cb7bfe8cSJoe Wallwork 
116720282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
1168*cb7bfe8cSJoe Wallwork @*/
116916522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut)
117020282da2SJoe Wallwork {
11713f07679eSJoe Wallwork   DM               dmDet;
117220282da2SJoe Wallwork   MPI_Comm         comm;
117316522936SJoe Wallwork   PetscBool        restrictAnisotropyFirst;
117420282da2SJoe Wallwork   PetscDS          ds;
117520282da2SJoe Wallwork   PetscErrorCode   ierr;
117620282da2SJoe Wallwork   PetscInt         dim, Nd, vStart, vEnd, v, i;
11773f07679eSJoe Wallwork   PetscScalar     *met, *det, integral, constants[1];
11783f07679eSJoe Wallwork   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, target, realIntegral;
11793f07679eSJoe Wallwork   Vec              determinant;
118020282da2SJoe Wallwork 
118120282da2SJoe Wallwork   PetscFunctionBegin;
118220282da2SJoe Wallwork 
118320282da2SJoe Wallwork   /* Extract metadata from dm */
118420282da2SJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
118520282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
118620282da2SJoe Wallwork   Nd = dim*dim;
118720282da2SJoe Wallwork 
118820282da2SJoe Wallwork   /* Set up metric and ensure it is SPD */
118916522936SJoe Wallwork   ierr = DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst);CHKERRQ(ierr);
11903f07679eSJoe Wallwork   ierr = DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant);CHKERRQ(ierr);
119120282da2SJoe Wallwork 
119220282da2SJoe Wallwork   /* Compute global normalization factor */
119316522936SJoe Wallwork   ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
119416522936SJoe Wallwork   ierr = DMPlexMetricGetNormalizationOrder(dm, &p);CHKERRQ(ierr);
119516522936SJoe Wallwork   constants[0] = p;
11963f07679eSJoe Wallwork   ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr);
11973f07679eSJoe Wallwork   ierr = DMGetDS(dmDet, &ds);CHKERRQ(ierr);
119820282da2SJoe Wallwork   ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
119920282da2SJoe Wallwork   ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
12003f07679eSJoe Wallwork   ierr = DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL);CHKERRQ(ierr);
12013f07679eSJoe Wallwork   realIntegral = PetscRealPart(integral);
12023f07679eSJoe Wallwork   if (realIntegral < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global metric normalization factor should be strictly positive, not %.4e Is the input metric positive-definite?", realIntegral);
12033f07679eSJoe Wallwork   factGlob = PetscPowReal(target/realIntegral, 2.0/dim);
120420282da2SJoe Wallwork 
120520282da2SJoe Wallwork   /* Apply local scaling */
120620282da2SJoe Wallwork   if (restrictSizes) {
120716522936SJoe Wallwork     ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr);
120816522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr);
120916522936SJoe Wallwork     h_min = PetscMax(h_min, 1.0e-30);
121016522936SJoe Wallwork     h_max = PetscMin(h_max, 1.0e+30);
121116522936SJoe Wallwork     if (h_min >= h_max) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max);
121216522936SJoe Wallwork   }
121316522936SJoe Wallwork   if (restrictAnisotropy && !restrictAnisotropyFirst) {
121416522936SJoe Wallwork     ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr);
121516522936SJoe Wallwork     a_max = PetscMin(a_max, 1.0e+30);
121620282da2SJoe Wallwork   }
121720282da2SJoe Wallwork   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
12183f07679eSJoe Wallwork   ierr = VecGetArray(determinant, &det);CHKERRQ(ierr);
121916522936SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
122020282da2SJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
12213f07679eSJoe Wallwork     PetscScalar *Mp, *detM;
12223f07679eSJoe Wallwork     PetscReal    fact;
122320282da2SJoe Wallwork 
122420282da2SJoe Wallwork     ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr);
12253f07679eSJoe Wallwork     ierr = DMPlexPointLocalRef(dmDet, v, det, &detM);CHKERRQ(ierr);
12263f07679eSJoe Wallwork     fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim));
122720282da2SJoe Wallwork     for (i = 0; i < Nd; ++i) Mp[i] *= fact;
12283f07679eSJoe Wallwork     if (restrictSizes) { ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, Mp, detM);CHKERRQ(ierr); }
122920282da2SJoe Wallwork   }
12303f07679eSJoe Wallwork   ierr = VecRestoreArray(determinant, &det);CHKERRQ(ierr);
123120282da2SJoe Wallwork   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
12323f07679eSJoe Wallwork   ierr = VecDestroy(&determinant);CHKERRQ(ierr);
12333f07679eSJoe Wallwork   ierr = DMDestroy(&dmDet);CHKERRQ(ierr);
123420282da2SJoe Wallwork 
123520282da2SJoe Wallwork   PetscFunctionReturn(0);
123620282da2SJoe Wallwork }
123720282da2SJoe Wallwork 
1238*cb7bfe8cSJoe Wallwork /*@
123920282da2SJoe Wallwork   DMPlexMetricAverage - Compute the average of a list of metrics
124020282da2SJoe Wallwork 
124120282da2SJoe Wallwork   Input Parameter:
124220282da2SJoe Wallwork + dm         - The DM
124320282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged
124420282da2SJoe Wallwork . weights    - Weights for the average
124520282da2SJoe Wallwork - metrics    - The metrics to be averaged
124620282da2SJoe Wallwork 
124720282da2SJoe Wallwork   Output Parameter:
124820282da2SJoe Wallwork . metricAvg  - The averaged metric
124920282da2SJoe Wallwork 
125020282da2SJoe Wallwork   Level: beginner
125120282da2SJoe Wallwork 
125220282da2SJoe Wallwork   Notes:
125320282da2SJoe Wallwork   The weights should sum to unity.
125420282da2SJoe Wallwork 
125520282da2SJoe Wallwork   If weights are not provided then an unweighted average is used.
125620282da2SJoe Wallwork 
125720282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
1258*cb7bfe8cSJoe Wallwork @*/
125920282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
126020282da2SJoe Wallwork {
126120282da2SJoe Wallwork   PetscBool      haveWeights = PETSC_TRUE;
126220282da2SJoe Wallwork   PetscErrorCode ierr;
126320282da2SJoe Wallwork   PetscInt       i;
126420282da2SJoe Wallwork   PetscReal      sum = 0.0, tol = 1.0e-10;
126520282da2SJoe Wallwork 
126620282da2SJoe Wallwork   PetscFunctionBegin;
126720282da2SJoe Wallwork   if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); }
126820282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr);
126920282da2SJoe Wallwork   ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr);
127020282da2SJoe Wallwork 
127120282da2SJoe Wallwork   /* Default to the unweighted case */
127220282da2SJoe Wallwork   if (!weights) {
127320282da2SJoe Wallwork     ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr);
127420282da2SJoe Wallwork     haveWeights = PETSC_FALSE;
127520282da2SJoe Wallwork     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
127620282da2SJoe Wallwork   }
127720282da2SJoe Wallwork 
127820282da2SJoe Wallwork   /* Check weights sum to unity */
127920282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { sum += weights[i]; }
128020282da2SJoe Wallwork   if (PetscAbsReal(sum - 1) > tol) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); }
128120282da2SJoe Wallwork 
128220282da2SJoe Wallwork   /* Compute metric average */
128320282da2SJoe Wallwork   for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); }
128420282da2SJoe Wallwork   if (!haveWeights) {ierr = PetscFree(weights); }
128520282da2SJoe Wallwork   PetscFunctionReturn(0);
128620282da2SJoe Wallwork }
128720282da2SJoe Wallwork 
1288*cb7bfe8cSJoe Wallwork /*@
128920282da2SJoe Wallwork   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
129020282da2SJoe Wallwork 
129120282da2SJoe Wallwork   Input Parameter:
129220282da2SJoe Wallwork + dm         - The DM
129320282da2SJoe Wallwork . metric1    - The first metric to be averaged
129420282da2SJoe Wallwork - metric2    - The second metric to be averaged
129520282da2SJoe Wallwork 
129620282da2SJoe Wallwork   Output Parameter:
129720282da2SJoe Wallwork . metricAvg  - The averaged metric
129820282da2SJoe Wallwork 
129920282da2SJoe Wallwork   Level: beginner
130020282da2SJoe Wallwork 
130120282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
1302*cb7bfe8cSJoe Wallwork @*/
130320282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
130420282da2SJoe Wallwork {
130520282da2SJoe Wallwork   PetscErrorCode ierr;
130620282da2SJoe Wallwork   PetscReal      weights[2] = {0.5, 0.5};
130720282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
130820282da2SJoe Wallwork 
130920282da2SJoe Wallwork   PetscFunctionBegin;
131020282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr);
131120282da2SJoe Wallwork   PetscFunctionReturn(0);
131220282da2SJoe Wallwork }
131320282da2SJoe Wallwork 
1314*cb7bfe8cSJoe Wallwork /*@
131520282da2SJoe Wallwork   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
131620282da2SJoe Wallwork 
131720282da2SJoe Wallwork   Input Parameter:
131820282da2SJoe Wallwork + dm         - The DM
131920282da2SJoe Wallwork . metric1    - The first metric to be averaged
132020282da2SJoe Wallwork . metric2    - The second metric to be averaged
132120282da2SJoe Wallwork - metric3    - The third metric to be averaged
132220282da2SJoe Wallwork 
132320282da2SJoe Wallwork   Output Parameter:
132420282da2SJoe Wallwork . metricAvg  - The averaged metric
132520282da2SJoe Wallwork 
132620282da2SJoe Wallwork   Level: beginner
132720282da2SJoe Wallwork 
132820282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
1329*cb7bfe8cSJoe Wallwork @*/
133020282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
133120282da2SJoe Wallwork {
133220282da2SJoe Wallwork   PetscErrorCode ierr;
133320282da2SJoe Wallwork   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
133420282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
133520282da2SJoe Wallwork 
133620282da2SJoe Wallwork   PetscFunctionBegin;
133720282da2SJoe Wallwork   ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr);
133820282da2SJoe Wallwork   PetscFunctionReturn(0);
133920282da2SJoe Wallwork }
134020282da2SJoe Wallwork 
134120282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
134220282da2SJoe Wallwork {
134320282da2SJoe Wallwork   PetscErrorCode ierr;
134420282da2SJoe Wallwork   PetscInt       i, j, k, l, m;
134520282da2SJoe Wallwork   PetscReal     *evals, *evals1;
134620282da2SJoe Wallwork   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
134720282da2SJoe Wallwork 
134820282da2SJoe Wallwork   PetscFunctionBegin;
134920282da2SJoe Wallwork   ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr);
135020282da2SJoe Wallwork   for (i = 0; i < dim; ++i) {
135120282da2SJoe Wallwork     for (j = 0; j < dim; ++j) {
135220282da2SJoe Wallwork       evecs[i*dim+j] = M1[i*dim+j];
135320282da2SJoe Wallwork     }
135420282da2SJoe Wallwork   }
135520282da2SJoe Wallwork   {
135620282da2SJoe Wallwork     PetscScalar *work;
135720282da2SJoe Wallwork     PetscBLASInt lwork;
135820282da2SJoe Wallwork 
135920282da2SJoe Wallwork     lwork = 5*dim;
136020282da2SJoe Wallwork     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
136120282da2SJoe Wallwork     {
136220282da2SJoe Wallwork       PetscBLASInt lierr, nb;
136320282da2SJoe Wallwork       PetscReal    sqrtk;
136420282da2SJoe Wallwork 
136520282da2SJoe Wallwork       /* Compute eigendecomposition of M1 */
136620282da2SJoe Wallwork       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
136720282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
136820282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
136920282da2SJoe Wallwork       {
137020282da2SJoe Wallwork         PetscReal *rwork;
137120282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
137220282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
137320282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
137420282da2SJoe Wallwork       }
137520282da2SJoe Wallwork #else
137620282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
137720282da2SJoe Wallwork #endif
137882490253SJoe Wallwork       if (lierr) {
137982490253SJoe Wallwork         LAPACKsyevFail(dim, M1);
138082490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
138182490253SJoe Wallwork       }
138220282da2SJoe Wallwork       ierr = PetscFPTrapPop();
138320282da2SJoe Wallwork 
138420282da2SJoe Wallwork       /* Compute square root and reciprocal */
138520282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
138620282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
138720282da2SJoe Wallwork           sqrtM1[i*dim+j] = 0.0;
138820282da2SJoe Wallwork           isqrtM1[i*dim+j] = 0.0;
138920282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
139020282da2SJoe Wallwork             sqrtk = PetscSqrtReal(evals1[k]);
139120282da2SJoe Wallwork             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
139220282da2SJoe Wallwork             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
139320282da2SJoe Wallwork           }
139420282da2SJoe Wallwork         }
139520282da2SJoe Wallwork       }
139620282da2SJoe Wallwork 
139720282da2SJoe Wallwork       /* Map into the space spanned by the eigenvectors of M1 */
139820282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
139920282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
140020282da2SJoe Wallwork           evecs[i*dim+j] = 0.0;
140120282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
140220282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
140320282da2SJoe Wallwork               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
140420282da2SJoe Wallwork             }
140520282da2SJoe Wallwork           }
140620282da2SJoe Wallwork         }
140720282da2SJoe Wallwork       }
140820282da2SJoe Wallwork 
140920282da2SJoe Wallwork       /* Compute eigendecomposition */
141020282da2SJoe Wallwork       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
141120282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX)
141220282da2SJoe Wallwork       {
141320282da2SJoe Wallwork         PetscReal *rwork;
141420282da2SJoe Wallwork         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
141520282da2SJoe Wallwork         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
141620282da2SJoe Wallwork         ierr = PetscFree(rwork);CHKERRQ(ierr);
141720282da2SJoe Wallwork       }
141820282da2SJoe Wallwork #else
141920282da2SJoe Wallwork       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
142020282da2SJoe Wallwork #endif
142182490253SJoe Wallwork       if (lierr) {
142282490253SJoe Wallwork         for (i = 0; i < dim; ++i) {
142382490253SJoe Wallwork           for (j = 0; j < dim; ++j) {
142482490253SJoe Wallwork             evecs[i*dim+j] = 0.0;
142582490253SJoe Wallwork             for (k = 0; k < dim; ++k) {
142682490253SJoe Wallwork               for (l = 0; l < dim; ++l) {
142782490253SJoe Wallwork                 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
142882490253SJoe Wallwork               }
142982490253SJoe Wallwork             }
143082490253SJoe Wallwork           }
143182490253SJoe Wallwork         }
143282490253SJoe Wallwork         LAPACKsyevFail(dim, evecs);
143382490253SJoe Wallwork         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
143482490253SJoe Wallwork       }
143520282da2SJoe Wallwork       ierr = PetscFPTrapPop();
143620282da2SJoe Wallwork 
143720282da2SJoe Wallwork       /* Modify eigenvalues */
143820282da2SJoe Wallwork       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
143920282da2SJoe Wallwork 
144020282da2SJoe Wallwork       /* Map back to get the intersection */
144120282da2SJoe Wallwork       for (i = 0; i < dim; ++i) {
144220282da2SJoe Wallwork         for (j = 0; j < dim; ++j) {
144320282da2SJoe Wallwork           M2[i*dim+j] = 0.0;
144420282da2SJoe Wallwork           for (k = 0; k < dim; ++k) {
144520282da2SJoe Wallwork             for (l = 0; l < dim; ++l) {
144620282da2SJoe Wallwork               for (m = 0; m < dim; ++m) {
144720282da2SJoe Wallwork                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
144820282da2SJoe Wallwork               }
144920282da2SJoe Wallwork             }
145020282da2SJoe Wallwork           }
145120282da2SJoe Wallwork         }
145220282da2SJoe Wallwork       }
145320282da2SJoe Wallwork     }
145420282da2SJoe Wallwork     ierr = PetscFree(work);CHKERRQ(ierr);
145520282da2SJoe Wallwork   }
145620282da2SJoe Wallwork   ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr);
145720282da2SJoe Wallwork   PetscFunctionReturn(0);
145820282da2SJoe Wallwork }
145920282da2SJoe Wallwork 
1460*cb7bfe8cSJoe Wallwork /*@
146120282da2SJoe Wallwork   DMPlexMetricIntersection - Compute the intersection of a list of metrics
146220282da2SJoe Wallwork 
146320282da2SJoe Wallwork   Input Parameter:
146420282da2SJoe Wallwork + dm         - The DM
146520282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected
146620282da2SJoe Wallwork - metrics    - The metrics to be intersected
146720282da2SJoe Wallwork 
146820282da2SJoe Wallwork   Output Parameter:
146920282da2SJoe Wallwork . metricInt  - The intersected metric
147020282da2SJoe Wallwork 
147120282da2SJoe Wallwork   Level: beginner
147220282da2SJoe Wallwork 
147320282da2SJoe Wallwork   Notes:
147420282da2SJoe Wallwork 
147520282da2SJoe Wallwork   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
147620282da2SJoe Wallwork 
147720282da2SJoe Wallwork   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
147820282da2SJoe Wallwork 
147920282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
1480*cb7bfe8cSJoe Wallwork @*/
148120282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
148220282da2SJoe Wallwork {
148320282da2SJoe Wallwork   PetscErrorCode ierr;
148420282da2SJoe Wallwork   PetscInt       dim, vStart, vEnd, v, i;
148520282da2SJoe Wallwork   PetscScalar   *met, *meti, *M, *Mi;
148620282da2SJoe Wallwork 
148720282da2SJoe Wallwork   PetscFunctionBegin;
148820282da2SJoe Wallwork   if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); }
148920282da2SJoe Wallwork 
149020282da2SJoe Wallwork   /* Extract metadata from dm */
149120282da2SJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
149220282da2SJoe Wallwork   ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr);
149320282da2SJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
149420282da2SJoe Wallwork 
149520282da2SJoe Wallwork   /* Copy over the first metric */
149620282da2SJoe Wallwork   ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr);
149720282da2SJoe Wallwork 
149820282da2SJoe Wallwork   /* Intersect subsequent metrics in turn */
149920282da2SJoe Wallwork   if (numMetrics > 1) {
150020282da2SJoe Wallwork     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
150120282da2SJoe Wallwork     for (i = 1; i < numMetrics; ++i) {
150220282da2SJoe Wallwork       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
150320282da2SJoe Wallwork       for (v = vStart; v < vEnd; ++v) {
150420282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr);
150520282da2SJoe Wallwork         ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr);
150620282da2SJoe Wallwork         ierr = DMPlexMetricIntersection_Private(dim, Mi, M);CHKERRQ(ierr);
150720282da2SJoe Wallwork       }
150820282da2SJoe Wallwork       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
150920282da2SJoe Wallwork     }
151020282da2SJoe Wallwork     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
151120282da2SJoe Wallwork   }
151220282da2SJoe Wallwork 
151320282da2SJoe Wallwork   PetscFunctionReturn(0);
151420282da2SJoe Wallwork }
151520282da2SJoe Wallwork 
1516*cb7bfe8cSJoe Wallwork /*@
151720282da2SJoe Wallwork   DMPlexMetricIntersection2 - Compute the intersection of two metrics
151820282da2SJoe Wallwork 
151920282da2SJoe Wallwork   Input Parameter:
152020282da2SJoe Wallwork + dm        - The DM
152120282da2SJoe Wallwork . metric1   - The first metric to be intersected
152220282da2SJoe Wallwork - metric2   - The second metric to be intersected
152320282da2SJoe Wallwork 
152420282da2SJoe Wallwork   Output Parameter:
152520282da2SJoe Wallwork . metricInt - The intersected metric
152620282da2SJoe Wallwork 
152720282da2SJoe Wallwork   Level: beginner
152820282da2SJoe Wallwork 
152920282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
1530*cb7bfe8cSJoe Wallwork @*/
153120282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
153220282da2SJoe Wallwork {
153320282da2SJoe Wallwork   PetscErrorCode ierr;
153420282da2SJoe Wallwork   Vec            metrics[2] = {metric1, metric2};
153520282da2SJoe Wallwork 
153620282da2SJoe Wallwork   PetscFunctionBegin;
153720282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr);
153820282da2SJoe Wallwork   PetscFunctionReturn(0);
153920282da2SJoe Wallwork }
154020282da2SJoe Wallwork 
1541*cb7bfe8cSJoe Wallwork /*@
154220282da2SJoe Wallwork   DMPlexMetricIntersection3 - Compute the intersection of three metrics
154320282da2SJoe Wallwork 
154420282da2SJoe Wallwork   Input Parameter:
154520282da2SJoe Wallwork + dm        - The DM
154620282da2SJoe Wallwork . metric1   - The first metric to be intersected
154720282da2SJoe Wallwork . metric2   - The second metric to be intersected
154820282da2SJoe Wallwork - metric3   - The third metric to be intersected
154920282da2SJoe Wallwork 
155020282da2SJoe Wallwork   Output Parameter:
155120282da2SJoe Wallwork . metricInt - The intersected metric
155220282da2SJoe Wallwork 
155320282da2SJoe Wallwork   Level: beginner
155420282da2SJoe Wallwork 
155520282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
1556*cb7bfe8cSJoe Wallwork @*/
155720282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
155820282da2SJoe Wallwork {
155920282da2SJoe Wallwork   PetscErrorCode ierr;
156020282da2SJoe Wallwork   Vec            metrics[3] = {metric1, metric2, metric3};
156120282da2SJoe Wallwork 
156220282da2SJoe Wallwork   PetscFunctionBegin;
156320282da2SJoe Wallwork   ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr);
156420282da2SJoe Wallwork   PetscFunctionReturn(0);
156520282da2SJoe Wallwork }
1566