120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 220282da2SJoe Wallwork #include <petscblaslapack.h> 320282da2SJoe Wallwork 4fe902aceSJoe Wallwork PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection; 5fe902aceSJoe Wallwork 6*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetFromOptions(DM dm) { 7bc00d7afSJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 831b70465SJoe Wallwork MPI_Comm comm; 993520041SJoe Wallwork PetscBool isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE; 1076f360caSJoe Wallwork PetscBool noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE; 11cc2eee55SJoe Wallwork PetscInt verbosity = -1, numIter = 3; 12ae8b063eSJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0, beta = 1.3, hausd = 0.01; 1331b70465SJoe Wallwork 1431b70465SJoe Wallwork PetscFunctionBegin; 15bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(PetscNew(&plex->metricCtx)); 169566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 17d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric"); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL)); 199566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetIsotropic(dm, isotropic)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL)); 219566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetUniform(dm, uniform)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL)); 239566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL)); 259566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetNoInsertion(dm, noInsert)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL)); 279566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetNoSwapping(dm, noSwap)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL)); 299566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetNoMovement(dm, noMove)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL)); 319566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetNoSurf(dm, noSurf)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0)); 339566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetNumIterations(dm, numIter)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10)); 359566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetVerbosity(dm, verbosity)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL)); 379566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetMinimumMagnitude(dm, h_min)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL)); 399566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetMaximumMagnitude(dm, h_max)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL)); 419566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetMaximumAnisotropy(dm, a_max)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL)); 439566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetNormalizationOrder(dm, p)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL)); 459566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetTargetComplexity(dm, target)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL)); 479566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetGradationFactor(dm, beta)); 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL)); 499566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetHausdorffNumber(dm, hausd)); 50d0609cedSBarry Smith PetscOptionsEnd(); 5131b70465SJoe Wallwork PetscFunctionReturn(0); 5231b70465SJoe Wallwork } 5331b70465SJoe Wallwork 54cb7bfe8cSJoe Wallwork /*@ 5531b70465SJoe Wallwork DMPlexMetricSetIsotropic - Record whether a metric is isotropic 5631b70465SJoe Wallwork 5731b70465SJoe Wallwork Input parameters: 5831b70465SJoe Wallwork + dm - The DM 5931b70465SJoe Wallwork - isotropic - Is the metric isotropic? 6031b70465SJoe Wallwork 6131b70465SJoe Wallwork Level: beginner 6231b70465SJoe Wallwork 63db781477SPatrick Sanan .seealso: `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetUniform()`, `DMPlexMetricSetRestrictAnisotropyFirst()` 64cb7bfe8cSJoe Wallwork @*/ 65*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic) { 6631b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 6731b70465SJoe Wallwork 6831b70465SJoe Wallwork PetscFunctionBegin; 69bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 7031b70465SJoe Wallwork plex->metricCtx->isotropic = isotropic; 7131b70465SJoe Wallwork PetscFunctionReturn(0); 7231b70465SJoe Wallwork } 7331b70465SJoe Wallwork 74cb7bfe8cSJoe Wallwork /*@ 7593520041SJoe Wallwork DMPlexMetricIsIsotropic - Is a metric isotropic? 7631b70465SJoe Wallwork 7731b70465SJoe Wallwork Input parameters: 7831b70465SJoe Wallwork . dm - The DM 7931b70465SJoe Wallwork 8031b70465SJoe Wallwork Output parameters: 8131b70465SJoe Wallwork . isotropic - Is the metric isotropic? 8231b70465SJoe Wallwork 8331b70465SJoe Wallwork Level: beginner 8431b70465SJoe Wallwork 85db781477SPatrick Sanan .seealso: `DMPlexMetricSetIsotropic()`, `DMPlexMetricIsUniform()`, `DMPlexMetricRestrictAnisotropyFirst()` 86cb7bfe8cSJoe Wallwork @*/ 87*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic) { 8831b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 8931b70465SJoe Wallwork 9031b70465SJoe Wallwork PetscFunctionBegin; 91bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 9231b70465SJoe Wallwork *isotropic = plex->metricCtx->isotropic; 9331b70465SJoe Wallwork PetscFunctionReturn(0); 9431b70465SJoe Wallwork } 9531b70465SJoe Wallwork 96cb7bfe8cSJoe Wallwork /*@ 9793520041SJoe Wallwork DMPlexMetricSetUniform - Record whether a metric is uniform 9893520041SJoe Wallwork 9993520041SJoe Wallwork Input parameters: 10093520041SJoe Wallwork + dm - The DM 10193520041SJoe Wallwork - uniform - Is the metric uniform? 10293520041SJoe Wallwork 10393520041SJoe Wallwork Level: beginner 10493520041SJoe Wallwork 10593520041SJoe Wallwork Notes: 10693520041SJoe Wallwork 10793520041SJoe Wallwork If the metric is specified as uniform then it is assumed to be isotropic, too. 10893520041SJoe Wallwork 109db781477SPatrick Sanan .seealso: `DMPlexMetricIsUniform()`, `DMPlexMetricSetIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()` 11093520041SJoe Wallwork @*/ 111*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform) { 11293520041SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 11393520041SJoe Wallwork 11493520041SJoe Wallwork PetscFunctionBegin; 115bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 11693520041SJoe Wallwork plex->metricCtx->uniform = uniform; 11793520041SJoe Wallwork if (uniform) plex->metricCtx->isotropic = uniform; 11893520041SJoe Wallwork PetscFunctionReturn(0); 11993520041SJoe Wallwork } 12093520041SJoe Wallwork 12193520041SJoe Wallwork /*@ 12293520041SJoe Wallwork DMPlexMetricIsUniform - Is a metric uniform? 12393520041SJoe Wallwork 12493520041SJoe Wallwork Input parameters: 12593520041SJoe Wallwork . dm - The DM 12693520041SJoe Wallwork 12793520041SJoe Wallwork Output parameters: 12893520041SJoe Wallwork . uniform - Is the metric uniform? 12993520041SJoe Wallwork 13093520041SJoe Wallwork Level: beginner 13193520041SJoe Wallwork 132db781477SPatrick Sanan .seealso: `DMPlexMetricSetUniform()`, `DMPlexMetricIsIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()` 13393520041SJoe Wallwork @*/ 134*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform) { 13593520041SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 13693520041SJoe Wallwork 13793520041SJoe Wallwork PetscFunctionBegin; 138bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 13993520041SJoe Wallwork *uniform = plex->metricCtx->uniform; 14093520041SJoe Wallwork PetscFunctionReturn(0); 14193520041SJoe Wallwork } 14293520041SJoe Wallwork 14393520041SJoe Wallwork /*@ 14431b70465SJoe Wallwork DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization 14531b70465SJoe Wallwork 14631b70465SJoe Wallwork Input parameters: 14731b70465SJoe Wallwork + dm - The DM 14831b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first? 14931b70465SJoe Wallwork 15031b70465SJoe Wallwork Level: beginner 15131b70465SJoe Wallwork 152db781477SPatrick Sanan .seealso: `DMPlexMetricSetIsotropic()`, `DMPlexMetricRestrictAnisotropyFirst()` 153cb7bfe8cSJoe Wallwork @*/ 154*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst) { 15531b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 15631b70465SJoe Wallwork 15731b70465SJoe Wallwork PetscFunctionBegin; 158bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 15931b70465SJoe Wallwork plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst; 16031b70465SJoe Wallwork PetscFunctionReturn(0); 16131b70465SJoe Wallwork } 16231b70465SJoe Wallwork 163cb7bfe8cSJoe Wallwork /*@ 16431b70465SJoe Wallwork DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after? 16531b70465SJoe Wallwork 16631b70465SJoe Wallwork Input parameters: 16731b70465SJoe Wallwork . dm - The DM 16831b70465SJoe Wallwork 16931b70465SJoe Wallwork Output parameters: 17031b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first? 17131b70465SJoe Wallwork 17231b70465SJoe Wallwork Level: beginner 17331b70465SJoe Wallwork 174db781477SPatrick Sanan .seealso: `DMPlexMetricIsIsotropic()`, `DMPlexMetricSetRestrictAnisotropyFirst()` 175cb7bfe8cSJoe Wallwork @*/ 176*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst) { 17731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 17831b70465SJoe Wallwork 17931b70465SJoe Wallwork PetscFunctionBegin; 180bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 18131b70465SJoe Wallwork *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst; 18231b70465SJoe Wallwork PetscFunctionReturn(0); 18331b70465SJoe Wallwork } 18431b70465SJoe Wallwork 185cb7bfe8cSJoe Wallwork /*@ 186cc2eee55SJoe Wallwork DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off? 187cc2eee55SJoe Wallwork 188cc2eee55SJoe Wallwork Input parameters: 189cc2eee55SJoe Wallwork + dm - The DM 190cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off? 191cc2eee55SJoe Wallwork 192cc2eee55SJoe Wallwork Level: beginner 193cc2eee55SJoe Wallwork 194cb7bfe8cSJoe Wallwork Notes: 195cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 196cb7bfe8cSJoe Wallwork 197db781477SPatrick Sanan .seealso: `DMPlexMetricNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()` 198cb7bfe8cSJoe Wallwork @*/ 199*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert) { 200cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 201cc2eee55SJoe Wallwork 202cc2eee55SJoe Wallwork PetscFunctionBegin; 203bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 204cc2eee55SJoe Wallwork plex->metricCtx->noInsert = noInsert; 205cc2eee55SJoe Wallwork PetscFunctionReturn(0); 206cc2eee55SJoe Wallwork } 207cc2eee55SJoe Wallwork 208cb7bfe8cSJoe Wallwork /*@ 209cc2eee55SJoe Wallwork DMPlexMetricNoInsertion - Are node insertion and deletion turned off? 210cc2eee55SJoe Wallwork 211cc2eee55SJoe Wallwork Input parameters: 212cc2eee55SJoe Wallwork . dm - The DM 213cc2eee55SJoe Wallwork 214cc2eee55SJoe Wallwork Output parameters: 215cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off? 216cc2eee55SJoe Wallwork 217cc2eee55SJoe Wallwork Level: beginner 218cc2eee55SJoe Wallwork 219cb7bfe8cSJoe Wallwork Notes: 220cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 221cb7bfe8cSJoe Wallwork 222db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()` 223cb7bfe8cSJoe Wallwork @*/ 224*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert) { 225cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 226cc2eee55SJoe Wallwork 227cc2eee55SJoe Wallwork PetscFunctionBegin; 228bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 229cc2eee55SJoe Wallwork *noInsert = plex->metricCtx->noInsert; 230cc2eee55SJoe Wallwork PetscFunctionReturn(0); 231cc2eee55SJoe Wallwork } 232cc2eee55SJoe Wallwork 233cb7bfe8cSJoe Wallwork /*@ 234cc2eee55SJoe Wallwork DMPlexMetricSetNoSwapping - Should facet swapping be turned off? 235cc2eee55SJoe Wallwork 236cc2eee55SJoe Wallwork Input parameters: 237cc2eee55SJoe Wallwork + dm - The DM 238cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off? 239cc2eee55SJoe Wallwork 240cc2eee55SJoe Wallwork Level: beginner 241cc2eee55SJoe Wallwork 242cb7bfe8cSJoe Wallwork Notes: 243cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 244cb7bfe8cSJoe Wallwork 245db781477SPatrick Sanan .seealso: `DMPlexMetricNoSwapping()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoSurf()` 246cb7bfe8cSJoe Wallwork @*/ 247*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap) { 248cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 249cc2eee55SJoe Wallwork 250cc2eee55SJoe Wallwork PetscFunctionBegin; 251bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 252cc2eee55SJoe Wallwork plex->metricCtx->noSwap = noSwap; 253cc2eee55SJoe Wallwork PetscFunctionReturn(0); 254cc2eee55SJoe Wallwork } 255cc2eee55SJoe Wallwork 256cb7bfe8cSJoe Wallwork /*@ 257cc2eee55SJoe Wallwork DMPlexMetricNoSwapping - Is facet swapping turned off? 258cc2eee55SJoe Wallwork 259cc2eee55SJoe Wallwork Input parameters: 260cc2eee55SJoe Wallwork . dm - The DM 261cc2eee55SJoe Wallwork 262cc2eee55SJoe Wallwork Output parameters: 263cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off? 264cc2eee55SJoe Wallwork 265cc2eee55SJoe Wallwork Level: beginner 266cc2eee55SJoe Wallwork 267cb7bfe8cSJoe Wallwork Notes: 268cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 269cb7bfe8cSJoe Wallwork 270db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoSwapping()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoSurf()` 271cb7bfe8cSJoe Wallwork @*/ 272*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap) { 273cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 274cc2eee55SJoe Wallwork 275cc2eee55SJoe Wallwork PetscFunctionBegin; 276bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 277cc2eee55SJoe Wallwork *noSwap = plex->metricCtx->noSwap; 278cc2eee55SJoe Wallwork PetscFunctionReturn(0); 279cc2eee55SJoe Wallwork } 280cc2eee55SJoe Wallwork 281cb7bfe8cSJoe Wallwork /*@ 282cc2eee55SJoe Wallwork DMPlexMetricSetNoMovement - Should node movement be turned off? 283cc2eee55SJoe Wallwork 284cc2eee55SJoe Wallwork Input parameters: 285cc2eee55SJoe Wallwork + dm - The DM 286cc2eee55SJoe Wallwork - noMove - Should node movement be turned off? 287cc2eee55SJoe Wallwork 288cc2eee55SJoe Wallwork Level: beginner 289cc2eee55SJoe Wallwork 290cb7bfe8cSJoe Wallwork Notes: 291cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 292cb7bfe8cSJoe Wallwork 293db781477SPatrick Sanan .seealso: `DMPlexMetricNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()`, `DMPlexMetricSetNoSurf()` 294cb7bfe8cSJoe Wallwork @*/ 295*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove) { 296cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 297cc2eee55SJoe Wallwork 298cc2eee55SJoe Wallwork PetscFunctionBegin; 299bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 300cc2eee55SJoe Wallwork plex->metricCtx->noMove = noMove; 301cc2eee55SJoe Wallwork PetscFunctionReturn(0); 302cc2eee55SJoe Wallwork } 303cc2eee55SJoe Wallwork 304cb7bfe8cSJoe Wallwork /*@ 305cc2eee55SJoe Wallwork DMPlexMetricNoMovement - Is node movement turned off? 306cc2eee55SJoe Wallwork 307cc2eee55SJoe Wallwork Input parameters: 308cc2eee55SJoe Wallwork . dm - The DM 309cc2eee55SJoe Wallwork 310cc2eee55SJoe Wallwork Output parameters: 311cc2eee55SJoe Wallwork . noMove - Is node movement turned off? 312cc2eee55SJoe Wallwork 313cc2eee55SJoe Wallwork Level: beginner 314cc2eee55SJoe Wallwork 315cb7bfe8cSJoe Wallwork Notes: 316cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 317cb7bfe8cSJoe Wallwork 318db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()`, `DMPlexMetricNoSurf()` 319cb7bfe8cSJoe Wallwork @*/ 320*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove) { 321cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 322cc2eee55SJoe Wallwork 323cc2eee55SJoe Wallwork PetscFunctionBegin; 324bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 325cc2eee55SJoe Wallwork *noMove = plex->metricCtx->noMove; 326cc2eee55SJoe Wallwork PetscFunctionReturn(0); 327cc2eee55SJoe Wallwork } 328cc2eee55SJoe Wallwork 329cb7bfe8cSJoe Wallwork /*@ 33076f360caSJoe Wallwork DMPlexMetricSetNoSurf - Should surface modification be turned off? 33176f360caSJoe Wallwork 33276f360caSJoe Wallwork Input parameters: 33376f360caSJoe Wallwork + dm - The DM 33476f360caSJoe Wallwork - noSurf - Should surface modification be turned off? 33576f360caSJoe Wallwork 33676f360caSJoe Wallwork Level: beginner 33776f360caSJoe Wallwork 33876f360caSJoe Wallwork Notes: 33976f360caSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 34076f360caSJoe Wallwork 341db781477SPatrick Sanan .seealso: `DMPlexMetricNoSurf()`, `DMPlexMetricSetNoMovement()`, `DMPlexMetricSetNoInsertion()`, `DMPlexMetricSetNoSwapping()` 34276f360caSJoe Wallwork @*/ 343*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf) { 34476f360caSJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 34576f360caSJoe Wallwork 34676f360caSJoe Wallwork PetscFunctionBegin; 347bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 34876f360caSJoe Wallwork plex->metricCtx->noSurf = noSurf; 34976f360caSJoe Wallwork PetscFunctionReturn(0); 35076f360caSJoe Wallwork } 35176f360caSJoe Wallwork 35276f360caSJoe Wallwork /*@ 35376f360caSJoe Wallwork DMPlexMetricNoSurf - Is surface modification turned off? 35476f360caSJoe Wallwork 35576f360caSJoe Wallwork Input parameters: 35676f360caSJoe Wallwork . dm - The DM 35776f360caSJoe Wallwork 35876f360caSJoe Wallwork Output parameters: 35976f360caSJoe Wallwork . noSurf - Is surface modification turned off? 36076f360caSJoe Wallwork 36176f360caSJoe Wallwork Level: beginner 36276f360caSJoe Wallwork 36376f360caSJoe Wallwork Notes: 36476f360caSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 36576f360caSJoe Wallwork 366db781477SPatrick Sanan .seealso: `DMPlexMetricSetNoSurf()`, `DMPlexMetricNoMovement()`, `DMPlexMetricNoInsertion()`, `DMPlexMetricNoSwapping()` 36776f360caSJoe Wallwork @*/ 368*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf) { 36976f360caSJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 37076f360caSJoe Wallwork 37176f360caSJoe Wallwork PetscFunctionBegin; 372bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 37376f360caSJoe Wallwork *noSurf = plex->metricCtx->noSurf; 37476f360caSJoe Wallwork PetscFunctionReturn(0); 37576f360caSJoe Wallwork } 37676f360caSJoe Wallwork 37776f360caSJoe Wallwork /*@ 37831b70465SJoe Wallwork DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude 37931b70465SJoe Wallwork 38031b70465SJoe Wallwork Input parameters: 38131b70465SJoe Wallwork + dm - The DM 38231b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude 38331b70465SJoe Wallwork 38431b70465SJoe Wallwork Level: beginner 38531b70465SJoe Wallwork 386db781477SPatrick Sanan .seealso: `DMPlexMetricGetMinimumMagnitude()`, `DMPlexMetricSetMaximumMagnitude()` 387cb7bfe8cSJoe Wallwork @*/ 388*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min) { 38931b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 39031b70465SJoe Wallwork 39131b70465SJoe Wallwork PetscFunctionBegin; 392bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 393bc00d7afSJoe Wallwork PetscCheck(h_min > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)"); 39431b70465SJoe Wallwork plex->metricCtx->h_min = h_min; 39531b70465SJoe Wallwork PetscFunctionReturn(0); 39631b70465SJoe Wallwork } 39731b70465SJoe Wallwork 398cb7bfe8cSJoe Wallwork /*@ 39931b70465SJoe Wallwork DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude 40031b70465SJoe Wallwork 40131b70465SJoe Wallwork Input parameters: 40231b70465SJoe Wallwork . dm - The DM 40331b70465SJoe Wallwork 404cc2eee55SJoe Wallwork Output parameters: 40531b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude 40631b70465SJoe Wallwork 40731b70465SJoe Wallwork Level: beginner 40831b70465SJoe Wallwork 409db781477SPatrick Sanan .seealso: `DMPlexMetricSetMinimumMagnitude()`, `DMPlexMetricGetMaximumMagnitude()` 410cb7bfe8cSJoe Wallwork @*/ 411*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min) { 41231b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 41331b70465SJoe Wallwork 41431b70465SJoe Wallwork PetscFunctionBegin; 415bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 41631b70465SJoe Wallwork *h_min = plex->metricCtx->h_min; 41731b70465SJoe Wallwork PetscFunctionReturn(0); 41831b70465SJoe Wallwork } 41931b70465SJoe Wallwork 420cb7bfe8cSJoe Wallwork /*@ 42131b70465SJoe Wallwork DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude 42231b70465SJoe Wallwork 42331b70465SJoe Wallwork Input parameters: 42431b70465SJoe Wallwork + dm - The DM 42531b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude 42631b70465SJoe Wallwork 42731b70465SJoe Wallwork Level: beginner 42831b70465SJoe Wallwork 429db781477SPatrick Sanan .seealso: `DMPlexMetricGetMaximumMagnitude()`, `DMPlexMetricSetMinimumMagnitude()` 430cb7bfe8cSJoe Wallwork @*/ 431*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max) { 43231b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 43331b70465SJoe Wallwork 43431b70465SJoe Wallwork PetscFunctionBegin; 435bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 436bc00d7afSJoe Wallwork PetscCheck(h_max > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric magnitudes must be in (0, inf)"); 43731b70465SJoe Wallwork plex->metricCtx->h_max = h_max; 43831b70465SJoe Wallwork PetscFunctionReturn(0); 43931b70465SJoe Wallwork } 44031b70465SJoe Wallwork 441cb7bfe8cSJoe Wallwork /*@ 44231b70465SJoe Wallwork DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude 44331b70465SJoe Wallwork 44431b70465SJoe Wallwork Input parameters: 44531b70465SJoe Wallwork . dm - The DM 44631b70465SJoe Wallwork 447cc2eee55SJoe Wallwork Output parameters: 44831b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude 44931b70465SJoe Wallwork 45031b70465SJoe Wallwork Level: beginner 45131b70465SJoe Wallwork 452db781477SPatrick Sanan .seealso: `DMPlexMetricSetMaximumMagnitude()`, `DMPlexMetricGetMinimumMagnitude()` 453cb7bfe8cSJoe Wallwork @*/ 454*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max) { 45531b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 45631b70465SJoe Wallwork 45731b70465SJoe Wallwork PetscFunctionBegin; 458bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 45931b70465SJoe Wallwork *h_max = plex->metricCtx->h_max; 46031b70465SJoe Wallwork PetscFunctionReturn(0); 46131b70465SJoe Wallwork } 46231b70465SJoe Wallwork 463cb7bfe8cSJoe Wallwork /*@ 46431b70465SJoe Wallwork DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy 46531b70465SJoe Wallwork 46631b70465SJoe Wallwork Input parameters: 46731b70465SJoe Wallwork + dm - The DM 46831b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy 46931b70465SJoe Wallwork 47031b70465SJoe Wallwork Level: beginner 47131b70465SJoe Wallwork 47231b70465SJoe Wallwork Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one. 47331b70465SJoe Wallwork 474db781477SPatrick Sanan .seealso: `DMPlexMetricGetMaximumAnisotropy()`, `DMPlexMetricSetMaximumMagnitude()` 475cb7bfe8cSJoe Wallwork @*/ 476*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max) { 47731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 47831b70465SJoe Wallwork 47931b70465SJoe Wallwork PetscFunctionBegin; 480bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 481bc00d7afSJoe Wallwork PetscCheck(a_max >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Anisotropy must be in [1, inf)"); 48231b70465SJoe Wallwork plex->metricCtx->a_max = a_max; 48331b70465SJoe Wallwork PetscFunctionReturn(0); 48431b70465SJoe Wallwork } 48531b70465SJoe Wallwork 486cb7bfe8cSJoe Wallwork /*@ 48731b70465SJoe Wallwork DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy 48831b70465SJoe Wallwork 48931b70465SJoe Wallwork Input parameters: 49031b70465SJoe Wallwork . dm - The DM 49131b70465SJoe Wallwork 492cc2eee55SJoe Wallwork Output parameters: 49331b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy 49431b70465SJoe Wallwork 49531b70465SJoe Wallwork Level: beginner 49631b70465SJoe Wallwork 497db781477SPatrick Sanan .seealso: `DMPlexMetricSetMaximumAnisotropy()`, `DMPlexMetricGetMaximumMagnitude()` 498cb7bfe8cSJoe Wallwork @*/ 499*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max) { 50031b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 50131b70465SJoe Wallwork 50231b70465SJoe Wallwork PetscFunctionBegin; 503bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 50431b70465SJoe Wallwork *a_max = plex->metricCtx->a_max; 50531b70465SJoe Wallwork PetscFunctionReturn(0); 50631b70465SJoe Wallwork } 50731b70465SJoe Wallwork 508cb7bfe8cSJoe Wallwork /*@ 50931b70465SJoe Wallwork DMPlexMetricSetTargetComplexity - Set the target metric complexity 51031b70465SJoe Wallwork 51131b70465SJoe Wallwork Input parameters: 51231b70465SJoe Wallwork + dm - The DM 51331b70465SJoe Wallwork - targetComplexity - The target metric complexity 51431b70465SJoe Wallwork 51531b70465SJoe Wallwork Level: beginner 51631b70465SJoe Wallwork 517db781477SPatrick Sanan .seealso: `DMPlexMetricGetTargetComplexity()`, `DMPlexMetricSetNormalizationOrder()` 518cb7bfe8cSJoe Wallwork @*/ 519*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity) { 52031b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 52131b70465SJoe Wallwork 52231b70465SJoe Wallwork PetscFunctionBegin; 523bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 524bc00d7afSJoe Wallwork PetscCheck(targetComplexity > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric complexity must be in (0, inf)"); 52531b70465SJoe Wallwork plex->metricCtx->targetComplexity = targetComplexity; 52631b70465SJoe Wallwork PetscFunctionReturn(0); 52731b70465SJoe Wallwork } 52831b70465SJoe Wallwork 529cb7bfe8cSJoe Wallwork /*@ 53031b70465SJoe Wallwork DMPlexMetricGetTargetComplexity - Get the target metric complexity 53131b70465SJoe Wallwork 53231b70465SJoe Wallwork Input parameters: 53331b70465SJoe Wallwork . dm - The DM 53431b70465SJoe Wallwork 535cc2eee55SJoe Wallwork Output parameters: 53631b70465SJoe Wallwork . targetComplexity - The target metric complexity 53731b70465SJoe Wallwork 53831b70465SJoe Wallwork Level: beginner 53931b70465SJoe Wallwork 540db781477SPatrick Sanan .seealso: `DMPlexMetricSetTargetComplexity()`, `DMPlexMetricGetNormalizationOrder()` 541cb7bfe8cSJoe Wallwork @*/ 542*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity) { 54331b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 54431b70465SJoe Wallwork 54531b70465SJoe Wallwork PetscFunctionBegin; 546bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 54731b70465SJoe Wallwork *targetComplexity = plex->metricCtx->targetComplexity; 54831b70465SJoe Wallwork PetscFunctionReturn(0); 54931b70465SJoe Wallwork } 55031b70465SJoe Wallwork 551cb7bfe8cSJoe Wallwork /*@ 55231b70465SJoe Wallwork DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization 55331b70465SJoe Wallwork 55431b70465SJoe Wallwork Input parameters: 55531b70465SJoe Wallwork + dm - The DM 55631b70465SJoe Wallwork - p - The normalization order 55731b70465SJoe Wallwork 55831b70465SJoe Wallwork Level: beginner 55931b70465SJoe Wallwork 560db781477SPatrick Sanan .seealso: `DMPlexMetricGetNormalizationOrder()`, `DMPlexMetricSetTargetComplexity()` 561cb7bfe8cSJoe Wallwork @*/ 562*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p) { 56331b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 56431b70465SJoe Wallwork 56531b70465SJoe Wallwork PetscFunctionBegin; 566bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 567bc00d7afSJoe Wallwork PetscCheck(p >= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Normalization order must be in [1, inf)"); 56831b70465SJoe Wallwork plex->metricCtx->p = p; 56931b70465SJoe Wallwork PetscFunctionReturn(0); 57031b70465SJoe Wallwork } 57131b70465SJoe Wallwork 572cb7bfe8cSJoe Wallwork /*@ 57331b70465SJoe Wallwork DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization 57431b70465SJoe Wallwork 57531b70465SJoe Wallwork Input parameters: 57631b70465SJoe Wallwork . dm - The DM 57731b70465SJoe Wallwork 578cc2eee55SJoe Wallwork Output parameters: 57931b70465SJoe Wallwork . p - The normalization order 58031b70465SJoe Wallwork 58131b70465SJoe Wallwork Level: beginner 58231b70465SJoe Wallwork 583db781477SPatrick Sanan .seealso: `DMPlexMetricSetNormalizationOrder()`, `DMPlexMetricGetTargetComplexity()` 584cb7bfe8cSJoe Wallwork @*/ 585*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p) { 58631b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 58731b70465SJoe Wallwork 58831b70465SJoe Wallwork PetscFunctionBegin; 589bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 59031b70465SJoe Wallwork *p = plex->metricCtx->p; 59131b70465SJoe Wallwork PetscFunctionReturn(0); 59231b70465SJoe Wallwork } 59331b70465SJoe Wallwork 594cb7bfe8cSJoe Wallwork /*@ 595cc2eee55SJoe Wallwork DMPlexMetricSetGradationFactor - Set the metric gradation factor 596cc2eee55SJoe Wallwork 597cc2eee55SJoe Wallwork Input parameters: 598cc2eee55SJoe Wallwork + dm - The DM 599cc2eee55SJoe Wallwork - beta - The metric gradation factor 600cc2eee55SJoe Wallwork 601cc2eee55SJoe Wallwork Level: beginner 602cc2eee55SJoe Wallwork 603cc2eee55SJoe Wallwork Notes: 604cc2eee55SJoe Wallwork 605cc2eee55SJoe Wallwork The gradation factor is the maximum tolerated length ratio between adjacent edges. 606cc2eee55SJoe Wallwork 607cc2eee55SJoe Wallwork Turn off gradation by passing the value -1. Otherwise, pass a positive value. 608cc2eee55SJoe Wallwork 609cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 610cb7bfe8cSJoe Wallwork 611db781477SPatrick Sanan .seealso: `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()` 612cb7bfe8cSJoe Wallwork @*/ 613*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta) { 614cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 615cc2eee55SJoe Wallwork 616cc2eee55SJoe Wallwork PetscFunctionBegin; 617bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 618bc00d7afSJoe Wallwork PetscCheck(beta > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric gradation factor must be in (0, inf)"); 619cc2eee55SJoe Wallwork plex->metricCtx->gradationFactor = beta; 620cc2eee55SJoe Wallwork PetscFunctionReturn(0); 621cc2eee55SJoe Wallwork } 622cc2eee55SJoe Wallwork 623cb7bfe8cSJoe Wallwork /*@ 624cc2eee55SJoe Wallwork DMPlexMetricGetGradationFactor - Get the metric gradation factor 625cc2eee55SJoe Wallwork 626cc2eee55SJoe Wallwork Input parameters: 627cc2eee55SJoe Wallwork . dm - The DM 628cc2eee55SJoe Wallwork 629cc2eee55SJoe Wallwork Output parameters: 630cc2eee55SJoe Wallwork . beta - The metric gradation factor 631cc2eee55SJoe Wallwork 632cc2eee55SJoe Wallwork Level: beginner 633cc2eee55SJoe Wallwork 634cb7bfe8cSJoe Wallwork Notes: 635cb7bfe8cSJoe Wallwork 636cb7bfe8cSJoe Wallwork The gradation factor is the maximum tolerated length ratio between adjacent edges. 637cb7bfe8cSJoe Wallwork 638cb7bfe8cSJoe Wallwork The value -1 implies that gradation is turned off. 639cb7bfe8cSJoe Wallwork 640cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 641cc2eee55SJoe Wallwork 642db781477SPatrick Sanan .seealso: `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()` 643cb7bfe8cSJoe Wallwork @*/ 644*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta) { 645cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 646cc2eee55SJoe Wallwork 647cc2eee55SJoe Wallwork PetscFunctionBegin; 648bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 649cc2eee55SJoe Wallwork *beta = plex->metricCtx->gradationFactor; 650cc2eee55SJoe Wallwork PetscFunctionReturn(0); 651cc2eee55SJoe Wallwork } 652cc2eee55SJoe Wallwork 653cb7bfe8cSJoe Wallwork /*@ 654ae8b063eSJoe Wallwork DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number 655ae8b063eSJoe Wallwork 656ae8b063eSJoe Wallwork Input parameters: 657ae8b063eSJoe Wallwork + dm - The DM 658ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number 659ae8b063eSJoe Wallwork 660ae8b063eSJoe Wallwork Level: beginner 661ae8b063eSJoe Wallwork 662ae8b063eSJoe Wallwork Notes: 663ae8b063eSJoe Wallwork 664ae8b063eSJoe Wallwork The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the 665ae8b063eSJoe Wallwork boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the 666ae8b063eSJoe Wallwork high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for 667ae8b063eSJoe Wallwork an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease 668ae8b063eSJoe Wallwork (resp. increase) the Hausdorff parameter. (Taken from 669ae8b063eSJoe Wallwork https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd). 670ae8b063eSJoe Wallwork 671ae8b063eSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 672ae8b063eSJoe Wallwork 673db781477SPatrick Sanan .seealso: `DMPlexMetricSetGradationFactor()`, `DMPlexMetricGetHausdorffNumber()` 674ae8b063eSJoe Wallwork @*/ 675*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd) { 676ae8b063eSJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 677ae8b063eSJoe Wallwork 678ae8b063eSJoe Wallwork PetscFunctionBegin; 679bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 680bc00d7afSJoe Wallwork PetscCheck(hausd > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Metric Hausdorff number must be in (0, inf)"); 681ae8b063eSJoe Wallwork plex->metricCtx->hausdorffNumber = hausd; 682ae8b063eSJoe Wallwork PetscFunctionReturn(0); 683ae8b063eSJoe Wallwork } 684ae8b063eSJoe Wallwork 685ae8b063eSJoe Wallwork /*@ 686ae8b063eSJoe Wallwork DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number 687ae8b063eSJoe Wallwork 688ae8b063eSJoe Wallwork Input parameters: 689ae8b063eSJoe Wallwork . dm - The DM 690ae8b063eSJoe Wallwork 691ae8b063eSJoe Wallwork Output parameters: 692ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number 693ae8b063eSJoe Wallwork 694ae8b063eSJoe Wallwork Level: beginner 695ae8b063eSJoe Wallwork 696ae8b063eSJoe Wallwork Notes: 697ae8b063eSJoe Wallwork 698ae8b063eSJoe Wallwork The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the 699ae8b063eSJoe Wallwork boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the 700ae8b063eSJoe Wallwork high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for 701ae8b063eSJoe Wallwork an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease 702ae8b063eSJoe Wallwork (resp. increase) the Hausdorff parameter. (Taken from 703ae8b063eSJoe Wallwork https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd). 704ae8b063eSJoe Wallwork 705ae8b063eSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 706ae8b063eSJoe Wallwork 707db781477SPatrick Sanan .seealso: `DMPlexMetricGetGradationFactor()`, `DMPlexMetricSetHausdorffNumber()` 708ae8b063eSJoe Wallwork @*/ 709*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd) { 710ae8b063eSJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 711ae8b063eSJoe Wallwork 712ae8b063eSJoe Wallwork PetscFunctionBegin; 713bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 714ae8b063eSJoe Wallwork *hausd = plex->metricCtx->hausdorffNumber; 715ae8b063eSJoe Wallwork PetscFunctionReturn(0); 716ae8b063eSJoe Wallwork } 717ae8b063eSJoe Wallwork 718ae8b063eSJoe Wallwork /*@ 719cc2eee55SJoe Wallwork DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package 720cc2eee55SJoe Wallwork 721cc2eee55SJoe Wallwork Input parameters: 722cc2eee55SJoe Wallwork + dm - The DM 723cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum 724cc2eee55SJoe Wallwork 725cb7bfe8cSJoe Wallwork Level: beginner 726cb7bfe8cSJoe Wallwork 727cb7bfe8cSJoe Wallwork Notes: 728cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 729cb7bfe8cSJoe Wallwork 730db781477SPatrick Sanan .seealso: `DMPlexMetricGetVerbosity()`, `DMPlexMetricSetNumIterations()` 731cb7bfe8cSJoe Wallwork @*/ 732*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity) { 733cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 734cc2eee55SJoe Wallwork 735cc2eee55SJoe Wallwork PetscFunctionBegin; 736bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 737cc2eee55SJoe Wallwork plex->metricCtx->verbosity = verbosity; 738cc2eee55SJoe Wallwork PetscFunctionReturn(0); 739cc2eee55SJoe Wallwork } 740cc2eee55SJoe Wallwork 741cb7bfe8cSJoe Wallwork /*@ 742cc2eee55SJoe Wallwork DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package 743cc2eee55SJoe Wallwork 744cc2eee55SJoe Wallwork Input parameters: 745cc2eee55SJoe Wallwork . dm - The DM 746cc2eee55SJoe Wallwork 747cc2eee55SJoe Wallwork Output parameters: 748cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum 749cc2eee55SJoe Wallwork 750cb7bfe8cSJoe Wallwork Level: beginner 751cb7bfe8cSJoe Wallwork 752cb7bfe8cSJoe Wallwork Notes: 753cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 754cb7bfe8cSJoe Wallwork 755db781477SPatrick Sanan .seealso: `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()` 756cb7bfe8cSJoe Wallwork @*/ 757*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity) { 758cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 759cc2eee55SJoe Wallwork 760cc2eee55SJoe Wallwork PetscFunctionBegin; 761bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 762cc2eee55SJoe Wallwork *verbosity = plex->metricCtx->verbosity; 763cc2eee55SJoe Wallwork PetscFunctionReturn(0); 764cc2eee55SJoe Wallwork } 765cc2eee55SJoe Wallwork 766cb7bfe8cSJoe Wallwork /*@ 767cc2eee55SJoe Wallwork DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations 768cc2eee55SJoe Wallwork 769cc2eee55SJoe Wallwork Input parameters: 770cc2eee55SJoe Wallwork + dm - The DM 771cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations 772cc2eee55SJoe Wallwork 773cb7bfe8cSJoe Wallwork Level: beginner 774cb7bfe8cSJoe Wallwork 775cb7bfe8cSJoe Wallwork Notes: 776cb7bfe8cSJoe Wallwork This is only used by ParMmg (not Pragmatic or Mmg). 777cc2eee55SJoe Wallwork 778db781477SPatrick Sanan .seealso: `DMPlexMetricSetVerbosity()`, `DMPlexMetricGetNumIterations()` 779cb7bfe8cSJoe Wallwork @*/ 780*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter) { 781cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 782cc2eee55SJoe Wallwork 783cc2eee55SJoe Wallwork PetscFunctionBegin; 784bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 785cc2eee55SJoe Wallwork plex->metricCtx->numIter = numIter; 786cc2eee55SJoe Wallwork PetscFunctionReturn(0); 787cc2eee55SJoe Wallwork } 788cc2eee55SJoe Wallwork 789cb7bfe8cSJoe Wallwork /*@ 790cc2eee55SJoe Wallwork DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations 791cc2eee55SJoe Wallwork 792cc2eee55SJoe Wallwork Input parameters: 793cc2eee55SJoe Wallwork . dm - The DM 794cc2eee55SJoe Wallwork 795cc2eee55SJoe Wallwork Output parameters: 796cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations 797cc2eee55SJoe Wallwork 798cb7bfe8cSJoe Wallwork Level: beginner 799cb7bfe8cSJoe Wallwork 800cb7bfe8cSJoe Wallwork Notes: 801cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic or Mmg). 802cc2eee55SJoe Wallwork 803db781477SPatrick Sanan .seealso: `DMPlexMetricSetNumIterations()`, `DMPlexMetricGetVerbosity()` 804cb7bfe8cSJoe Wallwork @*/ 805*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter) { 806cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *)dm->data; 807cc2eee55SJoe Wallwork 808cc2eee55SJoe Wallwork PetscFunctionBegin; 809bc00d7afSJoe Wallwork if (!plex->metricCtx) PetscCall(DMPlexMetricSetFromOptions(dm)); 810cc2eee55SJoe Wallwork *numIter = plex->metricCtx->numIter; 811cc2eee55SJoe Wallwork PetscFunctionReturn(0); 812cc2eee55SJoe Wallwork } 813cc2eee55SJoe Wallwork 814*9371c9d4SSatish Balay PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric) { 81520282da2SJoe Wallwork MPI_Comm comm; 81620282da2SJoe Wallwork PetscFE fe; 81720282da2SJoe Wallwork PetscInt dim; 81820282da2SJoe Wallwork 81920282da2SJoe Wallwork PetscFunctionBegin; 82020282da2SJoe Wallwork 82120282da2SJoe Wallwork /* Extract metadata from dm */ 8229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 8239566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 82420282da2SJoe Wallwork 82520282da2SJoe Wallwork /* Create a P1 field of the requested size */ 8269566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 8279566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe)); 8289566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 8299566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 8309566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, metric)); 83120282da2SJoe Wallwork 83220282da2SJoe Wallwork PetscFunctionReturn(0); 83320282da2SJoe Wallwork } 83420282da2SJoe Wallwork 835cb7bfe8cSJoe Wallwork /*@ 83620282da2SJoe Wallwork DMPlexMetricCreate - Create a Riemannian metric field 83720282da2SJoe Wallwork 83820282da2SJoe Wallwork Input parameters: 83920282da2SJoe Wallwork + dm - The DM 84020282da2SJoe Wallwork - f - The field number to use 84120282da2SJoe Wallwork 84220282da2SJoe Wallwork Output parameter: 84320282da2SJoe Wallwork . metric - The metric 84420282da2SJoe Wallwork 84520282da2SJoe Wallwork Level: beginner 84620282da2SJoe Wallwork 84731b70465SJoe Wallwork Notes: 84831b70465SJoe Wallwork 84931b70465SJoe Wallwork It is assumed that the DM is comprised of simplices. 85031b70465SJoe Wallwork 85131b70465SJoe Wallwork Command line options for Riemannian metrics: 85231b70465SJoe Wallwork 853cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic? 85493520041SJoe Wallwork . -dm_plex_metric_uniform - Is the metric uniform? 855cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 856cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 857cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 858cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max - Maximum tolerated anisotropy 859cb7bfe8cSJoe Wallwork . -dm_plex_metric_p - L-p normalization order 86067b8a455SSatish Balay - -dm_plex_metric_target_complexity - Target metric complexity 861cb7bfe8cSJoe Wallwork 862cb7bfe8cSJoe Wallwork Switching between remeshers can be achieved using 863cb7bfe8cSJoe Wallwork 86467b8a455SSatish Balay . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use 865cb7bfe8cSJoe Wallwork 866cb7bfe8cSJoe Wallwork Further options that are only relevant to Mmg and ParMmg: 867cb7bfe8cSJoe Wallwork 868cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation 869cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations - Number of parallel mesh adaptation iterations for ParMmg 870cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert - Should node insertion/deletion be turned off? 871cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap - Should facet swapping be turned off? 872cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move - Should node movement be turned off? 873cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity - Choose a verbosity level from -1 (silent) to 10 (maximum). 87420282da2SJoe Wallwork 875db781477SPatrick Sanan .seealso: `DMPlexMetricCreateUniform()`, `DMPlexMetricCreateIsotropic()` 876cb7bfe8cSJoe Wallwork @*/ 877*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric) { 87893520041SJoe Wallwork PetscBool isotropic, uniform; 87920282da2SJoe Wallwork PetscInt coordDim, Nd; 88020282da2SJoe Wallwork 88120282da2SJoe Wallwork PetscFunctionBegin; 8829566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &coordDim)); 88320282da2SJoe Wallwork Nd = coordDim * coordDim; 8849566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 8859566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 88693520041SJoe Wallwork if (uniform) { 88793520041SJoe Wallwork MPI_Comm comm; 88893520041SJoe Wallwork 8899566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 890bc00d7afSJoe Wallwork PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported"); 8919566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, metric)); 8929566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE)); 8939566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(*metric)); 894f6db075bSJoe Wallwork } else if (isotropic) PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric)); 895f6db075bSJoe Wallwork else PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric)); 89620282da2SJoe Wallwork PetscFunctionReturn(0); 89720282da2SJoe Wallwork } 89820282da2SJoe Wallwork 899cb7bfe8cSJoe Wallwork /*@ 90020282da2SJoe Wallwork DMPlexMetricCreateUniform - Construct a uniform isotropic metric 90120282da2SJoe Wallwork 90220282da2SJoe Wallwork Input parameters: 90320282da2SJoe Wallwork + dm - The DM 90420282da2SJoe Wallwork . f - The field number to use 90520282da2SJoe Wallwork - alpha - Scaling parameter for the diagonal 90620282da2SJoe Wallwork 90720282da2SJoe Wallwork Output parameter: 90820282da2SJoe Wallwork . metric - The uniform metric 90920282da2SJoe Wallwork 91020282da2SJoe Wallwork Level: beginner 91120282da2SJoe Wallwork 91293520041SJoe Wallwork Note: In this case, the metric is constant in space and so is only specified for a single vertex. 91320282da2SJoe Wallwork 914db781477SPatrick Sanan .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateIsotropic()` 915cb7bfe8cSJoe Wallwork @*/ 916*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric) { 91720282da2SJoe Wallwork PetscFunctionBegin; 9189566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE)); 9199566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, f, metric)); 9202c71b3e2SJacob Faibussowitsch PetscCheck(alpha, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined"); 921bc00d7afSJoe Wallwork PetscCheck(alpha > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Uniform metric scaling should be in (0, inf)"); 9229566063dSJacob Faibussowitsch PetscCall(VecSet(*metric, alpha)); 9239566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(*metric)); 9249566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(*metric)); 92520282da2SJoe Wallwork PetscFunctionReturn(0); 92620282da2SJoe Wallwork } 92720282da2SJoe Wallwork 928*9371c9d4SSatish Balay 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[]) { 92993520041SJoe Wallwork f0[0] = u[0]; 93093520041SJoe Wallwork } 93193520041SJoe Wallwork 932cb7bfe8cSJoe Wallwork /*@ 93320282da2SJoe Wallwork DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator 93420282da2SJoe Wallwork 93520282da2SJoe Wallwork Input parameters: 93620282da2SJoe Wallwork + dm - The DM 93720282da2SJoe Wallwork . f - The field number to use 93820282da2SJoe Wallwork - indicator - The error indicator 93920282da2SJoe Wallwork 94020282da2SJoe Wallwork Output parameter: 94120282da2SJoe Wallwork . metric - The isotropic metric 94220282da2SJoe Wallwork 94320282da2SJoe Wallwork Level: beginner 94420282da2SJoe Wallwork 94520282da2SJoe Wallwork Notes: 94620282da2SJoe Wallwork 94720282da2SJoe Wallwork It is assumed that the DM is comprised of simplices. 94820282da2SJoe Wallwork 94993520041SJoe Wallwork The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately. 95020282da2SJoe Wallwork 951db781477SPatrick Sanan .seealso: `DMPlexMetricCreate()`, `DMPlexMetricCreateUniform()` 952cb7bfe8cSJoe Wallwork @*/ 953*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric) { 95493520041SJoe Wallwork PetscInt m, n; 95520282da2SJoe Wallwork 95620282da2SJoe Wallwork PetscFunctionBegin; 9579566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE)); 9589566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, f, metric)); 9599566063dSJacob Faibussowitsch PetscCall(VecGetSize(indicator, &m)); 9609566063dSJacob Faibussowitsch PetscCall(VecGetSize(*metric, &n)); 9619566063dSJacob Faibussowitsch if (m == n) PetscCall(VecCopy(indicator, *metric)); 96293520041SJoe Wallwork else { 96393520041SJoe Wallwork DM dmIndi; 964*9371c9d4SSatish Balay void (*funcs[1])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 96593520041SJoe Wallwork 9669566063dSJacob Faibussowitsch PetscCall(VecGetDM(indicator, &dmIndi)); 96793520041SJoe Wallwork funcs[0] = identity; 9689566063dSJacob Faibussowitsch PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric)); 96920282da2SJoe Wallwork } 97020282da2SJoe Wallwork PetscFunctionReturn(0); 97120282da2SJoe Wallwork } 97220282da2SJoe Wallwork 973f6db075bSJoe Wallwork /*@ 974f6db075bSJoe Wallwork DMPlexMetricDeterminantCreate - Create the determinant field for a Riemannian metric 975f6db075bSJoe Wallwork 976f6db075bSJoe Wallwork Input parameters: 977f6db075bSJoe Wallwork + dm - The DM of the metric field 978f6db075bSJoe Wallwork - f - The field number to use 979f6db075bSJoe Wallwork 980f6db075bSJoe Wallwork Output parameter: 981f6db075bSJoe Wallwork + determinant - The determinant field 982f6db075bSJoe Wallwork - dmDet - The corresponding DM 983f6db075bSJoe Wallwork 984f6db075bSJoe Wallwork Level: beginner 985f6db075bSJoe Wallwork 986f6db075bSJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic(), DMPlexMetricCreate() 987f6db075bSJoe Wallwork @*/ 988*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricDeterminantCreate(DM dm, PetscInt f, Vec *determinant, DM *dmDet) { 989f6db075bSJoe Wallwork PetscBool uniform; 990f6db075bSJoe Wallwork 991f6db075bSJoe Wallwork PetscFunctionBegin; 992f6db075bSJoe Wallwork PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 993f6db075bSJoe Wallwork PetscCall(DMClone(dm, dmDet)); 994f6db075bSJoe Wallwork if (uniform) { 995f6db075bSJoe Wallwork MPI_Comm comm; 996f6db075bSJoe Wallwork 997f6db075bSJoe Wallwork PetscCall(PetscObjectGetComm((PetscObject)*dmDet, &comm)); 998f6db075bSJoe Wallwork PetscCall(VecCreate(comm, determinant)); 999f6db075bSJoe Wallwork PetscCall(VecSetSizes(*determinant, 1, PETSC_DECIDE)); 1000f6db075bSJoe Wallwork PetscCall(VecSetFromOptions(*determinant)); 1001f6db075bSJoe Wallwork } else PetscCall(DMPlexP1FieldCreate_Private(*dmDet, f, 1, determinant)); 1002f6db075bSJoe Wallwork PetscFunctionReturn(0); 1003f6db075bSJoe Wallwork } 1004f6db075bSJoe Wallwork 1005*9371c9d4SSatish Balay static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[]) { 100682490253SJoe Wallwork PetscInt i, j; 100782490253SJoe Wallwork 100882490253SJoe Wallwork PetscFunctionBegin; 100982490253SJoe Wallwork PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"); 101082490253SJoe Wallwork for (i = 0; i < dim; ++i) { 101182490253SJoe Wallwork if (i == 0) PetscPrintf(PETSC_COMM_SELF, " [["); 101282490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, " ["); 101382490253SJoe Wallwork for (j = 0; j < dim; ++j) { 101463a3b9bcSJacob Faibussowitsch if (j < dim - 1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", (double)PetscAbsScalar(Mpos[i * dim + j])); 101563a3b9bcSJacob Faibussowitsch else PetscPrintf(PETSC_COMM_SELF, "%15.8e", (double)PetscAbsScalar(Mpos[i * dim + j])); 101682490253SJoe Wallwork } 101782490253SJoe Wallwork if (i < dim - 1) PetscPrintf(PETSC_COMM_SELF, "]\n"); 101882490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, "]]\n"); 101982490253SJoe Wallwork } 102082490253SJoe Wallwork PetscFunctionReturn(0); 102182490253SJoe Wallwork } 102282490253SJoe Wallwork 1023*9371c9d4SSatish Balay static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp) { 102420282da2SJoe Wallwork PetscInt i, j, k; 102520282da2SJoe Wallwork PetscReal *eigs, max_eig, l_min = 1.0 / (h_max * h_max), l_max = 1.0 / (h_min * h_min), la_min = 1.0 / (a_max * a_max); 102620282da2SJoe Wallwork PetscScalar *Mpos; 102720282da2SJoe Wallwork 102820282da2SJoe Wallwork PetscFunctionBegin; 10299566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dim * dim, &Mpos, dim, &eigs)); 103020282da2SJoe Wallwork 103120282da2SJoe Wallwork /* Symmetrize */ 103220282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 103320282da2SJoe Wallwork Mpos[i * dim + i] = Mp[i * dim + i]; 103420282da2SJoe Wallwork for (j = i + 1; j < dim; ++j) { 103520282da2SJoe Wallwork Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]); 103620282da2SJoe Wallwork Mpos[j * dim + i] = Mpos[i * dim + j]; 103720282da2SJoe Wallwork } 103820282da2SJoe Wallwork } 103920282da2SJoe Wallwork 104020282da2SJoe Wallwork /* Compute eigendecomposition */ 104193520041SJoe Wallwork if (dim == 1) { 104293520041SJoe Wallwork /* Isotropic case */ 104393520041SJoe Wallwork eigs[0] = PetscRealPart(Mpos[0]); 104493520041SJoe Wallwork Mpos[0] = 1.0; 104593520041SJoe Wallwork } else { 104693520041SJoe Wallwork /* Anisotropic case */ 104720282da2SJoe Wallwork PetscScalar *work; 104820282da2SJoe Wallwork PetscBLASInt lwork; 104920282da2SJoe Wallwork 105020282da2SJoe Wallwork lwork = 5 * dim; 10519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * dim, &work)); 105220282da2SJoe Wallwork { 105320282da2SJoe Wallwork PetscBLASInt lierr; 105420282da2SJoe Wallwork PetscBLASInt nb; 105520282da2SJoe Wallwork 10569566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(dim, &nb)); 10579566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 105820282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 105920282da2SJoe Wallwork { 106020282da2SJoe Wallwork PetscReal *rwork; 10619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3 * dim, &rwork)); 1062792fecdfSBarry Smith PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, rwork, &lierr)); 10639566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 106420282da2SJoe Wallwork } 106520282da2SJoe Wallwork #else 1066792fecdfSBarry Smith PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, Mpos, &nb, eigs, work, &lwork, &lierr)); 106720282da2SJoe Wallwork #endif 106882490253SJoe Wallwork if (lierr) { 106982490253SJoe Wallwork for (i = 0; i < dim; ++i) { 107082490253SJoe Wallwork Mpos[i * dim + i] = Mp[i * dim + i]; 107182490253SJoe Wallwork for (j = i + 1; j < dim; ++j) { 107282490253SJoe Wallwork Mpos[i * dim + j] = 0.5 * (Mp[i * dim + j] + Mp[j * dim + i]); 107382490253SJoe Wallwork Mpos[j * dim + i] = Mpos[i * dim + j]; 107482490253SJoe Wallwork } 107582490253SJoe Wallwork } 107682490253SJoe Wallwork LAPACKsyevFail(dim, Mpos); 107798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr); 107882490253SJoe Wallwork } 10799566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 108020282da2SJoe Wallwork } 10819566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 108220282da2SJoe Wallwork } 108320282da2SJoe Wallwork 108420282da2SJoe Wallwork /* Reflect to positive orthant and enforce maximum and minimum size */ 108520282da2SJoe Wallwork max_eig = 0.0; 108620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 108720282da2SJoe Wallwork eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i]))); 108820282da2SJoe Wallwork max_eig = PetscMax(eigs[i], max_eig); 108920282da2SJoe Wallwork } 109020282da2SJoe Wallwork 10913f07679eSJoe Wallwork /* Enforce maximum anisotropy and compute determinant */ 10923f07679eSJoe Wallwork *detMp = 1.0; 109320282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 109420282da2SJoe Wallwork if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig * la_min); 10953f07679eSJoe Wallwork *detMp *= eigs[i]; 109620282da2SJoe Wallwork } 109720282da2SJoe Wallwork 109820282da2SJoe Wallwork /* Reconstruct Hessian */ 109920282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 110020282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 110120282da2SJoe Wallwork Mp[i * dim + j] = 0.0; 1102*9371c9d4SSatish Balay for (k = 0; k < dim; ++k) { Mp[i * dim + j] += Mpos[k * dim + i] * eigs[k] * Mpos[k * dim + j]; } 110320282da2SJoe Wallwork } 110420282da2SJoe Wallwork } 11059566063dSJacob Faibussowitsch PetscCall(PetscFree2(Mpos, eigs)); 110620282da2SJoe Wallwork 110720282da2SJoe Wallwork PetscFunctionReturn(0); 110820282da2SJoe Wallwork } 110920282da2SJoe Wallwork 1110cb7bfe8cSJoe Wallwork /*@ 111120282da2SJoe Wallwork DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric 111220282da2SJoe Wallwork 111320282da2SJoe Wallwork Input parameters: 111420282da2SJoe Wallwork + dm - The DM 11153f07679eSJoe Wallwork . metricIn - The metric 111699abec2bSJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 11173f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced? 111820282da2SJoe Wallwork 111920282da2SJoe Wallwork Output parameter: 11203f07679eSJoe Wallwork + metricOut - The metric 11213f07679eSJoe Wallwork - determinant - Its determinant 112220282da2SJoe Wallwork 112320282da2SJoe Wallwork Level: beginner 112420282da2SJoe Wallwork 1125cb7bfe8cSJoe Wallwork Notes: 1126cb7bfe8cSJoe Wallwork 1127cb7bfe8cSJoe Wallwork Relevant command line options: 1128cb7bfe8cSJoe Wallwork 112993520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic? 113093520041SJoe Wallwork . -dm_plex_metric_uniform - Is the metric uniform? 113193520041SJoe Wallwork . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 1132cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 1133cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max - Maximum tolerated anisotropy 1134cb7bfe8cSJoe Wallwork 1135db781477SPatrick Sanan .seealso: `DMPlexMetricNormalize()`, `DMPlexMetricIntersection()` 1136cb7bfe8cSJoe Wallwork @*/ 1137*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant) { 11383f07679eSJoe Wallwork DM dmDet; 113993520041SJoe Wallwork PetscBool isotropic, uniform; 114020282da2SJoe Wallwork PetscInt dim, vStart, vEnd, v; 11413f07679eSJoe Wallwork PetscScalar *met, *det; 114220282da2SJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0; 114320282da2SJoe Wallwork 114420282da2SJoe Wallwork PetscFunctionBegin; 11459566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0)); 114620282da2SJoe Wallwork 114720282da2SJoe Wallwork /* Extract metadata from dm */ 11489566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 114920282da2SJoe Wallwork if (restrictSizes) { 11509566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min)); 11519566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max)); 115299abec2bSJoe Wallwork h_min = PetscMax(h_min, 1.0e-30); 115399abec2bSJoe Wallwork h_max = PetscMin(h_max, 1.0e+30); 1154bc00d7afSJoe Wallwork PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude"); 115599abec2bSJoe Wallwork } 115699abec2bSJoe Wallwork if (restrictAnisotropy) { 11579566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max)); 115899abec2bSJoe Wallwork a_max = PetscMin(a_max, 1.0e+30); 115920282da2SJoe Wallwork } 116020282da2SJoe Wallwork 116193520041SJoe Wallwork /* Setup output metric */ 11625241a91fSJoe Wallwork PetscCall(VecCopy(metricIn, metricOut)); 116393520041SJoe Wallwork 116493520041SJoe Wallwork /* Enforce SPD and extract determinant */ 11655241a91fSJoe Wallwork PetscCall(VecGetArray(metricOut, &met)); 11669566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 11679566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 116893520041SJoe Wallwork if (uniform) { 1169d1a5b989SMatthew G. Knepley PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist"); 117093520041SJoe Wallwork 117193520041SJoe Wallwork /* Uniform case */ 11725241a91fSJoe Wallwork PetscCall(VecGetArray(determinant, &det)); 11739566063dSJacob Faibussowitsch PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det)); 11745241a91fSJoe Wallwork PetscCall(VecRestoreArray(determinant, &det)); 117593520041SJoe Wallwork } else { 117693520041SJoe Wallwork /* Spatially varying case */ 117793520041SJoe Wallwork PetscInt nrow; 117893520041SJoe Wallwork 117993520041SJoe Wallwork if (isotropic) nrow = 1; 118093520041SJoe Wallwork else nrow = dim; 11815241a91fSJoe Wallwork PetscCall(VecGetDM(determinant, &dmDet)); 11829566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 11835241a91fSJoe Wallwork PetscCall(VecGetArray(determinant, &det)); 118420282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 11853f07679eSJoe Wallwork PetscScalar *vmet, *vdet; 11869566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet)); 11879566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet)); 11889566063dSJacob Faibussowitsch PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet)); 118920282da2SJoe Wallwork } 11905241a91fSJoe Wallwork PetscCall(VecRestoreArray(determinant, &det)); 119193520041SJoe Wallwork } 11925241a91fSJoe Wallwork PetscCall(VecRestoreArray(metricOut, &met)); 1193fe902aceSJoe Wallwork 11949566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD, 0, 0, 0, 0)); 119520282da2SJoe Wallwork PetscFunctionReturn(0); 119620282da2SJoe Wallwork } 119720282da2SJoe Wallwork 1198*9371c9d4SSatish Balay 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[]) { 119920282da2SJoe Wallwork const PetscScalar p = constants[0]; 120020282da2SJoe Wallwork 1201985ad86eSJose E. Roman f0[0] = PetscPowScalar(u[0], p / (2.0 * p + dim)); 120220282da2SJoe Wallwork } 120320282da2SJoe Wallwork 1204cb7bfe8cSJoe Wallwork /*@ 120520282da2SJoe Wallwork DMPlexMetricNormalize - Apply L-p normalization to a metric 120620282da2SJoe Wallwork 120720282da2SJoe Wallwork Input parameters: 120820282da2SJoe Wallwork + dm - The DM 120920282da2SJoe Wallwork . metricIn - The unnormalized metric 121016522936SJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 121116522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced? 121220282da2SJoe Wallwork 121320282da2SJoe Wallwork Output parameter: 121420282da2SJoe Wallwork . metricOut - The normalized metric 121520282da2SJoe Wallwork 121620282da2SJoe Wallwork Level: beginner 121720282da2SJoe Wallwork 1218cb7bfe8cSJoe Wallwork Notes: 1219cb7bfe8cSJoe Wallwork 1220cb7bfe8cSJoe Wallwork Relevant command line options: 1221cb7bfe8cSJoe Wallwork 122293520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic? 122393520041SJoe Wallwork . -dm_plex_metric_uniform - Is the metric uniform? 122493520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 1225cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 1226cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 1227cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max - Maximum tolerated anisotropy 1228cb7bfe8cSJoe Wallwork . -dm_plex_metric_p - L-p normalization order 1229cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity - Target metric complexity 1230cb7bfe8cSJoe Wallwork 1231db781477SPatrick Sanan .seealso: `DMPlexMetricEnforceSPD()`, `DMPlexMetricIntersection()` 1232cb7bfe8cSJoe Wallwork @*/ 1233*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant) { 12343f07679eSJoe Wallwork DM dmDet; 123520282da2SJoe Wallwork MPI_Comm comm; 123693520041SJoe Wallwork PetscBool restrictAnisotropyFirst, isotropic, uniform; 123720282da2SJoe Wallwork PetscDS ds; 123820282da2SJoe Wallwork PetscInt dim, Nd, vStart, vEnd, v, i; 12393f07679eSJoe Wallwork PetscScalar *met, *det, integral, constants[1]; 124093520041SJoe Wallwork PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral; 124120282da2SJoe Wallwork 124220282da2SJoe Wallwork PetscFunctionBegin; 12439566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize, 0, 0, 0, 0)); 124420282da2SJoe Wallwork 124520282da2SJoe Wallwork /* Extract metadata from dm */ 12469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 12479566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 12489566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 12499566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 125093520041SJoe Wallwork if (isotropic) Nd = 1; 125193520041SJoe Wallwork else Nd = dim * dim; 125220282da2SJoe Wallwork 125320282da2SJoe Wallwork /* Set up metric and ensure it is SPD */ 12549566063dSJacob Faibussowitsch PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst)); 12555241a91fSJoe Wallwork PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, determinant)); 125620282da2SJoe Wallwork 125720282da2SJoe Wallwork /* Compute global normalization factor */ 12589566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetTargetComplexity(dm, &target)); 12599566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p)); 126016522936SJoe Wallwork constants[0] = p; 126193520041SJoe Wallwork if (uniform) { 1262bc00d7afSJoe Wallwork PetscCheck(isotropic, comm, PETSC_ERR_SUP, "Uniform anisotropic metrics not supported"); 126393520041SJoe Wallwork DM dmTmp; 126493520041SJoe Wallwork Vec tmp; 126593520041SJoe Wallwork 12669566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmTmp)); 12679566063dSJacob Faibussowitsch PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp)); 12689566063dSJacob Faibussowitsch PetscCall(VecGetArray(determinant, &det)); 12699566063dSJacob Faibussowitsch PetscCall(VecSet(tmp, det[0])); 12709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(determinant, &det)); 12719566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmTmp, &ds)); 12729566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants)); 12739566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, detMFunc)); 12749566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL)); 12759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 12769566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmTmp)); 127793520041SJoe Wallwork } else { 12789566063dSJacob Faibussowitsch PetscCall(VecGetDM(determinant, &dmDet)); 12799566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmDet, &ds)); 12809566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants)); 12819566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, detMFunc)); 12829566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL)); 128393520041SJoe Wallwork } 12843f07679eSJoe Wallwork realIntegral = PetscRealPart(integral); 1285bc00d7afSJoe Wallwork PetscCheck(realIntegral > 1.0e-30, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Global metric normalization factor must be in (0, inf). Is the input metric positive-definite?"); 12863f07679eSJoe Wallwork factGlob = PetscPowReal(target / realIntegral, 2.0 / dim); 128720282da2SJoe Wallwork 128820282da2SJoe Wallwork /* Apply local scaling */ 128920282da2SJoe Wallwork if (restrictSizes) { 12909566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min)); 12919566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max)); 129216522936SJoe Wallwork h_min = PetscMax(h_min, 1.0e-30); 129316522936SJoe Wallwork h_max = PetscMin(h_max, 1.0e+30); 1294bc00d7afSJoe Wallwork PetscCheck(h_min < h_max, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum metric magnitude should be smaller than maximum metric magnitude"); 129516522936SJoe Wallwork } 129616522936SJoe Wallwork if (restrictAnisotropy && !restrictAnisotropyFirst) { 12979566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max)); 129816522936SJoe Wallwork a_max = PetscMin(a_max, 1.0e+30); 129920282da2SJoe Wallwork } 13005241a91fSJoe Wallwork PetscCall(VecGetArray(metricOut, &met)); 13019566063dSJacob Faibussowitsch PetscCall(VecGetArray(determinant, &det)); 130293520041SJoe Wallwork if (uniform) { 130393520041SJoe Wallwork /* Uniform case */ 130493520041SJoe Wallwork met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0 / (2 * p + dim)); 13059566063dSJacob Faibussowitsch if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det)); 130693520041SJoe Wallwork } else { 130793520041SJoe Wallwork /* Spatially varying case */ 130893520041SJoe Wallwork PetscInt nrow; 130993520041SJoe Wallwork 131093520041SJoe Wallwork if (isotropic) nrow = 1; 131193520041SJoe Wallwork else nrow = dim; 13129566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 13135241a91fSJoe Wallwork PetscCall(VecGetDM(determinant, &dmDet)); 131420282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 13153f07679eSJoe Wallwork PetscScalar *Mp, *detM; 131620282da2SJoe Wallwork 13179566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp)); 13189566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM)); 13193f07679eSJoe Wallwork fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0 / (2 * p + dim)); 132020282da2SJoe Wallwork for (i = 0; i < Nd; ++i) Mp[i] *= fact; 13219566063dSJacob Faibussowitsch if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM)); 132293520041SJoe Wallwork } 132320282da2SJoe Wallwork } 13249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(determinant, &det)); 13255241a91fSJoe Wallwork PetscCall(VecRestoreArray(metricOut, &met)); 132620282da2SJoe Wallwork 13279566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize, 0, 0, 0, 0)); 132820282da2SJoe Wallwork PetscFunctionReturn(0); 132920282da2SJoe Wallwork } 133020282da2SJoe Wallwork 1331cb7bfe8cSJoe Wallwork /*@ 133220282da2SJoe Wallwork DMPlexMetricAverage - Compute the average of a list of metrics 133320282da2SJoe Wallwork 1334f1a722f8SMatthew G. Knepley Input Parameters: 133520282da2SJoe Wallwork + dm - The DM 133620282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged 133720282da2SJoe Wallwork . weights - Weights for the average 133820282da2SJoe Wallwork - metrics - The metrics to be averaged 133920282da2SJoe Wallwork 134020282da2SJoe Wallwork Output Parameter: 134120282da2SJoe Wallwork . metricAvg - The averaged metric 134220282da2SJoe Wallwork 134320282da2SJoe Wallwork Level: beginner 134420282da2SJoe Wallwork 134520282da2SJoe Wallwork Notes: 134620282da2SJoe Wallwork The weights should sum to unity. 134720282da2SJoe Wallwork 134820282da2SJoe Wallwork If weights are not provided then an unweighted average is used. 134920282da2SJoe Wallwork 1350db781477SPatrick Sanan .seealso: `DMPlexMetricAverage2()`, `DMPlexMetricAverage3()`, `DMPlexMetricIntersection()` 1351cb7bfe8cSJoe Wallwork @*/ 1352*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg) { 135320282da2SJoe Wallwork PetscBool haveWeights = PETSC_TRUE; 135493520041SJoe Wallwork PetscInt i, m, n; 135520282da2SJoe Wallwork PetscReal sum = 0.0, tol = 1.0e-10; 135620282da2SJoe Wallwork 135720282da2SJoe Wallwork PetscFunctionBegin; 13589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage, 0, 0, 0, 0)); 135963a3b9bcSJacob Faibussowitsch PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %" PetscInt_FMT " < 1 metrics", numMetrics); 13605241a91fSJoe Wallwork PetscCall(VecSet(metricAvg, 0.0)); 13615241a91fSJoe Wallwork PetscCall(VecGetSize(metricAvg, &m)); 136293520041SJoe Wallwork for (i = 0; i < numMetrics; ++i) { 13639566063dSJacob Faibussowitsch PetscCall(VecGetSize(metrics[i], &n)); 13645f80ce2aSJacob Faibussowitsch PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented"); 136593520041SJoe Wallwork } 136620282da2SJoe Wallwork 136720282da2SJoe Wallwork /* Default to the unweighted case */ 136820282da2SJoe Wallwork if (!weights) { 13699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numMetrics, &weights)); 137020282da2SJoe Wallwork haveWeights = PETSC_FALSE; 137120282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) { weights[i] = 1.0 / numMetrics; } 137220282da2SJoe Wallwork } 137320282da2SJoe Wallwork 137420282da2SJoe Wallwork /* Check weights sum to unity */ 137593520041SJoe Wallwork for (i = 0; i < numMetrics; ++i) sum += weights[i]; 13765f80ce2aSJacob Faibussowitsch PetscCheck(PetscAbsReal(sum - 1) <= tol, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); 137720282da2SJoe Wallwork 137820282da2SJoe Wallwork /* Compute metric average */ 13795241a91fSJoe Wallwork for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(metricAvg, weights[i], metrics[i])); 13809566063dSJacob Faibussowitsch if (!haveWeights) PetscCall(PetscFree(weights)); 1381fe902aceSJoe Wallwork 13829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage, 0, 0, 0, 0)); 138320282da2SJoe Wallwork PetscFunctionReturn(0); 138420282da2SJoe Wallwork } 138520282da2SJoe Wallwork 1386cb7bfe8cSJoe Wallwork /*@ 138720282da2SJoe Wallwork DMPlexMetricAverage2 - Compute the unweighted average of two metrics 138820282da2SJoe Wallwork 1389f1a722f8SMatthew G. Knepley Input Parameters: 139020282da2SJoe Wallwork + dm - The DM 139120282da2SJoe Wallwork . metric1 - The first metric to be averaged 139220282da2SJoe Wallwork - metric2 - The second metric to be averaged 139320282da2SJoe Wallwork 139420282da2SJoe Wallwork Output Parameter: 139520282da2SJoe Wallwork . metricAvg - The averaged metric 139620282da2SJoe Wallwork 139720282da2SJoe Wallwork Level: beginner 139820282da2SJoe Wallwork 1399db781477SPatrick Sanan .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage3()` 1400cb7bfe8cSJoe Wallwork @*/ 1401*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec metricAvg) { 140220282da2SJoe Wallwork PetscReal weights[2] = {0.5, 0.5}; 140320282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 140420282da2SJoe Wallwork 140520282da2SJoe Wallwork PetscFunctionBegin; 14069566063dSJacob Faibussowitsch PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg)); 140720282da2SJoe Wallwork PetscFunctionReturn(0); 140820282da2SJoe Wallwork } 140920282da2SJoe Wallwork 1410cb7bfe8cSJoe Wallwork /*@ 141120282da2SJoe Wallwork DMPlexMetricAverage3 - Compute the unweighted average of three metrics 141220282da2SJoe Wallwork 1413f1a722f8SMatthew G. Knepley Input Parameters: 141420282da2SJoe Wallwork + dm - The DM 141520282da2SJoe Wallwork . metric1 - The first metric to be averaged 141620282da2SJoe Wallwork . metric2 - The second metric to be averaged 141720282da2SJoe Wallwork - metric3 - The third metric to be averaged 141820282da2SJoe Wallwork 141920282da2SJoe Wallwork Output Parameter: 142020282da2SJoe Wallwork . metricAvg - The averaged metric 142120282da2SJoe Wallwork 142220282da2SJoe Wallwork Level: beginner 142320282da2SJoe Wallwork 1424db781477SPatrick Sanan .seealso: `DMPlexMetricAverage()`, `DMPlexMetricAverage2()` 1425cb7bfe8cSJoe Wallwork @*/ 1426*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricAvg) { 142720282da2SJoe Wallwork PetscReal weights[3] = {1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0}; 142820282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 142920282da2SJoe Wallwork 143020282da2SJoe Wallwork PetscFunctionBegin; 14319566063dSJacob Faibussowitsch PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg)); 143220282da2SJoe Wallwork PetscFunctionReturn(0); 143320282da2SJoe Wallwork } 143420282da2SJoe Wallwork 1435*9371c9d4SSatish Balay static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[]) { 143620282da2SJoe Wallwork PetscInt i, j, k, l, m; 1437e2606525SJoe Wallwork PetscReal *evals; 143820282da2SJoe Wallwork PetscScalar *evecs, *sqrtM1, *isqrtM1; 143920282da2SJoe Wallwork 144020282da2SJoe Wallwork PetscFunctionBegin; 144193520041SJoe Wallwork 144293520041SJoe Wallwork /* Isotropic case */ 144393520041SJoe Wallwork if (dim == 1) { 14446f71502aSJoe Wallwork M2[0] = (PetscScalar)PetscMax(PetscRealPart(M1[0]), PetscRealPart(M2[0])); 144593520041SJoe Wallwork PetscFunctionReturn(0); 144693520041SJoe Wallwork } 144793520041SJoe Wallwork 144893520041SJoe Wallwork /* Anisotropic case */ 1449e2606525SJoe Wallwork PetscCall(PetscMalloc4(dim * dim, &evecs, dim * dim, &sqrtM1, dim * dim, &isqrtM1, dim, &evals)); 145020282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 1451*9371c9d4SSatish Balay for (j = 0; j < dim; ++j) { evecs[i * dim + j] = M1[i * dim + j]; } 145220282da2SJoe Wallwork } 145320282da2SJoe Wallwork { 145420282da2SJoe Wallwork PetscScalar *work; 145520282da2SJoe Wallwork PetscBLASInt lwork; 145620282da2SJoe Wallwork 145720282da2SJoe Wallwork lwork = 5 * dim; 14589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * dim, &work)); 145920282da2SJoe Wallwork { 146020282da2SJoe Wallwork PetscBLASInt lierr, nb; 1461e2606525SJoe Wallwork PetscReal sqrtj; 146220282da2SJoe Wallwork 146320282da2SJoe Wallwork /* Compute eigendecomposition of M1 */ 14649566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(dim, &nb)); 14659566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 146620282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 146720282da2SJoe Wallwork { 146820282da2SJoe Wallwork PetscReal *rwork; 14699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3 * dim, &rwork)); 1470792fecdfSBarry Smith PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 14719566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 147220282da2SJoe Wallwork } 147320282da2SJoe Wallwork #else 1474792fecdfSBarry Smith PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 147520282da2SJoe Wallwork #endif 147682490253SJoe Wallwork if (lierr) { 147782490253SJoe Wallwork LAPACKsyevFail(dim, M1); 147898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr); 147982490253SJoe Wallwork } 14809566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 148120282da2SJoe Wallwork 1482e2606525SJoe Wallwork /* Compute square root and the reciprocal thereof */ 148320282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 148420282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 1485e2606525SJoe Wallwork sqrtM1[i * dim + k] = 0.0; 1486e2606525SJoe Wallwork isqrtM1[i * dim + k] = 0.0; 1487e2606525SJoe Wallwork for (j = 0; j < dim; ++j) { 1488e2606525SJoe Wallwork sqrtj = PetscSqrtReal(evals[j]); 1489e2606525SJoe Wallwork sqrtM1[i * dim + k] += evecs[j * dim + i] * sqrtj * evecs[j * dim + k]; 1490e2606525SJoe Wallwork isqrtM1[i * dim + k] += evecs[j * dim + i] * (1.0 / sqrtj) * evecs[j * dim + k]; 149120282da2SJoe Wallwork } 149220282da2SJoe Wallwork } 149320282da2SJoe Wallwork } 149420282da2SJoe Wallwork 1495e2606525SJoe Wallwork /* Map M2 into the space spanned by the eigenvectors of M1 */ 149620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 149720282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 1498e2606525SJoe Wallwork evecs[i * dim + l] = 0.0; 1499e2606525SJoe Wallwork for (j = 0; j < dim; ++j) { 1500*9371c9d4SSatish Balay for (k = 0; k < dim; ++k) { evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l]; } 150120282da2SJoe Wallwork } 150220282da2SJoe Wallwork } 150320282da2SJoe Wallwork } 150420282da2SJoe Wallwork 150520282da2SJoe Wallwork /* Compute eigendecomposition */ 15069566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 150720282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 150820282da2SJoe Wallwork { 150920282da2SJoe Wallwork PetscReal *rwork; 15109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3 * dim, &rwork)); 1511792fecdfSBarry Smith PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 15129566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 151320282da2SJoe Wallwork } 151420282da2SJoe Wallwork #else 1515792fecdfSBarry Smith PetscCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 151620282da2SJoe Wallwork #endif 151782490253SJoe Wallwork if (lierr) { 151882490253SJoe Wallwork for (i = 0; i < dim; ++i) { 151982490253SJoe Wallwork for (l = 0; l < dim; ++l) { 1520e2606525SJoe Wallwork evecs[i * dim + l] = 0.0; 1521e2606525SJoe Wallwork for (j = 0; j < dim; ++j) { 1522*9371c9d4SSatish Balay for (k = 0; k < dim; ++k) { evecs[i * dim + l] += isqrtM1[j * dim + i] * M2[j * dim + k] * isqrtM1[k * dim + l]; } 152382490253SJoe Wallwork } 152482490253SJoe Wallwork } 152582490253SJoe Wallwork } 152682490253SJoe Wallwork LAPACKsyevFail(dim, evecs); 152798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr); 152882490253SJoe Wallwork } 15299566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 153020282da2SJoe Wallwork 153120282da2SJoe Wallwork /* Modify eigenvalues */ 1532e2606525SJoe Wallwork for (i = 0; i < dim; ++i) evals[i] = PetscMax(evals[i], 1.0); 153320282da2SJoe Wallwork 153420282da2SJoe Wallwork /* Map back to get the intersection */ 153520282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 1536e2606525SJoe Wallwork for (m = 0; m < dim; ++m) { 1537e2606525SJoe Wallwork M2[i * dim + m] = 0.0; 153820282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 153920282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 1540*9371c9d4SSatish Balay 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]; } 154120282da2SJoe Wallwork } 154220282da2SJoe Wallwork } 154320282da2SJoe Wallwork } 154420282da2SJoe Wallwork } 154520282da2SJoe Wallwork } 15469566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 154720282da2SJoe Wallwork } 1548e2606525SJoe Wallwork PetscCall(PetscFree4(evecs, sqrtM1, isqrtM1, evals)); 154920282da2SJoe Wallwork PetscFunctionReturn(0); 155020282da2SJoe Wallwork } 155120282da2SJoe Wallwork 1552cb7bfe8cSJoe Wallwork /*@ 155320282da2SJoe Wallwork DMPlexMetricIntersection - Compute the intersection of a list of metrics 155420282da2SJoe Wallwork 1555f1a722f8SMatthew G. Knepley Input Parameters: 155620282da2SJoe Wallwork + dm - The DM 155720282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected 155820282da2SJoe Wallwork - metrics - The metrics to be intersected 155920282da2SJoe Wallwork 156020282da2SJoe Wallwork Output Parameter: 156120282da2SJoe Wallwork . metricInt - The intersected metric 156220282da2SJoe Wallwork 156320282da2SJoe Wallwork Level: beginner 156420282da2SJoe Wallwork 156520282da2SJoe Wallwork Notes: 156620282da2SJoe Wallwork 1567e2606525SJoe Wallwork The intersection of a list of metrics has the minimal ellipsoid which fits within the ellipsoids of the component metrics. 156820282da2SJoe Wallwork 1569e2606525SJoe Wallwork The implementation used here is only consistent with the minimal ellipsoid definition in the case numMetrics = 2. 157020282da2SJoe Wallwork 1571db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection2()`, `DMPlexMetricIntersection3()`, `DMPlexMetricAverage()` 1572cb7bfe8cSJoe Wallwork @*/ 1573*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt) { 157493520041SJoe Wallwork PetscBool isotropic, uniform; 157593520041SJoe Wallwork PetscInt v, i, m, n; 157693520041SJoe Wallwork PetscScalar *met, *meti; 157720282da2SJoe Wallwork 157820282da2SJoe Wallwork PetscFunctionBegin; 15799566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection, 0, 0, 0, 0)); 158063a3b9bcSJacob Faibussowitsch PetscCheck(numMetrics >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %" PetscInt_FMT " < 1 metrics", numMetrics); 158120282da2SJoe Wallwork 158220282da2SJoe Wallwork /* Copy over the first metric */ 15835241a91fSJoe Wallwork PetscCall(VecCopy(metrics[0], metricInt)); 158493520041SJoe Wallwork if (numMetrics == 1) PetscFunctionReturn(0); 15855241a91fSJoe Wallwork PetscCall(VecGetSize(metricInt, &m)); 158693520041SJoe Wallwork for (i = 0; i < numMetrics; ++i) { 15879566063dSJacob Faibussowitsch PetscCall(VecGetSize(metrics[i], &n)); 158808401ef6SPierre Jolivet PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented"); 158993520041SJoe Wallwork } 159020282da2SJoe Wallwork 159120282da2SJoe Wallwork /* Intersect subsequent metrics in turn */ 15929566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 15939566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 159493520041SJoe Wallwork if (uniform) { 159593520041SJoe Wallwork /* Uniform case */ 15965241a91fSJoe Wallwork PetscCall(VecGetArray(metricInt, &met)); 159793520041SJoe Wallwork for (i = 1; i < numMetrics; ++i) { 15989566063dSJacob Faibussowitsch PetscCall(VecGetArray(metrics[i], &meti)); 15999566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection_Private(1, meti, met)); 16009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(metrics[i], &meti)); 160193520041SJoe Wallwork } 16025241a91fSJoe Wallwork PetscCall(VecRestoreArray(metricInt, &met)); 160393520041SJoe Wallwork } else { 160493520041SJoe Wallwork /* Spatially varying case */ 160593520041SJoe Wallwork PetscInt dim, vStart, vEnd, nrow; 160693520041SJoe Wallwork PetscScalar *M, *Mi; 160793520041SJoe Wallwork 16089566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 160993520041SJoe Wallwork if (isotropic) nrow = 1; 161093520041SJoe Wallwork else nrow = dim; 16119566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 16125241a91fSJoe Wallwork PetscCall(VecGetArray(metricInt, &met)); 161320282da2SJoe Wallwork for (i = 1; i < numMetrics; ++i) { 16149566063dSJacob Faibussowitsch PetscCall(VecGetArray(metrics[i], &meti)); 161520282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 16169566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &M)); 16179566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi)); 16189566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M)); 161920282da2SJoe Wallwork } 16209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(metrics[i], &meti)); 162120282da2SJoe Wallwork } 16225241a91fSJoe Wallwork PetscCall(VecRestoreArray(metricInt, &met)); 162320282da2SJoe Wallwork } 1624fe902aceSJoe Wallwork 16259566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection, 0, 0, 0, 0)); 162620282da2SJoe Wallwork PetscFunctionReturn(0); 162720282da2SJoe Wallwork } 162820282da2SJoe Wallwork 1629cb7bfe8cSJoe Wallwork /*@ 163020282da2SJoe Wallwork DMPlexMetricIntersection2 - Compute the intersection of two metrics 163120282da2SJoe Wallwork 1632f1a722f8SMatthew G. Knepley Input Parameters: 163320282da2SJoe Wallwork + dm - The DM 163420282da2SJoe Wallwork . metric1 - The first metric to be intersected 163520282da2SJoe Wallwork - metric2 - The second metric to be intersected 163620282da2SJoe Wallwork 163720282da2SJoe Wallwork Output Parameter: 163820282da2SJoe Wallwork . metricInt - The intersected metric 163920282da2SJoe Wallwork 164020282da2SJoe Wallwork Level: beginner 164120282da2SJoe Wallwork 1642db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection3()` 1643cb7bfe8cSJoe Wallwork @*/ 1644*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec metricInt) { 164520282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 164620282da2SJoe Wallwork 164720282da2SJoe Wallwork PetscFunctionBegin; 16489566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt)); 164920282da2SJoe Wallwork PetscFunctionReturn(0); 165020282da2SJoe Wallwork } 165120282da2SJoe Wallwork 1652cb7bfe8cSJoe Wallwork /*@ 165320282da2SJoe Wallwork DMPlexMetricIntersection3 - Compute the intersection of three metrics 165420282da2SJoe Wallwork 1655f1a722f8SMatthew G. Knepley Input Parameters: 165620282da2SJoe Wallwork + dm - The DM 165720282da2SJoe Wallwork . metric1 - The first metric to be intersected 165820282da2SJoe Wallwork . metric2 - The second metric to be intersected 165920282da2SJoe Wallwork - metric3 - The third metric to be intersected 166020282da2SJoe Wallwork 166120282da2SJoe Wallwork Output Parameter: 166220282da2SJoe Wallwork . metricInt - The intersected metric 166320282da2SJoe Wallwork 166420282da2SJoe Wallwork Level: beginner 166520282da2SJoe Wallwork 1666db781477SPatrick Sanan .seealso: `DMPlexMetricIntersection()`, `DMPlexMetricIntersection2()` 1667cb7bfe8cSJoe Wallwork @*/ 1668*9371c9d4SSatish Balay PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec metricInt) { 166920282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 167020282da2SJoe Wallwork 167120282da2SJoe Wallwork PetscFunctionBegin; 16729566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt)); 167320282da2SJoe Wallwork PetscFunctionReturn(0); 167420282da2SJoe Wallwork } 1675