Actual source code: plexmetric.c
  1: #include <petsc/private/dmpleximpl.h>
  2: #include <petscblaslapack.h>
  4: PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection;
  6: PetscErrorCode DMPlexMetricSetFromOptions(DM dm)
  7: {
  8:   DM_Plex  *plex = (DM_Plex *)dm->data;
  9:   MPI_Comm  comm;
 10:   PetscBool isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE;
 11:   PetscBool noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE;
 12:   PetscInt  verbosity = -1, numIter = 3;
 13:   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;
 15:   PetscFunctionBegin;
 16:   if (!plex->metricCtx) PetscCall(PetscNew(&plex->metricCtx));
 17:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
 18:   PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");
 19:   PetscCall(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL));
 20:   PetscCall(DMPlexMetricSetIsotropic(dm, isotropic));
 21:   PetscCall(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL));
 22:   PetscCall(DMPlexMetricSetUniform(dm, uniform));
 23:   PetscCall(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL));
 24:   PetscCall(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst));
 25:   PetscCall(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL));
 26:   PetscCall(DMPlexMetricSetNoInsertion(dm, noInsert));
 27:   PetscCall(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL));
 28:   PetscCall(DMPlexMetricSetNoSwapping(dm, noSwap));
 29:   PetscCall(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL));
 30:   PetscCall(DMPlexMetricSetNoMovement(dm, noMove));
 31:   PetscCall(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL));
 32:   PetscCall(DMPlexMetricSetNoSurf(dm, noSurf));
 33:   PetscCall(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0));
 34:   PetscCall(DMPlexMetricSetNumIterations(dm, numIter));
 35:   PetscCall(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10));
 36:   PetscCall(DMPlexMetricSetVerbosity(dm, verbosity));
 37:   PetscCall(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL));
 38:   PetscCall(DMPlexMetricSetMinimumMagnitude(dm, h_min));
 39:   PetscCall(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL));
 40:   PetscCall(DMPlexMetricSetMaximumMagnitude(dm, h_max));
 41:   PetscCall(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL));
 42:   PetscCall(DMPlexMetricSetMaximumAnisotropy(dm, a_max));
 43:   PetscCall(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL));
 44:   PetscCall(DMPlexMetricSetNormalizationOrder(dm, p));
 45:   PetscCall(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL));
 46:   PetscCall(DMPlexMetricSetTargetComplexity(dm, target));
 47:   PetscCall(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL));
 48:   PetscCall(DMPlexMetricSetGradationFactor(dm, beta));
 49:   PetscCall(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL));
 50:   PetscCall(DMPlexMetricSetHausdorffNumber(dm, hausd));
 51:   PetscOptionsEnd();
 52:   PetscFunctionReturn(PETSC_SUCCESS);
 53: }
 55: /*@
 56:   DMPlexMetricSetIsotropic - Record whether a metric is isotropic
 58:   Input Parameters:
 59: + dm        - The `DM`
 60: - isotropic - Is the metric isotropic?
 62:   Level: beginner
 64: .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetUniform()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
 65: @*/
 66: PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic)
 67: {
 68:   DM_Plex *plex = (DM_Plex *)dm->data;
 70:   PetscFunctionBegin;
 71:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
 72:   plex->metricCtx->isotropic = isotropic;
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }
 76: /*@
 77:   DMPlexMetricIsIsotropic - Is a metric isotropic?
 79:   Input Parameters:
 80: . dm - The `DM`
 82:   Output Parameters:
 83: . isotropic - Is the metric isotropic?
 85:   Level: beginner
 87: .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricIsUniform()`, `DMPlexMetricRestrictAnisotropyFirst()`
 88: @*/
 89: PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic)
 90: {
 91:   DM_Plex *plex = (DM_Plex *)dm->data;
 93:   PetscFunctionBegin;
 94:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
 95:   *isotropic = plex->metricCtx->isotropic;
 96:   PetscFunctionReturn(PETSC_SUCCESS);
 97: }
 99: /*@
100:   DMPlexMetricSetUniform - Record whether a metric is uniform
102:   Input Parameters:
103: + dm      - The `DM`
104: - uniform - Is the metric uniform?
106:   Level: beginner
108:   Note:
109:   If the metric is specified as uniform then it is assumed to be isotropic, too.
111: .seealso: `DMPLEX`, `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
112: @*/
113: PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
114: {
115:   DM_Plex *plex = (DM_Plex *)dm->data;
117:   PetscFunctionBegin;
118:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
119:   plex->metricCtx->uniform = uniform;
120:   if (uniform) plex->metricCtx->isotropic = uniform;
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }
124: /*@
125:   DMPlexMetricIsUniform - Is a metric uniform?
127:   Input Parameters:
128: . dm - The `DM`
130:   Output Parameters:
131: . uniform - Is the metric uniform?
133:   Level: beginner
135: .seealso: `DMPLEX`, `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
136: @*/
137: PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
138: {
139:   DM_Plex *plex = (DM_Plex *)dm->data;
141:   PetscFunctionBegin;
142:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
143:   *uniform = plex->metricCtx->uniform;
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }
147: /*@
148:   DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
150:   Input Parameters:
151: + dm                      - The `DM`
152: - restrictAnisotropyFirst - Should anisotropy be normalized first?
154:   Level: beginner
156: .seealso: `DMPLEX`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
157: @*/
158: PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
159: {
160:   DM_Plex *plex = (DM_Plex *)dm->data;
162:   PetscFunctionBegin;
163:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
164:   plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
165:   PetscFunctionReturn(PETSC_SUCCESS);
166: }
168: /*@
169:   DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
171:   Input Parameters:
172: . dm - The `DM`
174:   Output Parameters:
175: . restrictAnisotropyFirst - Is anisotropy be normalized first?
177:   Level: beginner
179: .seealso: `DMPLEX`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
180: @*/
181: PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
182: {
183:   DM_Plex *plex = (DM_Plex *)dm->data;
185:   PetscFunctionBegin;
186:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
187:   *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
188:   PetscFunctionReturn(PETSC_SUCCESS);
189: }
191: /*@
192:   DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
194:   Input Parameters:
195: + dm       - The `DM`
196: - noInsert - Should node insertion and deletion be turned off?
198:   Level: beginner
200:   Note:
201:   This is only used by Mmg and ParMmg (not Pragmatic).
203: .seealso: `DMPLEX`, `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
204: @*/
205: PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
206: {
207:   DM_Plex *plex = (DM_Plex *)dm->data;
209:   PetscFunctionBegin;
210:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
211:   plex->metricCtx->noInsert = noInsert;
212:   PetscFunctionReturn(PETSC_SUCCESS);
213: }
215: /*@
216:   DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
218:   Input Parameters:
219: . dm - The `DM`
221:   Output Parameters:
222: . noInsert - Are node insertion and deletion turned off?
224:   Level: beginner
226:   Note:
227:   This is only used by Mmg and ParMmg (not Pragmatic).
229: .seealso: `DMPLEX`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
230: @*/
231: PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
232: {
233:   DM_Plex *plex = (DM_Plex *)dm->data;
235:   PetscFunctionBegin;
236:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
237:   *noInsert = plex->metricCtx->noInsert;
238:   PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: /*@
242:   DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
244:   Input Parameters:
245: + dm     - The `DM`
246: - noSwap - Should facet swapping be turned off?
248:   Level: beginner
250:   Note:
251:   This is only used by Mmg and ParMmg (not Pragmatic).
253: .seealso: `DMPLEX`, `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
254: @*/
255: PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
256: {
257:   DM_Plex *plex = (DM_Plex *)dm->data;
259:   PetscFunctionBegin;
260:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
261:   plex->metricCtx->noSwap = noSwap;
262:   PetscFunctionReturn(PETSC_SUCCESS);
263: }
265: /*@
266:   DMPlexMetricNoSwapping - Is facet swapping turned off?
268:   Input Parameters:
269: . dm - The `DM`
271:   Output Parameters:
272: . noSwap - Is facet swapping turned off?
274:   Level: beginner
276:   Note:
277:   This is only used by Mmg and ParMmg (not Pragmatic).
279: .seealso: `DMPLEX`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
280: @*/
281: PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
282: {
283:   DM_Plex *plex = (DM_Plex *)dm->data;
285:   PetscFunctionBegin;
286:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
287:   *noSwap = plex->metricCtx->noSwap;
288:   PetscFunctionReturn(PETSC_SUCCESS);
289: }
291: /*@
292:   DMPlexMetricSetNoMovement - Should node movement be turned off?
294:   Input Parameters:
295: + dm     - The `DM`
296: - noMove - Should node movement be turned off?
298:   Level: beginner
300:   Note:
301:   This is only used by Mmg and ParMmg (not Pragmatic).
303: .seealso: `DMPLEX`, `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()`
304: @*/
305: PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
306: {
307:   DM_Plex *plex = (DM_Plex *)dm->data;
309:   PetscFunctionBegin;
310:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
311:   plex->metricCtx->noMove = noMove;
312:   PetscFunctionReturn(PETSC_SUCCESS);
313: }
315: /*@
316:   DMPlexMetricNoMovement - Is node movement turned off?
318:   Input Parameters:
319: . dm - The `DM`
321:   Output Parameters:
322: . noMove - Is node movement turned off?
324:   Level: beginner
326:   Note:
327:   This is only used by Mmg and ParMmg (not Pragmatic).
329: .seealso: `DMPLEX`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()`
330: @*/
331: PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
332: {
333:   DM_Plex *plex = (DM_Plex *)dm->data;
335:   PetscFunctionBegin;
336:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
337:   *noMove = plex->metricCtx->noMove;
338:   PetscFunctionReturn(PETSC_SUCCESS);
339: }
341: /*@
342:   DMPlexMetricSetNoSurf - Should surface modification be turned off?
344:   Input Parameters:
345: + dm     - The `DM`
346: - noSurf - Should surface modification be turned off?
348:   Level: beginner
350:   Note:
351:   This is only used by Mmg and ParMmg (not Pragmatic).
353: .seealso: `DMPLEX`, `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`
354: @*/
355: PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
356: {
357:   DM_Plex *plex = (DM_Plex *)dm->data;
359:   PetscFunctionBegin;
360:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
361:   plex->metricCtx->noSurf = noSurf;
362:   PetscFunctionReturn(PETSC_SUCCESS);
363: }
365: /*@
366:   DMPlexMetricNoSurf - Is surface modification turned off?
368:   Input Parameters:
369: . dm - The `DM`
371:   Output Parameters:
372: . noSurf - Is surface modification turned off?
374:   Level: beginner
376:   Note:
377:   This is only used by Mmg and ParMmg (not Pragmatic).
379: .seealso: `DMPLEX`, `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`
380: @*/
381: PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
382: {
383:   DM_Plex *plex = (DM_Plex *)dm->data;
385:   PetscFunctionBegin;
386:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
387:   *noSurf = plex->metricCtx->noSurf;
388:   PetscFunctionReturn(PETSC_SUCCESS);
389: }
391: /*@
392:   DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
394:   Input Parameters:
395: + dm    - The `DM`
396: - h_min - The minimum tolerated metric magnitude
398:   Level: beginner
400: .seealso: `DMPLEX`, `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()`
401: @*/
402: PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
403: {
404:   DM_Plex *plex = (DM_Plex *)dm->data;
406:   PetscFunctionBegin;
407:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
408:   PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
409:   plex->metricCtx->h_min = h_min;
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }
413: /*@
414:   DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
416:   Input Parameters:
417: . dm - The `DM`
419:   Output Parameters:
420: . h_min - The minimum tolerated metric magnitude
422:   Level: beginner
424: .seealso: `DMPLEX`, `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()`
425: @*/
426: PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
427: {
428:   DM_Plex *plex = (DM_Plex *)dm->data;
430:   PetscFunctionBegin;
431:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
432:   *h_min = plex->metricCtx->h_min;
433:   PetscFunctionReturn(PETSC_SUCCESS);
434: }
436: /*@
437:   DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
439:   Input Parameters:
440: + dm    - The `DM`
441: - h_max - The maximum tolerated metric magnitude
443:   Level: beginner
445: .seealso: `DMPLEX`, `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()`
446: @*/
447: PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
448: {
449:   DM_Plex *plex = (DM_Plex *)dm->data;
451:   PetscFunctionBegin;
452:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
453:   PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
454:   plex->metricCtx->h_max = h_max;
455:   PetscFunctionReturn(PETSC_SUCCESS);
456: }
458: /*@
459:   DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
461:   Input Parameters:
462: . dm - The `DM`
464:   Output Parameters:
465: . h_max - The maximum tolerated metric magnitude
467:   Level: beginner
469: .seealso: `DMPLEX`, `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()`
470: @*/
471: PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
472: {
473:   DM_Plex *plex = (DM_Plex *)dm->data;
475:   PetscFunctionBegin;
476:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
477:   *h_max = plex->metricCtx->h_max;
478:   PetscFunctionReturn(PETSC_SUCCESS);
479: }
481: /*@
482:   DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
484:   Input Parameters:
485: + dm    - The `DM`
486: - a_max - The maximum tolerated metric anisotropy
488:   Level: beginner
490:   Note:
491:   If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
493: .seealso: `DMPLEX`, `DMPlexMetricGetMaximumAnisotropy()`, `DMPlexMetricSetMaximumMagnitude()`
494: @*/
495: PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max)
496: {
497:   DM_Plex *plex = (DM_Plex *)dm->data;
499:   PetscFunctionBegin;
500:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
501:   PetscCheck(a_max >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Anisotropy must be in [1, inf)");
502:   plex->metricCtx->a_max = a_max;
503:   PetscFunctionReturn(PETSC_SUCCESS);
504: }
506: /*@
507:   DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy
509:   Input Parameters:
510: . dm - The `DM`
512:   Output Parameters:
513: . a_max - The maximum tolerated metric anisotropy
515:   Level: beginner
517: .seealso: `DMPLEX`, `DMPlexMetricSetMaximumAnisotropy()`, `DMPlexMetricGetMaximumMagnitude()`
518: @*/
519: PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max)
520: {
521:   DM_Plex *plex = (DM_Plex *)dm->data;
523:   PetscFunctionBegin;
524:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
525:   *a_max = plex->metricCtx->a_max;
526:   PetscFunctionReturn(PETSC_SUCCESS);
527: }
529: /*@
530:   DMPlexMetricSetTargetComplexity - Set the target metric complexity
532:   Input Parameters:
533: + dm               - The `DM`
534: - targetComplexity - The target metric complexity
536:   Level: beginner
538: .seealso: `DMPLEX`, `DMPlexMetricGetTargetComplexity()`, `DMPlexMetricSetNormalizationOrder()`
539: @*/
540: PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity)
541: {
542:   DM_Plex *plex = (DM_Plex *)dm->data;
544:   PetscFunctionBegin;
545:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
546:   PetscCheck(targetComplexity > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric complexity must be in (0, inf)");
547:   plex->metricCtx->targetComplexity = targetComplexity;
548:   PetscFunctionReturn(PETSC_SUCCESS);
549: }
551: /*@
552:   DMPlexMetricGetTargetComplexity - Get the target metric complexity
554:   Input Parameters:
555: . dm - The `DM`
557:   Output Parameters:
558: . targetComplexity - The target metric complexity
560:   Level: beginner
562: .seealso: `DMPLEX`, `DMPlexMetricSetTargetComplexity()`, `DMPlexMetricGetNormalizationOrder()`
563: @*/
564: PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity)
565: {
566:   DM_Plex *plex = (DM_Plex *)dm->data;
568:   PetscFunctionBegin;
569:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
570:   *targetComplexity = plex->metricCtx->targetComplexity;
571:   PetscFunctionReturn(PETSC_SUCCESS);
572: }
574: /*@
575:   DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization
577:   Input Parameters:
578: + dm - The `DM`
579: - p  - The normalization order
581:   Level: beginner
583: .seealso: `DMPLEX`, `DMPlexMetricGetNormalizationOrder()`, `DMPlexMetricSetTargetComplexity()`
584: @*/
585: PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p)
586: {
587:   DM_Plex *plex = (DM_Plex *)dm->data;
589:   PetscFunctionBegin;
590:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
591:   PetscCheck(p >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Normalization order must be in [1, inf)");
592:   plex->metricCtx->p = p;
593:   PetscFunctionReturn(PETSC_SUCCESS);
594: }
596: /*@
597:   DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization
599:   Input Parameters:
600: . dm - The `DM`
602:   Output Parameters:
603: . p - The normalization order
605:   Level: beginner
607: .seealso: `DMPLEX`, `DMPlexMetricSetNormalizationOrder()`, `DMPlexMetricGetTargetComplexity()`
608: @*/
609: PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p)
610: {
611:   DM_Plex *plex = (DM_Plex *)dm->data;
613:   PetscFunctionBegin;
614:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
615:   *p = plex->metricCtx->p;
616:   PetscFunctionReturn(PETSC_SUCCESS);
617: }
619: /*@
620:   DMPlexMetricSetGradationFactor - Set the metric gradation factor
622:   Input Parameters:
623: + dm   - The `DM`
624: - beta - The metric gradation factor
626:   Level: beginner
628:   Notes:
629:   The gradation factor is the maximum tolerated length ratio between adjacent edges.
631:   Turn off gradation by passing the value -1. Otherwise, pass a positive value.
633:   This is only used by Mmg and ParMmg (not Pragmatic).
635: .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
636: @*/
637: PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
638: {
639:   DM_Plex *plex = (DM_Plex *)dm->data;
641:   PetscFunctionBegin;
642:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
643:   PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)");
644:   plex->metricCtx->gradationFactor = beta;
645:   PetscFunctionReturn(PETSC_SUCCESS);
646: }
648: /*@
649:   DMPlexMetricGetGradationFactor - Get the metric gradation factor
651:   Input Parameters:
652: . dm - The `DM`
654:   Output Parameters:
655: . beta - The metric gradation factor
657:   Level: beginner
659:   Notes:
660:   The gradation factor is the maximum tolerated length ratio between adjacent edges.
662:   The value -1 implies that gradation is turned off.
664:   This is only used by Mmg and ParMmg (not Pragmatic).
666: .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
667: @*/
668: PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
669: {
670:   DM_Plex *plex = (DM_Plex *)dm->data;
672:   PetscFunctionBegin;
673:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
674:   *beta = plex->metricCtx->gradationFactor;
675:   PetscFunctionReturn(PETSC_SUCCESS);
676: }
678: /*@
679:   DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
681:   Input Parameters:
682: + dm    - The `DM`
683: - hausd - The metric Hausdorff number
685:   Level: beginner
687:   Notes:
688:   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
689:   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
690:   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
691:   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
692:   (resp. increase) the Hausdorff parameter. (Taken from
693:   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
695:   This is only used by Mmg and ParMmg (not Pragmatic).
697: .seealso: `DMPLEX`, `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
698: @*/
699: PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
700: {
701:   DM_Plex *plex = (DM_Plex *)dm->data;
703:   PetscFunctionBegin;
704:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
705:   PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)");
706:   plex->metricCtx->hausdorffNumber = hausd;
707:   PetscFunctionReturn(PETSC_SUCCESS);
708: }
710: /*@
711:   DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
713:   Input Parameters:
714: . dm - The `DM`
716:   Output Parameters:
717: . hausd - The metric Hausdorff number
719:   Level: beginner
721:   Notes:
722:   The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
723:   boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
724:   high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
725:   an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
726:   (resp. increase) the Hausdorff parameter. (Taken from
727:   https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
729:   This is only used by Mmg and ParMmg (not Pragmatic).
731: .seealso: `DMPLEX`, `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
732: @*/
733: PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
734: {
735:   DM_Plex *plex = (DM_Plex *)dm->data;
737:   PetscFunctionBegin;
738:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
739:   *hausd = plex->metricCtx->hausdorffNumber;
740:   PetscFunctionReturn(PETSC_SUCCESS);
741: }
743: /*@
744:   DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
746:   Input Parameters:
747: + dm        - The `DM`
748: - verbosity - The verbosity, where -1 is silent and 10 is maximum
750:   Level: beginner
752:   Note:
753:   This is only used by Mmg and ParMmg (not Pragmatic).
755: .seealso: `DMPLEX`, `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()`
756: @*/
757: PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
758: {
759:   DM_Plex *plex = (DM_Plex *)dm->data;
761:   PetscFunctionBegin;
762:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
763:   plex->metricCtx->verbosity = verbosity;
764:   PetscFunctionReturn(PETSC_SUCCESS);
765: }
767: /*@
768:   DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
770:   Input Parameters:
771: . dm - The `DM`
773:   Output Parameters:
774: . verbosity - The verbosity, where -1 is silent and 10 is maximum
776:   Level: beginner
778:   Note:
779:   This is only used by Mmg and ParMmg (not Pragmatic).
781: .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
782: @*/
783: PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
784: {
785:   DM_Plex *plex = (DM_Plex *)dm->data;
787:   PetscFunctionBegin;
788:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
789:   *verbosity = plex->metricCtx->verbosity;
790:   PetscFunctionReturn(PETSC_SUCCESS);
791: }
793: /*@
794:   DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
796:   Input Parameters:
797: + dm      - The `DM`
798: - numIter - the number of parallel adaptation iterations
800:   Level: beginner
802:   Note:
803:   This is only used by ParMmg (not Pragmatic or Mmg).
805: .seealso: `DMPLEX`, `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
806: @*/
807: PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
808: {
809:   DM_Plex *plex = (DM_Plex *)dm->data;
811:   PetscFunctionBegin;
812:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
813:   plex->metricCtx->numIter = numIter;
814:   PetscFunctionReturn(PETSC_SUCCESS);
815: }
817: /*@
818:   DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
820:   Input Parameters:
821: . dm - The `DM`
823:   Output Parameters:
824: . numIter - the number of parallel adaptation iterations
826:   Level: beginner
828:   Note:
829:   This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
831: .seealso: `DMPLEX`, `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()`
832: @*/
833: PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
834: {
835:   DM_Plex *plex = (DM_Plex *)dm->data;
837:   PetscFunctionBegin;
838:   if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
839:   *numIter = plex->metricCtx->numIter;
840:   PetscFunctionReturn(PETSC_SUCCESS);
841: }
843: static PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
844: {
845:   MPI_Comm comm;
846:   PetscFE  fe;
847:   PetscInt dim;
849:   PetscFunctionBegin;
850:   /* Extract metadata from dm */
851:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
852:   PetscCall(DMGetDimension(dm, &dim));
854:   /* Create a P1 field of the requested size */
855:   PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
856:   PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
857:   PetscCall(DMCreateDS(dm));
858:   PetscCall(PetscFEDestroy(&fe));
859:   PetscCall(DMCreateLocalVector(dm, metric));
860:   PetscFunctionReturn(PETSC_SUCCESS);
861: }
863: /*@
864:   DMPlexMetricCreate - Create a Riemannian metric field
866:   Input Parameters:
867: + dm - The `DM`
868: - f  - The field number to use
870:   Output Parameter:
871: . metric - The metric
873:   Options Database Key:
874: . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use
876:   Options Database Keys for Mmg and ParMmg:
877: + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation
878: . -dm_plex_metric_num_iterations   - Number of parallel mesh adaptation iterations for ParMmg
879: . -dm_plex_metric_no_insert        - Should node insertion/deletion be turned off?
880: . -dm_plex_metric_no_swap          - Should facet swapping be turned off?
881: . -dm_plex_metric_no_move          - Should node movement be turned off?
882: - -dm_plex_metric_verbosity        - Choose a verbosity level from -1 (silent) to 10 (maximum).
884:   Options Database Keys for Riemannian metrics:
885: + -dm_plex_metric_isotropic                 - Is the metric isotropic?
886: . -dm_plex_metric_uniform                   - Is the metric uniform?
887: . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
888: . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
889: . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
890: . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
891: . -dm_plex_metric_p                         - L-p normalization order
892: - -dm_plex_metric_target_complexity         - Target metric complexity
894:   Level: beginner
896:   Note:
897:   It is assumed that the `DM` is comprised of simplices.
899: .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`
900: @*/
901: PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
902: {
903:   PetscBool isotropic, uniform;
904:   PetscInt  coordDim, Nd;
906:   PetscFunctionBegin;
907:   PetscCall(DMGetCoordinateDim(dm, &coordDim));
908:   Nd = coordDim * coordDim;
909:   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
910:   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
911:   if (uniform) {
912:     MPI_Comm comm;
914:     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
915:     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
916:     PetscCall(VecCreate(comm, metric));
917:     PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
918:     PetscCall(VecSetFromOptions(*metric));
919:   } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
920:   else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
921:   PetscFunctionReturn(PETSC_SUCCESS);
922: }
924: /*@
925:   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
927:   Input Parameters:
928: + dm    - The `DM`
929: . f     - The field number to use
930: - alpha - Scaling parameter for the diagonal
932:   Output Parameter:
933: . metric - The uniform metric
935:   Level: beginner
937:   Note:
938:   In this case, the metric is constant in space and so is only specified for a single vertex.
940: .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
941: @*/
942: PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
943: {
944:   PetscFunctionBegin;
945:   PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
946:   PetscCall(DMPlexMetricCreate(dm, f, metric));
947:   PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
948:   PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
949:   PetscCall(VecSet(*metric, alpha));
950:   PetscCall(VecAssemblyBegin(*metric));
951:   PetscCall(VecAssemblyEnd(*metric));
952:   PetscFunctionReturn(PETSC_SUCCESS);
953: }
955: 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[])
956: {
957:   f0[0] = u[0];
958: }
960: /*@
961:   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
963:   Input Parameters:
964: + dm        - The `DM`
965: . f         - The field number to use
966: - indicator - The error indicator
968:   Output Parameter:
969: . metric - The isotropic metric
971:   Level: beginner
973:   Notes:
974:   It is assumed that the `DM` is comprised of simplices.
976:   The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
978: .seealso: `DMPLEX`, `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
979: @*/
980: PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
981: {
982:   PetscInt m, n;
984:   PetscFunctionBegin;
985:   PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
986:   PetscCall(DMPlexMetricCreate(dm, f, metric));
987:   PetscCall(VecGetSize(indicator, &m));
988:   PetscCall(VecGetSize(*metric, &n));
989:   if (m == n) PetscCall(VecCopy(indicator, *metric));
990:   else {
991:     DM dmIndi;
992:     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[]);
994:     PetscCall(VecGetDM(indicator, &dmIndi));
995:     funcs[0] = identity;
996:     PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
997:   }
998:   PetscFunctionReturn(PETSC_SUCCESS);
999: }
1001: /*@
1002:   DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1004:   Input Parameters:
1005: + dm - The `DM` of the metric field
1006: - f  - The field number to use
1008:   Output Parameters:
1009: + determinant - The determinant field
1010: - dmDet       - The corresponding `DM`
1012:   Level: beginner
1014: .seealso: `DMPLEX`, `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`, `DMPlexMetricCreate()`
1015: @*/
1016: PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1017: {
1018:   PetscBool uniform;
1020:   PetscFunctionBegin;
1021:   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1022:   PetscCall(DMClone(dm, dmDet));
1023:   if (uniform) {
1024:     MPI_Comm comm;
1026:     PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm));
1027:     PetscCall(VecCreate(comm, determinant));
1028:     PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1029:     PetscCall(VecSetFromOptions(*determinant));
1030:   } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
1031:   PetscFunctionReturn(PETSC_SUCCESS);
1032: }
1034: static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1035: {
1036:   PetscInt i, j;
1038:   PetscFunctionBegin;
1039:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"));
1040:   for (i = 0; i < dim; ++i) {
1041:     if (i == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    [["));
1042:     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "     ["));
1043:     for (j = 0; j < dim; ++j) {
1044:       if (j < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])));
1045:       else PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])));
1046:     }
1047:     if (i < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
1048:     else PetscCall(PetscPrintf(PETSC_COMM_SELF, "]]\n"));
1049:   }
1050:   PetscFunctionReturn(PETSC_SUCCESS);
1051: }
1053: static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1054: {
1055:   PetscInt     i, j, k;
1056:   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);
1057:   PetscScalar *Mpos;
1059:   PetscFunctionBegin;
1060:   PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs));
1062:   /* Symmetrize */
1063:   for (i = 0; i < dim; ++i) {
1064:     Mpos[i * dim + i] = Mp[i * dim + i];
1065:     for (j = i + 1; j < dim; ++j) {
1066:       Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
1067:       Mpos[j * dim + i] = Mpos[i * dim + j];
1068:     }
1069:   }
1071:   /* Compute eigendecomposition */
1072:   if (dim == 1) {
1073:     /* Isotropic case */
1074:     eigs[0] = PetscRealPart(Mpos[0]);
1075:     Mpos[0] = 1.0;
1076:   } else {
1077:     /* Anisotropic case */
1078:     PetscScalar *work;
1079:     PetscBLASInt lwork = (PetscBLASInt)(5 * dim);
1081:     PetscCall(PetscMalloc1(5 * dim, &work));
1082:     {
1083:       PetscBLASInt lierr;
1084:       PetscBLASInt nb;
1086:       PetscCall(PetscBLASIntCast(dim, &nb));
1087:       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1088: #if defined(PETSC_USE_COMPLEX)
1089:       {
1090:         PetscReal *rwork;
1091:         PetscCall(PetscMalloc1(3 * dim, &rwork));
1092:         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr));
1093:         PetscCall(PetscFree(rwork));
1094:       }
1095: #else
1096:       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr));
1097: #endif
1098:       if (lierr) {
1099:         for (i = 0; i < dim; ++i) {
1100:           Mpos[i * dim + i] = Mp[i * dim + i];
1101:           for (j = i + 1; j < dim; ++j) {
1102:             Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
1103:             Mpos[j * dim + i] = Mpos[i * dim + j];
1104:           }
1105:         }
1106:         PetscCall(LAPACKsyevFail(dim, Mpos));
1107:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
1108:       }
1109:       PetscCall(PetscFPTrapPop());
1110:     }
1111:     PetscCall(PetscFree(work));
1112:   }
1114:   /* Reflect to positive orthant and enforce maximum and minimum size */
1115:   max_eig = 0.0;
1116:   for (i = 0; i < dim; ++i) {
1117:     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
1118:     max_eig = PetscMax(eigs[i], max_eig);
1119:   }
1121:   /* Enforce maximum anisotropy and compute determinant */
1122:   *detMp = 1.0;
1123:   for (i = 0; i < dim; ++i) {
1124:     if (a_max >= 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min);
1125:     *detMp *= eigs[i];
1126:   }
1128:   /* Reconstruct Hessian */
1129:   for (i = 0; i < dim; ++i) {
1130:     for (j = 0; j < dim; ++j) {
1131:       Mp[i * dim + j] = 0.0;
1132:       for (k = 0; k < dim; ++k) Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j];
1133:     }
1134:   }
1135:   PetscCall(PetscFree2(Mpos, eigs));
1136:   PetscFunctionReturn(PETSC_SUCCESS);
1137: }
1139: /*@
1140:   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
1142:   Input Parameters:
1143: + dm                 - The `DM`
1144: . metricIn           - The metric
1145: . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
1146: - restrictAnisotropy - Should maximum anisotropy be enforced?
1148:   Output Parameters:
1149: + metricOut   - The metric
1150: - determinant - Its determinant
1152:   Options Database Keys:
1153: + -dm_plex_metric_isotropic - Is the metric isotropic?
1154: . -dm_plex_metric_uniform   - Is the metric uniform?
1155: . -dm_plex_metric_h_min     - Minimum tolerated metric magnitude
1156: . -dm_plex_metric_h_max     - Maximum tolerated metric magnitude
1157: - -dm_plex_metric_a_max     - Maximum tolerated anisotropy
1159:   Level: beginner
1161: .seealso: `DMPLEX`, `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1162: @*/
1163: PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1164: {
1165:   DM           dmDet;
1166:   PetscBool    isotropic, uniform;
1167:   PetscInt     dim, vStart, vEnd, v;
1168:   PetscScalar *met, *det;
1169:   PetscReal    h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+15;
1171:   PetscFunctionBegin;
1172:   PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1174:   /* Extract metadata from dm */
1175:   PetscCall(DMGetDimension(dm, &dim));
1176:   if (restrictSizes) {
1177:     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1178:     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1179:     h_min = PetscMax(h_min, 1.0e-30);
1180:     h_max = PetscMin(h_max, 1.0e+30);
1181:     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1182:   }
1183:   if (restrictAnisotropy) {
1184:     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1185:     a_max = PetscMin(a_max, 1.0e+30);
1186:   }
1188:   /* Setup output metric */
1189:   PetscCall(VecCopy(metricIn, metricOut));
1191:   /* Enforce SPD and extract determinant */
1192:   PetscCall(VecGetArray(metricOut, &met));
1193:   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1194:   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1195:   if (uniform) {
1196:     PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
1198:     /* Uniform case */
1199:     PetscCall(VecGetArray(determinant, &det));
1200:     PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1201:     PetscCall(VecRestoreArray(determinant, &det));
1202:   } else {
1203:     /* Spatially varying case */
1204:     PetscInt nrow;
1206:     if (isotropic) nrow = 1;
1207:     else nrow = dim;
1208:     PetscCall(VecGetDM(determinant, &dmDet));
1209:     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1210:     PetscCall(VecGetArray(determinant, &det));
1211:     for (v = vStart; v < vEnd; ++v) {
1212:       PetscScalar *vmet, *vdet;
1213:       PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
1214:       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
1215:       PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
1216:     }
1217:     PetscCall(VecRestoreArray(determinant, &det));
1218:   }
1219:   PetscCall(VecRestoreArray(metricOut, &met));
1221:   PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1222:   PetscFunctionReturn(PETSC_SUCCESS);
1223: }
1225: 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[])
1226: {
1227:   const PetscScalar p = constants[0];
1229:   f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim));
1230: }
1232: /*@
1233:   DMPlexMetricNormalize - Apply L-p normalization to a metric
1235:   Input Parameters:
1236: + dm                 - The `DM`
1237: . metricIn           - The unnormalized metric
1238: . restrictSizes      - Should maximum/minimum metric magnitudes be enforced?
1239: - restrictAnisotropy - Should maximum metric anisotropy be enforced?
1241:   Output Parameters:
1242: + metricOut   - The normalized metric
1243: - determinant - computed determinant
1245:   Options Database Keys:
1246: + -dm_plex_metric_isotropic                 - Is the metric isotropic?
1247: . -dm_plex_metric_uniform                   - Is the metric uniform?
1248: . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1249: . -dm_plex_metric_h_min                     - Minimum tolerated metric magnitude
1250: . -dm_plex_metric_h_max                     - Maximum tolerated metric magnitude
1251: . -dm_plex_metric_a_max                     - Maximum tolerated anisotropy
1252: . -dm_plex_metric_p                         - L-p normalization order
1253: - -dm_plex_metric_target_complexity         - Target metric complexity
1255:   Level: beginner
1257: .seealso: `DMPLEX`, `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1258: @*/
1259: PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1260: {
1261:   DM           dmDet;
1262:   MPI_Comm     comm;
1263:   PetscBool    restrictAnisotropyFirst, isotropic, uniform;
1264:   PetscDS      ds;
1265:   PetscInt     dim, Nd, vStart, vEnd, v, i;
1266:   PetscScalar *met, *det, integral, constants[1];
1267:   PetscReal    p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
1269:   PetscFunctionBegin;
1270:   PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1272:   /* Extract metadata from dm */
1273:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1274:   PetscCall(DMGetDimension(dm, &dim));
1275:   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1276:   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1277:   if (isotropic) Nd = 1;
1278:   else Nd = dim * dim;
1280:   /* Set up metric and ensure it is SPD */
1281:   PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
1282:   PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
1284:   /* Compute global normalization factor */
1285:   PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
1286:   PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
1287:   constants[0] = p;
1288:   if (uniform) {
1289:     PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
1290:     DM  dmTmp;
1291:     Vec tmp;
1293:     PetscCall(DMClone(dm, &dmTmp));
1294:     PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
1295:     PetscCall(VecGetArray(determinant, &det));
1296:     PetscCall(VecSet(tmp, det[0]));
1297:     PetscCall(VecRestoreArray(determinant, &det));
1298:     PetscCall(DMGetDS(dmTmp, &ds));
1299:     PetscCall(PetscDSSetConstants(ds, 1, constants));
1300:     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1301:     PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
1302:     PetscCall(VecDestroy(&tmp));
1303:     PetscCall(DMDestroy(&dmTmp));
1304:   } else {
1305:     PetscCall(VecGetDM(determinant, &dmDet));
1306:     PetscCall(DMGetDS(dmDet, &ds));
1307:     PetscCall(PetscDSSetConstants(ds, 1, constants));
1308:     PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1309:     PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
1310:   }
1311:   realIntegral = PetscRealPart(integral);
1312:   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?");
1313:   factGlob = PetscPowReal(target / realIntegral, 2.0 / dim);
1315:   /* Apply local scaling */
1316:   if (restrictSizes) {
1317:     PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1318:     PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1319:     h_min = PetscMax(h_min, 1.0e-30);
1320:     h_max = PetscMin(h_max, 1.0e+30);
1321:     PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1322:   }
1323:   if (restrictAnisotropy && !restrictAnisotropyFirst) {
1324:     PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1325:     a_max = PetscMin(a_max, 1.0e+30);
1326:   }
1327:   PetscCall(VecGetArray(metricOut, &met));
1328:   PetscCall(VecGetArray(determinant, &det));
1329:   if (uniform) {
1330:     /* Uniform case */
1331:     met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim));
1332:     if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1333:   } else {
1334:     /* Spatially varying case */
1335:     PetscInt nrow;
1337:     if (isotropic) nrow = 1;
1338:     else nrow = dim;
1339:     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1340:     PetscCall(VecGetDM(determinant, &dmDet));
1341:     for (v = vStart; v < vEnd; ++v) {
1342:       PetscScalar *Mp, *detM;
1344:       PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
1345:       PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
1346:       fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim));
1347:       for (i = 0; i < Nd; ++i) Mp[i] *= fact;
1348:       if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
1349:     }
1350:   }
1351:   PetscCall(VecRestoreArray(determinant, &det));
1352:   PetscCall(VecRestoreArray(metricOut, &met));
1354:   PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1355:   PetscFunctionReturn(PETSC_SUCCESS);
1356: }
1358: /*@
1359:   DMPlexMetricAverage - Compute the average of a list of metrics
1361:   Input Parameters:
1362: + dm         - The `DM`
1363: . numMetrics - The number of metrics to be averaged
1364: . weights    - Weights for the average
1365: - metrics    - The metrics to be averaged
1367:   Output Parameter:
1368: . metricAvg - The averaged metric
1370:   Level: beginner
1372:   Notes:
1373:   The weights should sum to unity.
1375:   If weights are not provided then an unweighted average is used.
1377: .seealso: `DMPLEX`, `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1378: @*/
1379: PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1380: {
1381:   PetscBool haveWeights = PETSC_TRUE;
1382:   PetscInt  i, m, n;
1383:   PetscReal sum = 0.0, tol = 1.0e-10;
1385:   PetscFunctionBegin;
1386:   PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0));
1387:   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
1388:   PetscCall(VecSet(metricAvg, 0.0));
1389:   PetscCall(VecGetSize(metricAvg, &m));
1390:   for (i = 0; i < numMetrics; ++i) {
1391:     PetscCall(VecGetSize(metrics[i], &n));
1392:     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
1393:   }
1395:   /* Default to the unweighted case */
1396:   if (!weights) {
1397:     PetscCall(PetscMalloc1(numMetrics, &weights));
1398:     haveWeights = PETSC_FALSE;
1399:     for (i = 0; i < numMetrics; ++i) weights[i] = 1.0 / numMetrics;
1400:   }
1402:   /* Check weights sum to unity */
1403:   for (i = 0; i < numMetrics; ++i) sum += weights[i];
1404:   PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
1406:   /* Compute metric average */
1407:   for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
1408:   if (!haveWeights) PetscCall(PetscFree(weights));
1410:   PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0));
1411:   PetscFunctionReturn(PETSC_SUCCESS);
1412: }
1414: /*@
1415:   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
1417:   Input Parameters:
1418: + dm      - The `DM`
1419: . metric1 - The first metric to be averaged
1420: - metric2 - The second metric to be averaged
1422:   Output Parameter:
1423: . metricAvg - The averaged metric
1425:   Level: beginner
1427: .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1428: @*/
1429: PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1430: {
1431:   PetscReal weights[2] = {0.5, 0.5};
1432:   Vec       metrics[2] = {metric1, metric2};
1434:   PetscFunctionBegin;
1435:   PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
1436:   PetscFunctionReturn(PETSC_SUCCESS);
1437: }
1439: /*@
1440:   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
1442:   Input Parameters:
1443: + dm      - The `DM`
1444: . metric1 - The first metric to be averaged
1445: . metric2 - The second metric to be averaged
1446: - metric3 - The third metric to be averaged
1448:   Output Parameter:
1449: . metricAvg - The averaged metric
1451:   Level: beginner
1453: .seealso: `DMPLEX`, `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1454: @*/
1455: PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1456: {
1457:   PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0};
1458:   Vec       metrics[3] = {metric1, metric2, metric3};
1460:   PetscFunctionBegin;
1461:   PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
1462:   PetscFunctionReturn(PETSC_SUCCESS);
1463: }
1465: static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1466: {
1467:   PetscInt     i, j, k, l, m;
1468:   PetscReal   *evals;
1469:   PetscScalar *evecs, *sqrtM1, *isqrtM1;
1471:   PetscFunctionBegin;
1472:   /* Isotropic case */
1473:   if (dim == 1) {
1474:     M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
1475:     PetscFunctionReturn(PETSC_SUCCESS);
1476:   }
1478:   /* Anisotropic case */
1479:   PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals));
1480:   for (i = 0; i < dim; ++i) {
1481:     for (j = 0; j < dim; ++j) evecs[i * dim + j] = M1[i * dim + j];
1482:   }
1483:   {
1484:     PetscScalar *work;
1485:     PetscBLASInt lwork = (PetscBLASInt)(5 * dim);
1486:     PetscCall(PetscMalloc1(5 * dim, &work));
1487:     {
1488:       PetscBLASInt lierr, nb;
1489:       PetscReal    sqrtj;
1491:       /* Compute eigendecomposition of M1 */
1492:       PetscCall(PetscBLASIntCast(dim, &nb));
1493:       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1494: #if defined(PETSC_USE_COMPLEX)
1495:       {
1496:         PetscReal *rwork;
1497:         PetscCall(PetscMalloc1(3 * dim, &rwork));
1498:         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1499:         PetscCall(PetscFree(rwork));
1500:       }
1501: #else
1502:       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1503: #endif
1504:       if (lierr) {
1505:         PetscCall(LAPACKsyevFail(dim, M1));
1506:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
1507:       }
1508:       PetscCall(PetscFPTrapPop());
1510:       /* Compute square root and the reciprocal thereof */
1511:       for (i = 0; i < dim; ++i) {
1512:         for (k = 0; k < dim; ++k) {
1513:           sqrtM1[i * dim + k]  = 0.0;
1514:           isqrtM1[i * dim + k] = 0.0;
1515:           for (j = 0; j < dim; ++j) {
1516:             sqrtj = PetscSqrtReal(evals[j]);
1517:             sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k];
1518:             isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k];
1519:           }
1520:         }
1521:       }
1523:       /* Map M2 into the space spanned by the eigenvectors of M1 */
1524:       for (i = 0; i < dim; ++i) {
1525:         for (l = 0; l < dim; ++l) {
1526:           evecs[i * dim + l] = 0.0;
1527:           for (j = 0; j < dim; ++j) {
1528:             for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
1529:           }
1530:         }
1531:       }
1533:       /* Compute eigendecomposition */
1534:       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1535: #if defined(PETSC_USE_COMPLEX)
1536:       {
1537:         PetscReal *rwork;
1538:         PetscCall(PetscMalloc1(3 * dim, &rwork));
1539:         PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1540:         PetscCall(PetscFree(rwork));
1541:       }
1542: #else
1543:       PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1544: #endif
1545:       if (lierr) {
1546:         for (i = 0; i < dim; ++i) {
1547:           for (l = 0; l < dim; ++l) {
1548:             evecs[i * dim + l] = 0.0;
1549:             for (j = 0; j < dim; ++j) {
1550:               for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
1551:             }
1552:           }
1553:         }
1554:         PetscCall(LAPACKsyevFail(dim, evecs));
1555:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
1556:       }
1557:       PetscCall(PetscFPTrapPop());
1559:       /* Modify eigenvalues */
1560:       for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
1562:       /* Map back to get the intersection */
1563:       for (i = 0; i < dim; ++i) {
1564:         for (m = 0; m < dim; ++m) {
1565:           M2[i * dim + m] = 0.0;
1566:           for (j = 0; j < dim; ++j) {
1567:             for (k = 0; k < dim; ++k) {
1568:               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];
1569:             }
1570:           }
1571:         }
1572:       }
1573:     }
1574:     PetscCall(PetscFree(work));
1575:   }
1576:   PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
1577:   PetscFunctionReturn(PETSC_SUCCESS);
1578: }
1580: /*@
1581:   DMPlexMetricIntersection - Compute the intersection of a list of metrics
1583:   Input Parameters:
1584: + dm         - The `DM`
1585: . numMetrics - The number of metrics to be intersected
1586: - metrics    - The metrics to be intersected
1588:   Output Parameter:
1589: . metricInt - The intersected metric
1591:   Level: beginner
1593:   Notes:
1594:   The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
1596:   The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
1598: .seealso: `DMPLEX`, `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1599: @*/
1600: PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1601: {
1602:   PetscBool    isotropic, uniform;
1603:   PetscInt     v, i, m, n;
1604:   PetscScalar *met, *meti;
1606:   PetscFunctionBegin;
1607:   PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1608:   PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
1610:   /* Copy over the first metric */
1611:   PetscCall(VecCopy(metrics[0], metricInt));
1612:   if (numMetrics == 1) PetscFunctionReturn(PETSC_SUCCESS);
1613:   PetscCall(VecGetSize(metricInt, &m));
1614:   for (i = 0; i < numMetrics; ++i) {
1615:     PetscCall(VecGetSize(metrics[i], &n));
1616:     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
1617:   }
1619:   /* Intersect subsequent metrics in turn */
1620:   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1621:   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1622:   if (uniform) {
1623:     /* Uniform case */
1624:     PetscCall(VecGetArray(metricInt, &met));
1625:     for (i = 1; i < numMetrics; ++i) {
1626:       PetscCall(VecGetArray(metrics[i], &meti));
1627:       PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
1628:       PetscCall(VecRestoreArray(metrics[i], &meti));
1629:     }
1630:     PetscCall(VecRestoreArray(metricInt, &met));
1631:   } else {
1632:     /* Spatially varying case */
1633:     PetscInt     dim, vStart, vEnd, nrow;
1634:     PetscScalar *M, *Mi;
1636:     PetscCall(DMGetDimension(dm, &dim));
1637:     if (isotropic) nrow = 1;
1638:     else nrow = dim;
1639:     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1640:     PetscCall(VecGetArray(metricInt, &met));
1641:     for (i = 1; i < numMetrics; ++i) {
1642:       PetscCall(VecGetArray(metrics[i], &meti));
1643:       for (v = vStart; v < vEnd; ++v) {
1644:         PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
1645:         PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
1646:         PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
1647:       }
1648:       PetscCall(VecRestoreArray(metrics[i], &meti));
1649:     }
1650:     PetscCall(VecRestoreArray(metricInt, &met));
1651:   }
1653:   PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1654:   PetscFunctionReturn(PETSC_SUCCESS);
1655: }
1657: /*@
1658:   DMPlexMetricIntersection2 - Compute the intersection of two metrics
1660:   Input Parameters:
1661: + dm      - The `DM`
1662: . metric1 - The first metric to be intersected
1663: - metric2 - The second metric to be intersected
1665:   Output Parameter:
1666: . metricInt - The intersected metric
1668:   Level: beginner
1670: .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1671: @*/
1672: PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1673: {
1674:   Vec metrics[2] = {metric1, metric2};
1676:   PetscFunctionBegin;
1677:   PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
1678:   PetscFunctionReturn(PETSC_SUCCESS);
1679: }
1681: /*@
1682:   DMPlexMetricIntersection3 - Compute the intersection of three metrics
1684:   Input Parameters:
1685: + dm      - The `DM`
1686: . metric1 - The first metric to be intersected
1687: . metric2 - The second metric to be intersected
1688: - metric3 - The third metric to be intersected
1690:   Output Parameter:
1691: . metricInt - The intersected metric
1693:   Level: beginner
1695: .seealso: `DMPLEX`, `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1696: @*/
1697: PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1698: {
1699:   Vec metrics[3] = {metric1, metric2, metric3};
1701:   PetscFunctionBegin;
1702:   PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
1703:   PetscFunctionReturn(PETSC_SUCCESS);
1704: }