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 631b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm) 731b70465SJoe Wallwork { 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; 1131b70465SJoe Wallwork PetscErrorCode ierr; 12cc2eee55SJoe Wallwork PetscInt verbosity = -1, numIter = 3; 13ae8b063eSJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0, beta = 1.3, hausd = 0.01; 1431b70465SJoe Wallwork 1531b70465SJoe Wallwork PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 179566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");PetscCall(ierr); 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)); 509566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 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 6393520041SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst() 64cb7bfe8cSJoe Wallwork @*/ 6531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic) 6631b70465SJoe Wallwork { 6731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 6831b70465SJoe Wallwork 6931b70465SJoe Wallwork PetscFunctionBegin; 7031b70465SJoe Wallwork if (!plex->metricCtx) { 719566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 729566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 7331b70465SJoe Wallwork } 7431b70465SJoe Wallwork plex->metricCtx->isotropic = isotropic; 7531b70465SJoe Wallwork PetscFunctionReturn(0); 7631b70465SJoe Wallwork } 7731b70465SJoe Wallwork 78cb7bfe8cSJoe Wallwork /*@ 7993520041SJoe Wallwork DMPlexMetricIsIsotropic - Is a metric isotropic? 8031b70465SJoe Wallwork 8131b70465SJoe Wallwork Input parameters: 8231b70465SJoe Wallwork . dm - The DM 8331b70465SJoe Wallwork 8431b70465SJoe Wallwork Output parameters: 8531b70465SJoe Wallwork . isotropic - Is the metric isotropic? 8631b70465SJoe Wallwork 8731b70465SJoe Wallwork Level: beginner 8831b70465SJoe Wallwork 8993520041SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst() 90cb7bfe8cSJoe Wallwork @*/ 9131b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic) 9231b70465SJoe Wallwork { 9331b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 9431b70465SJoe Wallwork 9531b70465SJoe Wallwork PetscFunctionBegin; 9631b70465SJoe Wallwork if (!plex->metricCtx) { 979566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 989566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 9931b70465SJoe Wallwork } 10031b70465SJoe Wallwork *isotropic = plex->metricCtx->isotropic; 10131b70465SJoe Wallwork PetscFunctionReturn(0); 10231b70465SJoe Wallwork } 10331b70465SJoe Wallwork 104cb7bfe8cSJoe Wallwork /*@ 10593520041SJoe Wallwork DMPlexMetricSetUniform - Record whether a metric is uniform 10693520041SJoe Wallwork 10793520041SJoe Wallwork Input parameters: 10893520041SJoe Wallwork + dm - The DM 10993520041SJoe Wallwork - uniform - Is the metric uniform? 11093520041SJoe Wallwork 11193520041SJoe Wallwork Level: beginner 11293520041SJoe Wallwork 11393520041SJoe Wallwork Notes: 11493520041SJoe Wallwork 11593520041SJoe Wallwork If the metric is specified as uniform then it is assumed to be isotropic, too. 11693520041SJoe Wallwork 11793520041SJoe Wallwork .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 11893520041SJoe Wallwork @*/ 11993520041SJoe Wallwork PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform) 12093520041SJoe Wallwork { 12193520041SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 12293520041SJoe Wallwork 12393520041SJoe Wallwork PetscFunctionBegin; 12493520041SJoe Wallwork if (!plex->metricCtx) { 1259566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 1269566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 12793520041SJoe Wallwork } 12893520041SJoe Wallwork plex->metricCtx->uniform = uniform; 12993520041SJoe Wallwork if (uniform) plex->metricCtx->isotropic = uniform; 13093520041SJoe Wallwork PetscFunctionReturn(0); 13193520041SJoe Wallwork } 13293520041SJoe Wallwork 13393520041SJoe Wallwork /*@ 13493520041SJoe Wallwork DMPlexMetricIsUniform - Is a metric uniform? 13593520041SJoe Wallwork 13693520041SJoe Wallwork Input parameters: 13793520041SJoe Wallwork . dm - The DM 13893520041SJoe Wallwork 13993520041SJoe Wallwork Output parameters: 14093520041SJoe Wallwork . uniform - Is the metric uniform? 14193520041SJoe Wallwork 14293520041SJoe Wallwork Level: beginner 14393520041SJoe Wallwork 14493520041SJoe Wallwork .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 14593520041SJoe Wallwork @*/ 14693520041SJoe Wallwork PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform) 14793520041SJoe Wallwork { 14893520041SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 14993520041SJoe Wallwork 15093520041SJoe Wallwork PetscFunctionBegin; 15193520041SJoe Wallwork if (!plex->metricCtx) { 1529566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 1539566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 15493520041SJoe Wallwork } 15593520041SJoe Wallwork *uniform = plex->metricCtx->uniform; 15693520041SJoe Wallwork PetscFunctionReturn(0); 15793520041SJoe Wallwork } 15893520041SJoe Wallwork 15993520041SJoe Wallwork /*@ 16031b70465SJoe Wallwork DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization 16131b70465SJoe Wallwork 16231b70465SJoe Wallwork Input parameters: 16331b70465SJoe Wallwork + dm - The DM 16431b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first? 16531b70465SJoe Wallwork 16631b70465SJoe Wallwork Level: beginner 16731b70465SJoe Wallwork 16831b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 169cb7bfe8cSJoe Wallwork @*/ 17031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst) 17131b70465SJoe Wallwork { 17231b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 17331b70465SJoe Wallwork 17431b70465SJoe Wallwork PetscFunctionBegin; 17531b70465SJoe Wallwork if (!plex->metricCtx) { 1769566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 1779566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 17831b70465SJoe Wallwork } 17931b70465SJoe Wallwork plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst; 18031b70465SJoe Wallwork PetscFunctionReturn(0); 18131b70465SJoe Wallwork } 18231b70465SJoe Wallwork 183cb7bfe8cSJoe Wallwork /*@ 18431b70465SJoe Wallwork DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after? 18531b70465SJoe Wallwork 18631b70465SJoe Wallwork Input parameters: 18731b70465SJoe Wallwork . dm - The DM 18831b70465SJoe Wallwork 18931b70465SJoe Wallwork Output parameters: 19031b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first? 19131b70465SJoe Wallwork 19231b70465SJoe Wallwork Level: beginner 19331b70465SJoe Wallwork 19431b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 195cb7bfe8cSJoe Wallwork @*/ 19631b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst) 19731b70465SJoe Wallwork { 19831b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 19931b70465SJoe Wallwork 20031b70465SJoe Wallwork PetscFunctionBegin; 20131b70465SJoe Wallwork if (!plex->metricCtx) { 2029566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 2039566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 20431b70465SJoe Wallwork } 20531b70465SJoe Wallwork *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst; 20631b70465SJoe Wallwork PetscFunctionReturn(0); 20731b70465SJoe Wallwork } 20831b70465SJoe Wallwork 209cb7bfe8cSJoe Wallwork /*@ 210cc2eee55SJoe Wallwork DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off? 211cc2eee55SJoe Wallwork 212cc2eee55SJoe Wallwork Input parameters: 213cc2eee55SJoe Wallwork + dm - The DM 214cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off? 215cc2eee55SJoe Wallwork 216cc2eee55SJoe Wallwork Level: beginner 217cc2eee55SJoe Wallwork 218cb7bfe8cSJoe Wallwork Notes: 219cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 220cb7bfe8cSJoe Wallwork 22176f360caSJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf() 222cb7bfe8cSJoe Wallwork @*/ 223cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert) 224cc2eee55SJoe Wallwork { 225cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 226cc2eee55SJoe Wallwork 227cc2eee55SJoe Wallwork PetscFunctionBegin; 228cc2eee55SJoe Wallwork if (!plex->metricCtx) { 2299566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 2309566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 231cc2eee55SJoe Wallwork } 232cc2eee55SJoe Wallwork plex->metricCtx->noInsert = noInsert; 233cc2eee55SJoe Wallwork PetscFunctionReturn(0); 234cc2eee55SJoe Wallwork } 235cc2eee55SJoe Wallwork 236cb7bfe8cSJoe Wallwork /*@ 237cc2eee55SJoe Wallwork DMPlexMetricNoInsertion - Are node insertion and deletion turned off? 238cc2eee55SJoe Wallwork 239cc2eee55SJoe Wallwork Input parameters: 240cc2eee55SJoe Wallwork . dm - The DM 241cc2eee55SJoe Wallwork 242cc2eee55SJoe Wallwork Output parameters: 243cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off? 244cc2eee55SJoe Wallwork 245cc2eee55SJoe Wallwork Level: beginner 246cc2eee55SJoe Wallwork 247cb7bfe8cSJoe Wallwork Notes: 248cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 249cb7bfe8cSJoe Wallwork 25076f360caSJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf() 251cb7bfe8cSJoe Wallwork @*/ 252cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert) 253cc2eee55SJoe Wallwork { 254cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 255cc2eee55SJoe Wallwork 256cc2eee55SJoe Wallwork PetscFunctionBegin; 257cc2eee55SJoe Wallwork if (!plex->metricCtx) { 2589566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 2599566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 260cc2eee55SJoe Wallwork } 261cc2eee55SJoe Wallwork *noInsert = plex->metricCtx->noInsert; 262cc2eee55SJoe Wallwork PetscFunctionReturn(0); 263cc2eee55SJoe Wallwork } 264cc2eee55SJoe Wallwork 265cb7bfe8cSJoe Wallwork /*@ 266cc2eee55SJoe Wallwork DMPlexMetricSetNoSwapping - Should facet swapping be turned off? 267cc2eee55SJoe Wallwork 268cc2eee55SJoe Wallwork Input parameters: 269cc2eee55SJoe Wallwork + dm - The DM 270cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off? 271cc2eee55SJoe Wallwork 272cc2eee55SJoe Wallwork Level: beginner 273cc2eee55SJoe Wallwork 274cb7bfe8cSJoe Wallwork Notes: 275cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 276cb7bfe8cSJoe Wallwork 27776f360caSJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf() 278cb7bfe8cSJoe Wallwork @*/ 279cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap) 280cc2eee55SJoe Wallwork { 281cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 282cc2eee55SJoe Wallwork 283cc2eee55SJoe Wallwork PetscFunctionBegin; 284cc2eee55SJoe Wallwork if (!plex->metricCtx) { 2859566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 2869566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 287cc2eee55SJoe Wallwork } 288cc2eee55SJoe Wallwork plex->metricCtx->noSwap = noSwap; 289cc2eee55SJoe Wallwork PetscFunctionReturn(0); 290cc2eee55SJoe Wallwork } 291cc2eee55SJoe Wallwork 292cb7bfe8cSJoe Wallwork /*@ 293cc2eee55SJoe Wallwork DMPlexMetricNoSwapping - Is facet swapping turned off? 294cc2eee55SJoe Wallwork 295cc2eee55SJoe Wallwork Input parameters: 296cc2eee55SJoe Wallwork . dm - The DM 297cc2eee55SJoe Wallwork 298cc2eee55SJoe Wallwork Output parameters: 299cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off? 300cc2eee55SJoe Wallwork 301cc2eee55SJoe Wallwork Level: beginner 302cc2eee55SJoe Wallwork 303cb7bfe8cSJoe Wallwork Notes: 304cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 305cb7bfe8cSJoe Wallwork 30676f360caSJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf() 307cb7bfe8cSJoe Wallwork @*/ 308cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap) 309cc2eee55SJoe Wallwork { 310cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 311cc2eee55SJoe Wallwork 312cc2eee55SJoe Wallwork PetscFunctionBegin; 313cc2eee55SJoe Wallwork if (!plex->metricCtx) { 3149566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 3159566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 316cc2eee55SJoe Wallwork } 317cc2eee55SJoe Wallwork *noSwap = plex->metricCtx->noSwap; 318cc2eee55SJoe Wallwork PetscFunctionReturn(0); 319cc2eee55SJoe Wallwork } 320cc2eee55SJoe Wallwork 321cb7bfe8cSJoe Wallwork /*@ 322cc2eee55SJoe Wallwork DMPlexMetricSetNoMovement - Should node movement be turned off? 323cc2eee55SJoe Wallwork 324cc2eee55SJoe Wallwork Input parameters: 325cc2eee55SJoe Wallwork + dm - The DM 326cc2eee55SJoe Wallwork - noMove - Should node movement be turned off? 327cc2eee55SJoe Wallwork 328cc2eee55SJoe Wallwork Level: beginner 329cc2eee55SJoe Wallwork 330cb7bfe8cSJoe Wallwork Notes: 331cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 332cb7bfe8cSJoe Wallwork 33376f360caSJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoSurf() 334cb7bfe8cSJoe Wallwork @*/ 335cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove) 336cc2eee55SJoe Wallwork { 337cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 338cc2eee55SJoe Wallwork 339cc2eee55SJoe Wallwork PetscFunctionBegin; 340cc2eee55SJoe Wallwork if (!plex->metricCtx) { 3419566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 3429566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 343cc2eee55SJoe Wallwork } 344cc2eee55SJoe Wallwork plex->metricCtx->noMove = noMove; 345cc2eee55SJoe Wallwork PetscFunctionReturn(0); 346cc2eee55SJoe Wallwork } 347cc2eee55SJoe Wallwork 348cb7bfe8cSJoe Wallwork /*@ 349cc2eee55SJoe Wallwork DMPlexMetricNoMovement - Is node movement turned off? 350cc2eee55SJoe Wallwork 351cc2eee55SJoe Wallwork Input parameters: 352cc2eee55SJoe Wallwork . dm - The DM 353cc2eee55SJoe Wallwork 354cc2eee55SJoe Wallwork Output parameters: 355cc2eee55SJoe Wallwork . noMove - Is node movement turned off? 356cc2eee55SJoe Wallwork 357cc2eee55SJoe Wallwork Level: beginner 358cc2eee55SJoe Wallwork 359cb7bfe8cSJoe Wallwork Notes: 360cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 361cb7bfe8cSJoe Wallwork 36276f360caSJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoSurf() 363cb7bfe8cSJoe Wallwork @*/ 364cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove) 365cc2eee55SJoe Wallwork { 366cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 367cc2eee55SJoe Wallwork 368cc2eee55SJoe Wallwork PetscFunctionBegin; 369cc2eee55SJoe Wallwork if (!plex->metricCtx) { 3709566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 3719566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 372cc2eee55SJoe Wallwork } 373cc2eee55SJoe Wallwork *noMove = plex->metricCtx->noMove; 374cc2eee55SJoe Wallwork PetscFunctionReturn(0); 375cc2eee55SJoe Wallwork } 376cc2eee55SJoe Wallwork 377cb7bfe8cSJoe Wallwork /*@ 37876f360caSJoe Wallwork DMPlexMetricSetNoSurf - Should surface modification be turned off? 37976f360caSJoe Wallwork 38076f360caSJoe Wallwork Input parameters: 38176f360caSJoe Wallwork + dm - The DM 38276f360caSJoe Wallwork - noSurf - Should surface modification be turned off? 38376f360caSJoe Wallwork 38476f360caSJoe Wallwork Level: beginner 38576f360caSJoe Wallwork 38676f360caSJoe Wallwork Notes: 38776f360caSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 38876f360caSJoe Wallwork 38976f360caSJoe Wallwork .seealso: DMPlexMetricNoSurf(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping() 39076f360caSJoe Wallwork @*/ 39176f360caSJoe Wallwork PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf) 39276f360caSJoe Wallwork { 39376f360caSJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 39476f360caSJoe Wallwork 39576f360caSJoe Wallwork PetscFunctionBegin; 39676f360caSJoe Wallwork if (!plex->metricCtx) { 3979566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 3989566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 39976f360caSJoe Wallwork } 40076f360caSJoe Wallwork plex->metricCtx->noSurf = noSurf; 40176f360caSJoe Wallwork PetscFunctionReturn(0); 40276f360caSJoe Wallwork } 40376f360caSJoe Wallwork 40476f360caSJoe Wallwork /*@ 40576f360caSJoe Wallwork DMPlexMetricNoSurf - Is surface modification turned off? 40676f360caSJoe Wallwork 40776f360caSJoe Wallwork Input parameters: 40876f360caSJoe Wallwork . dm - The DM 40976f360caSJoe Wallwork 41076f360caSJoe Wallwork Output parameters: 41176f360caSJoe Wallwork . noSurf - Is surface modification turned off? 41276f360caSJoe Wallwork 41376f360caSJoe Wallwork Level: beginner 41476f360caSJoe Wallwork 41576f360caSJoe Wallwork Notes: 41676f360caSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 41776f360caSJoe Wallwork 41876f360caSJoe Wallwork .seealso: DMPlexMetricSetNoSurf(), DMPlexMetricNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping() 41976f360caSJoe Wallwork @*/ 42076f360caSJoe Wallwork PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf) 42176f360caSJoe Wallwork { 42276f360caSJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 42376f360caSJoe Wallwork 42476f360caSJoe Wallwork PetscFunctionBegin; 42576f360caSJoe Wallwork if (!plex->metricCtx) { 4269566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 4279566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 42876f360caSJoe Wallwork } 42976f360caSJoe Wallwork *noSurf = plex->metricCtx->noSurf; 43076f360caSJoe Wallwork PetscFunctionReturn(0); 43176f360caSJoe Wallwork } 43276f360caSJoe Wallwork 43376f360caSJoe Wallwork /*@ 43431b70465SJoe Wallwork DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude 43531b70465SJoe Wallwork 43631b70465SJoe Wallwork Input parameters: 43731b70465SJoe Wallwork + dm - The DM 43831b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude 43931b70465SJoe Wallwork 44031b70465SJoe Wallwork Level: beginner 44131b70465SJoe Wallwork 44231b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude() 443cb7bfe8cSJoe Wallwork @*/ 44431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min) 44531b70465SJoe Wallwork { 44631b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 44731b70465SJoe Wallwork 44831b70465SJoe Wallwork PetscFunctionBegin; 44931b70465SJoe Wallwork if (!plex->metricCtx) { 4509566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 4519566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 45231b70465SJoe Wallwork } 453*08401ef6SPierre Jolivet PetscCheck(h_min > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min); 45431b70465SJoe Wallwork plex->metricCtx->h_min = h_min; 45531b70465SJoe Wallwork PetscFunctionReturn(0); 45631b70465SJoe Wallwork } 45731b70465SJoe Wallwork 458cb7bfe8cSJoe Wallwork /*@ 45931b70465SJoe Wallwork DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude 46031b70465SJoe Wallwork 46131b70465SJoe Wallwork Input parameters: 46231b70465SJoe Wallwork . dm - The DM 46331b70465SJoe Wallwork 464cc2eee55SJoe Wallwork Output parameters: 46531b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude 46631b70465SJoe Wallwork 46731b70465SJoe Wallwork Level: beginner 46831b70465SJoe Wallwork 46931b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude() 470cb7bfe8cSJoe Wallwork @*/ 47131b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min) 47231b70465SJoe Wallwork { 47331b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 47431b70465SJoe Wallwork 47531b70465SJoe Wallwork PetscFunctionBegin; 47631b70465SJoe Wallwork if (!plex->metricCtx) { 4779566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 4789566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 47931b70465SJoe Wallwork } 48031b70465SJoe Wallwork *h_min = plex->metricCtx->h_min; 48131b70465SJoe Wallwork PetscFunctionReturn(0); 48231b70465SJoe Wallwork } 48331b70465SJoe Wallwork 484cb7bfe8cSJoe Wallwork /*@ 48531b70465SJoe Wallwork DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude 48631b70465SJoe Wallwork 48731b70465SJoe Wallwork Input parameters: 48831b70465SJoe Wallwork + dm - The DM 48931b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude 49031b70465SJoe Wallwork 49131b70465SJoe Wallwork Level: beginner 49231b70465SJoe Wallwork 49331b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude() 494cb7bfe8cSJoe Wallwork @*/ 49531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max) 49631b70465SJoe Wallwork { 49731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 49831b70465SJoe Wallwork 49931b70465SJoe Wallwork PetscFunctionBegin; 50031b70465SJoe Wallwork if (!plex->metricCtx) { 5019566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 5029566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 50331b70465SJoe Wallwork } 504*08401ef6SPierre Jolivet PetscCheck(h_max > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max); 50531b70465SJoe Wallwork plex->metricCtx->h_max = h_max; 50631b70465SJoe Wallwork PetscFunctionReturn(0); 50731b70465SJoe Wallwork } 50831b70465SJoe Wallwork 509cb7bfe8cSJoe Wallwork /*@ 51031b70465SJoe Wallwork DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude 51131b70465SJoe Wallwork 51231b70465SJoe Wallwork Input parameters: 51331b70465SJoe Wallwork . dm - The DM 51431b70465SJoe Wallwork 515cc2eee55SJoe Wallwork Output parameters: 51631b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude 51731b70465SJoe Wallwork 51831b70465SJoe Wallwork Level: beginner 51931b70465SJoe Wallwork 52031b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude() 521cb7bfe8cSJoe Wallwork @*/ 52231b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max) 52331b70465SJoe Wallwork { 52431b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 52531b70465SJoe Wallwork 52631b70465SJoe Wallwork PetscFunctionBegin; 52731b70465SJoe Wallwork if (!plex->metricCtx) { 5289566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 5299566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 53031b70465SJoe Wallwork } 53131b70465SJoe Wallwork *h_max = plex->metricCtx->h_max; 53231b70465SJoe Wallwork PetscFunctionReturn(0); 53331b70465SJoe Wallwork } 53431b70465SJoe Wallwork 535cb7bfe8cSJoe Wallwork /*@ 53631b70465SJoe Wallwork DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy 53731b70465SJoe Wallwork 53831b70465SJoe Wallwork Input parameters: 53931b70465SJoe Wallwork + dm - The DM 54031b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy 54131b70465SJoe Wallwork 54231b70465SJoe Wallwork Level: beginner 54331b70465SJoe Wallwork 54431b70465SJoe Wallwork Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one. 54531b70465SJoe Wallwork 54631b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude() 547cb7bfe8cSJoe Wallwork @*/ 54831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max) 54931b70465SJoe Wallwork { 55031b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 55131b70465SJoe Wallwork 55231b70465SJoe Wallwork PetscFunctionBegin; 55331b70465SJoe Wallwork if (!plex->metricCtx) { 5549566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 5559566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 55631b70465SJoe Wallwork } 557*08401ef6SPierre Jolivet PetscCheck(a_max >= 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max); 55831b70465SJoe Wallwork plex->metricCtx->a_max = a_max; 55931b70465SJoe Wallwork PetscFunctionReturn(0); 56031b70465SJoe Wallwork } 56131b70465SJoe Wallwork 562cb7bfe8cSJoe Wallwork /*@ 56331b70465SJoe Wallwork DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy 56431b70465SJoe Wallwork 56531b70465SJoe Wallwork Input parameters: 56631b70465SJoe Wallwork . dm - The DM 56731b70465SJoe Wallwork 568cc2eee55SJoe Wallwork Output parameters: 56931b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy 57031b70465SJoe Wallwork 57131b70465SJoe Wallwork Level: beginner 57231b70465SJoe Wallwork 57331b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude() 574cb7bfe8cSJoe Wallwork @*/ 57531b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max) 57631b70465SJoe Wallwork { 57731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 57831b70465SJoe Wallwork 57931b70465SJoe Wallwork PetscFunctionBegin; 58031b70465SJoe Wallwork if (!plex->metricCtx) { 5819566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 5829566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 58331b70465SJoe Wallwork } 58431b70465SJoe Wallwork *a_max = plex->metricCtx->a_max; 58531b70465SJoe Wallwork PetscFunctionReturn(0); 58631b70465SJoe Wallwork } 58731b70465SJoe Wallwork 588cb7bfe8cSJoe Wallwork /*@ 58931b70465SJoe Wallwork DMPlexMetricSetTargetComplexity - Set the target metric complexity 59031b70465SJoe Wallwork 59131b70465SJoe Wallwork Input parameters: 59231b70465SJoe Wallwork + dm - The DM 59331b70465SJoe Wallwork - targetComplexity - The target metric complexity 59431b70465SJoe Wallwork 59531b70465SJoe Wallwork Level: beginner 59631b70465SJoe Wallwork 59731b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder() 598cb7bfe8cSJoe Wallwork @*/ 59931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity) 60031b70465SJoe Wallwork { 60131b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 60231b70465SJoe Wallwork 60331b70465SJoe Wallwork PetscFunctionBegin; 60431b70465SJoe Wallwork if (!plex->metricCtx) { 6059566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 6069566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 60731b70465SJoe Wallwork } 608*08401ef6SPierre Jolivet PetscCheck(targetComplexity > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive"); 60931b70465SJoe Wallwork plex->metricCtx->targetComplexity = targetComplexity; 61031b70465SJoe Wallwork PetscFunctionReturn(0); 61131b70465SJoe Wallwork } 61231b70465SJoe Wallwork 613cb7bfe8cSJoe Wallwork /*@ 61431b70465SJoe Wallwork DMPlexMetricGetTargetComplexity - Get the target metric complexity 61531b70465SJoe Wallwork 61631b70465SJoe Wallwork Input parameters: 61731b70465SJoe Wallwork . dm - The DM 61831b70465SJoe Wallwork 619cc2eee55SJoe Wallwork Output parameters: 62031b70465SJoe Wallwork . targetComplexity - The target metric complexity 62131b70465SJoe Wallwork 62231b70465SJoe Wallwork Level: beginner 62331b70465SJoe Wallwork 62431b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder() 625cb7bfe8cSJoe Wallwork @*/ 62631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity) 62731b70465SJoe Wallwork { 62831b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 62931b70465SJoe Wallwork 63031b70465SJoe Wallwork PetscFunctionBegin; 63131b70465SJoe Wallwork if (!plex->metricCtx) { 6329566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 6339566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 63431b70465SJoe Wallwork } 63531b70465SJoe Wallwork *targetComplexity = plex->metricCtx->targetComplexity; 63631b70465SJoe Wallwork PetscFunctionReturn(0); 63731b70465SJoe Wallwork } 63831b70465SJoe Wallwork 639cb7bfe8cSJoe Wallwork /*@ 64031b70465SJoe Wallwork DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization 64131b70465SJoe Wallwork 64231b70465SJoe Wallwork Input parameters: 64331b70465SJoe Wallwork + dm - The DM 64431b70465SJoe Wallwork - p - The normalization order 64531b70465SJoe Wallwork 64631b70465SJoe Wallwork Level: beginner 64731b70465SJoe Wallwork 64831b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity() 649cb7bfe8cSJoe Wallwork @*/ 65031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p) 65131b70465SJoe Wallwork { 65231b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 65331b70465SJoe Wallwork 65431b70465SJoe Wallwork PetscFunctionBegin; 65531b70465SJoe Wallwork if (!plex->metricCtx) { 6569566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 6579566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 65831b70465SJoe Wallwork } 659*08401ef6SPierre Jolivet PetscCheck(p >= 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater"); 66031b70465SJoe Wallwork plex->metricCtx->p = p; 66131b70465SJoe Wallwork PetscFunctionReturn(0); 66231b70465SJoe Wallwork } 66331b70465SJoe Wallwork 664cb7bfe8cSJoe Wallwork /*@ 66531b70465SJoe Wallwork DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization 66631b70465SJoe Wallwork 66731b70465SJoe Wallwork Input parameters: 66831b70465SJoe Wallwork . dm - The DM 66931b70465SJoe Wallwork 670cc2eee55SJoe Wallwork Output parameters: 67131b70465SJoe Wallwork . p - The normalization order 67231b70465SJoe Wallwork 67331b70465SJoe Wallwork Level: beginner 67431b70465SJoe Wallwork 67531b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity() 676cb7bfe8cSJoe Wallwork @*/ 67731b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p) 67831b70465SJoe Wallwork { 67931b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 68031b70465SJoe Wallwork 68131b70465SJoe Wallwork PetscFunctionBegin; 68231b70465SJoe Wallwork if (!plex->metricCtx) { 6839566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 6849566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 68531b70465SJoe Wallwork } 68631b70465SJoe Wallwork *p = plex->metricCtx->p; 68731b70465SJoe Wallwork PetscFunctionReturn(0); 68831b70465SJoe Wallwork } 68931b70465SJoe Wallwork 690cb7bfe8cSJoe Wallwork /*@ 691cc2eee55SJoe Wallwork DMPlexMetricSetGradationFactor - Set the metric gradation factor 692cc2eee55SJoe Wallwork 693cc2eee55SJoe Wallwork Input parameters: 694cc2eee55SJoe Wallwork + dm - The DM 695cc2eee55SJoe Wallwork - beta - The metric gradation factor 696cc2eee55SJoe Wallwork 697cc2eee55SJoe Wallwork Level: beginner 698cc2eee55SJoe Wallwork 699cc2eee55SJoe Wallwork Notes: 700cc2eee55SJoe Wallwork 701cc2eee55SJoe Wallwork The gradation factor is the maximum tolerated length ratio between adjacent edges. 702cc2eee55SJoe Wallwork 703cc2eee55SJoe Wallwork Turn off gradation by passing the value -1. Otherwise, pass a positive value. 704cc2eee55SJoe Wallwork 705cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 706cb7bfe8cSJoe Wallwork 707ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber() 708cb7bfe8cSJoe Wallwork @*/ 709cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta) 710cc2eee55SJoe Wallwork { 711cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 712cc2eee55SJoe Wallwork 713cc2eee55SJoe Wallwork PetscFunctionBegin; 714cc2eee55SJoe Wallwork if (!plex->metricCtx) { 7159566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 7169566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 717cc2eee55SJoe Wallwork } 718cc2eee55SJoe Wallwork plex->metricCtx->gradationFactor = beta; 719cc2eee55SJoe Wallwork PetscFunctionReturn(0); 720cc2eee55SJoe Wallwork } 721cc2eee55SJoe Wallwork 722cb7bfe8cSJoe Wallwork /*@ 723cc2eee55SJoe Wallwork DMPlexMetricGetGradationFactor - Get the metric gradation factor 724cc2eee55SJoe Wallwork 725cc2eee55SJoe Wallwork Input parameters: 726cc2eee55SJoe Wallwork . dm - The DM 727cc2eee55SJoe Wallwork 728cc2eee55SJoe Wallwork Output parameters: 729cc2eee55SJoe Wallwork . beta - The metric gradation factor 730cc2eee55SJoe Wallwork 731cc2eee55SJoe Wallwork Level: beginner 732cc2eee55SJoe Wallwork 733cb7bfe8cSJoe Wallwork Notes: 734cb7bfe8cSJoe Wallwork 735cb7bfe8cSJoe Wallwork The gradation factor is the maximum tolerated length ratio between adjacent edges. 736cb7bfe8cSJoe Wallwork 737cb7bfe8cSJoe Wallwork The value -1 implies that gradation is turned off. 738cb7bfe8cSJoe Wallwork 739cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 740cc2eee55SJoe Wallwork 741ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber() 742cb7bfe8cSJoe Wallwork @*/ 743cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta) 744cc2eee55SJoe Wallwork { 745cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 746cc2eee55SJoe Wallwork 747cc2eee55SJoe Wallwork PetscFunctionBegin; 748cc2eee55SJoe Wallwork if (!plex->metricCtx) { 7499566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 7509566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 751cc2eee55SJoe Wallwork } 752cc2eee55SJoe Wallwork *beta = plex->metricCtx->gradationFactor; 753cc2eee55SJoe Wallwork PetscFunctionReturn(0); 754cc2eee55SJoe Wallwork } 755cc2eee55SJoe Wallwork 756cb7bfe8cSJoe Wallwork /*@ 757ae8b063eSJoe Wallwork DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number 758ae8b063eSJoe Wallwork 759ae8b063eSJoe Wallwork Input parameters: 760ae8b063eSJoe Wallwork + dm - The DM 761ae8b063eSJoe Wallwork - hausd - The metric Hausdorff number 762ae8b063eSJoe Wallwork 763ae8b063eSJoe Wallwork Level: beginner 764ae8b063eSJoe Wallwork 765ae8b063eSJoe Wallwork Notes: 766ae8b063eSJoe Wallwork 767ae8b063eSJoe Wallwork The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the 768ae8b063eSJoe Wallwork boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the 769ae8b063eSJoe Wallwork high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for 770ae8b063eSJoe Wallwork an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease 771ae8b063eSJoe Wallwork (resp. increase) the Hausdorff parameter. (Taken from 772ae8b063eSJoe Wallwork https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd). 773ae8b063eSJoe Wallwork 774ae8b063eSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 775ae8b063eSJoe Wallwork 776ae8b063eSJoe Wallwork .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber() 777ae8b063eSJoe Wallwork @*/ 778ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd) 779ae8b063eSJoe Wallwork { 780ae8b063eSJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 781ae8b063eSJoe Wallwork 782ae8b063eSJoe Wallwork PetscFunctionBegin; 783ae8b063eSJoe Wallwork if (!plex->metricCtx) { 7849566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 7859566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 786ae8b063eSJoe Wallwork } 787ae8b063eSJoe Wallwork plex->metricCtx->hausdorffNumber = hausd; 788ae8b063eSJoe Wallwork PetscFunctionReturn(0); 789ae8b063eSJoe Wallwork } 790ae8b063eSJoe Wallwork 791ae8b063eSJoe Wallwork /*@ 792ae8b063eSJoe Wallwork DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number 793ae8b063eSJoe Wallwork 794ae8b063eSJoe Wallwork Input parameters: 795ae8b063eSJoe Wallwork . dm - The DM 796ae8b063eSJoe Wallwork 797ae8b063eSJoe Wallwork Output parameters: 798ae8b063eSJoe Wallwork . hausd - The metric Hausdorff number 799ae8b063eSJoe Wallwork 800ae8b063eSJoe Wallwork Level: beginner 801ae8b063eSJoe Wallwork 802ae8b063eSJoe Wallwork Notes: 803ae8b063eSJoe Wallwork 804ae8b063eSJoe Wallwork The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the 805ae8b063eSJoe Wallwork boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the 806ae8b063eSJoe Wallwork high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for 807ae8b063eSJoe Wallwork an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease 808ae8b063eSJoe Wallwork (resp. increase) the Hausdorff parameter. (Taken from 809ae8b063eSJoe Wallwork https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd). 810ae8b063eSJoe Wallwork 811ae8b063eSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 812ae8b063eSJoe Wallwork 813ae8b063eSJoe Wallwork .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber() 814ae8b063eSJoe Wallwork @*/ 815ae8b063eSJoe Wallwork PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd) 816ae8b063eSJoe Wallwork { 817ae8b063eSJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 818ae8b063eSJoe Wallwork 819ae8b063eSJoe Wallwork PetscFunctionBegin; 820ae8b063eSJoe Wallwork if (!plex->metricCtx) { 8219566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 8229566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 823ae8b063eSJoe Wallwork } 824ae8b063eSJoe Wallwork *hausd = plex->metricCtx->hausdorffNumber; 825ae8b063eSJoe Wallwork PetscFunctionReturn(0); 826ae8b063eSJoe Wallwork } 827ae8b063eSJoe Wallwork 828ae8b063eSJoe Wallwork /*@ 829cc2eee55SJoe Wallwork DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package 830cc2eee55SJoe Wallwork 831cc2eee55SJoe Wallwork Input parameters: 832cc2eee55SJoe Wallwork + dm - The DM 833cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum 834cc2eee55SJoe Wallwork 835cb7bfe8cSJoe Wallwork Level: beginner 836cb7bfe8cSJoe Wallwork 837cb7bfe8cSJoe Wallwork Notes: 838cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 839cb7bfe8cSJoe Wallwork 840cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations() 841cb7bfe8cSJoe Wallwork @*/ 842cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity) 843cc2eee55SJoe Wallwork { 844cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 845cc2eee55SJoe Wallwork 846cc2eee55SJoe Wallwork PetscFunctionBegin; 847cc2eee55SJoe Wallwork if (!plex->metricCtx) { 8489566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 8499566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 850cc2eee55SJoe Wallwork } 851cc2eee55SJoe Wallwork plex->metricCtx->verbosity = verbosity; 852cc2eee55SJoe Wallwork PetscFunctionReturn(0); 853cc2eee55SJoe Wallwork } 854cc2eee55SJoe Wallwork 855cb7bfe8cSJoe Wallwork /*@ 856cc2eee55SJoe Wallwork DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package 857cc2eee55SJoe Wallwork 858cc2eee55SJoe Wallwork Input parameters: 859cc2eee55SJoe Wallwork . dm - The DM 860cc2eee55SJoe Wallwork 861cc2eee55SJoe Wallwork Output parameters: 862cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum 863cc2eee55SJoe Wallwork 864cb7bfe8cSJoe Wallwork Level: beginner 865cb7bfe8cSJoe Wallwork 866cb7bfe8cSJoe Wallwork Notes: 867cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic). 868cb7bfe8cSJoe Wallwork 869cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations() 870cb7bfe8cSJoe Wallwork @*/ 871cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity) 872cc2eee55SJoe Wallwork { 873cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 874cc2eee55SJoe Wallwork 875cc2eee55SJoe Wallwork PetscFunctionBegin; 876cc2eee55SJoe Wallwork if (!plex->metricCtx) { 8779566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 8789566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 879cc2eee55SJoe Wallwork } 880cc2eee55SJoe Wallwork *verbosity = plex->metricCtx->verbosity; 881cc2eee55SJoe Wallwork PetscFunctionReturn(0); 882cc2eee55SJoe Wallwork } 883cc2eee55SJoe Wallwork 884cb7bfe8cSJoe Wallwork /*@ 885cc2eee55SJoe Wallwork DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations 886cc2eee55SJoe Wallwork 887cc2eee55SJoe Wallwork Input parameters: 888cc2eee55SJoe Wallwork + dm - The DM 889cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations 890cc2eee55SJoe Wallwork 891cb7bfe8cSJoe Wallwork Level: beginner 892cb7bfe8cSJoe Wallwork 893cb7bfe8cSJoe Wallwork Notes: 894cb7bfe8cSJoe Wallwork This is only used by ParMmg (not Pragmatic or Mmg). 895cc2eee55SJoe Wallwork 896cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations() 897cb7bfe8cSJoe Wallwork @*/ 898cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter) 899cc2eee55SJoe Wallwork { 900cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 901cc2eee55SJoe Wallwork 902cc2eee55SJoe Wallwork PetscFunctionBegin; 903cc2eee55SJoe Wallwork if (!plex->metricCtx) { 9049566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 9059566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 906cc2eee55SJoe Wallwork } 907cc2eee55SJoe Wallwork plex->metricCtx->numIter = numIter; 908cc2eee55SJoe Wallwork PetscFunctionReturn(0); 909cc2eee55SJoe Wallwork } 910cc2eee55SJoe Wallwork 911cb7bfe8cSJoe Wallwork /*@ 912cc2eee55SJoe Wallwork DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations 913cc2eee55SJoe Wallwork 914cc2eee55SJoe Wallwork Input parameters: 915cc2eee55SJoe Wallwork . dm - The DM 916cc2eee55SJoe Wallwork 917cc2eee55SJoe Wallwork Output parameters: 918cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations 919cc2eee55SJoe Wallwork 920cb7bfe8cSJoe Wallwork Level: beginner 921cb7bfe8cSJoe Wallwork 922cb7bfe8cSJoe Wallwork Notes: 923cb7bfe8cSJoe Wallwork This is only used by Mmg and ParMmg (not Pragmatic or Mmg). 924cc2eee55SJoe Wallwork 925cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity() 926cb7bfe8cSJoe Wallwork @*/ 927cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter) 928cc2eee55SJoe Wallwork { 929cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 930cc2eee55SJoe Wallwork 931cc2eee55SJoe Wallwork PetscFunctionBegin; 932cc2eee55SJoe Wallwork if (!plex->metricCtx) { 9339566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 9349566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 935cc2eee55SJoe Wallwork } 936cc2eee55SJoe Wallwork *numIter = plex->metricCtx->numIter; 937cc2eee55SJoe Wallwork PetscFunctionReturn(0); 938cc2eee55SJoe Wallwork } 939cc2eee55SJoe Wallwork 94020282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric) 94120282da2SJoe Wallwork { 94220282da2SJoe Wallwork MPI_Comm comm; 94320282da2SJoe Wallwork PetscFE fe; 94420282da2SJoe Wallwork PetscInt dim; 94520282da2SJoe Wallwork 94620282da2SJoe Wallwork PetscFunctionBegin; 94720282da2SJoe Wallwork 94820282da2SJoe Wallwork /* Extract metadata from dm */ 9499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 9509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 95120282da2SJoe Wallwork 95220282da2SJoe Wallwork /* Create a P1 field of the requested size */ 9539566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 9549566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe)); 9559566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 9569566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 9579566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, metric)); 95820282da2SJoe Wallwork 95920282da2SJoe Wallwork PetscFunctionReturn(0); 96020282da2SJoe Wallwork } 96120282da2SJoe Wallwork 962cb7bfe8cSJoe Wallwork /*@ 96320282da2SJoe Wallwork DMPlexMetricCreate - Create a Riemannian metric field 96420282da2SJoe Wallwork 96520282da2SJoe Wallwork Input parameters: 96620282da2SJoe Wallwork + dm - The DM 96720282da2SJoe Wallwork - f - The field number to use 96820282da2SJoe Wallwork 96920282da2SJoe Wallwork Output parameter: 97020282da2SJoe Wallwork . metric - The metric 97120282da2SJoe Wallwork 97220282da2SJoe Wallwork Level: beginner 97320282da2SJoe Wallwork 97431b70465SJoe Wallwork Notes: 97531b70465SJoe Wallwork 97631b70465SJoe Wallwork It is assumed that the DM is comprised of simplices. 97731b70465SJoe Wallwork 97831b70465SJoe Wallwork Command line options for Riemannian metrics: 97931b70465SJoe Wallwork 980cb7bfe8cSJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic? 98193520041SJoe Wallwork . -dm_plex_metric_uniform - Is the metric uniform? 982cb7bfe8cSJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 983cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 984cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 985cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max - Maximum tolerated anisotropy 986cb7bfe8cSJoe Wallwork . -dm_plex_metric_p - L-p normalization order 98767b8a455SSatish Balay - -dm_plex_metric_target_complexity - Target metric complexity 988cb7bfe8cSJoe Wallwork 989cb7bfe8cSJoe Wallwork Switching between remeshers can be achieved using 990cb7bfe8cSJoe Wallwork 99167b8a455SSatish Balay . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use 992cb7bfe8cSJoe Wallwork 993cb7bfe8cSJoe Wallwork Further options that are only relevant to Mmg and ParMmg: 994cb7bfe8cSJoe Wallwork 995cb7bfe8cSJoe Wallwork + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation 996cb7bfe8cSJoe Wallwork . -dm_plex_metric_num_iterations - Number of parallel mesh adaptation iterations for ParMmg 997cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_insert - Should node insertion/deletion be turned off? 998cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_swap - Should facet swapping be turned off? 999cb7bfe8cSJoe Wallwork . -dm_plex_metric_no_move - Should node movement be turned off? 1000cb7bfe8cSJoe Wallwork - -dm_plex_metric_verbosity - Choose a verbosity level from -1 (silent) to 10 (maximum). 100120282da2SJoe Wallwork 100220282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic() 1003cb7bfe8cSJoe Wallwork @*/ 100420282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric) 100520282da2SJoe Wallwork { 100631b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 100793520041SJoe Wallwork PetscBool isotropic, uniform; 100820282da2SJoe Wallwork PetscInt coordDim, Nd; 100920282da2SJoe Wallwork 101020282da2SJoe Wallwork PetscFunctionBegin; 101131b70465SJoe Wallwork if (!plex->metricCtx) { 10129566063dSJacob Faibussowitsch PetscCall(PetscNew(&plex->metricCtx)); 10139566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 101431b70465SJoe Wallwork } 10159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &coordDim)); 101620282da2SJoe Wallwork Nd = coordDim*coordDim; 10179566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 10189566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 101993520041SJoe Wallwork if (uniform) { 102093520041SJoe Wallwork MPI_Comm comm; 102193520041SJoe Wallwork 10229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 102328b400f6SJacob Faibussowitsch PetscCheck(isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported"); 10249566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, metric)); 10259566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*metric, 1, PETSC_DECIDE)); 10269566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(*metric)); 102793520041SJoe Wallwork } else if (isotropic) { 10289566063dSJacob Faibussowitsch PetscCall(DMPlexP1FieldCreate_Private(dm, f, 1, metric)); 102993520041SJoe Wallwork } else { 10309566063dSJacob Faibussowitsch PetscCall(DMPlexP1FieldCreate_Private(dm, f, Nd, metric)); 103193520041SJoe Wallwork } 103220282da2SJoe Wallwork PetscFunctionReturn(0); 103320282da2SJoe Wallwork } 103420282da2SJoe Wallwork 1035cb7bfe8cSJoe Wallwork /*@ 103620282da2SJoe Wallwork DMPlexMetricCreateUniform - Construct a uniform isotropic metric 103720282da2SJoe Wallwork 103820282da2SJoe Wallwork Input parameters: 103920282da2SJoe Wallwork + dm - The DM 104020282da2SJoe Wallwork . f - The field number to use 104120282da2SJoe Wallwork - alpha - Scaling parameter for the diagonal 104220282da2SJoe Wallwork 104320282da2SJoe Wallwork Output parameter: 104420282da2SJoe Wallwork . metric - The uniform metric 104520282da2SJoe Wallwork 104620282da2SJoe Wallwork Level: beginner 104720282da2SJoe Wallwork 104893520041SJoe Wallwork Note: In this case, the metric is constant in space and so is only specified for a single vertex. 104920282da2SJoe Wallwork 105020282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic() 1051cb7bfe8cSJoe Wallwork @*/ 105220282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric) 105320282da2SJoe Wallwork { 105420282da2SJoe Wallwork PetscFunctionBegin; 10559566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetUniform(dm, PETSC_TRUE)); 10569566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, f, metric)); 10572c71b3e2SJacob Faibussowitsch PetscCheck(alpha,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined"); 10582c71b3e2SJacob Faibussowitsch PetscCheck(alpha >= 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha); 10599566063dSJacob Faibussowitsch PetscCall(VecSet(*metric, alpha)); 10609566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(*metric)); 10619566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(*metric)); 106220282da2SJoe Wallwork PetscFunctionReturn(0); 106320282da2SJoe Wallwork } 106420282da2SJoe Wallwork 106593520041SJoe Wallwork static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 106693520041SJoe Wallwork const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 106793520041SJoe Wallwork const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 106893520041SJoe Wallwork PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 106993520041SJoe Wallwork { 107093520041SJoe Wallwork f0[0] = u[0]; 107193520041SJoe Wallwork } 107293520041SJoe Wallwork 1073cb7bfe8cSJoe Wallwork /*@ 107420282da2SJoe Wallwork DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator 107520282da2SJoe Wallwork 107620282da2SJoe Wallwork Input parameters: 107720282da2SJoe Wallwork + dm - The DM 107820282da2SJoe Wallwork . f - The field number to use 107920282da2SJoe Wallwork - indicator - The error indicator 108020282da2SJoe Wallwork 108120282da2SJoe Wallwork Output parameter: 108220282da2SJoe Wallwork . metric - The isotropic metric 108320282da2SJoe Wallwork 108420282da2SJoe Wallwork Level: beginner 108520282da2SJoe Wallwork 108620282da2SJoe Wallwork Notes: 108720282da2SJoe Wallwork 108820282da2SJoe Wallwork It is assumed that the DM is comprised of simplices. 108920282da2SJoe Wallwork 109093520041SJoe Wallwork The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately. 109120282da2SJoe Wallwork 109220282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform() 1093cb7bfe8cSJoe Wallwork @*/ 109420282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric) 109520282da2SJoe Wallwork { 109693520041SJoe Wallwork PetscInt m, n; 109720282da2SJoe Wallwork 109820282da2SJoe Wallwork PetscFunctionBegin; 10999566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetIsotropic(dm, PETSC_TRUE)); 11009566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, f, metric)); 11019566063dSJacob Faibussowitsch PetscCall(VecGetSize(indicator, &m)); 11029566063dSJacob Faibussowitsch PetscCall(VecGetSize(*metric, &n)); 11039566063dSJacob Faibussowitsch if (m == n) PetscCall(VecCopy(indicator, *metric)); 110493520041SJoe Wallwork else { 110593520041SJoe Wallwork DM dmIndi; 110693520041SJoe Wallwork void (*funcs[1])(PetscInt, PetscInt, PetscInt, 110793520041SJoe Wallwork const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 110893520041SJoe Wallwork const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 110993520041SJoe Wallwork PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 111093520041SJoe Wallwork 11119566063dSJacob Faibussowitsch PetscCall(VecGetDM(indicator, &dmIndi)); 111293520041SJoe Wallwork funcs[0] = identity; 11139566063dSJacob Faibussowitsch PetscCall(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric)); 111420282da2SJoe Wallwork } 111520282da2SJoe Wallwork PetscFunctionReturn(0); 111620282da2SJoe Wallwork } 111720282da2SJoe Wallwork 111882490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[]) 111982490253SJoe Wallwork { 112082490253SJoe Wallwork PetscInt i, j; 112182490253SJoe Wallwork 112282490253SJoe Wallwork PetscFunctionBegin; 112382490253SJoe Wallwork PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"); 112482490253SJoe Wallwork for (i = 0; i < dim; ++i) { 112582490253SJoe Wallwork if (i == 0) PetscPrintf(PETSC_COMM_SELF, " [["); 112682490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, " ["); 112782490253SJoe Wallwork for (j = 0; j < dim; ++j) { 112882490253SJoe Wallwork if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]); 112982490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]); 113082490253SJoe Wallwork } 113182490253SJoe Wallwork if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n"); 113282490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, "]]\n"); 113382490253SJoe Wallwork } 113482490253SJoe Wallwork PetscFunctionReturn(0); 113582490253SJoe Wallwork } 113682490253SJoe Wallwork 11373f07679eSJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp) 113820282da2SJoe Wallwork { 113920282da2SJoe Wallwork PetscInt i, j, k; 114020282da2SJoe 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); 114120282da2SJoe Wallwork PetscScalar *Mpos; 114220282da2SJoe Wallwork 114320282da2SJoe Wallwork PetscFunctionBegin; 11449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dim*dim, &Mpos, dim, &eigs)); 114520282da2SJoe Wallwork 114620282da2SJoe Wallwork /* Symmetrize */ 114720282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 114820282da2SJoe Wallwork Mpos[i*dim+i] = Mp[i*dim+i]; 114920282da2SJoe Wallwork for (j = i+1; j < dim; ++j) { 115020282da2SJoe Wallwork Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 115120282da2SJoe Wallwork Mpos[j*dim+i] = Mpos[i*dim+j]; 115220282da2SJoe Wallwork } 115320282da2SJoe Wallwork } 115420282da2SJoe Wallwork 115520282da2SJoe Wallwork /* Compute eigendecomposition */ 115693520041SJoe Wallwork if (dim == 1) { 115793520041SJoe Wallwork 115893520041SJoe Wallwork /* Isotropic case */ 115993520041SJoe Wallwork eigs[0] = PetscRealPart(Mpos[0]); 116093520041SJoe Wallwork Mpos[0] = 1.0; 116193520041SJoe Wallwork } else { 116293520041SJoe Wallwork 116393520041SJoe Wallwork /* Anisotropic case */ 116420282da2SJoe Wallwork PetscScalar *work; 116520282da2SJoe Wallwork PetscBLASInt lwork; 116620282da2SJoe Wallwork 116720282da2SJoe Wallwork lwork = 5*dim; 11689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5*dim, &work)); 116920282da2SJoe Wallwork { 117020282da2SJoe Wallwork PetscBLASInt lierr; 117120282da2SJoe Wallwork PetscBLASInt nb; 117220282da2SJoe Wallwork 11739566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(dim, &nb)); 11749566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 117520282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 117620282da2SJoe Wallwork { 117720282da2SJoe Wallwork PetscReal *rwork; 11789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*dim, &rwork)); 117920282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr)); 11809566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 118120282da2SJoe Wallwork } 118220282da2SJoe Wallwork #else 118320282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr)); 118420282da2SJoe Wallwork #endif 118582490253SJoe Wallwork if (lierr) { 118682490253SJoe Wallwork for (i = 0; i < dim; ++i) { 118782490253SJoe Wallwork Mpos[i*dim+i] = Mp[i*dim+i]; 118882490253SJoe Wallwork for (j = i+1; j < dim; ++j) { 118982490253SJoe Wallwork Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 119082490253SJoe Wallwork Mpos[j*dim+i] = Mpos[i*dim+j]; 119182490253SJoe Wallwork } 119282490253SJoe Wallwork } 119382490253SJoe Wallwork LAPACKsyevFail(dim, Mpos); 119498921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 119582490253SJoe Wallwork } 11969566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 119720282da2SJoe Wallwork } 11989566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 119920282da2SJoe Wallwork } 120020282da2SJoe Wallwork 120120282da2SJoe Wallwork /* Reflect to positive orthant and enforce maximum and minimum size */ 120220282da2SJoe Wallwork max_eig = 0.0; 120320282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 120420282da2SJoe Wallwork eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i]))); 120520282da2SJoe Wallwork max_eig = PetscMax(eigs[i], max_eig); 120620282da2SJoe Wallwork } 120720282da2SJoe Wallwork 12083f07679eSJoe Wallwork /* Enforce maximum anisotropy and compute determinant */ 12093f07679eSJoe Wallwork *detMp = 1.0; 121020282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 121120282da2SJoe Wallwork if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min); 12123f07679eSJoe Wallwork *detMp *= eigs[i]; 121320282da2SJoe Wallwork } 121420282da2SJoe Wallwork 121520282da2SJoe Wallwork /* Reconstruct Hessian */ 121620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 121720282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 121820282da2SJoe Wallwork Mp[i*dim+j] = 0.0; 121920282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 122020282da2SJoe Wallwork Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j]; 122120282da2SJoe Wallwork } 122220282da2SJoe Wallwork } 122320282da2SJoe Wallwork } 12249566063dSJacob Faibussowitsch PetscCall(PetscFree2(Mpos, eigs)); 122520282da2SJoe Wallwork 122620282da2SJoe Wallwork PetscFunctionReturn(0); 122720282da2SJoe Wallwork } 122820282da2SJoe Wallwork 1229cb7bfe8cSJoe Wallwork /*@ 123020282da2SJoe Wallwork DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric 123120282da2SJoe Wallwork 123220282da2SJoe Wallwork Input parameters: 123320282da2SJoe Wallwork + dm - The DM 12343f07679eSJoe Wallwork . metricIn - The metric 123599abec2bSJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 12363f07679eSJoe Wallwork - restrictAnisotropy - Should maximum anisotropy be enforced? 123720282da2SJoe Wallwork 123820282da2SJoe Wallwork Output parameter: 12393f07679eSJoe Wallwork + metricOut - The metric 12403f07679eSJoe Wallwork - determinant - Its determinant 124120282da2SJoe Wallwork 124220282da2SJoe Wallwork Level: beginner 124320282da2SJoe Wallwork 1244cb7bfe8cSJoe Wallwork Notes: 1245cb7bfe8cSJoe Wallwork 1246cb7bfe8cSJoe Wallwork Relevant command line options: 1247cb7bfe8cSJoe Wallwork 124893520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic? 124993520041SJoe Wallwork . -dm_plex_metric_uniform - Is the metric uniform? 125093520041SJoe Wallwork . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 1251cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 1252cb7bfe8cSJoe Wallwork - -dm_plex_metric_a_max - Maximum tolerated anisotropy 1253cb7bfe8cSJoe Wallwork 125420282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection() 1255cb7bfe8cSJoe Wallwork @*/ 12563f07679eSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant) 125720282da2SJoe Wallwork { 12583f07679eSJoe Wallwork DM dmDet; 125993520041SJoe Wallwork PetscBool isotropic, uniform; 126020282da2SJoe Wallwork PetscInt dim, vStart, vEnd, v; 12613f07679eSJoe Wallwork PetscScalar *met, *det; 126220282da2SJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0; 126320282da2SJoe Wallwork 126420282da2SJoe Wallwork PetscFunctionBegin; 12659566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0)); 126620282da2SJoe Wallwork 126720282da2SJoe Wallwork /* Extract metadata from dm */ 12689566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 126920282da2SJoe Wallwork if (restrictSizes) { 12709566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min)); 12719566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max)); 127299abec2bSJoe Wallwork h_min = PetscMax(h_min, 1.0e-30); 127399abec2bSJoe Wallwork h_max = PetscMin(h_max, 1.0e+30); 127454c59aa7SJacob Faibussowitsch PetscCheck(h_min < h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max); 127599abec2bSJoe Wallwork } 127699abec2bSJoe Wallwork if (restrictAnisotropy) { 12779566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max)); 127899abec2bSJoe Wallwork a_max = PetscMin(a_max, 1.0e+30); 127920282da2SJoe Wallwork } 128020282da2SJoe Wallwork 128193520041SJoe Wallwork /* Setup output metric */ 12829566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, 0, metricOut)); 12839566063dSJacob Faibussowitsch PetscCall(VecCopy(metricIn, *metricOut)); 128493520041SJoe Wallwork 128593520041SJoe Wallwork /* Enforce SPD and extract determinant */ 12869566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metricOut, &met)); 12879566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 12889566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 128993520041SJoe Wallwork if (uniform) { 1290d1a5b989SMatthew G. Knepley PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist"); 129193520041SJoe Wallwork 129293520041SJoe Wallwork /* Uniform case */ 12939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(metricIn, determinant)); 12949566063dSJacob Faibussowitsch PetscCall(VecGetArray(*determinant, &det)); 12959566063dSJacob Faibussowitsch PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det)); 12969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*determinant, &det)); 129793520041SJoe Wallwork } else { 129893520041SJoe Wallwork 129993520041SJoe Wallwork /* Spatially varying case */ 130093520041SJoe Wallwork PetscInt nrow; 130193520041SJoe Wallwork 130293520041SJoe Wallwork if (isotropic) nrow = 1; 130393520041SJoe Wallwork else nrow = dim; 13049566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmDet)); 13059566063dSJacob Faibussowitsch PetscCall(DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant)); 13069566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 13079566063dSJacob Faibussowitsch PetscCall(VecGetArray(*determinant, &det)); 130820282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 13093f07679eSJoe Wallwork PetscScalar *vmet, *vdet; 13109566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &vmet)); 13119566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dmDet, v, det, &vdet)); 13129566063dSJacob Faibussowitsch PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet)); 131320282da2SJoe Wallwork } 13149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*determinant, &det)); 131593520041SJoe Wallwork } 13169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metricOut, &met)); 1317fe902aceSJoe Wallwork 13189566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0)); 131920282da2SJoe Wallwork PetscFunctionReturn(0); 132020282da2SJoe Wallwork } 132120282da2SJoe Wallwork 132220282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 132320282da2SJoe Wallwork const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 132420282da2SJoe Wallwork const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 132520282da2SJoe Wallwork PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 132620282da2SJoe Wallwork { 132720282da2SJoe Wallwork const PetscScalar p = constants[0]; 132820282da2SJoe Wallwork 1329985ad86eSJose E. Roman f0[0] = PetscPowScalar(u[0], p/(2.0*p + dim)); 133020282da2SJoe Wallwork } 133120282da2SJoe Wallwork 1332cb7bfe8cSJoe Wallwork /*@ 133320282da2SJoe Wallwork DMPlexMetricNormalize - Apply L-p normalization to a metric 133420282da2SJoe Wallwork 133520282da2SJoe Wallwork Input parameters: 133620282da2SJoe Wallwork + dm - The DM 133720282da2SJoe Wallwork . metricIn - The unnormalized metric 133816522936SJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 133916522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced? 134020282da2SJoe Wallwork 134120282da2SJoe Wallwork Output parameter: 134220282da2SJoe Wallwork . metricOut - The normalized metric 134320282da2SJoe Wallwork 134420282da2SJoe Wallwork Level: beginner 134520282da2SJoe Wallwork 1346cb7bfe8cSJoe Wallwork Notes: 1347cb7bfe8cSJoe Wallwork 1348cb7bfe8cSJoe Wallwork Relevant command line options: 1349cb7bfe8cSJoe Wallwork 135093520041SJoe Wallwork + -dm_plex_metric_isotropic - Is the metric isotropic? 135193520041SJoe Wallwork . -dm_plex_metric_uniform - Is the metric uniform? 135293520041SJoe Wallwork . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 1353cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 1354cb7bfe8cSJoe Wallwork . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 1355cb7bfe8cSJoe Wallwork . -dm_plex_metric_a_max - Maximum tolerated anisotropy 1356cb7bfe8cSJoe Wallwork . -dm_plex_metric_p - L-p normalization order 1357cb7bfe8cSJoe Wallwork - -dm_plex_metric_target_complexity - Target metric complexity 1358cb7bfe8cSJoe Wallwork 135920282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection() 1360cb7bfe8cSJoe Wallwork @*/ 136116522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut) 136220282da2SJoe Wallwork { 13633f07679eSJoe Wallwork DM dmDet; 136420282da2SJoe Wallwork MPI_Comm comm; 136593520041SJoe Wallwork PetscBool restrictAnisotropyFirst, isotropic, uniform; 136620282da2SJoe Wallwork PetscDS ds; 136720282da2SJoe Wallwork PetscInt dim, Nd, vStart, vEnd, v, i; 13683f07679eSJoe Wallwork PetscScalar *met, *det, integral, constants[1]; 136993520041SJoe Wallwork PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral; 13703f07679eSJoe Wallwork Vec determinant; 137120282da2SJoe Wallwork 137220282da2SJoe Wallwork PetscFunctionBegin; 13739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0)); 137420282da2SJoe Wallwork 137520282da2SJoe Wallwork /* Extract metadata from dm */ 13769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 13779566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 13789566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 13799566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 138093520041SJoe Wallwork if (isotropic) Nd = 1; 138193520041SJoe Wallwork else Nd = dim*dim; 138220282da2SJoe Wallwork 138320282da2SJoe Wallwork /* Set up metric and ensure it is SPD */ 13849566063dSJacob Faibussowitsch PetscCall(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst)); 13859566063dSJacob Faibussowitsch PetscCall(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant)); 138620282da2SJoe Wallwork 138720282da2SJoe Wallwork /* Compute global normalization factor */ 13889566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetTargetComplexity(dm, &target)); 13899566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetNormalizationOrder(dm, &p)); 139016522936SJoe Wallwork constants[0] = p; 139193520041SJoe Wallwork if (uniform) { 139228b400f6SJacob Faibussowitsch PetscCheck(isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported"); 139393520041SJoe Wallwork DM dmTmp; 139493520041SJoe Wallwork Vec tmp; 139593520041SJoe Wallwork 13969566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmTmp)); 13979566063dSJacob Faibussowitsch PetscCall(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp)); 13989566063dSJacob Faibussowitsch PetscCall(VecGetArray(determinant, &det)); 13999566063dSJacob Faibussowitsch PetscCall(VecSet(tmp, det[0])); 14009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(determinant, &det)); 14019566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmTmp, &ds)); 14029566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants)); 14039566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, detMFunc)); 14049566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL)); 14059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 14069566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmTmp)); 140793520041SJoe Wallwork } else { 14089566063dSJacob Faibussowitsch PetscCall(VecGetDM(determinant, &dmDet)); 14099566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmDet, &ds)); 14109566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants)); 14119566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, detMFunc)); 14129566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL)); 141393520041SJoe Wallwork } 14143f07679eSJoe Wallwork realIntegral = PetscRealPart(integral); 1415*08401ef6SPierre Jolivet PetscCheck(realIntegral >= 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global metric normalization factor should be strictly positive, not %.4e Is the input metric positive-definite?", realIntegral); 14163f07679eSJoe Wallwork factGlob = PetscPowReal(target/realIntegral, 2.0/dim); 141720282da2SJoe Wallwork 141820282da2SJoe Wallwork /* Apply local scaling */ 141920282da2SJoe Wallwork if (restrictSizes) { 14209566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMinimumMagnitude(dm, &h_min)); 14219566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumMagnitude(dm, &h_max)); 142216522936SJoe Wallwork h_min = PetscMax(h_min, 1.0e-30); 142316522936SJoe Wallwork h_max = PetscMin(h_max, 1.0e+30); 1424*08401ef6SPierre Jolivet PetscCheck(h_min < h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max); 142516522936SJoe Wallwork } 142616522936SJoe Wallwork if (restrictAnisotropy && !restrictAnisotropyFirst) { 14279566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetMaximumAnisotropy(dm, &a_max)); 142816522936SJoe Wallwork a_max = PetscMin(a_max, 1.0e+30); 142920282da2SJoe Wallwork } 14309566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metricOut, &met)); 14319566063dSJacob Faibussowitsch PetscCall(VecGetArray(determinant, &det)); 143293520041SJoe Wallwork if (uniform) { 143393520041SJoe Wallwork 143493520041SJoe Wallwork /* Uniform case */ 143593520041SJoe Wallwork met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim)); 14369566063dSJacob Faibussowitsch if (restrictSizes) PetscCall(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det)); 143793520041SJoe Wallwork } else { 143893520041SJoe Wallwork 143993520041SJoe Wallwork /* Spatially varying case */ 144093520041SJoe Wallwork PetscInt nrow; 144193520041SJoe Wallwork 144293520041SJoe Wallwork if (isotropic) nrow = 1; 144393520041SJoe Wallwork else nrow = dim; 14449566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 144520282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 14463f07679eSJoe Wallwork PetscScalar *Mp, *detM; 144720282da2SJoe Wallwork 14489566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &Mp)); 14499566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dmDet, v, det, &detM)); 14503f07679eSJoe Wallwork fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim)); 145120282da2SJoe Wallwork for (i = 0; i < Nd; ++i) Mp[i] *= fact; 14529566063dSJacob Faibussowitsch if (restrictSizes) PetscCall(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM)); 145393520041SJoe Wallwork } 145420282da2SJoe Wallwork } 14559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(determinant, &det)); 14569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metricOut, &met)); 14579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&determinant)); 14589566063dSJacob Faibussowitsch if (!uniform) PetscCall(DMDestroy(&dmDet)); 145920282da2SJoe Wallwork 14609566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0)); 146120282da2SJoe Wallwork PetscFunctionReturn(0); 146220282da2SJoe Wallwork } 146320282da2SJoe Wallwork 1464cb7bfe8cSJoe Wallwork /*@ 146520282da2SJoe Wallwork DMPlexMetricAverage - Compute the average of a list of metrics 146620282da2SJoe Wallwork 1467f1a722f8SMatthew G. Knepley Input Parameters: 146820282da2SJoe Wallwork + dm - The DM 146920282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged 147020282da2SJoe Wallwork . weights - Weights for the average 147120282da2SJoe Wallwork - metrics - The metrics to be averaged 147220282da2SJoe Wallwork 147320282da2SJoe Wallwork Output Parameter: 147420282da2SJoe Wallwork . metricAvg - The averaged metric 147520282da2SJoe Wallwork 147620282da2SJoe Wallwork Level: beginner 147720282da2SJoe Wallwork 147820282da2SJoe Wallwork Notes: 147920282da2SJoe Wallwork The weights should sum to unity. 148020282da2SJoe Wallwork 148120282da2SJoe Wallwork If weights are not provided then an unweighted average is used. 148220282da2SJoe Wallwork 148320282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection() 1484cb7bfe8cSJoe Wallwork @*/ 148520282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg) 148620282da2SJoe Wallwork { 148720282da2SJoe Wallwork PetscBool haveWeights = PETSC_TRUE; 148893520041SJoe Wallwork PetscInt i, m, n; 148920282da2SJoe Wallwork PetscReal sum = 0.0, tol = 1.0e-10; 149020282da2SJoe Wallwork 149120282da2SJoe Wallwork PetscFunctionBegin; 14929566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0)); 14932c71b3e2SJacob Faibussowitsch PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); 14949566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, 0, metricAvg)); 14959566063dSJacob Faibussowitsch PetscCall(VecSet(*metricAvg, 0.0)); 14969566063dSJacob Faibussowitsch PetscCall(VecGetSize(*metricAvg, &m)); 149793520041SJoe Wallwork for (i = 0; i < numMetrics; ++i) { 14989566063dSJacob Faibussowitsch PetscCall(VecGetSize(metrics[i], &n)); 14995f80ce2aSJacob Faibussowitsch PetscCheck(m == n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented"); 150093520041SJoe Wallwork } 150120282da2SJoe Wallwork 150220282da2SJoe Wallwork /* Default to the unweighted case */ 150320282da2SJoe Wallwork if (!weights) { 15049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numMetrics, &weights)); 150520282da2SJoe Wallwork haveWeights = PETSC_FALSE; 150620282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; } 150720282da2SJoe Wallwork } 150820282da2SJoe Wallwork 150920282da2SJoe Wallwork /* Check weights sum to unity */ 151093520041SJoe Wallwork for (i = 0; i < numMetrics; ++i) sum += weights[i]; 15115f80ce2aSJacob Faibussowitsch PetscCheck(PetscAbsReal(sum - 1) <= tol,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); 151220282da2SJoe Wallwork 151320282da2SJoe Wallwork /* Compute metric average */ 15149566063dSJacob Faibussowitsch for (i = 0; i < numMetrics; ++i) PetscCall(VecAXPY(*metricAvg, weights[i], metrics[i])); 15159566063dSJacob Faibussowitsch if (!haveWeights) PetscCall(PetscFree(weights)); 1516fe902aceSJoe Wallwork 15179566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0)); 151820282da2SJoe Wallwork PetscFunctionReturn(0); 151920282da2SJoe Wallwork } 152020282da2SJoe Wallwork 1521cb7bfe8cSJoe Wallwork /*@ 152220282da2SJoe Wallwork DMPlexMetricAverage2 - Compute the unweighted average of two metrics 152320282da2SJoe Wallwork 1524f1a722f8SMatthew G. Knepley Input Parameters: 152520282da2SJoe Wallwork + dm - The DM 152620282da2SJoe Wallwork . metric1 - The first metric to be averaged 152720282da2SJoe Wallwork - metric2 - The second metric to be averaged 152820282da2SJoe Wallwork 152920282da2SJoe Wallwork Output Parameter: 153020282da2SJoe Wallwork . metricAvg - The averaged metric 153120282da2SJoe Wallwork 153220282da2SJoe Wallwork Level: beginner 153320282da2SJoe Wallwork 153420282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3() 1535cb7bfe8cSJoe Wallwork @*/ 153620282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg) 153720282da2SJoe Wallwork { 153820282da2SJoe Wallwork PetscReal weights[2] = {0.5, 0.5}; 153920282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 154020282da2SJoe Wallwork 154120282da2SJoe Wallwork PetscFunctionBegin; 15429566063dSJacob Faibussowitsch PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg)); 154320282da2SJoe Wallwork PetscFunctionReturn(0); 154420282da2SJoe Wallwork } 154520282da2SJoe Wallwork 1546cb7bfe8cSJoe Wallwork /*@ 154720282da2SJoe Wallwork DMPlexMetricAverage3 - Compute the unweighted average of three metrics 154820282da2SJoe Wallwork 1549f1a722f8SMatthew G. Knepley Input Parameters: 155020282da2SJoe Wallwork + dm - The DM 155120282da2SJoe Wallwork . metric1 - The first metric to be averaged 155220282da2SJoe Wallwork . metric2 - The second metric to be averaged 155320282da2SJoe Wallwork - metric3 - The third metric to be averaged 155420282da2SJoe Wallwork 155520282da2SJoe Wallwork Output Parameter: 155620282da2SJoe Wallwork . metricAvg - The averaged metric 155720282da2SJoe Wallwork 155820282da2SJoe Wallwork Level: beginner 155920282da2SJoe Wallwork 156020282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2() 1561cb7bfe8cSJoe Wallwork @*/ 156220282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg) 156320282da2SJoe Wallwork { 156420282da2SJoe Wallwork PetscReal weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0}; 156520282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 156620282da2SJoe Wallwork 156720282da2SJoe Wallwork PetscFunctionBegin; 15689566063dSJacob Faibussowitsch PetscCall(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg)); 156920282da2SJoe Wallwork PetscFunctionReturn(0); 157020282da2SJoe Wallwork } 157120282da2SJoe Wallwork 157220282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[]) 157320282da2SJoe Wallwork { 157420282da2SJoe Wallwork PetscInt i, j, k, l, m; 157520282da2SJoe Wallwork PetscReal *evals, *evals1; 157620282da2SJoe Wallwork PetscScalar *evecs, *sqrtM1, *isqrtM1; 157720282da2SJoe Wallwork 157820282da2SJoe Wallwork PetscFunctionBegin; 157993520041SJoe Wallwork 158093520041SJoe Wallwork /* Isotropic case */ 158193520041SJoe Wallwork if (dim == 1) { 158293520041SJoe Wallwork M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0])); 158393520041SJoe Wallwork PetscFunctionReturn(0); 158493520041SJoe Wallwork } 158593520041SJoe Wallwork 158693520041SJoe Wallwork /* Anisotropic case */ 15879566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1)); 158820282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 158920282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 159020282da2SJoe Wallwork evecs[i*dim+j] = M1[i*dim+j]; 159120282da2SJoe Wallwork } 159220282da2SJoe Wallwork } 159320282da2SJoe Wallwork { 159420282da2SJoe Wallwork PetscScalar *work; 159520282da2SJoe Wallwork PetscBLASInt lwork; 159620282da2SJoe Wallwork 159720282da2SJoe Wallwork lwork = 5*dim; 15989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5*dim, &work)); 159920282da2SJoe Wallwork { 160020282da2SJoe Wallwork PetscBLASInt lierr, nb; 160120282da2SJoe Wallwork PetscReal sqrtk; 160220282da2SJoe Wallwork 160320282da2SJoe Wallwork /* Compute eigendecomposition of M1 */ 16049566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(dim, &nb)); 16059566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 160620282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 160720282da2SJoe Wallwork { 160820282da2SJoe Wallwork PetscReal *rwork; 16099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*dim, &rwork)); 161020282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr)); 16119566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 161220282da2SJoe Wallwork } 161320282da2SJoe Wallwork #else 161420282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr)); 161520282da2SJoe Wallwork #endif 161682490253SJoe Wallwork if (lierr) { 161782490253SJoe Wallwork LAPACKsyevFail(dim, M1); 161898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 161982490253SJoe Wallwork } 16209566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 162120282da2SJoe Wallwork 162220282da2SJoe Wallwork /* Compute square root and reciprocal */ 162320282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 162420282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 162520282da2SJoe Wallwork sqrtM1[i*dim+j] = 0.0; 162620282da2SJoe Wallwork isqrtM1[i*dim+j] = 0.0; 162720282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 162820282da2SJoe Wallwork sqrtk = PetscSqrtReal(evals1[k]); 162920282da2SJoe Wallwork sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j]; 163020282da2SJoe Wallwork isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j]; 163120282da2SJoe Wallwork } 163220282da2SJoe Wallwork } 163320282da2SJoe Wallwork } 163420282da2SJoe Wallwork 163520282da2SJoe Wallwork /* Map into the space spanned by the eigenvectors of M1 */ 163620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 163720282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 163820282da2SJoe Wallwork evecs[i*dim+j] = 0.0; 163920282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 164020282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 164120282da2SJoe Wallwork evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 164220282da2SJoe Wallwork } 164320282da2SJoe Wallwork } 164420282da2SJoe Wallwork } 164520282da2SJoe Wallwork } 164620282da2SJoe Wallwork 164720282da2SJoe Wallwork /* Compute eigendecomposition */ 16489566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 164920282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 165020282da2SJoe Wallwork { 165120282da2SJoe Wallwork PetscReal *rwork; 16529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*dim, &rwork)); 165320282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 16549566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 165520282da2SJoe Wallwork } 165620282da2SJoe Wallwork #else 165720282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 165820282da2SJoe Wallwork #endif 165982490253SJoe Wallwork if (lierr) { 166082490253SJoe Wallwork for (i = 0; i < dim; ++i) { 166182490253SJoe Wallwork for (j = 0; j < dim; ++j) { 166282490253SJoe Wallwork evecs[i*dim+j] = 0.0; 166382490253SJoe Wallwork for (k = 0; k < dim; ++k) { 166482490253SJoe Wallwork for (l = 0; l < dim; ++l) { 166582490253SJoe Wallwork evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 166682490253SJoe Wallwork } 166782490253SJoe Wallwork } 166882490253SJoe Wallwork } 166982490253SJoe Wallwork } 167082490253SJoe Wallwork LAPACKsyevFail(dim, evecs); 167198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 167282490253SJoe Wallwork } 16739566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 167420282da2SJoe Wallwork 167520282da2SJoe Wallwork /* Modify eigenvalues */ 167620282da2SJoe Wallwork for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]); 167720282da2SJoe Wallwork 167820282da2SJoe Wallwork /* Map back to get the intersection */ 167920282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 168020282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 168120282da2SJoe Wallwork M2[i*dim+j] = 0.0; 168220282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 168320282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 168420282da2SJoe Wallwork for (m = 0; m < dim; ++m) { 168520282da2SJoe Wallwork M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m]; 168620282da2SJoe Wallwork } 168720282da2SJoe Wallwork } 168820282da2SJoe Wallwork } 168920282da2SJoe Wallwork } 169020282da2SJoe Wallwork } 169120282da2SJoe Wallwork } 16929566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 169320282da2SJoe Wallwork } 16949566063dSJacob Faibussowitsch PetscCall(PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1)); 169520282da2SJoe Wallwork PetscFunctionReturn(0); 169620282da2SJoe Wallwork } 169720282da2SJoe Wallwork 1698cb7bfe8cSJoe Wallwork /*@ 169920282da2SJoe Wallwork DMPlexMetricIntersection - Compute the intersection of a list of metrics 170020282da2SJoe Wallwork 1701f1a722f8SMatthew G. Knepley Input Parameters: 170220282da2SJoe Wallwork + dm - The DM 170320282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected 170420282da2SJoe Wallwork - metrics - The metrics to be intersected 170520282da2SJoe Wallwork 170620282da2SJoe Wallwork Output Parameter: 170720282da2SJoe Wallwork . metricInt - The intersected metric 170820282da2SJoe Wallwork 170920282da2SJoe Wallwork Level: beginner 171020282da2SJoe Wallwork 171120282da2SJoe Wallwork Notes: 171220282da2SJoe Wallwork 171320282da2SJoe Wallwork The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics. 171420282da2SJoe Wallwork 171520282da2SJoe Wallwork The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2. 171620282da2SJoe Wallwork 171720282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage() 1718cb7bfe8cSJoe Wallwork @*/ 171920282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt) 172020282da2SJoe Wallwork { 172193520041SJoe Wallwork PetscBool isotropic, uniform; 172293520041SJoe Wallwork PetscInt v, i, m, n; 172393520041SJoe Wallwork PetscScalar *met, *meti; 172420282da2SJoe Wallwork 172520282da2SJoe Wallwork PetscFunctionBegin; 17269566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0)); 17272c71b3e2SJacob Faibussowitsch PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); 172820282da2SJoe Wallwork 172920282da2SJoe Wallwork /* Copy over the first metric */ 17309566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, 0, metricInt)); 17319566063dSJacob Faibussowitsch PetscCall(VecCopy(metrics[0], *metricInt)); 173293520041SJoe Wallwork if (numMetrics == 1) PetscFunctionReturn(0); 17339566063dSJacob Faibussowitsch PetscCall(VecGetSize(*metricInt, &m)); 173493520041SJoe Wallwork for (i = 0; i < numMetrics; ++i) { 17359566063dSJacob Faibussowitsch PetscCall(VecGetSize(metrics[i], &n)); 1736*08401ef6SPierre Jolivet PetscCheck(m == n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented"); 173793520041SJoe Wallwork } 173820282da2SJoe Wallwork 173920282da2SJoe Wallwork /* Intersect subsequent metrics in turn */ 17409566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 17419566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 174293520041SJoe Wallwork if (uniform) { 174393520041SJoe Wallwork 174493520041SJoe Wallwork /* Uniform case */ 17459566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metricInt, &met)); 174693520041SJoe Wallwork for (i = 1; i < numMetrics; ++i) { 17479566063dSJacob Faibussowitsch PetscCall(VecGetArray(metrics[i], &meti)); 17489566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection_Private(1, meti, met)); 17499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(metrics[i], &meti)); 175093520041SJoe Wallwork } 17519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metricInt, &met)); 175293520041SJoe Wallwork } else { 175393520041SJoe Wallwork 175493520041SJoe Wallwork /* Spatially varying case */ 175593520041SJoe Wallwork PetscInt dim, vStart, vEnd, nrow; 175693520041SJoe Wallwork PetscScalar *M, *Mi; 175793520041SJoe Wallwork 17589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 175993520041SJoe Wallwork if (isotropic) nrow = 1; 176093520041SJoe Wallwork else nrow = dim; 17619566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 17629566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metricInt, &met)); 176320282da2SJoe Wallwork for (i = 1; i < numMetrics; ++i) { 17649566063dSJacob Faibussowitsch PetscCall(VecGetArray(metrics[i], &meti)); 176520282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 17669566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &M)); 17679566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, meti, &Mi)); 17689566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection_Private(nrow, Mi, M)); 176920282da2SJoe Wallwork } 17709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(metrics[i], &meti)); 177120282da2SJoe Wallwork } 17729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metricInt, &met)); 177320282da2SJoe Wallwork } 1774fe902aceSJoe Wallwork 17759566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0)); 177620282da2SJoe Wallwork PetscFunctionReturn(0); 177720282da2SJoe Wallwork } 177820282da2SJoe Wallwork 1779cb7bfe8cSJoe Wallwork /*@ 178020282da2SJoe Wallwork DMPlexMetricIntersection2 - Compute the intersection of two metrics 178120282da2SJoe Wallwork 1782f1a722f8SMatthew G. Knepley Input Parameters: 178320282da2SJoe Wallwork + dm - The DM 178420282da2SJoe Wallwork . metric1 - The first metric to be intersected 178520282da2SJoe Wallwork - metric2 - The second metric to be intersected 178620282da2SJoe Wallwork 178720282da2SJoe Wallwork Output Parameter: 178820282da2SJoe Wallwork . metricInt - The intersected metric 178920282da2SJoe Wallwork 179020282da2SJoe Wallwork Level: beginner 179120282da2SJoe Wallwork 179220282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3() 1793cb7bfe8cSJoe Wallwork @*/ 179420282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt) 179520282da2SJoe Wallwork { 179620282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 179720282da2SJoe Wallwork 179820282da2SJoe Wallwork PetscFunctionBegin; 17999566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricInt)); 180020282da2SJoe Wallwork PetscFunctionReturn(0); 180120282da2SJoe Wallwork } 180220282da2SJoe Wallwork 1803cb7bfe8cSJoe Wallwork /*@ 180420282da2SJoe Wallwork DMPlexMetricIntersection3 - Compute the intersection of three metrics 180520282da2SJoe Wallwork 1806f1a722f8SMatthew G. Knepley Input Parameters: 180720282da2SJoe Wallwork + dm - The DM 180820282da2SJoe Wallwork . metric1 - The first metric to be intersected 180920282da2SJoe Wallwork . metric2 - The second metric to be intersected 181020282da2SJoe Wallwork - metric3 - The third metric to be intersected 181120282da2SJoe Wallwork 181220282da2SJoe Wallwork Output Parameter: 181320282da2SJoe Wallwork . metricInt - The intersected metric 181420282da2SJoe Wallwork 181520282da2SJoe Wallwork Level: beginner 181620282da2SJoe Wallwork 181720282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2() 1818cb7bfe8cSJoe Wallwork @*/ 181920282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt) 182020282da2SJoe Wallwork { 182120282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 182220282da2SJoe Wallwork 182320282da2SJoe Wallwork PetscFunctionBegin; 18249566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIntersection(dm, 3, metrics, metricInt)); 182520282da2SJoe Wallwork PetscFunctionReturn(0); 182620282da2SJoe Wallwork } 1827