120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 220282da2SJoe Wallwork #include <petscblaslapack.h> 320282da2SJoe Wallwork 431b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm) 531b70465SJoe Wallwork { 631b70465SJoe Wallwork MPI_Comm comm; 731b70465SJoe Wallwork PetscBool isotropic = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE; 8cc2eee55SJoe Wallwork PetscBool noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE; 931b70465SJoe Wallwork PetscErrorCode ierr; 10cc2eee55SJoe Wallwork PetscInt verbosity = -1, numIter = 3; 11cc2eee55SJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0, beta = 1.3; 1231b70465SJoe Wallwork 1331b70465SJoe Wallwork PetscFunctionBegin; 1431b70465SJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 1531b70465SJoe Wallwork ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr); 1631b70465SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr); 17cc2eee55SJoe Wallwork ierr = DMPlexMetricSetIsotropic(dm, isotropic);CHKERRQ(ierr); 1831b70465SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr); 19cc2eee55SJoe Wallwork ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst);CHKERRQ(ierr); 20cc2eee55SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL);CHKERRQ(ierr); 21cc2eee55SJoe Wallwork ierr = DMPlexMetricSetNoInsertion(dm, noInsert);CHKERRQ(ierr); 22cc2eee55SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL);CHKERRQ(ierr); 23cc2eee55SJoe Wallwork ierr = DMPlexMetricSetNoSwapping(dm, noSwap);CHKERRQ(ierr); 24cc2eee55SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL);CHKERRQ(ierr); 25cc2eee55SJoe Wallwork ierr = DMPlexMetricSetNoMovement(dm, noMove);CHKERRQ(ierr); 26cc2eee55SJoe Wallwork ierr = PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0);CHKERRQ(ierr); 27cc2eee55SJoe Wallwork ierr = DMPlexMetricSetNumIterations(dm, numIter);CHKERRQ(ierr); 28cc2eee55SJoe Wallwork ierr = PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10);CHKERRQ(ierr); 29cc2eee55SJoe Wallwork ierr = DMPlexMetricSetVerbosity(dm, verbosity);CHKERRQ(ierr); 3031b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr); 3131b70465SJoe Wallwork ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr); 3231b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr); 3331b70465SJoe Wallwork ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr); 3431b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr); 3531b70465SJoe Wallwork ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr); 3631b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr); 3731b70465SJoe Wallwork ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr); 3831b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr); 3931b70465SJoe Wallwork ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr); 40cc2eee55SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL);CHKERRQ(ierr); 41cc2eee55SJoe Wallwork ierr = DMPlexMetricSetGradationFactor(dm, beta);CHKERRQ(ierr); 4231b70465SJoe Wallwork ierr = PetscOptionsEnd();CHKERRQ(ierr); 4331b70465SJoe Wallwork PetscFunctionReturn(0); 4431b70465SJoe Wallwork } 4531b70465SJoe Wallwork 4631b70465SJoe Wallwork /* 4731b70465SJoe Wallwork DMPlexMetricSetIsotropic - Record whether a metric is isotropic 4831b70465SJoe Wallwork 4931b70465SJoe Wallwork Input parameters: 5031b70465SJoe Wallwork + dm - The DM 5131b70465SJoe Wallwork - isotropic - Is the metric isotropic? 5231b70465SJoe Wallwork 5331b70465SJoe Wallwork Level: beginner 5431b70465SJoe Wallwork 5531b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 5631b70465SJoe Wallwork */ 5731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic) 5831b70465SJoe Wallwork { 5931b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 6031b70465SJoe Wallwork PetscErrorCode ierr; 6131b70465SJoe Wallwork 6231b70465SJoe Wallwork PetscFunctionBegin; 6331b70465SJoe Wallwork if (!plex->metricCtx) { 6431b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 6531b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 6631b70465SJoe Wallwork } 6731b70465SJoe Wallwork plex->metricCtx->isotropic = isotropic; 6831b70465SJoe Wallwork PetscFunctionReturn(0); 6931b70465SJoe Wallwork } 7031b70465SJoe Wallwork 7131b70465SJoe Wallwork /* 7231b70465SJoe Wallwork DMPlexMetricIsIsotropic - Is a metric is isotropic? 7331b70465SJoe Wallwork 7431b70465SJoe Wallwork Input parameters: 7531b70465SJoe Wallwork . dm - The DM 7631b70465SJoe Wallwork 7731b70465SJoe Wallwork Output parameters: 7831b70465SJoe Wallwork . isotropic - Is the metric isotropic? 7931b70465SJoe Wallwork 8031b70465SJoe Wallwork Level: beginner 8131b70465SJoe Wallwork 8231b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 8331b70465SJoe Wallwork */ 8431b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic) 8531b70465SJoe Wallwork { 8631b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 8731b70465SJoe Wallwork PetscErrorCode ierr; 8831b70465SJoe Wallwork 8931b70465SJoe Wallwork PetscFunctionBegin; 9031b70465SJoe Wallwork if (!plex->metricCtx) { 9131b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 9231b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 9331b70465SJoe Wallwork } 9431b70465SJoe Wallwork *isotropic = plex->metricCtx->isotropic; 9531b70465SJoe Wallwork PetscFunctionReturn(0); 9631b70465SJoe Wallwork } 9731b70465SJoe Wallwork 9831b70465SJoe Wallwork /* 9931b70465SJoe Wallwork DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization 10031b70465SJoe Wallwork 10131b70465SJoe Wallwork Input parameters: 10231b70465SJoe Wallwork + dm - The DM 10331b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first? 10431b70465SJoe Wallwork 10531b70465SJoe Wallwork Level: beginner 10631b70465SJoe Wallwork 10731b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 10831b70465SJoe Wallwork */ 10931b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst) 11031b70465SJoe Wallwork { 11131b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 11231b70465SJoe Wallwork PetscErrorCode ierr; 11331b70465SJoe Wallwork 11431b70465SJoe Wallwork PetscFunctionBegin; 11531b70465SJoe Wallwork if (!plex->metricCtx) { 11631b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 11731b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 11831b70465SJoe Wallwork } 11931b70465SJoe Wallwork plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst; 12031b70465SJoe Wallwork PetscFunctionReturn(0); 12131b70465SJoe Wallwork } 12231b70465SJoe Wallwork 12331b70465SJoe Wallwork /* 12431b70465SJoe Wallwork DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after? 12531b70465SJoe Wallwork 12631b70465SJoe Wallwork Input parameters: 12731b70465SJoe Wallwork . dm - The DM 12831b70465SJoe Wallwork 12931b70465SJoe Wallwork Output parameters: 13031b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first? 13131b70465SJoe Wallwork 13231b70465SJoe Wallwork Level: beginner 13331b70465SJoe Wallwork 13431b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 13531b70465SJoe Wallwork */ 13631b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst) 13731b70465SJoe Wallwork { 13831b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 13931b70465SJoe Wallwork PetscErrorCode ierr; 14031b70465SJoe Wallwork 14131b70465SJoe Wallwork PetscFunctionBegin; 14231b70465SJoe Wallwork if (!plex->metricCtx) { 14331b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 14431b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 14531b70465SJoe Wallwork } 14631b70465SJoe Wallwork *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst; 14731b70465SJoe Wallwork PetscFunctionReturn(0); 14831b70465SJoe Wallwork } 14931b70465SJoe Wallwork 15031b70465SJoe Wallwork /* 151cc2eee55SJoe Wallwork DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off? 152cc2eee55SJoe Wallwork 153cc2eee55SJoe Wallwork Input parameters: 154cc2eee55SJoe Wallwork + dm - The DM 155cc2eee55SJoe Wallwork - noInsert - Should node insertion and deletion be turned off? 156cc2eee55SJoe Wallwork 157cc2eee55SJoe Wallwork Level: beginner 158cc2eee55SJoe Wallwork 159cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement() 160cc2eee55SJoe Wallwork */ 161cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert) 162cc2eee55SJoe Wallwork { 163cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 164cc2eee55SJoe Wallwork PetscErrorCode ierr; 165cc2eee55SJoe Wallwork 166cc2eee55SJoe Wallwork PetscFunctionBegin; 167cc2eee55SJoe Wallwork if (!plex->metricCtx) { 168cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 169cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 170cc2eee55SJoe Wallwork } 171cc2eee55SJoe Wallwork plex->metricCtx->noInsert = noInsert; 172cc2eee55SJoe Wallwork PetscFunctionReturn(0); 173cc2eee55SJoe Wallwork } 174cc2eee55SJoe Wallwork 175cc2eee55SJoe Wallwork /* 176cc2eee55SJoe Wallwork DMPlexMetricNoInsertion - Are node insertion and deletion turned off? 177cc2eee55SJoe Wallwork 178cc2eee55SJoe Wallwork Input parameters: 179cc2eee55SJoe Wallwork . dm - The DM 180cc2eee55SJoe Wallwork 181cc2eee55SJoe Wallwork Output parameters: 182cc2eee55SJoe Wallwork . noInsert - Are node insertion and deletion turned off? 183cc2eee55SJoe Wallwork 184cc2eee55SJoe Wallwork Level: beginner 185cc2eee55SJoe Wallwork 186cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement() 187cc2eee55SJoe Wallwork */ 188cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert) 189cc2eee55SJoe Wallwork { 190cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 191cc2eee55SJoe Wallwork PetscErrorCode ierr; 192cc2eee55SJoe Wallwork 193cc2eee55SJoe Wallwork PetscFunctionBegin; 194cc2eee55SJoe Wallwork if (!plex->metricCtx) { 195cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 196cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 197cc2eee55SJoe Wallwork } 198cc2eee55SJoe Wallwork *noInsert = plex->metricCtx->noInsert; 199cc2eee55SJoe Wallwork PetscFunctionReturn(0); 200cc2eee55SJoe Wallwork } 201cc2eee55SJoe Wallwork 202cc2eee55SJoe Wallwork /* 203cc2eee55SJoe Wallwork DMPlexMetricSetNoSwapping - Should facet swapping be turned off? 204cc2eee55SJoe Wallwork 205cc2eee55SJoe Wallwork Input parameters: 206cc2eee55SJoe Wallwork + dm - The DM 207cc2eee55SJoe Wallwork - noSwap - Should facet swapping be turned off? 208cc2eee55SJoe Wallwork 209cc2eee55SJoe Wallwork Level: beginner 210cc2eee55SJoe Wallwork 211cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement() 212cc2eee55SJoe Wallwork */ 213cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap) 214cc2eee55SJoe Wallwork { 215cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 216cc2eee55SJoe Wallwork PetscErrorCode ierr; 217cc2eee55SJoe Wallwork 218cc2eee55SJoe Wallwork PetscFunctionBegin; 219cc2eee55SJoe Wallwork if (!plex->metricCtx) { 220cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 221cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 222cc2eee55SJoe Wallwork } 223cc2eee55SJoe Wallwork plex->metricCtx->noSwap = noSwap; 224cc2eee55SJoe Wallwork PetscFunctionReturn(0); 225cc2eee55SJoe Wallwork } 226cc2eee55SJoe Wallwork 227cc2eee55SJoe Wallwork /* 228cc2eee55SJoe Wallwork DMPlexMetricNoSwapping - Is facet swapping turned off? 229cc2eee55SJoe Wallwork 230cc2eee55SJoe Wallwork Input parameters: 231cc2eee55SJoe Wallwork . dm - The DM 232cc2eee55SJoe Wallwork 233cc2eee55SJoe Wallwork Output parameters: 234cc2eee55SJoe Wallwork . noSwap - Is facet swapping turned off? 235cc2eee55SJoe Wallwork 236cc2eee55SJoe Wallwork Level: beginner 237cc2eee55SJoe Wallwork 238cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement() 239cc2eee55SJoe Wallwork */ 240cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap) 241cc2eee55SJoe Wallwork { 242cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 243cc2eee55SJoe Wallwork PetscErrorCode ierr; 244cc2eee55SJoe Wallwork 245cc2eee55SJoe Wallwork PetscFunctionBegin; 246cc2eee55SJoe Wallwork if (!plex->metricCtx) { 247cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 248cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 249cc2eee55SJoe Wallwork } 250cc2eee55SJoe Wallwork *noSwap = plex->metricCtx->noSwap; 251cc2eee55SJoe Wallwork PetscFunctionReturn(0); 252cc2eee55SJoe Wallwork } 253cc2eee55SJoe Wallwork 254cc2eee55SJoe Wallwork /* 255cc2eee55SJoe Wallwork DMPlexMetricSetNoMovement - Should node movement be turned off? 256cc2eee55SJoe Wallwork 257cc2eee55SJoe Wallwork Input parameters: 258cc2eee55SJoe Wallwork + dm - The DM 259cc2eee55SJoe Wallwork - noMove - Should node movement be turned off? 260cc2eee55SJoe Wallwork 261cc2eee55SJoe Wallwork Level: beginner 262cc2eee55SJoe Wallwork 263cc2eee55SJoe Wallwork .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping() 264cc2eee55SJoe Wallwork */ 265cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove) 266cc2eee55SJoe Wallwork { 267cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 268cc2eee55SJoe Wallwork PetscErrorCode ierr; 269cc2eee55SJoe Wallwork 270cc2eee55SJoe Wallwork PetscFunctionBegin; 271cc2eee55SJoe Wallwork if (!plex->metricCtx) { 272cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 273cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 274cc2eee55SJoe Wallwork } 275cc2eee55SJoe Wallwork plex->metricCtx->noMove = noMove; 276cc2eee55SJoe Wallwork PetscFunctionReturn(0); 277cc2eee55SJoe Wallwork } 278cc2eee55SJoe Wallwork 279cc2eee55SJoe Wallwork /* 280cc2eee55SJoe Wallwork DMPlexMetricNoMovement - Is node movement turned off? 281cc2eee55SJoe Wallwork 282cc2eee55SJoe Wallwork Input parameters: 283cc2eee55SJoe Wallwork . dm - The DM 284cc2eee55SJoe Wallwork 285cc2eee55SJoe Wallwork Output parameters: 286cc2eee55SJoe Wallwork . noMove - Is node movement turned off? 287cc2eee55SJoe Wallwork 288cc2eee55SJoe Wallwork Level: beginner 289cc2eee55SJoe Wallwork 290cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping() 291cc2eee55SJoe Wallwork */ 292cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove) 293cc2eee55SJoe Wallwork { 294cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 295cc2eee55SJoe Wallwork PetscErrorCode ierr; 296cc2eee55SJoe Wallwork 297cc2eee55SJoe Wallwork PetscFunctionBegin; 298cc2eee55SJoe Wallwork if (!plex->metricCtx) { 299cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 300cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 301cc2eee55SJoe Wallwork } 302cc2eee55SJoe Wallwork *noMove = plex->metricCtx->noMove; 303cc2eee55SJoe Wallwork PetscFunctionReturn(0); 304cc2eee55SJoe Wallwork } 305cc2eee55SJoe Wallwork 306cc2eee55SJoe Wallwork /* 30731b70465SJoe Wallwork DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude 30831b70465SJoe Wallwork 30931b70465SJoe Wallwork Input parameters: 31031b70465SJoe Wallwork + dm - The DM 31131b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude 31231b70465SJoe Wallwork 31331b70465SJoe Wallwork Level: beginner 31431b70465SJoe Wallwork 31531b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude() 31631b70465SJoe Wallwork */ 31731b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min) 31831b70465SJoe Wallwork { 31931b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 32031b70465SJoe Wallwork PetscErrorCode ierr; 32131b70465SJoe Wallwork 32231b70465SJoe Wallwork PetscFunctionBegin; 32331b70465SJoe Wallwork if (!plex->metricCtx) { 32431b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 32531b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 32631b70465SJoe Wallwork } 32731b70465SJoe Wallwork if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min); 32831b70465SJoe Wallwork plex->metricCtx->h_min = h_min; 32931b70465SJoe Wallwork PetscFunctionReturn(0); 33031b70465SJoe Wallwork } 33131b70465SJoe Wallwork 33231b70465SJoe Wallwork /* 33331b70465SJoe Wallwork DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude 33431b70465SJoe Wallwork 33531b70465SJoe Wallwork Input parameters: 33631b70465SJoe Wallwork . dm - The DM 33731b70465SJoe Wallwork 338cc2eee55SJoe Wallwork Output parameters: 33931b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude 34031b70465SJoe Wallwork 34131b70465SJoe Wallwork Level: beginner 34231b70465SJoe Wallwork 34331b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude() 34431b70465SJoe Wallwork */ 34531b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min) 34631b70465SJoe Wallwork { 34731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 34831b70465SJoe Wallwork PetscErrorCode ierr; 34931b70465SJoe Wallwork 35031b70465SJoe Wallwork PetscFunctionBegin; 35131b70465SJoe Wallwork if (!plex->metricCtx) { 35231b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 35331b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 35431b70465SJoe Wallwork } 35531b70465SJoe Wallwork *h_min = plex->metricCtx->h_min; 35631b70465SJoe Wallwork PetscFunctionReturn(0); 35731b70465SJoe Wallwork } 35831b70465SJoe Wallwork 35931b70465SJoe Wallwork /* 36031b70465SJoe Wallwork DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude 36131b70465SJoe Wallwork 36231b70465SJoe Wallwork Input parameters: 36331b70465SJoe Wallwork + dm - The DM 36431b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude 36531b70465SJoe Wallwork 36631b70465SJoe Wallwork Level: beginner 36731b70465SJoe Wallwork 36831b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude() 36931b70465SJoe Wallwork */ 37031b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max) 37131b70465SJoe Wallwork { 37231b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 37331b70465SJoe Wallwork PetscErrorCode ierr; 37431b70465SJoe Wallwork 37531b70465SJoe Wallwork PetscFunctionBegin; 37631b70465SJoe Wallwork if (!plex->metricCtx) { 37731b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 37831b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 37931b70465SJoe Wallwork } 38031b70465SJoe Wallwork if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max); 38131b70465SJoe Wallwork plex->metricCtx->h_max = h_max; 38231b70465SJoe Wallwork PetscFunctionReturn(0); 38331b70465SJoe Wallwork } 38431b70465SJoe Wallwork 38531b70465SJoe Wallwork /* 38631b70465SJoe Wallwork DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude 38731b70465SJoe Wallwork 38831b70465SJoe Wallwork Input parameters: 38931b70465SJoe Wallwork . dm - The DM 39031b70465SJoe Wallwork 391cc2eee55SJoe Wallwork Output parameters: 39231b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude 39331b70465SJoe Wallwork 39431b70465SJoe Wallwork Level: beginner 39531b70465SJoe Wallwork 39631b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude() 39731b70465SJoe Wallwork */ 39831b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max) 39931b70465SJoe Wallwork { 40031b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 40131b70465SJoe Wallwork PetscErrorCode ierr; 40231b70465SJoe Wallwork 40331b70465SJoe Wallwork PetscFunctionBegin; 40431b70465SJoe Wallwork if (!plex->metricCtx) { 40531b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 40631b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 40731b70465SJoe Wallwork } 40831b70465SJoe Wallwork *h_max = plex->metricCtx->h_max; 40931b70465SJoe Wallwork PetscFunctionReturn(0); 41031b70465SJoe Wallwork } 41131b70465SJoe Wallwork 41231b70465SJoe Wallwork /* 41331b70465SJoe Wallwork DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy 41431b70465SJoe Wallwork 41531b70465SJoe Wallwork Input parameters: 41631b70465SJoe Wallwork + dm - The DM 41731b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy 41831b70465SJoe Wallwork 41931b70465SJoe Wallwork Level: beginner 42031b70465SJoe Wallwork 42131b70465SJoe Wallwork Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one. 42231b70465SJoe Wallwork 42331b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude() 42431b70465SJoe Wallwork */ 42531b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max) 42631b70465SJoe Wallwork { 42731b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 42831b70465SJoe Wallwork PetscErrorCode ierr; 42931b70465SJoe Wallwork 43031b70465SJoe Wallwork PetscFunctionBegin; 43131b70465SJoe Wallwork if (!plex->metricCtx) { 43231b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 43331b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 43431b70465SJoe Wallwork } 43531b70465SJoe Wallwork if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max); 43631b70465SJoe Wallwork plex->metricCtx->a_max = a_max; 43731b70465SJoe Wallwork PetscFunctionReturn(0); 43831b70465SJoe Wallwork } 43931b70465SJoe Wallwork 44031b70465SJoe Wallwork /* 44131b70465SJoe Wallwork DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy 44231b70465SJoe Wallwork 44331b70465SJoe Wallwork Input parameters: 44431b70465SJoe Wallwork . dm - The DM 44531b70465SJoe Wallwork 446cc2eee55SJoe Wallwork Output parameters: 44731b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy 44831b70465SJoe Wallwork 44931b70465SJoe Wallwork Level: beginner 45031b70465SJoe Wallwork 45131b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude() 45231b70465SJoe Wallwork */ 45331b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max) 45431b70465SJoe Wallwork { 45531b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 45631b70465SJoe Wallwork PetscErrorCode ierr; 45731b70465SJoe Wallwork 45831b70465SJoe Wallwork PetscFunctionBegin; 45931b70465SJoe Wallwork if (!plex->metricCtx) { 46031b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 46131b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 46231b70465SJoe Wallwork } 46331b70465SJoe Wallwork *a_max = plex->metricCtx->a_max; 46431b70465SJoe Wallwork PetscFunctionReturn(0); 46531b70465SJoe Wallwork } 46631b70465SJoe Wallwork 46731b70465SJoe Wallwork /* 46831b70465SJoe Wallwork DMPlexMetricSetTargetComplexity - Set the target metric complexity 46931b70465SJoe Wallwork 47031b70465SJoe Wallwork Input parameters: 47131b70465SJoe Wallwork + dm - The DM 47231b70465SJoe Wallwork - targetComplexity - The target metric complexity 47331b70465SJoe Wallwork 47431b70465SJoe Wallwork Level: beginner 47531b70465SJoe Wallwork 47631b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder() 47731b70465SJoe Wallwork */ 47831b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity) 47931b70465SJoe Wallwork { 48031b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 48131b70465SJoe Wallwork PetscErrorCode ierr; 48231b70465SJoe Wallwork 48331b70465SJoe Wallwork PetscFunctionBegin; 48431b70465SJoe Wallwork if (!plex->metricCtx) { 48531b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 48631b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 48731b70465SJoe Wallwork } 48831b70465SJoe Wallwork if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive"); 48931b70465SJoe Wallwork plex->metricCtx->targetComplexity = targetComplexity; 49031b70465SJoe Wallwork PetscFunctionReturn(0); 49131b70465SJoe Wallwork } 49231b70465SJoe Wallwork 49331b70465SJoe Wallwork /* 49431b70465SJoe Wallwork DMPlexMetricGetTargetComplexity - Get the target metric complexity 49531b70465SJoe Wallwork 49631b70465SJoe Wallwork Input parameters: 49731b70465SJoe Wallwork . dm - The DM 49831b70465SJoe Wallwork 499cc2eee55SJoe Wallwork Output parameters: 50031b70465SJoe Wallwork . targetComplexity - The target metric complexity 50131b70465SJoe Wallwork 50231b70465SJoe Wallwork Level: beginner 50331b70465SJoe Wallwork 50431b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder() 50531b70465SJoe Wallwork */ 50631b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity) 50731b70465SJoe Wallwork { 50831b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 50931b70465SJoe Wallwork PetscErrorCode ierr; 51031b70465SJoe Wallwork 51131b70465SJoe Wallwork PetscFunctionBegin; 51231b70465SJoe Wallwork if (!plex->metricCtx) { 51331b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 51431b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 51531b70465SJoe Wallwork } 51631b70465SJoe Wallwork *targetComplexity = plex->metricCtx->targetComplexity; 51731b70465SJoe Wallwork PetscFunctionReturn(0); 51831b70465SJoe Wallwork } 51931b70465SJoe Wallwork 52031b70465SJoe Wallwork /* 52131b70465SJoe Wallwork DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization 52231b70465SJoe Wallwork 52331b70465SJoe Wallwork Input parameters: 52431b70465SJoe Wallwork + dm - The DM 52531b70465SJoe Wallwork - p - The normalization order 52631b70465SJoe Wallwork 52731b70465SJoe Wallwork Level: beginner 52831b70465SJoe Wallwork 52931b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity() 53031b70465SJoe Wallwork */ 53131b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p) 53231b70465SJoe Wallwork { 53331b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 53431b70465SJoe Wallwork PetscErrorCode ierr; 53531b70465SJoe Wallwork 53631b70465SJoe Wallwork PetscFunctionBegin; 53731b70465SJoe Wallwork if (!plex->metricCtx) { 53831b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 53931b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 54031b70465SJoe Wallwork } 54131b70465SJoe Wallwork if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater"); 54231b70465SJoe Wallwork plex->metricCtx->p = p; 54331b70465SJoe Wallwork PetscFunctionReturn(0); 54431b70465SJoe Wallwork } 54531b70465SJoe Wallwork 54631b70465SJoe Wallwork /* 54731b70465SJoe Wallwork DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization 54831b70465SJoe Wallwork 54931b70465SJoe Wallwork Input parameters: 55031b70465SJoe Wallwork . dm - The DM 55131b70465SJoe Wallwork 552cc2eee55SJoe Wallwork Output parameters: 55331b70465SJoe Wallwork . p - The normalization order 55431b70465SJoe Wallwork 55531b70465SJoe Wallwork Level: beginner 55631b70465SJoe Wallwork 55731b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity() 55831b70465SJoe Wallwork */ 55931b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p) 56031b70465SJoe Wallwork { 56131b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 56231b70465SJoe Wallwork PetscErrorCode ierr; 56331b70465SJoe Wallwork 56431b70465SJoe Wallwork PetscFunctionBegin; 56531b70465SJoe Wallwork if (!plex->metricCtx) { 56631b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 56731b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 56831b70465SJoe Wallwork } 56931b70465SJoe Wallwork *p = plex->metricCtx->p; 57031b70465SJoe Wallwork PetscFunctionReturn(0); 57131b70465SJoe Wallwork } 57231b70465SJoe Wallwork 573cc2eee55SJoe Wallwork /* 574cc2eee55SJoe Wallwork DMPlexMetricSetGradationFactor - Set the metric gradation factor 575cc2eee55SJoe Wallwork 576cc2eee55SJoe Wallwork Input parameters: 577cc2eee55SJoe Wallwork + dm - The DM 578cc2eee55SJoe Wallwork - beta - The metric gradation factor 579cc2eee55SJoe Wallwork 580cc2eee55SJoe Wallwork Level: beginner 581cc2eee55SJoe Wallwork 582cc2eee55SJoe Wallwork Notes: 583cc2eee55SJoe Wallwork 584cc2eee55SJoe Wallwork The gradation factor is the maximum tolerated length ratio between adjacent edges. 585cc2eee55SJoe Wallwork 586cc2eee55SJoe Wallwork Turn off gradation by passing the value -1. Otherwise, pass a positive value. 587cc2eee55SJoe Wallwork 588cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetGradationFactor() 589cc2eee55SJoe Wallwork */ 590cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta) 591cc2eee55SJoe Wallwork { 592cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 593cc2eee55SJoe Wallwork PetscErrorCode ierr; 594cc2eee55SJoe Wallwork 595cc2eee55SJoe Wallwork PetscFunctionBegin; 596cc2eee55SJoe Wallwork if (!plex->metricCtx) { 597cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 598cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 599cc2eee55SJoe Wallwork } 600cc2eee55SJoe Wallwork plex->metricCtx->gradationFactor = beta; 601cc2eee55SJoe Wallwork PetscFunctionReturn(0); 602cc2eee55SJoe Wallwork } 603cc2eee55SJoe Wallwork 604cc2eee55SJoe Wallwork /* 605cc2eee55SJoe Wallwork DMPlexMetricGetGradationFactor - Get the metric gradation factor 606cc2eee55SJoe Wallwork 607cc2eee55SJoe Wallwork Input parameters: 608cc2eee55SJoe Wallwork . dm - The DM 609cc2eee55SJoe Wallwork 610cc2eee55SJoe Wallwork Output parameters: 611cc2eee55SJoe Wallwork . beta - The metric gradation factor 612cc2eee55SJoe Wallwork 613cc2eee55SJoe Wallwork Level: beginner 614cc2eee55SJoe Wallwork 615cc2eee55SJoe Wallwork Note: The gradation factor is the maximum tolerated length ratio between adjacent edges. 616cc2eee55SJoe Wallwork 617cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetGradationFactor() 618cc2eee55SJoe Wallwork */ 619cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta) 620cc2eee55SJoe Wallwork { 621cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 622cc2eee55SJoe Wallwork PetscErrorCode ierr; 623cc2eee55SJoe Wallwork 624cc2eee55SJoe Wallwork PetscFunctionBegin; 625cc2eee55SJoe Wallwork if (!plex->metricCtx) { 626cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 627cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 628cc2eee55SJoe Wallwork } 629cc2eee55SJoe Wallwork *beta = plex->metricCtx->gradationFactor; 630cc2eee55SJoe Wallwork PetscFunctionReturn(0); 631cc2eee55SJoe Wallwork } 632cc2eee55SJoe Wallwork 633cc2eee55SJoe Wallwork /* 634cc2eee55SJoe Wallwork DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package 635cc2eee55SJoe Wallwork 636cc2eee55SJoe Wallwork Input parameters: 637cc2eee55SJoe Wallwork + dm - The DM 638cc2eee55SJoe Wallwork - verbosity - The verbosity, where -1 is silent and 10 is maximum 639cc2eee55SJoe Wallwork 640cc2eee55SJoe Wallwork .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations() 641cc2eee55SJoe Wallwork */ 642cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity) 643cc2eee55SJoe Wallwork { 644cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 645cc2eee55SJoe Wallwork PetscErrorCode ierr; 646cc2eee55SJoe Wallwork 647cc2eee55SJoe Wallwork PetscFunctionBegin; 648cc2eee55SJoe Wallwork if (!plex->metricCtx) { 649cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 650cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 651cc2eee55SJoe Wallwork } 652cc2eee55SJoe Wallwork plex->metricCtx->verbosity = verbosity; 653cc2eee55SJoe Wallwork PetscFunctionReturn(0); 654cc2eee55SJoe Wallwork } 655cc2eee55SJoe Wallwork 656cc2eee55SJoe Wallwork /* 657cc2eee55SJoe Wallwork DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package 658cc2eee55SJoe Wallwork 659cc2eee55SJoe Wallwork Input parameters: 660cc2eee55SJoe Wallwork . dm - The DM 661cc2eee55SJoe Wallwork 662cc2eee55SJoe Wallwork Output parameters: 663cc2eee55SJoe Wallwork . verbosity - The verbosity, where -1 is silent and 10 is maximum 664cc2eee55SJoe Wallwork 665cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations() 666cc2eee55SJoe Wallwork */ 667cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity) 668cc2eee55SJoe Wallwork { 669cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 670cc2eee55SJoe Wallwork PetscErrorCode ierr; 671cc2eee55SJoe Wallwork 672cc2eee55SJoe Wallwork PetscFunctionBegin; 673cc2eee55SJoe Wallwork if (!plex->metricCtx) { 674cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 675cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 676cc2eee55SJoe Wallwork } 677cc2eee55SJoe Wallwork *verbosity = plex->metricCtx->verbosity; 678cc2eee55SJoe Wallwork PetscFunctionReturn(0); 679cc2eee55SJoe Wallwork } 680cc2eee55SJoe Wallwork 681cc2eee55SJoe Wallwork /* 682cc2eee55SJoe Wallwork DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations 683cc2eee55SJoe Wallwork 684cc2eee55SJoe Wallwork Input parameters: 685cc2eee55SJoe Wallwork + dm - The DM 686cc2eee55SJoe Wallwork - numIter - the number of parallel adaptation iterations 687cc2eee55SJoe Wallwork 688cc2eee55SJoe Wallwork Note: This option is only used by ParMmg, not Mmg or Pragmatic. 689cc2eee55SJoe Wallwork 690cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations() 691cc2eee55SJoe Wallwork */ 692cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter) 693cc2eee55SJoe Wallwork { 694cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 695cc2eee55SJoe Wallwork PetscErrorCode ierr; 696cc2eee55SJoe Wallwork 697cc2eee55SJoe Wallwork PetscFunctionBegin; 698cc2eee55SJoe Wallwork if (!plex->metricCtx) { 699cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 700cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 701cc2eee55SJoe Wallwork } 702cc2eee55SJoe Wallwork plex->metricCtx->numIter = numIter; 703cc2eee55SJoe Wallwork PetscFunctionReturn(0); 704cc2eee55SJoe Wallwork } 705cc2eee55SJoe Wallwork 706cc2eee55SJoe Wallwork /* 707cc2eee55SJoe Wallwork DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations 708cc2eee55SJoe Wallwork 709cc2eee55SJoe Wallwork Input parameters: 710cc2eee55SJoe Wallwork . dm - The DM 711cc2eee55SJoe Wallwork 712cc2eee55SJoe Wallwork Output parameters: 713cc2eee55SJoe Wallwork . numIter - the number of parallel adaptation iterations 714cc2eee55SJoe Wallwork 715cc2eee55SJoe Wallwork Note: This option is only used by ParMmg, not Mmg or Pragmatic. 716cc2eee55SJoe Wallwork 717cc2eee55SJoe Wallwork .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity() 718cc2eee55SJoe Wallwork */ 719cc2eee55SJoe Wallwork PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter) 720cc2eee55SJoe Wallwork { 721cc2eee55SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 722cc2eee55SJoe Wallwork PetscErrorCode ierr; 723cc2eee55SJoe Wallwork 724cc2eee55SJoe Wallwork PetscFunctionBegin; 725cc2eee55SJoe Wallwork if (!plex->metricCtx) { 726cc2eee55SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 727cc2eee55SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 728cc2eee55SJoe Wallwork } 729cc2eee55SJoe Wallwork *numIter = plex->metricCtx->numIter; 730cc2eee55SJoe Wallwork PetscFunctionReturn(0); 731cc2eee55SJoe Wallwork } 732cc2eee55SJoe Wallwork 73320282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric) 73420282da2SJoe Wallwork { 73520282da2SJoe Wallwork MPI_Comm comm; 73620282da2SJoe Wallwork PetscErrorCode ierr; 73720282da2SJoe Wallwork PetscFE fe; 73820282da2SJoe Wallwork PetscInt dim; 73920282da2SJoe Wallwork 74020282da2SJoe Wallwork PetscFunctionBegin; 74120282da2SJoe Wallwork 74220282da2SJoe Wallwork /* Extract metadata from dm */ 74320282da2SJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 74420282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 74520282da2SJoe Wallwork 74620282da2SJoe Wallwork /* Create a P1 field of the requested size */ 74720282da2SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 74820282da2SJoe Wallwork ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr); 74920282da2SJoe Wallwork ierr = DMCreateDS(dm);CHKERRQ(ierr); 75020282da2SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 75120282da2SJoe Wallwork ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr); 75220282da2SJoe Wallwork 75320282da2SJoe Wallwork PetscFunctionReturn(0); 75420282da2SJoe Wallwork } 75520282da2SJoe Wallwork 75620282da2SJoe Wallwork /* 75720282da2SJoe Wallwork DMPlexMetricCreate - Create a Riemannian metric field 75820282da2SJoe Wallwork 75920282da2SJoe Wallwork Input parameters: 76020282da2SJoe Wallwork + dm - The DM 76120282da2SJoe Wallwork - f - The field number to use 76220282da2SJoe Wallwork 76320282da2SJoe Wallwork Output parameter: 76420282da2SJoe Wallwork . metric - The metric 76520282da2SJoe Wallwork 76620282da2SJoe Wallwork Level: beginner 76720282da2SJoe Wallwork 76831b70465SJoe Wallwork Notes: 76931b70465SJoe Wallwork 77031b70465SJoe Wallwork It is assumed that the DM is comprised of simplices. 77131b70465SJoe Wallwork 77231b70465SJoe Wallwork Command line options for Riemannian metrics: 77331b70465SJoe Wallwork 77431b70465SJoe Wallwork -dm_plex_metric_isotropic - Is the metric isotropic? 77531b70465SJoe Wallwork -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 77631b70465SJoe Wallwork -dm_plex_metric_h_min - Minimum tolerated metric magnitude 77731b70465SJoe Wallwork -dm_plex_metric_h_max - Maximum tolerated metric magnitude 77831b70465SJoe Wallwork -dm_plex_metric_a_max - Maximum tolerated anisotropy 77931b70465SJoe Wallwork -dm_plex_metric_p - L-p normalization order 78031b70465SJoe Wallwork -dm_plex_metric_target_complexity - Target metric complexity 78120282da2SJoe Wallwork 78220282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic() 78320282da2SJoe Wallwork */ 78420282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric) 78520282da2SJoe Wallwork { 78631b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 78720282da2SJoe Wallwork PetscErrorCode ierr; 78820282da2SJoe Wallwork PetscInt coordDim, Nd; 78920282da2SJoe Wallwork 79020282da2SJoe Wallwork PetscFunctionBegin; 79131b70465SJoe Wallwork if (!plex->metricCtx) { 79231b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 79331b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 79431b70465SJoe Wallwork } 79520282da2SJoe Wallwork ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 79620282da2SJoe Wallwork Nd = coordDim*coordDim; 79720282da2SJoe Wallwork ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr); 79820282da2SJoe Wallwork PetscFunctionReturn(0); 79920282da2SJoe Wallwork } 80020282da2SJoe Wallwork 80120282da2SJoe Wallwork typedef struct { 80220282da2SJoe Wallwork PetscReal scaling; /* Scaling for uniform metric diagonal */ 80320282da2SJoe Wallwork } DMPlexMetricUniformCtx; 80420282da2SJoe Wallwork 80520282da2SJoe Wallwork static PetscErrorCode diagonal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 80620282da2SJoe Wallwork { 80720282da2SJoe Wallwork DMPlexMetricUniformCtx *user = (DMPlexMetricUniformCtx*)ctx; 80820282da2SJoe Wallwork PetscInt i, j; 80920282da2SJoe Wallwork 81020282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 81120282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 81220282da2SJoe Wallwork if (i == j) u[i+dim*j] = user->scaling; 81320282da2SJoe Wallwork else u[i+dim*j] = 0.0; 81420282da2SJoe Wallwork } 81520282da2SJoe Wallwork } 81620282da2SJoe Wallwork return 0; 81720282da2SJoe Wallwork } 81820282da2SJoe Wallwork 81920282da2SJoe Wallwork /* 82020282da2SJoe Wallwork DMPlexMetricCreateUniform - Construct a uniform isotropic metric 82120282da2SJoe Wallwork 82220282da2SJoe Wallwork Input parameters: 82320282da2SJoe Wallwork + dm - The DM 82420282da2SJoe Wallwork . f - The field number to use 82520282da2SJoe Wallwork - alpha - Scaling parameter for the diagonal 82620282da2SJoe Wallwork 82720282da2SJoe Wallwork Output parameter: 82820282da2SJoe Wallwork . metric - The uniform metric 82920282da2SJoe Wallwork 83020282da2SJoe Wallwork Level: beginner 83120282da2SJoe Wallwork 83220282da2SJoe Wallwork Note: It is assumed that the DM is comprised of simplices. 83320282da2SJoe Wallwork 83420282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic() 83520282da2SJoe Wallwork */ 83620282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric) 83720282da2SJoe Wallwork { 83820282da2SJoe Wallwork DMPlexMetricUniformCtx user; 83920282da2SJoe Wallwork PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 84020282da2SJoe Wallwork PetscErrorCode ierr; 84120282da2SJoe Wallwork void *ctxs[1]; 84220282da2SJoe Wallwork 84320282da2SJoe Wallwork PetscFunctionBegin; 84420282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr); 84520282da2SJoe Wallwork if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined"); 84620282da2SJoe Wallwork if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha); 84720282da2SJoe Wallwork else user.scaling = alpha; 84820282da2SJoe Wallwork funcs[0] = diagonal; 84920282da2SJoe Wallwork ctxs[0] = &user; 85020282da2SJoe Wallwork ierr = DMProjectFunctionLocal(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, *metric);CHKERRQ(ierr); 85120282da2SJoe Wallwork PetscFunctionReturn(0); 85220282da2SJoe Wallwork } 85320282da2SJoe Wallwork 85420282da2SJoe Wallwork /* 85520282da2SJoe Wallwork DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator 85620282da2SJoe Wallwork 85720282da2SJoe Wallwork Input parameters: 85820282da2SJoe Wallwork + dm - The DM 85920282da2SJoe Wallwork . f - The field number to use 86020282da2SJoe Wallwork - indicator - The error indicator 86120282da2SJoe Wallwork 86220282da2SJoe Wallwork Output parameter: 86320282da2SJoe Wallwork . metric - The isotropic metric 86420282da2SJoe Wallwork 86520282da2SJoe Wallwork Level: beginner 86620282da2SJoe Wallwork 86720282da2SJoe Wallwork Notes: 86820282da2SJoe Wallwork 86920282da2SJoe Wallwork It is assumed that the DM is comprised of simplices. 87020282da2SJoe Wallwork 87120282da2SJoe Wallwork The indicator needs to be a scalar field defined at *vertices*. 87220282da2SJoe Wallwork 87320282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform() 87420282da2SJoe Wallwork */ 87520282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric) 87620282da2SJoe Wallwork { 87720282da2SJoe Wallwork DM dmIndi; 87820282da2SJoe Wallwork PetscErrorCode ierr; 87920282da2SJoe Wallwork PetscInt dim, d, vStart, vEnd, v; 88020282da2SJoe Wallwork const PetscScalar *indi; 88120282da2SJoe Wallwork PetscScalar *met; 88220282da2SJoe Wallwork 88320282da2SJoe Wallwork PetscFunctionBegin; 88420282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 88520282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr); 88620282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 88720282da2SJoe Wallwork ierr = VecGetArrayRead(indicator, &indi);CHKERRQ(ierr); 88820282da2SJoe Wallwork ierr = VecGetArrayWrite(*metric, &met);CHKERRQ(ierr); 88920282da2SJoe Wallwork ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr); 89020282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 89120282da2SJoe Wallwork PetscScalar *vindi, *vmet; 89220282da2SJoe Wallwork ierr = DMPlexPointLocalRead(dmIndi, v, indi, &vindi);CHKERRQ(ierr); 89320282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr); 89420282da2SJoe Wallwork for (d = 0; d < dim; ++d) vmet[d*(dim+1)] = vindi[0]; 89520282da2SJoe Wallwork } 89620282da2SJoe Wallwork ierr = VecRestoreArrayWrite(*metric, &met);CHKERRQ(ierr); 89720282da2SJoe Wallwork ierr = VecRestoreArrayRead(indicator, &indi);CHKERRQ(ierr); 89820282da2SJoe Wallwork PetscFunctionReturn(0); 89920282da2SJoe Wallwork } 90020282da2SJoe Wallwork 901*82490253SJoe Wallwork static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[]) 902*82490253SJoe Wallwork { 903*82490253SJoe Wallwork PetscInt i, j; 904*82490253SJoe Wallwork 905*82490253SJoe Wallwork PetscFunctionBegin; 906*82490253SJoe Wallwork PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"); 907*82490253SJoe Wallwork for (i = 0; i < dim; ++i) { 908*82490253SJoe Wallwork if (i == 0) PetscPrintf(PETSC_COMM_SELF, " [["); 909*82490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, " ["); 910*82490253SJoe Wallwork for (j = 0; j < dim; ++j) { 911*82490253SJoe Wallwork if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]); 912*82490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]); 913*82490253SJoe Wallwork } 914*82490253SJoe Wallwork if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n"); 915*82490253SJoe Wallwork else PetscPrintf(PETSC_COMM_SELF, "]]\n"); 916*82490253SJoe Wallwork } 917*82490253SJoe Wallwork PetscFunctionReturn(0); 918*82490253SJoe Wallwork } 919*82490253SJoe Wallwork 92020282da2SJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[]) 92120282da2SJoe Wallwork { 92220282da2SJoe Wallwork PetscErrorCode ierr; 92320282da2SJoe Wallwork PetscInt i, j, k; 92420282da2SJoe 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); 92520282da2SJoe Wallwork PetscScalar *Mpos; 92620282da2SJoe Wallwork 92720282da2SJoe Wallwork PetscFunctionBegin; 92820282da2SJoe Wallwork ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr); 92920282da2SJoe Wallwork 93020282da2SJoe Wallwork /* Symmetrize */ 93120282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 93220282da2SJoe Wallwork Mpos[i*dim+i] = Mp[i*dim+i]; 93320282da2SJoe Wallwork for (j = i+1; j < dim; ++j) { 93420282da2SJoe Wallwork Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 93520282da2SJoe Wallwork Mpos[j*dim+i] = Mpos[i*dim+j]; 93620282da2SJoe Wallwork } 93720282da2SJoe Wallwork } 93820282da2SJoe Wallwork 93920282da2SJoe Wallwork /* Compute eigendecomposition */ 94020282da2SJoe Wallwork { 94120282da2SJoe Wallwork PetscScalar *work; 94220282da2SJoe Wallwork PetscBLASInt lwork; 94320282da2SJoe Wallwork 94420282da2SJoe Wallwork lwork = 5*dim; 94520282da2SJoe Wallwork ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr); 94620282da2SJoe Wallwork { 94720282da2SJoe Wallwork PetscBLASInt lierr; 94820282da2SJoe Wallwork PetscBLASInt nb; 94920282da2SJoe Wallwork 95020282da2SJoe Wallwork ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr); 95120282da2SJoe Wallwork ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 95220282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 95320282da2SJoe Wallwork { 95420282da2SJoe Wallwork PetscReal *rwork; 95520282da2SJoe Wallwork ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 95620282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr)); 95720282da2SJoe Wallwork ierr = PetscFree(rwork);CHKERRQ(ierr); 95820282da2SJoe Wallwork } 95920282da2SJoe Wallwork #else 96020282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr)); 96120282da2SJoe Wallwork #endif 962*82490253SJoe Wallwork if (lierr) { 963*82490253SJoe Wallwork for (i = 0; i < dim; ++i) { 964*82490253SJoe Wallwork Mpos[i*dim+i] = Mp[i*dim+i]; 965*82490253SJoe Wallwork for (j = i+1; j < dim; ++j) { 966*82490253SJoe Wallwork Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 967*82490253SJoe Wallwork Mpos[j*dim+i] = Mpos[i*dim+j]; 968*82490253SJoe Wallwork } 969*82490253SJoe Wallwork } 970*82490253SJoe Wallwork LAPACKsyevFail(dim, Mpos); 971*82490253SJoe Wallwork SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 972*82490253SJoe Wallwork } 97320282da2SJoe Wallwork ierr = PetscFPTrapPop();CHKERRQ(ierr); 97420282da2SJoe Wallwork } 97520282da2SJoe Wallwork ierr = PetscFree(work);CHKERRQ(ierr); 97620282da2SJoe Wallwork } 97720282da2SJoe Wallwork 97820282da2SJoe Wallwork /* Reflect to positive orthant and enforce maximum and minimum size */ 97920282da2SJoe Wallwork max_eig = 0.0; 98020282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 98120282da2SJoe Wallwork eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i]))); 98220282da2SJoe Wallwork max_eig = PetscMax(eigs[i], max_eig); 98320282da2SJoe Wallwork } 98420282da2SJoe Wallwork 98520282da2SJoe Wallwork /* Enforce maximum anisotropy */ 98620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 98720282da2SJoe Wallwork if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min); 98820282da2SJoe Wallwork } 98920282da2SJoe Wallwork 99020282da2SJoe Wallwork /* Reconstruct Hessian */ 99120282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 99220282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 99320282da2SJoe Wallwork Mp[i*dim+j] = 0.0; 99420282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 99520282da2SJoe Wallwork Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j]; 99620282da2SJoe Wallwork } 99720282da2SJoe Wallwork } 99820282da2SJoe Wallwork } 99920282da2SJoe Wallwork ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr); 100020282da2SJoe Wallwork 100120282da2SJoe Wallwork PetscFunctionReturn(0); 100220282da2SJoe Wallwork } 100320282da2SJoe Wallwork 100420282da2SJoe Wallwork /* 100520282da2SJoe Wallwork DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric 100620282da2SJoe Wallwork 100720282da2SJoe Wallwork Input parameters: 100820282da2SJoe Wallwork + dm - The DM 100999abec2bSJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 101099abec2bSJoe Wallwork . restrictAnisotropy - Should maximum anisotropy be enforced? 101120282da2SJoe Wallwork - metric - The metric 101220282da2SJoe Wallwork 101320282da2SJoe Wallwork Output parameter: 101420282da2SJoe Wallwork . metric - The metric 101520282da2SJoe Wallwork 101620282da2SJoe Wallwork Level: beginner 101720282da2SJoe Wallwork 101820282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection() 101920282da2SJoe Wallwork */ 102099abec2bSJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metric) 102120282da2SJoe Wallwork { 102220282da2SJoe Wallwork PetscErrorCode ierr; 102320282da2SJoe Wallwork PetscInt dim, vStart, vEnd, v; 102420282da2SJoe Wallwork PetscScalar *met; 102520282da2SJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0; 102620282da2SJoe Wallwork 102720282da2SJoe Wallwork PetscFunctionBegin; 102820282da2SJoe Wallwork 102920282da2SJoe Wallwork /* Extract metadata from dm */ 103020282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 103120282da2SJoe Wallwork if (restrictSizes) { 103299abec2bSJoe Wallwork ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr); 103399abec2bSJoe Wallwork ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr); 103499abec2bSJoe Wallwork h_min = PetscMax(h_min, 1.0e-30); 103599abec2bSJoe Wallwork h_max = PetscMin(h_max, 1.0e+30); 103699abec2bSJoe Wallwork if (h_min >= h_max) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max); 103799abec2bSJoe Wallwork } 103899abec2bSJoe Wallwork if (restrictAnisotropy) { 103999abec2bSJoe Wallwork ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr); 104099abec2bSJoe Wallwork a_max = PetscMin(a_max, 1.0e+30); 104120282da2SJoe Wallwork } 104220282da2SJoe Wallwork 104320282da2SJoe Wallwork /* Enforce SPD */ 104420282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 104520282da2SJoe Wallwork ierr = VecGetArray(metric, &met);CHKERRQ(ierr); 104620282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 104720282da2SJoe Wallwork PetscScalar *vmet; 104820282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr); 104920282da2SJoe Wallwork ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, vmet);CHKERRQ(ierr); 105020282da2SJoe Wallwork } 105120282da2SJoe Wallwork ierr = VecRestoreArray(metric, &met);CHKERRQ(ierr); 105220282da2SJoe Wallwork 105320282da2SJoe Wallwork PetscFunctionReturn(0); 105420282da2SJoe Wallwork } 105520282da2SJoe Wallwork 105620282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 105720282da2SJoe Wallwork const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 105820282da2SJoe Wallwork const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 105920282da2SJoe Wallwork PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 106020282da2SJoe Wallwork { 106120282da2SJoe Wallwork const PetscScalar p = constants[0]; 106220282da2SJoe Wallwork PetscReal detH = 0.0; 106320282da2SJoe Wallwork 106420282da2SJoe Wallwork if (dim == 2) DMPlex_Det2D_Scalar_Internal(&detH, u); 106520282da2SJoe Wallwork else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detH, u); 106620282da2SJoe Wallwork f0[0] = PetscPowReal(detH, p/(2.0*p + dim)); 106720282da2SJoe Wallwork } 106820282da2SJoe Wallwork 106920282da2SJoe Wallwork /* 107020282da2SJoe Wallwork DMPlexMetricNormalize - Apply L-p normalization to a metric 107120282da2SJoe Wallwork 107220282da2SJoe Wallwork Input parameters: 107320282da2SJoe Wallwork + dm - The DM 107420282da2SJoe Wallwork . metricIn - The unnormalized metric 107516522936SJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 107616522936SJoe Wallwork - restrictAnisotropy - Should maximum metric anisotropy be enforced? 107720282da2SJoe Wallwork 107820282da2SJoe Wallwork Output parameter: 107920282da2SJoe Wallwork . metricOut - The normalized metric 108020282da2SJoe Wallwork 108120282da2SJoe Wallwork Level: beginner 108220282da2SJoe Wallwork 108320282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection() 108420282da2SJoe Wallwork */ 108516522936SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut) 108620282da2SJoe Wallwork { 108720282da2SJoe Wallwork MPI_Comm comm; 108816522936SJoe Wallwork PetscBool restrictAnisotropyFirst; 108920282da2SJoe Wallwork PetscDS ds; 109020282da2SJoe Wallwork PetscErrorCode ierr; 109120282da2SJoe Wallwork PetscInt dim, Nd, vStart, vEnd, v, i; 109220282da2SJoe Wallwork PetscScalar *met, integral, constants[1]; 109320282da2SJoe Wallwork PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, target; 109420282da2SJoe Wallwork 109520282da2SJoe Wallwork PetscFunctionBegin; 109620282da2SJoe Wallwork 109720282da2SJoe Wallwork /* Extract metadata from dm */ 109820282da2SJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 109920282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 110020282da2SJoe Wallwork Nd = dim*dim; 110120282da2SJoe Wallwork 110220282da2SJoe Wallwork /* Set up metric and ensure it is SPD */ 110320282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr); 110420282da2SJoe Wallwork ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr); 110516522936SJoe Wallwork ierr = DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst);CHKERRQ(ierr); 110616522936SJoe Wallwork ierr = DMPlexMetricEnforceSPD(dm, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), *metricOut);CHKERRQ(ierr); 110720282da2SJoe Wallwork 110820282da2SJoe Wallwork /* Compute global normalization factor */ 110916522936SJoe Wallwork ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr); 111016522936SJoe Wallwork ierr = DMPlexMetricGetNormalizationOrder(dm, &p);CHKERRQ(ierr); 111116522936SJoe Wallwork constants[0] = p; 111220282da2SJoe Wallwork ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 111320282da2SJoe Wallwork ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr); 111420282da2SJoe Wallwork ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr); 111520282da2SJoe Wallwork ierr = DMPlexComputeIntegralFEM(dm, *metricOut, &integral, NULL);CHKERRQ(ierr); 111620282da2SJoe Wallwork factGlob = PetscPowReal(target/PetscRealPart(integral), 2.0/dim); 111720282da2SJoe Wallwork 111820282da2SJoe Wallwork /* Apply local scaling */ 111920282da2SJoe Wallwork if (restrictSizes) { 112016522936SJoe Wallwork ierr = DMPlexMetricGetMinimumMagnitude(dm, &h_min);CHKERRQ(ierr); 112116522936SJoe Wallwork ierr = DMPlexMetricGetMaximumMagnitude(dm, &h_max);CHKERRQ(ierr); 112216522936SJoe Wallwork h_min = PetscMax(h_min, 1.0e-30); 112316522936SJoe Wallwork h_max = PetscMin(h_max, 1.0e+30); 112416522936SJoe Wallwork if (h_min >= h_max) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max); 112516522936SJoe Wallwork } 112616522936SJoe Wallwork if (restrictAnisotropy && !restrictAnisotropyFirst) { 112716522936SJoe Wallwork ierr = DMPlexMetricGetMaximumAnisotropy(dm, &a_max);CHKERRQ(ierr); 112816522936SJoe Wallwork a_max = PetscMin(a_max, 1.0e+30); 112920282da2SJoe Wallwork } 113020282da2SJoe Wallwork ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr); 113116522936SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 113220282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 113320282da2SJoe Wallwork PetscScalar *Mp; 113420282da2SJoe Wallwork PetscReal detM, fact; 113520282da2SJoe Wallwork 113620282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr); 113720282da2SJoe Wallwork if (dim == 2) DMPlex_Det2D_Scalar_Internal(&detM, Mp); 113820282da2SJoe Wallwork else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detM, Mp); 113920282da2SJoe Wallwork else SETERRQ1(comm, PETSC_ERR_SUP, "Dimension %d not supported", dim); 114020282da2SJoe Wallwork fact = factGlob * PetscPowReal(PetscAbsReal(detM), -1.0/(2*p+dim)); 114120282da2SJoe Wallwork for (i = 0; i < Nd; ++i) Mp[i] *= fact; 114220282da2SJoe Wallwork if (restrictSizes) { ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, Mp);CHKERRQ(ierr); } 114320282da2SJoe Wallwork } 114420282da2SJoe Wallwork ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr); 114520282da2SJoe Wallwork 114620282da2SJoe Wallwork PetscFunctionReturn(0); 114720282da2SJoe Wallwork } 114820282da2SJoe Wallwork 114920282da2SJoe Wallwork /* 115020282da2SJoe Wallwork DMPlexMetricAverage - Compute the average of a list of metrics 115120282da2SJoe Wallwork 115220282da2SJoe Wallwork Input Parameter: 115320282da2SJoe Wallwork + dm - The DM 115420282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged 115520282da2SJoe Wallwork . weights - Weights for the average 115620282da2SJoe Wallwork - metrics - The metrics to be averaged 115720282da2SJoe Wallwork 115820282da2SJoe Wallwork Output Parameter: 115920282da2SJoe Wallwork . metricAvg - The averaged metric 116020282da2SJoe Wallwork 116120282da2SJoe Wallwork Level: beginner 116220282da2SJoe Wallwork 116320282da2SJoe Wallwork Notes: 116420282da2SJoe Wallwork The weights should sum to unity. 116520282da2SJoe Wallwork 116620282da2SJoe Wallwork If weights are not provided then an unweighted average is used. 116720282da2SJoe Wallwork 116820282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection() 116920282da2SJoe Wallwork */ 117020282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg) 117120282da2SJoe Wallwork { 117220282da2SJoe Wallwork PetscBool haveWeights = PETSC_TRUE; 117320282da2SJoe Wallwork PetscErrorCode ierr; 117420282da2SJoe Wallwork PetscInt i; 117520282da2SJoe Wallwork PetscReal sum = 0.0, tol = 1.0e-10; 117620282da2SJoe Wallwork 117720282da2SJoe Wallwork PetscFunctionBegin; 117820282da2SJoe Wallwork if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); } 117920282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr); 118020282da2SJoe Wallwork ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr); 118120282da2SJoe Wallwork 118220282da2SJoe Wallwork /* Default to the unweighted case */ 118320282da2SJoe Wallwork if (!weights) { 118420282da2SJoe Wallwork ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr); 118520282da2SJoe Wallwork haveWeights = PETSC_FALSE; 118620282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; } 118720282da2SJoe Wallwork } 118820282da2SJoe Wallwork 118920282da2SJoe Wallwork /* Check weights sum to unity */ 119020282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) { sum += weights[i]; } 119120282da2SJoe Wallwork if (PetscAbsReal(sum - 1) > tol) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); } 119220282da2SJoe Wallwork 119320282da2SJoe Wallwork /* Compute metric average */ 119420282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); } 119520282da2SJoe Wallwork if (!haveWeights) {ierr = PetscFree(weights); } 119620282da2SJoe Wallwork PetscFunctionReturn(0); 119720282da2SJoe Wallwork } 119820282da2SJoe Wallwork 119920282da2SJoe Wallwork /* 120020282da2SJoe Wallwork DMPlexMetricAverage2 - Compute the unweighted average of two metrics 120120282da2SJoe Wallwork 120220282da2SJoe Wallwork Input Parameter: 120320282da2SJoe Wallwork + dm - The DM 120420282da2SJoe Wallwork . metric1 - The first metric to be averaged 120520282da2SJoe Wallwork - metric2 - The second metric to be averaged 120620282da2SJoe Wallwork 120720282da2SJoe Wallwork Output Parameter: 120820282da2SJoe Wallwork . metricAvg - The averaged metric 120920282da2SJoe Wallwork 121020282da2SJoe Wallwork Level: beginner 121120282da2SJoe Wallwork 121220282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3() 121320282da2SJoe Wallwork */ 121420282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg) 121520282da2SJoe Wallwork { 121620282da2SJoe Wallwork PetscErrorCode ierr; 121720282da2SJoe Wallwork PetscReal weights[2] = {0.5, 0.5}; 121820282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 121920282da2SJoe Wallwork 122020282da2SJoe Wallwork PetscFunctionBegin; 122120282da2SJoe Wallwork ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr); 122220282da2SJoe Wallwork PetscFunctionReturn(0); 122320282da2SJoe Wallwork } 122420282da2SJoe Wallwork 122520282da2SJoe Wallwork /* 122620282da2SJoe Wallwork DMPlexMetricAverage3 - Compute the unweighted average of three metrics 122720282da2SJoe Wallwork 122820282da2SJoe Wallwork Input Parameter: 122920282da2SJoe Wallwork + dm - The DM 123020282da2SJoe Wallwork . metric1 - The first metric to be averaged 123120282da2SJoe Wallwork . metric2 - The second metric to be averaged 123220282da2SJoe Wallwork - metric3 - The third metric to be averaged 123320282da2SJoe Wallwork 123420282da2SJoe Wallwork Output Parameter: 123520282da2SJoe Wallwork . metricAvg - The averaged metric 123620282da2SJoe Wallwork 123720282da2SJoe Wallwork Level: beginner 123820282da2SJoe Wallwork 123920282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2() 124020282da2SJoe Wallwork */ 124120282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg) 124220282da2SJoe Wallwork { 124320282da2SJoe Wallwork PetscErrorCode ierr; 124420282da2SJoe Wallwork PetscReal weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0}; 124520282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 124620282da2SJoe Wallwork 124720282da2SJoe Wallwork PetscFunctionBegin; 124820282da2SJoe Wallwork ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr); 124920282da2SJoe Wallwork PetscFunctionReturn(0); 125020282da2SJoe Wallwork } 125120282da2SJoe Wallwork 125220282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[]) 125320282da2SJoe Wallwork { 125420282da2SJoe Wallwork PetscErrorCode ierr; 125520282da2SJoe Wallwork PetscInt i, j, k, l, m; 125620282da2SJoe Wallwork PetscReal *evals, *evals1; 125720282da2SJoe Wallwork PetscScalar *evecs, *sqrtM1, *isqrtM1; 125820282da2SJoe Wallwork 125920282da2SJoe Wallwork PetscFunctionBegin; 126020282da2SJoe Wallwork ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr); 126120282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 126220282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 126320282da2SJoe Wallwork evecs[i*dim+j] = M1[i*dim+j]; 126420282da2SJoe Wallwork } 126520282da2SJoe Wallwork } 126620282da2SJoe Wallwork { 126720282da2SJoe Wallwork PetscScalar *work; 126820282da2SJoe Wallwork PetscBLASInt lwork; 126920282da2SJoe Wallwork 127020282da2SJoe Wallwork lwork = 5*dim; 127120282da2SJoe Wallwork ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr); 127220282da2SJoe Wallwork { 127320282da2SJoe Wallwork PetscBLASInt lierr, nb; 127420282da2SJoe Wallwork PetscReal sqrtk; 127520282da2SJoe Wallwork 127620282da2SJoe Wallwork /* Compute eigendecomposition of M1 */ 127720282da2SJoe Wallwork ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr); 127820282da2SJoe Wallwork ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 127920282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 128020282da2SJoe Wallwork { 128120282da2SJoe Wallwork PetscReal *rwork; 128220282da2SJoe Wallwork ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 128320282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr)); 128420282da2SJoe Wallwork ierr = PetscFree(rwork);CHKERRQ(ierr); 128520282da2SJoe Wallwork } 128620282da2SJoe Wallwork #else 128720282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr)); 128820282da2SJoe Wallwork #endif 1289*82490253SJoe Wallwork if (lierr) { 1290*82490253SJoe Wallwork LAPACKsyevFail(dim, M1); 1291*82490253SJoe Wallwork SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 1292*82490253SJoe Wallwork } 129320282da2SJoe Wallwork ierr = PetscFPTrapPop(); 129420282da2SJoe Wallwork 129520282da2SJoe Wallwork /* Compute square root and reciprocal */ 129620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 129720282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 129820282da2SJoe Wallwork sqrtM1[i*dim+j] = 0.0; 129920282da2SJoe Wallwork isqrtM1[i*dim+j] = 0.0; 130020282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 130120282da2SJoe Wallwork sqrtk = PetscSqrtReal(evals1[k]); 130220282da2SJoe Wallwork sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j]; 130320282da2SJoe Wallwork isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j]; 130420282da2SJoe Wallwork } 130520282da2SJoe Wallwork } 130620282da2SJoe Wallwork } 130720282da2SJoe Wallwork 130820282da2SJoe Wallwork /* Map into the space spanned by the eigenvectors of M1 */ 130920282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 131020282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 131120282da2SJoe Wallwork evecs[i*dim+j] = 0.0; 131220282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 131320282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 131420282da2SJoe Wallwork evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 131520282da2SJoe Wallwork } 131620282da2SJoe Wallwork } 131720282da2SJoe Wallwork } 131820282da2SJoe Wallwork } 131920282da2SJoe Wallwork 132020282da2SJoe Wallwork /* Compute eigendecomposition */ 132120282da2SJoe Wallwork ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 132220282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 132320282da2SJoe Wallwork { 132420282da2SJoe Wallwork PetscReal *rwork; 132520282da2SJoe Wallwork ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 132620282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 132720282da2SJoe Wallwork ierr = PetscFree(rwork);CHKERRQ(ierr); 132820282da2SJoe Wallwork } 132920282da2SJoe Wallwork #else 133020282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 133120282da2SJoe Wallwork #endif 1332*82490253SJoe Wallwork if (lierr) { 1333*82490253SJoe Wallwork for (i = 0; i < dim; ++i) { 1334*82490253SJoe Wallwork for (j = 0; j < dim; ++j) { 1335*82490253SJoe Wallwork evecs[i*dim+j] = 0.0; 1336*82490253SJoe Wallwork for (k = 0; k < dim; ++k) { 1337*82490253SJoe Wallwork for (l = 0; l < dim; ++l) { 1338*82490253SJoe Wallwork evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 1339*82490253SJoe Wallwork } 1340*82490253SJoe Wallwork } 1341*82490253SJoe Wallwork } 1342*82490253SJoe Wallwork } 1343*82490253SJoe Wallwork LAPACKsyevFail(dim, evecs); 1344*82490253SJoe Wallwork SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 1345*82490253SJoe Wallwork } 134620282da2SJoe Wallwork ierr = PetscFPTrapPop(); 134720282da2SJoe Wallwork 134820282da2SJoe Wallwork /* Modify eigenvalues */ 134920282da2SJoe Wallwork for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]); 135020282da2SJoe Wallwork 135120282da2SJoe Wallwork /* Map back to get the intersection */ 135220282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 135320282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 135420282da2SJoe Wallwork M2[i*dim+j] = 0.0; 135520282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 135620282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 135720282da2SJoe Wallwork for (m = 0; m < dim; ++m) { 135820282da2SJoe Wallwork M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m]; 135920282da2SJoe Wallwork } 136020282da2SJoe Wallwork } 136120282da2SJoe Wallwork } 136220282da2SJoe Wallwork } 136320282da2SJoe Wallwork } 136420282da2SJoe Wallwork } 136520282da2SJoe Wallwork ierr = PetscFree(work);CHKERRQ(ierr); 136620282da2SJoe Wallwork } 136720282da2SJoe Wallwork ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr); 136820282da2SJoe Wallwork PetscFunctionReturn(0); 136920282da2SJoe Wallwork } 137020282da2SJoe Wallwork 137120282da2SJoe Wallwork /* 137220282da2SJoe Wallwork DMPlexMetricIntersection - Compute the intersection of a list of metrics 137320282da2SJoe Wallwork 137420282da2SJoe Wallwork Input Parameter: 137520282da2SJoe Wallwork + dm - The DM 137620282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected 137720282da2SJoe Wallwork - metrics - The metrics to be intersected 137820282da2SJoe Wallwork 137920282da2SJoe Wallwork Output Parameter: 138020282da2SJoe Wallwork . metricInt - The intersected metric 138120282da2SJoe Wallwork 138220282da2SJoe Wallwork Level: beginner 138320282da2SJoe Wallwork 138420282da2SJoe Wallwork Notes: 138520282da2SJoe Wallwork 138620282da2SJoe Wallwork The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics. 138720282da2SJoe Wallwork 138820282da2SJoe Wallwork The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2. 138920282da2SJoe Wallwork 139020282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage() 139120282da2SJoe Wallwork */ 139220282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt) 139320282da2SJoe Wallwork { 139420282da2SJoe Wallwork PetscErrorCode ierr; 139520282da2SJoe Wallwork PetscInt dim, vStart, vEnd, v, i; 139620282da2SJoe Wallwork PetscScalar *met, *meti, *M, *Mi; 139720282da2SJoe Wallwork 139820282da2SJoe Wallwork PetscFunctionBegin; 139920282da2SJoe Wallwork if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); } 140020282da2SJoe Wallwork 140120282da2SJoe Wallwork /* Extract metadata from dm */ 140220282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 140320282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr); 140420282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 140520282da2SJoe Wallwork 140620282da2SJoe Wallwork /* Copy over the first metric */ 140720282da2SJoe Wallwork ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr); 140820282da2SJoe Wallwork 140920282da2SJoe Wallwork /* Intersect subsequent metrics in turn */ 141020282da2SJoe Wallwork if (numMetrics > 1) { 141120282da2SJoe Wallwork ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr); 141220282da2SJoe Wallwork for (i = 1; i < numMetrics; ++i) { 141320282da2SJoe Wallwork ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr); 141420282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 141520282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr); 141620282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr); 141720282da2SJoe Wallwork ierr = DMPlexMetricIntersection_Private(dim, Mi, M);CHKERRQ(ierr); 141820282da2SJoe Wallwork } 141920282da2SJoe Wallwork ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr); 142020282da2SJoe Wallwork } 142120282da2SJoe Wallwork ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr); 142220282da2SJoe Wallwork } 142320282da2SJoe Wallwork 142420282da2SJoe Wallwork PetscFunctionReturn(0); 142520282da2SJoe Wallwork } 142620282da2SJoe Wallwork 142720282da2SJoe Wallwork /* 142820282da2SJoe Wallwork DMPlexMetricIntersection2 - Compute the intersection of two metrics 142920282da2SJoe Wallwork 143020282da2SJoe Wallwork Input Parameter: 143120282da2SJoe Wallwork + dm - The DM 143220282da2SJoe Wallwork . metric1 - The first metric to be intersected 143320282da2SJoe Wallwork - metric2 - The second metric to be intersected 143420282da2SJoe Wallwork 143520282da2SJoe Wallwork Output Parameter: 143620282da2SJoe Wallwork . metricInt - The intersected metric 143720282da2SJoe Wallwork 143820282da2SJoe Wallwork Level: beginner 143920282da2SJoe Wallwork 144020282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3() 144120282da2SJoe Wallwork */ 144220282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt) 144320282da2SJoe Wallwork { 144420282da2SJoe Wallwork PetscErrorCode ierr; 144520282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 144620282da2SJoe Wallwork 144720282da2SJoe Wallwork PetscFunctionBegin; 144820282da2SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr); 144920282da2SJoe Wallwork PetscFunctionReturn(0); 145020282da2SJoe Wallwork } 145120282da2SJoe Wallwork 145220282da2SJoe Wallwork /* 145320282da2SJoe Wallwork DMPlexMetricIntersection3 - Compute the intersection of three metrics 145420282da2SJoe Wallwork 145520282da2SJoe Wallwork Input Parameter: 145620282da2SJoe Wallwork + dm - The DM 145720282da2SJoe Wallwork . metric1 - The first metric to be intersected 145820282da2SJoe Wallwork . metric2 - The second metric to be intersected 145920282da2SJoe Wallwork - metric3 - The third metric to be intersected 146020282da2SJoe Wallwork 146120282da2SJoe Wallwork Output Parameter: 146220282da2SJoe Wallwork . metricInt - The intersected metric 146320282da2SJoe Wallwork 146420282da2SJoe Wallwork Level: beginner 146520282da2SJoe Wallwork 146620282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2() 146720282da2SJoe Wallwork */ 146820282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt) 146920282da2SJoe Wallwork { 147020282da2SJoe Wallwork PetscErrorCode ierr; 147120282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 147220282da2SJoe Wallwork 147320282da2SJoe Wallwork PetscFunctionBegin; 147420282da2SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr); 147520282da2SJoe Wallwork PetscFunctionReturn(0); 147620282da2SJoe Wallwork } 1477