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: `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: `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: Notes:
110: If the metric is specified as uniform then it is assumed to be isotropic, too.
112: .seealso: `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
113: @*/
114: PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform)
115: {
116: DM_Plex *plex = (DM_Plex *)dm->data;
118: PetscFunctionBegin;
119: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
120: plex->metricCtx->uniform = uniform;
121: if (uniform) plex->metricCtx->isotropic = uniform;
122: PetscFunctionReturn(PETSC_SUCCESS);
123: }
125: /*@
126: DMPlexMetricIsUniform - Is a metric uniform?
128: Input parameters:
129: . dm - The DM
131: Output parameters:
132: . uniform - Is the metric uniform?
134: Level: beginner
136: .seealso: `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
137: @*/
138: PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform)
139: {
140: DM_Plex *plex = (DM_Plex *)dm->data;
142: PetscFunctionBegin;
143: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
144: *uniform = plex->metricCtx->uniform;
145: PetscFunctionReturn(PETSC_SUCCESS);
146: }
148: /*@
149: DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization
151: Input parameters:
152: + dm - The DM
153: - restrictAnisotropyFirst - Should anisotropy be normalized first?
155: Level: beginner
157: .seealso: `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()`
158: @*/
159: PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst)
160: {
161: DM_Plex *plex = (DM_Plex *)dm->data;
163: PetscFunctionBegin;
164: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
165: plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst;
166: PetscFunctionReturn(PETSC_SUCCESS);
167: }
169: /*@
170: DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after?
172: Input parameters:
173: . dm - The DM
175: Output parameters:
176: . restrictAnisotropyFirst - Is anisotropy be normalized first?
178: Level: beginner
180: .seealso: `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()`
181: @*/
182: PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst)
183: {
184: DM_Plex *plex = (DM_Plex *)dm->data;
186: PetscFunctionBegin;
187: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
188: *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst;
189: PetscFunctionReturn(PETSC_SUCCESS);
190: }
192: /*@
193: DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off?
195: Input parameters:
196: + dm - The DM
197: - noInsert - Should node insertion and deletion be turned off?
199: Level: beginner
201: Notes:
202: This is only used by Mmg and ParMmg (not Pragmatic).
204: .seealso: `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
205: @*/
206: PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert)
207: {
208: DM_Plex *plex = (DM_Plex *)dm->data;
210: PetscFunctionBegin;
211: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
212: plex->metricCtx->noInsert = noInsert;
213: PetscFunctionReturn(PETSC_SUCCESS);
214: }
216: /*@
217: DMPlexMetricNoInsertion - Are node insertion and deletion turned off?
219: Input parameters:
220: . dm - The DM
222: Output parameters:
223: . noInsert - Are node insertion and deletion turned off?
225: Level: beginner
227: Notes:
228: This is only used by Mmg and ParMmg (not Pragmatic).
230: .seealso: `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
231: @*/
232: PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert)
233: {
234: DM_Plex *plex = (DM_Plex *)dm->data;
236: PetscFunctionBegin;
237: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
238: *noInsert = plex->metricCtx->noInsert;
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }
242: /*@
243: DMPlexMetricSetNoSwapping - Should facet swapping be turned off?
245: Input parameters:
246: + dm - The DM
247: - noSwap - Should facet swapping be turned off?
249: Level: beginner
251: Notes:
252: This is only used by Mmg and ParMmg (not Pragmatic).
254: .seealso: `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()`
255: @*/
256: PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap)
257: {
258: DM_Plex *plex = (DM_Plex *)dm->data;
260: PetscFunctionBegin;
261: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
262: plex->metricCtx->noSwap = noSwap;
263: PetscFunctionReturn(PETSC_SUCCESS);
264: }
266: /*@
267: DMPlexMetricNoSwapping - Is facet swapping turned off?
269: Input parameters:
270: . dm - The DM
272: Output parameters:
273: . noSwap - Is facet swapping turned off?
275: Level: beginner
277: Notes:
278: This is only used by Mmg and ParMmg (not Pragmatic).
280: .seealso: `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()`
281: @*/
282: PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap)
283: {
284: DM_Plex *plex = (DM_Plex *)dm->data;
286: PetscFunctionBegin;
287: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
288: *noSwap = plex->metricCtx->noSwap;
289: PetscFunctionReturn(PETSC_SUCCESS);
290: }
292: /*@
293: DMPlexMetricSetNoMovement - Should node movement be turned off?
295: Input parameters:
296: + dm - The DM
297: - noMove - Should node movement be turned off?
299: Level: beginner
301: Notes:
302: This is only used by Mmg and ParMmg (not Pragmatic).
304: .seealso: `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()`
305: @*/
306: PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove)
307: {
308: DM_Plex *plex = (DM_Plex *)dm->data;
310: PetscFunctionBegin;
311: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
312: plex->metricCtx->noMove = noMove;
313: PetscFunctionReturn(PETSC_SUCCESS);
314: }
316: /*@
317: DMPlexMetricNoMovement - Is node movement turned off?
319: Input parameters:
320: . dm - The DM
322: Output parameters:
323: . noMove - Is node movement turned off?
325: Level: beginner
327: Notes:
328: This is only used by Mmg and ParMmg (not Pragmatic).
330: .seealso: `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()`
331: @*/
332: PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove)
333: {
334: DM_Plex *plex = (DM_Plex *)dm->data;
336: PetscFunctionBegin;
337: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
338: *noMove = plex->metricCtx->noMove;
339: PetscFunctionReturn(PETSC_SUCCESS);
340: }
342: /*@
343: DMPlexMetricSetNoSurf - Should surface modification be turned off?
345: Input parameters:
346: + dm - The DM
347: - noSurf - Should surface modification be turned off?
349: Level: beginner
351: Notes:
352: This is only used by Mmg and ParMmg (not Pragmatic).
354: .seealso: `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`
355: @*/
356: PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf)
357: {
358: DM_Plex *plex = (DM_Plex *)dm->data;
360: PetscFunctionBegin;
361: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
362: plex->metricCtx->noSurf = noSurf;
363: PetscFunctionReturn(PETSC_SUCCESS);
364: }
366: /*@
367: DMPlexMetricNoSurf - Is surface modification turned off?
369: Input parameters:
370: . dm - The DM
372: Output parameters:
373: . noSurf - Is surface modification turned off?
375: Level: beginner
377: Notes:
378: This is only used by Mmg and ParMmg (not Pragmatic).
380: .seealso: `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`
381: @*/
382: PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf)
383: {
384: DM_Plex *plex = (DM_Plex *)dm->data;
386: PetscFunctionBegin;
387: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
388: *noSurf = plex->metricCtx->noSurf;
389: PetscFunctionReturn(PETSC_SUCCESS);
390: }
392: /*@
393: DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude
395: Input parameters:
396: + dm - The DM
397: - h_min - The minimum tolerated metric magnitude
399: Level: beginner
401: .seealso: `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()`
402: @*/
403: PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min)
404: {
405: DM_Plex *plex = (DM_Plex *)dm->data;
407: PetscFunctionBegin;
408: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
409: PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
410: plex->metricCtx->h_min = h_min;
411: PetscFunctionReturn(PETSC_SUCCESS);
412: }
414: /*@
415: DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude
417: Input parameters:
418: . dm - The DM
420: Output parameters:
421: . h_min - The minimum tolerated metric magnitude
423: Level: beginner
425: .seealso: `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()`
426: @*/
427: PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min)
428: {
429: DM_Plex *plex = (DM_Plex *)dm->data;
431: PetscFunctionBegin;
432: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
433: *h_min = plex->metricCtx->h_min;
434: PetscFunctionReturn(PETSC_SUCCESS);
435: }
437: /*@
438: DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude
440: Input parameters:
441: + dm - The DM
442: - h_max - The maximum tolerated metric magnitude
444: Level: beginner
446: .seealso: `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()`
447: @*/
448: PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max)
449: {
450: DM_Plex *plex = (DM_Plex *)dm->data;
452: PetscFunctionBegin;
453: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
454: PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)");
455: plex->metricCtx->h_max = h_max;
456: PetscFunctionReturn(PETSC_SUCCESS);
457: }
459: /*@
460: DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude
462: Input parameters:
463: . dm - The DM
465: Output parameters:
466: . h_max - The maximum tolerated metric magnitude
468: Level: beginner
470: .seealso: `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()`
471: @*/
472: PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max)
473: {
474: DM_Plex *plex = (DM_Plex *)dm->data;
476: PetscFunctionBegin;
477: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
478: *h_max = plex->metricCtx->h_max;
479: PetscFunctionReturn(PETSC_SUCCESS);
480: }
482: /*@
483: DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy
485: Input parameters:
486: + dm - The DM
487: - a_max - The maximum tolerated metric anisotropy
489: Level: beginner
491: Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one.
493: .seealso: `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: `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: `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: `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: `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: `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:
630: The gradation factor is the maximum tolerated length ratio between adjacent edges.
632: Turn off gradation by passing the value -1. Otherwise, pass a positive value.
634: This is only used by Mmg and ParMmg (not Pragmatic).
636: .seealso: `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
637: @*/
638: PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta)
639: {
640: DM_Plex *plex = (DM_Plex *)dm->data;
642: PetscFunctionBegin;
643: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
644: PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)");
645: plex->metricCtx->gradationFactor = beta;
646: PetscFunctionReturn(PETSC_SUCCESS);
647: }
649: /*@
650: DMPlexMetricGetGradationFactor - Get the metric gradation factor
652: Input parameters:
653: . dm - The DM
655: Output parameters:
656: . beta - The metric gradation factor
658: Level: beginner
660: Notes:
662: The gradation factor is the maximum tolerated length ratio between adjacent edges.
664: The value -1 implies that gradation is turned off.
666: This is only used by Mmg and ParMmg (not Pragmatic).
668: .seealso: `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
669: @*/
670: PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta)
671: {
672: DM_Plex *plex = (DM_Plex *)dm->data;
674: PetscFunctionBegin;
675: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
676: *beta = plex->metricCtx->gradationFactor;
677: PetscFunctionReturn(PETSC_SUCCESS);
678: }
680: /*@
681: DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number
683: Input parameters:
684: + dm - The DM
685: - hausd - The metric Hausdorff number
687: Level: beginner
689: Notes:
691: The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
692: boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
693: high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
694: an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
695: (resp. increase) the Hausdorff parameter. (Taken from
696: https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
698: This is only used by Mmg and ParMmg (not Pragmatic).
700: .seealso: `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()`
701: @*/
702: PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd)
703: {
704: DM_Plex *plex = (DM_Plex *)dm->data;
706: PetscFunctionBegin;
707: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
708: PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)");
709: plex->metricCtx->hausdorffNumber = hausd;
710: PetscFunctionReturn(PETSC_SUCCESS);
711: }
713: /*@
714: DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number
716: Input parameters:
717: . dm - The DM
719: Output parameters:
720: . hausd - The metric Hausdorff number
722: Level: beginner
724: Notes:
726: The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the
727: boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the
728: high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for
729: an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease
730: (resp. increase) the Hausdorff parameter. (Taken from
731: https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd).
733: This is only used by Mmg and ParMmg (not Pragmatic).
735: .seealso: `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()`
736: @*/
737: PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd)
738: {
739: DM_Plex *plex = (DM_Plex *)dm->data;
741: PetscFunctionBegin;
742: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
743: *hausd = plex->metricCtx->hausdorffNumber;
744: PetscFunctionReturn(PETSC_SUCCESS);
745: }
747: /*@
748: DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package
750: Input parameters:
751: + dm - The DM
752: - verbosity - The verbosity, where -1 is silent and 10 is maximum
754: Level: beginner
756: Notes:
757: This is only used by Mmg and ParMmg (not Pragmatic).
759: .seealso: `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()`
760: @*/
761: PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity)
762: {
763: DM_Plex *plex = (DM_Plex *)dm->data;
765: PetscFunctionBegin;
766: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
767: plex->metricCtx->verbosity = verbosity;
768: PetscFunctionReturn(PETSC_SUCCESS);
769: }
771: /*@
772: DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package
774: Input parameters:
775: . dm - The DM
777: Output parameters:
778: . verbosity - The verbosity, where -1 is silent and 10 is maximum
780: Level: beginner
782: Notes:
783: This is only used by Mmg and ParMmg (not Pragmatic).
785: .seealso: `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
786: @*/
787: PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity)
788: {
789: DM_Plex *plex = (DM_Plex *)dm->data;
791: PetscFunctionBegin;
792: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
793: *verbosity = plex->metricCtx->verbosity;
794: PetscFunctionReturn(PETSC_SUCCESS);
795: }
797: /*@
798: DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations
800: Input parameters:
801: + dm - The DM
802: - numIter - the number of parallel adaptation iterations
804: Level: beginner
806: Notes:
807: This is only used by ParMmg (not Pragmatic or Mmg).
809: .seealso: `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()`
810: @*/
811: PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter)
812: {
813: DM_Plex *plex = (DM_Plex *)dm->data;
815: PetscFunctionBegin;
816: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
817: plex->metricCtx->numIter = numIter;
818: PetscFunctionReturn(PETSC_SUCCESS);
819: }
821: /*@
822: DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations
824: Input parameters:
825: . dm - The DM
827: Output parameters:
828: . numIter - the number of parallel adaptation iterations
830: Level: beginner
832: Notes:
833: This is only used by Mmg and ParMmg (not Pragmatic or Mmg).
835: .seealso: `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()`
836: @*/
837: PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter)
838: {
839: DM_Plex *plex = (DM_Plex *)dm->data;
841: PetscFunctionBegin;
842: if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm));
843: *numIter = plex->metricCtx->numIter;
844: PetscFunctionReturn(PETSC_SUCCESS);
845: }
847: PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
848: {
849: MPI_Comm comm;
850: PetscFE fe;
851: PetscInt dim;
853: PetscFunctionBegin;
855: /* Extract metadata from dm */
856: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
857: PetscCall(DMGetDimension(dm, &dim));
859: /* Create a P1 field of the requested size */
860: PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
861: PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
862: PetscCall(DMCreateDS(dm));
863: PetscCall(PetscFEDestroy(&fe));
864: PetscCall(DMCreateLocalVector(dm, metric));
866: PetscFunctionReturn(PETSC_SUCCESS);
867: }
869: /*@
870: DMPlexMetricCreate - Create a Riemannian metric field
872: Input parameters:
873: + dm - The DM
874: - f - The field number to use
876: Output parameter:
877: . metric - The metric
879: Level: beginner
881: Notes:
883: It is assumed that the DM is comprised of simplices.
885: Command line options for Riemannian metrics:
887: + -dm_plex_metric_isotropic - Is the metric isotropic?
888: . -dm_plex_metric_uniform - Is the metric uniform?
889: . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
890: . -dm_plex_metric_h_min - Minimum tolerated metric magnitude
891: . -dm_plex_metric_h_max - Maximum tolerated metric magnitude
892: . -dm_plex_metric_a_max - Maximum tolerated anisotropy
893: . -dm_plex_metric_p - L-p normalization order
894: - -dm_plex_metric_target_complexity - Target metric complexity
896: Switching between remeshers can be achieved using
898: . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use
900: Further options that are only relevant to Mmg and ParMmg:
902: + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation
903: . -dm_plex_metric_num_iterations - Number of parallel mesh adaptation iterations for ParMmg
904: . -dm_plex_metric_no_insert - Should node insertion/deletion be turned off?
905: . -dm_plex_metric_no_swap - Should facet swapping be turned off?
906: . -dm_plex_metric_no_move - Should node movement be turned off?
907: - -dm_plex_metric_verbosity - Choose a verbosity level from -1 (silent) to 10 (maximum).
909: .seealso: `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()`
910: @*/
911: PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
912: {
913: PetscBool isotropic, uniform;
914: PetscInt coordDim, Nd;
916: PetscFunctionBegin;
917: PetscCall(DMGetCoordinateDim(dm, &coordDim));
918: Nd = coordDim * coordDim;
919: PetscCall(DMPlexMetricIsUniform(dm, &uniform));
920: PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
921: if (uniform) {
922: MPI_Comm comm;
924: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
925: PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
926: PetscCall(VecCreate(comm, metric));
927: PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE));
928: PetscCall(VecSetFromOptions(*metric));
929: } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric));
930: else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric));
931: PetscFunctionReturn(PETSC_SUCCESS);
932: }
934: /*@
935: DMPlexMetricCreateUniform - Construct a uniform isotropic metric
937: Input parameters:
938: + dm - The DM
939: . f - The field number to use
940: - alpha - Scaling parameter for the diagonal
942: Output parameter:
943: . metric - The uniform metric
945: Level: beginner
947: Note: In this case, the metric is constant in space and so is only specified for a single vertex.
949: .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()`
950: @*/
951: PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
952: {
953: PetscFunctionBegin;
954: PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE));
955: PetscCall(DMPlexMetricCreate(dm, f, metric));
956: PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
957: PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)");
958: PetscCall(VecSet(*metric, alpha));
959: PetscCall(VecAssemblyBegin(*metric));
960: PetscCall(VecAssemblyEnd(*metric));
961: PetscFunctionReturn(PETSC_SUCCESS);
962: }
964: 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[])
965: {
966: f0[0] = u[0];
967: }
969: /*@
970: DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
972: Input parameters:
973: + dm - The DM
974: . f - The field number to use
975: - indicator - The error indicator
977: Output parameter:
978: . metric - The isotropic metric
980: Level: beginner
982: Notes:
984: It is assumed that the DM is comprised of simplices.
986: The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately.
988: .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()`
989: @*/
990: PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
991: {
992: PetscInt m, n;
994: PetscFunctionBegin;
995: PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE));
996: PetscCall(DMPlexMetricCreate(dm, f, metric));
997: PetscCall(VecGetSize(indicator, &m));
998: PetscCall(VecGetSize(*metric, &n));
999: if (m == n) PetscCall(VecCopy(indicator, *metric));
1000: else {
1001: DM dmIndi;
1002: 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[]);
1004: PetscCall(VecGetDM(indicator, &dmIndi));
1005: funcs[0] = identity;
1006: PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric));
1007: }
1008: PetscFunctionReturn(PETSC_SUCCESS);
1009: }
1011: /*@
1012: DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric
1014: Input parameters:
1015: + dm - The DM of the metric field
1016: - f - The field number to use
1018: Output parameter:
1019: + determinant - The determinant field
1020: - dmDet - The corresponding DM
1022: Level: beginner
1024: .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic(), DMPlexMetricCreate()
1025: @*/
1026: PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet)
1027: {
1028: PetscBool uniform;
1030: PetscFunctionBegin;
1031: PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1032: PetscCall(DMClone(dm, dmDet));
1033: if (uniform) {
1034: MPI_Comm comm;
1036: PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm));
1037: PetscCall(VecCreate(comm, determinant));
1038: PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE));
1039: PetscCall(VecSetFromOptions(*determinant));
1040: } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant));
1041: PetscFunctionReturn(PETSC_SUCCESS);
1042: }
1044: static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[])
1045: {
1046: PetscInt i, j;
1048: PetscFunctionBegin;
1049: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"));
1050: for (i = 0; i < dim; ++i) {
1051: if (i == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " [["));
1052: else PetscCall(PetscPrintf(PETSC_COMM_SELF, " ["));
1053: for (j = 0; j < dim; ++j) {
1054: if (j < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])));
1055: else PetscCall(PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])));
1056: }
1057: if (i < dim - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
1058: else PetscCall(PetscPrintf(PETSC_COMM_SELF, "]]\n"));
1059: }
1060: PetscFunctionReturn(PETSC_SUCCESS);
1061: }
1063: static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp)
1064: {
1065: PetscInt i, j, k;
1066: 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);
1067: PetscScalar *Mpos;
1069: PetscFunctionBegin;
1070: PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs));
1072: /* Symmetrize */
1073: for (i = 0; i < dim; ++i) {
1074: Mpos[i * dim + i] = Mp[i * dim + i];
1075: for (j = i + 1; j < dim; ++j) {
1076: Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
1077: Mpos[j * dim + i] = Mpos[i * dim + j];
1078: }
1079: }
1081: /* Compute eigendecomposition */
1082: if (dim == 1) {
1083: /* Isotropic case */
1084: eigs[0] = PetscRealPart(Mpos[0]);
1085: Mpos[0] = 1.0;
1086: } else {
1087: /* Anisotropic case */
1088: PetscScalar *work;
1089: PetscBLASInt lwork;
1091: lwork = 5 * dim;
1092: PetscCall(PetscMalloc1(5 * dim, &work));
1093: {
1094: PetscBLASInt lierr;
1095: PetscBLASInt nb;
1097: PetscCall(PetscBLASIntCast(dim, &nb));
1098: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1099: #if defined(PETSC_USE_COMPLEX)
1100: {
1101: PetscReal *rwork;
1102: PetscCall(PetscMalloc1(3 * dim, &rwork));
1103: PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr));
1104: PetscCall(PetscFree(rwork));
1105: }
1106: #else
1107: PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr));
1108: #endif
1109: if (lierr) {
1110: for (i = 0; i < dim; ++i) {
1111: Mpos[i * dim + i] = Mp[i * dim + i];
1112: for (j = i + 1; j < dim; ++j) {
1113: Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]);
1114: Mpos[j * dim + i] = Mpos[i * dim + j];
1115: }
1116: }
1117: PetscCall(LAPACKsyevFail(dim, Mpos));
1118: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
1119: }
1120: PetscCall(PetscFPTrapPop());
1121: }
1122: PetscCall(PetscFree(work));
1123: }
1125: /* Reflect to positive orthant and enforce maximum and minimum size */
1126: max_eig = 0.0;
1127: for (i = 0; i < dim; ++i) {
1128: eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
1129: max_eig = PetscMax(eigs[i], max_eig);
1130: }
1132: /* Enforce maximum anisotropy and compute determinant */
1133: *detMp = 1.0;
1134: for (i = 0; i < dim; ++i) {
1135: if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min);
1136: *detMp *= eigs[i];
1137: }
1139: /* Reconstruct Hessian */
1140: for (i = 0; i < dim; ++i) {
1141: for (j = 0; j < dim; ++j) {
1142: Mp[i * dim + j] = 0.0;
1143: for (k = 0; k < dim; ++k) Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j];
1144: }
1145: }
1146: PetscCall(PetscFree2(Mpos, eigs));
1148: PetscFunctionReturn(PETSC_SUCCESS);
1149: }
1151: /*@
1152: DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
1154: Input parameters:
1155: + dm - The DM
1156: . metricIn - The metric
1157: . restrictSizes - Should maximum/minimum metric magnitudes be enforced?
1158: - restrictAnisotropy - Should maximum anisotropy be enforced?
1160: Output parameter:
1161: + metricOut - The metric
1162: - determinant - Its determinant
1164: Level: beginner
1166: Notes:
1168: Relevant command line options:
1170: + -dm_plex_metric_isotropic - Is the metric isotropic?
1171: . -dm_plex_metric_uniform - Is the metric uniform?
1172: . -dm_plex_metric_h_min - Minimum tolerated metric magnitude
1173: . -dm_plex_metric_h_max - Maximum tolerated metric magnitude
1174: - -dm_plex_metric_a_max - Maximum tolerated anisotropy
1176: .seealso: `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()`
1177: @*/
1178: PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1179: {
1180: DM dmDet;
1181: PetscBool isotropic, uniform;
1182: PetscInt dim, vStart, vEnd, v;
1183: PetscScalar *met, *det;
1184: PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
1186: PetscFunctionBegin;
1187: PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1189: /* Extract metadata from dm */
1190: PetscCall(DMGetDimension(dm, &dim));
1191: if (restrictSizes) {
1192: PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1193: PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1194: h_min = PetscMax(h_min, 1.0e-30);
1195: h_max = PetscMin(h_max, 1.0e+30);
1196: PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1197: }
1198: if (restrictAnisotropy) {
1199: PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1200: a_max = PetscMin(a_max, 1.0e+30);
1201: }
1203: /* Setup output metric */
1204: PetscCall(VecCopy(metricIn, metricOut));
1206: /* Enforce SPD and extract determinant */
1207: PetscCall(VecGetArray(metricOut, &met));
1208: PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1209: PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1210: if (uniform) {
1211: PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist");
1213: /* Uniform case */
1214: PetscCall(VecGetArray(determinant, &det));
1215: PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1216: PetscCall(VecRestoreArray(determinant, &det));
1217: } else {
1218: /* Spatially varying case */
1219: PetscInt nrow;
1221: if (isotropic) nrow = 1;
1222: else nrow = dim;
1223: PetscCall(VecGetDM(determinant, &dmDet));
1224: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1225: PetscCall(VecGetArray(determinant, &det));
1226: for (v = vStart; v < vEnd; ++v) {
1227: PetscScalar *vmet, *vdet;
1228: PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet));
1229: PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet));
1230: PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet));
1231: }
1232: PetscCall(VecRestoreArray(determinant, &det));
1233: }
1234: PetscCall(VecRestoreArray(metricOut, &met));
1236: PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0));
1237: PetscFunctionReturn(PETSC_SUCCESS);
1238: }
1240: 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[])
1241: {
1242: const PetscScalar p = constants[0];
1244: f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim));
1245: }
1247: /*@
1248: DMPlexMetricNormalize - Apply L-p normalization to a metric
1250: Input parameters:
1251: + dm - The DM
1252: . metricIn - The unnormalized metric
1253: . restrictSizes - Should maximum/minimum metric magnitudes be enforced?
1254: - restrictAnisotropy - Should maximum metric anisotropy be enforced?
1256: Output parameter:
1257: . metricOut - The normalized metric
1259: Level: beginner
1261: Notes:
1263: Relevant command line options:
1265: + -dm_plex_metric_isotropic - Is the metric isotropic?
1266: . -dm_plex_metric_uniform - Is the metric uniform?
1267: . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization?
1268: . -dm_plex_metric_h_min - Minimum tolerated metric magnitude
1269: . -dm_plex_metric_h_max - Maximum tolerated metric magnitude
1270: . -dm_plex_metric_a_max - Maximum tolerated anisotropy
1271: . -dm_plex_metric_p - L-p normalization order
1272: - -dm_plex_metric_target_complexity - Target metric complexity
1274: .seealso: `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()`
1275: @*/
1276: PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant)
1277: {
1278: DM dmDet;
1279: MPI_Comm comm;
1280: PetscBool restrictAnisotropyFirst, isotropic, uniform;
1281: PetscDS ds;
1282: PetscInt dim, Nd, vStart, vEnd, v, i;
1283: PetscScalar *met, *det, integral, constants[1];
1284: PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral;
1286: PetscFunctionBegin;
1287: PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1289: /* Extract metadata from dm */
1290: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1291: PetscCall(DMGetDimension(dm, &dim));
1292: PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1293: PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1294: if (isotropic) Nd = 1;
1295: else Nd = dim * dim;
1297: /* Set up metric and ensure it is SPD */
1298: PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst));
1299: PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant));
1301: /* Compute global normalization factor */
1302: PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
1303: PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p));
1304: constants[0] = p;
1305: if (uniform) {
1306: PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported");
1307: DM dmTmp;
1308: Vec tmp;
1310: PetscCall(DMClone(dm, &dmTmp));
1311: PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp));
1312: PetscCall(VecGetArray(determinant, &det));
1313: PetscCall(VecSet(tmp, det[0]));
1314: PetscCall(VecRestoreArray(determinant, &det));
1315: PetscCall(DMGetDS(dmTmp, &ds));
1316: PetscCall(PetscDSSetConstants(ds, 1, constants));
1317: PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1318: PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL));
1319: PetscCall(VecDestroy(&tmp));
1320: PetscCall(DMDestroy(&dmTmp));
1321: } else {
1322: PetscCall(VecGetDM(determinant, &dmDet));
1323: PetscCall(DMGetDS(dmDet, &ds));
1324: PetscCall(PetscDSSetConstants(ds, 1, constants));
1325: PetscCall(PetscDSSetObjective(ds, 0, detMFunc));
1326: PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL));
1327: }
1328: realIntegral = PetscRealPart(integral);
1329: 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?");
1330: factGlob = PetscPowReal(target / realIntegral, 2.0 / dim);
1332: /* Apply local scaling */
1333: if (restrictSizes) {
1334: PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min));
1335: PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max));
1336: h_min = PetscMax(h_min, 1.0e-30);
1337: h_max = PetscMin(h_max, 1.0e+30);
1338: PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude");
1339: }
1340: if (restrictAnisotropy && !restrictAnisotropyFirst) {
1341: PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max));
1342: a_max = PetscMin(a_max, 1.0e+30);
1343: }
1344: PetscCall(VecGetArray(metricOut, &met));
1345: PetscCall(VecGetArray(determinant, &det));
1346: if (uniform) {
1347: /* Uniform case */
1348: met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim));
1349: if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det));
1350: } else {
1351: /* Spatially varying case */
1352: PetscInt nrow;
1354: if (isotropic) nrow = 1;
1355: else nrow = dim;
1356: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1357: PetscCall(VecGetDM(determinant, &dmDet));
1358: for (v = vStart; v < vEnd; ++v) {
1359: PetscScalar *Mp, *detM;
1361: PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp));
1362: PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM));
1363: fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim));
1364: for (i = 0; i < Nd; ++i) Mp[i] *= fact;
1365: if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM));
1366: }
1367: }
1368: PetscCall(VecRestoreArray(determinant, &det));
1369: PetscCall(VecRestoreArray(metricOut, &met));
1371: PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0));
1372: PetscFunctionReturn(PETSC_SUCCESS);
1373: }
1375: /*@
1376: DMPlexMetricAverage - Compute the average of a list of metrics
1378: Input Parameters:
1379: + dm - The DM
1380: . numMetrics - The number of metrics to be averaged
1381: . weights - Weights for the average
1382: - metrics - The metrics to be averaged
1384: Output Parameter:
1385: . metricAvg - The averaged metric
1387: Level: beginner
1389: Notes:
1390: The weights should sum to unity.
1392: If weights are not provided then an unweighted average is used.
1394: .seealso: `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()`
1395: @*/
1396: PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg)
1397: {
1398: PetscBool haveWeights = PETSC_TRUE;
1399: PetscInt i, m, n;
1400: PetscReal sum = 0.0, tol = 1.0e-10;
1402: PetscFunctionBegin;
1403: PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0));
1404: PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics);
1405: PetscCall(VecSet(metricAvg, 0.0));
1406: PetscCall(VecGetSize(metricAvg, &m));
1407: for (i = 0; i < numMetrics; ++i) {
1408: PetscCall(VecGetSize(metrics[i], &n));
1409: PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented");
1410: }
1412: /* Default to the unweighted case */
1413: if (!weights) {
1414: PetscCall(PetscMalloc1(numMetrics, &weights));
1415: haveWeights = PETSC_FALSE;
1416: for (i = 0; i < numMetrics; ++i) weights[i] = 1.0 / numMetrics;
1417: }
1419: /* Check weights sum to unity */
1420: for (i = 0; i < numMetrics; ++i) sum += weights[i];
1421: PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity");
1423: /* Compute metric average */
1424: for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i]));
1425: if (!haveWeights) PetscCall(PetscFree(weights));
1427: PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0));
1428: PetscFunctionReturn(PETSC_SUCCESS);
1429: }
1431: /*@
1432: DMPlexMetricAverage2 - Compute the unweighted average of two metrics
1434: Input Parameters:
1435: + dm - The DM
1436: . metric1 - The first metric to be averaged
1437: - metric2 - The second metric to be averaged
1439: Output Parameter:
1440: . metricAvg - The averaged metric
1442: Level: beginner
1444: .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage3()`
1445: @*/
1446: PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg)
1447: {
1448: PetscReal weights[2] = {0.5, 0.5};
1449: Vec metrics[2] = {metric1, metric2};
1451: PetscFunctionBegin;
1452: PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg));
1453: PetscFunctionReturn(PETSC_SUCCESS);
1454: }
1456: /*@
1457: DMPlexMetricAverage3 - Compute the unweighted average of three metrics
1459: Input Parameters:
1460: + dm - The DM
1461: . metric1 - The first metric to be averaged
1462: . metric2 - The second metric to be averaged
1463: - metric3 - The third metric to be averaged
1465: Output Parameter:
1466: . metricAvg - The averaged metric
1468: Level: beginner
1470: .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage2()`
1471: @*/
1472: PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg)
1473: {
1474: PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0};
1475: Vec metrics[3] = {metric1, metric2, metric3};
1477: PetscFunctionBegin;
1478: PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg));
1479: PetscFunctionReturn(PETSC_SUCCESS);
1480: }
1482: static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
1483: {
1484: PetscInt i, j, k, l, m;
1485: PetscReal *evals;
1486: PetscScalar *evecs, *sqrtM1, *isqrtM1;
1488: PetscFunctionBegin;
1490: /* Isotropic case */
1491: if (dim == 1) {
1492: M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0]));
1493: PetscFunctionReturn(PETSC_SUCCESS);
1494: }
1496: /* Anisotropic case */
1497: PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals));
1498: for (i = 0; i < dim; ++i) {
1499: for (j = 0; j < dim; ++j) evecs[i * dim + j] = M1[i * dim + j];
1500: }
1501: {
1502: PetscScalar *work;
1503: PetscBLASInt lwork;
1505: lwork = 5 * dim;
1506: PetscCall(PetscMalloc1(5 * dim, &work));
1507: {
1508: PetscBLASInt lierr, nb;
1509: PetscReal sqrtj;
1511: /* Compute eigendecomposition of M1 */
1512: PetscCall(PetscBLASIntCast(dim, &nb));
1513: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1514: #if defined(PETSC_USE_COMPLEX)
1515: {
1516: PetscReal *rwork;
1517: PetscCall(PetscMalloc1(3 * dim, &rwork));
1518: PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1519: PetscCall(PetscFree(rwork));
1520: }
1521: #else
1522: PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1523: #endif
1524: if (lierr) {
1525: PetscCall(LAPACKsyevFail(dim, M1));
1526: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
1527: }
1528: PetscCall(PetscFPTrapPop());
1530: /* Compute square root and the reciprocal thereof */
1531: for (i = 0; i < dim; ++i) {
1532: for (k = 0; k < dim; ++k) {
1533: sqrtM1[i * dim + k] = 0.0;
1534: isqrtM1[i * dim + k] = 0.0;
1535: for (j = 0; j < dim; ++j) {
1536: sqrtj = PetscSqrtReal(evals[j]);
1537: sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k];
1538: isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k];
1539: }
1540: }
1541: }
1543: /* Map M2 into the space spanned by the eigenvectors of M1 */
1544: for (i = 0; i < dim; ++i) {
1545: for (l = 0; l < dim; ++l) {
1546: evecs[i * dim + l] = 0.0;
1547: for (j = 0; j < dim; ++j) {
1548: for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
1549: }
1550: }
1551: }
1553: /* Compute eigendecomposition */
1554: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1555: #if defined(PETSC_USE_COMPLEX)
1556: {
1557: PetscReal *rwork;
1558: PetscCall(PetscMalloc1(3 * dim, &rwork));
1559: PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
1560: PetscCall(PetscFree(rwork));
1561: }
1562: #else
1563: PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
1564: #endif
1565: if (lierr) {
1566: for (i = 0; i < dim; ++i) {
1567: for (l = 0; l < dim; ++l) {
1568: evecs[i * dim + l] = 0.0;
1569: for (j = 0; j < dim; ++j) {
1570: for (k = 0; k < dim; ++k) evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l];
1571: }
1572: }
1573: }
1574: PetscCall(LAPACKsyevFail(dim, evecs));
1575: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
1576: }
1577: PetscCall(PetscFPTrapPop());
1579: /* Modify eigenvalues */
1580: for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0);
1582: /* Map back to get the intersection */
1583: for (i = 0; i < dim; ++i) {
1584: for (m = 0; m < dim; ++m) {
1585: M2[i * dim + m] = 0.0;
1586: for (j = 0; j < dim; ++j) {
1587: for (k = 0; k < dim; ++k) {
1588: 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];
1589: }
1590: }
1591: }
1592: }
1593: }
1594: PetscCall(PetscFree(work));
1595: }
1596: PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals));
1597: PetscFunctionReturn(PETSC_SUCCESS);
1598: }
1600: /*@
1601: DMPlexMetricIntersection - Compute the intersection of a list of metrics
1603: Input Parameters:
1604: + dm - The DM
1605: . numMetrics - The number of metrics to be intersected
1606: - metrics - The metrics to be intersected
1608: Output Parameter:
1609: . metricInt - The intersected metric
1611: Level: beginner
1613: Notes:
1615: The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics.
1617: The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2.
1619: .seealso: `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()`
1620: @*/
1621: PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt)
1622: {
1623: PetscBool isotropic, uniform;
1624: PetscInt v, i, m, n;
1625: PetscScalar *met, *meti;
1627: PetscFunctionBegin;
1628: PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1629: PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics);
1631: /* Copy over the first metric */
1632: PetscCall(VecCopy(metrics[0], metricInt));
1633: if (numMetrics == 1) PetscFunctionReturn(PETSC_SUCCESS);
1634: PetscCall(VecGetSize(metricInt, &m));
1635: for (i = 0; i < numMetrics; ++i) {
1636: PetscCall(VecGetSize(metrics[i], &n));
1637: PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented");
1638: }
1640: /* Intersect subsequent metrics in turn */
1641: PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1642: PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1643: if (uniform) {
1644: /* Uniform case */
1645: PetscCall(VecGetArray(metricInt, &met));
1646: for (i = 1; i < numMetrics; ++i) {
1647: PetscCall(VecGetArray(metrics[i], &meti));
1648: PetscCall(DMPlexMetricIntersection_Private(1, meti, met));
1649: PetscCall(VecRestoreArray(metrics[i], &meti));
1650: }
1651: PetscCall(VecRestoreArray(metricInt, &met));
1652: } else {
1653: /* Spatially varying case */
1654: PetscInt dim, vStart, vEnd, nrow;
1655: PetscScalar *M, *Mi;
1657: PetscCall(DMGetDimension(dm, &dim));
1658: if (isotropic) nrow = 1;
1659: else nrow = dim;
1660: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1661: PetscCall(VecGetArray(metricInt, &met));
1662: for (i = 1; i < numMetrics; ++i) {
1663: PetscCall(VecGetArray(metrics[i], &meti));
1664: for (v = vStart; v < vEnd; ++v) {
1665: PetscCall(DMPlexPointLocalRef(dm, v, met, &M));
1666: PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi));
1667: PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M));
1668: }
1669: PetscCall(VecRestoreArray(metrics[i], &meti));
1670: }
1671: PetscCall(VecRestoreArray(metricInt, &met));
1672: }
1674: PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0));
1675: PetscFunctionReturn(PETSC_SUCCESS);
1676: }
1678: /*@
1679: DMPlexMetricIntersection2 - Compute the intersection of two metrics
1681: Input Parameters:
1682: + dm - The DM
1683: . metric1 - The first metric to be intersected
1684: - metric2 - The second metric to be intersected
1686: Output Parameter:
1687: . metricInt - The intersected metric
1689: Level: beginner
1691: .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()`
1692: @*/
1693: PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt)
1694: {
1695: Vec metrics[2] = {metric1, metric2};
1697: PetscFunctionBegin;
1698: PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt));
1699: PetscFunctionReturn(PETSC_SUCCESS);
1700: }
1702: /*@
1703: DMPlexMetricIntersection3 - Compute the intersection of three metrics
1705: Input Parameters:
1706: + dm - The DM
1707: . metric1 - The first metric to be intersected
1708: . metric2 - The second metric to be intersected
1709: - metric3 - The third metric to be intersected
1711: Output Parameter:
1712: . metricInt - The intersected metric
1714: Level: beginner
1716: .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()`
1717: @*/
1718: PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt)
1719: {
1720: Vec metrics[3] = {metric1, metric2, metric3};
1722: PetscFunctionBegin;
1723: PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt));
1724: PetscFunctionReturn(PETSC_SUCCESS);
1725: }