120282da2SJoe Wallwork #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 220282da2SJoe Wallwork #include <petscblaslapack.h> 320282da2SJoe Wallwork 4*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetFromOptions(DM dm) 5*31b70465SJoe Wallwork { 6*31b70465SJoe Wallwork MPI_Comm comm; 7*31b70465SJoe Wallwork PetscBool isotropic = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE; 8*31b70465SJoe Wallwork PetscErrorCode ierr; 9*31b70465SJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0; 10*31b70465SJoe Wallwork 11*31b70465SJoe Wallwork PetscFunctionBegin; 12*31b70465SJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 13*31b70465SJoe Wallwork ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr); 14*31b70465SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr); 15*31b70465SJoe Wallwork ierr = DMPlexMetricSetIsotropic(dm, isotropic); 16*31b70465SJoe Wallwork ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr); 17*31b70465SJoe Wallwork ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst); 18*31b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr); 19*31b70465SJoe Wallwork ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr); 20*31b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr); 21*31b70465SJoe Wallwork ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr); 22*31b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr); 23*31b70465SJoe Wallwork ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr); 24*31b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr); 25*31b70465SJoe Wallwork ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr); 26*31b70465SJoe Wallwork ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr); 27*31b70465SJoe Wallwork ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr); 28*31b70465SJoe Wallwork ierr = PetscOptionsEnd();CHKERRQ(ierr); 29*31b70465SJoe Wallwork PetscFunctionReturn(0); 30*31b70465SJoe Wallwork } 31*31b70465SJoe Wallwork 32*31b70465SJoe Wallwork /* 33*31b70465SJoe Wallwork DMPlexMetricSetIsotropic - Record whether a metric is isotropic 34*31b70465SJoe Wallwork 35*31b70465SJoe Wallwork Input parameters: 36*31b70465SJoe Wallwork + dm - The DM 37*31b70465SJoe Wallwork - isotropic - Is the metric isotropic? 38*31b70465SJoe Wallwork 39*31b70465SJoe Wallwork Level: beginner 40*31b70465SJoe Wallwork 41*31b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 42*31b70465SJoe Wallwork */ 43*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic) 44*31b70465SJoe Wallwork { 45*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 46*31b70465SJoe Wallwork PetscErrorCode ierr; 47*31b70465SJoe Wallwork 48*31b70465SJoe Wallwork PetscFunctionBegin; 49*31b70465SJoe Wallwork if (!plex->metricCtx) { 50*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 51*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 52*31b70465SJoe Wallwork } 53*31b70465SJoe Wallwork plex->metricCtx->isotropic = isotropic; 54*31b70465SJoe Wallwork PetscFunctionReturn(0); 55*31b70465SJoe Wallwork } 56*31b70465SJoe Wallwork 57*31b70465SJoe Wallwork /* 58*31b70465SJoe Wallwork DMPlexMetricIsIsotropic - Is a metric is isotropic? 59*31b70465SJoe Wallwork 60*31b70465SJoe Wallwork Input parameters: 61*31b70465SJoe Wallwork . dm - The DM 62*31b70465SJoe Wallwork 63*31b70465SJoe Wallwork Output parameters: 64*31b70465SJoe Wallwork . isotropic - Is the metric isotropic? 65*31b70465SJoe Wallwork 66*31b70465SJoe Wallwork Level: beginner 67*31b70465SJoe Wallwork 68*31b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 69*31b70465SJoe Wallwork */ 70*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic) 71*31b70465SJoe Wallwork { 72*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 73*31b70465SJoe Wallwork PetscErrorCode ierr; 74*31b70465SJoe Wallwork 75*31b70465SJoe Wallwork PetscFunctionBegin; 76*31b70465SJoe Wallwork if (!plex->metricCtx) { 77*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 78*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 79*31b70465SJoe Wallwork } 80*31b70465SJoe Wallwork *isotropic = plex->metricCtx->isotropic; 81*31b70465SJoe Wallwork PetscFunctionReturn(0); 82*31b70465SJoe Wallwork } 83*31b70465SJoe Wallwork 84*31b70465SJoe Wallwork /* 85*31b70465SJoe Wallwork DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization 86*31b70465SJoe Wallwork 87*31b70465SJoe Wallwork Input parameters: 88*31b70465SJoe Wallwork + dm - The DM 89*31b70465SJoe Wallwork - restrictAnisotropyFirst - Should anisotropy be normalized first? 90*31b70465SJoe Wallwork 91*31b70465SJoe Wallwork Level: beginner 92*31b70465SJoe Wallwork 93*31b70465SJoe Wallwork .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 94*31b70465SJoe Wallwork */ 95*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst) 96*31b70465SJoe Wallwork { 97*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 98*31b70465SJoe Wallwork PetscErrorCode ierr; 99*31b70465SJoe Wallwork 100*31b70465SJoe Wallwork PetscFunctionBegin; 101*31b70465SJoe Wallwork if (!plex->metricCtx) { 102*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 103*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 104*31b70465SJoe Wallwork } 105*31b70465SJoe Wallwork plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst; 106*31b70465SJoe Wallwork PetscFunctionReturn(0); 107*31b70465SJoe Wallwork } 108*31b70465SJoe Wallwork 109*31b70465SJoe Wallwork /* 110*31b70465SJoe Wallwork DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after? 111*31b70465SJoe Wallwork 112*31b70465SJoe Wallwork Input parameters: 113*31b70465SJoe Wallwork . dm - The DM 114*31b70465SJoe Wallwork 115*31b70465SJoe Wallwork Output parameters: 116*31b70465SJoe Wallwork . restrictAnisotropyFirst - Is anisotropy be normalized first? 117*31b70465SJoe Wallwork 118*31b70465SJoe Wallwork Level: beginner 119*31b70465SJoe Wallwork 120*31b70465SJoe Wallwork .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 121*31b70465SJoe Wallwork */ 122*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst) 123*31b70465SJoe Wallwork { 124*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 125*31b70465SJoe Wallwork PetscErrorCode ierr; 126*31b70465SJoe Wallwork 127*31b70465SJoe Wallwork PetscFunctionBegin; 128*31b70465SJoe Wallwork if (!plex->metricCtx) { 129*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 130*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 131*31b70465SJoe Wallwork } 132*31b70465SJoe Wallwork *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst; 133*31b70465SJoe Wallwork PetscFunctionReturn(0); 134*31b70465SJoe Wallwork } 135*31b70465SJoe Wallwork 136*31b70465SJoe Wallwork /* 137*31b70465SJoe Wallwork DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude 138*31b70465SJoe Wallwork 139*31b70465SJoe Wallwork Input parameters: 140*31b70465SJoe Wallwork + dm - The DM 141*31b70465SJoe Wallwork - h_min - The minimum tolerated metric magnitude 142*31b70465SJoe Wallwork 143*31b70465SJoe Wallwork Level: beginner 144*31b70465SJoe Wallwork 145*31b70465SJoe Wallwork .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude() 146*31b70465SJoe Wallwork */ 147*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min) 148*31b70465SJoe Wallwork { 149*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 150*31b70465SJoe Wallwork PetscErrorCode ierr; 151*31b70465SJoe Wallwork 152*31b70465SJoe Wallwork PetscFunctionBegin; 153*31b70465SJoe Wallwork if (!plex->metricCtx) { 154*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 155*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 156*31b70465SJoe Wallwork } 157*31b70465SJoe Wallwork if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min); 158*31b70465SJoe Wallwork plex->metricCtx->h_min = h_min; 159*31b70465SJoe Wallwork PetscFunctionReturn(0); 160*31b70465SJoe Wallwork } 161*31b70465SJoe Wallwork 162*31b70465SJoe Wallwork /* 163*31b70465SJoe Wallwork DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude 164*31b70465SJoe Wallwork 165*31b70465SJoe Wallwork Input parameters: 166*31b70465SJoe Wallwork . dm - The DM 167*31b70465SJoe Wallwork 168*31b70465SJoe Wallwork Input parameters: 169*31b70465SJoe Wallwork . h_min - The minimum tolerated metric magnitude 170*31b70465SJoe Wallwork 171*31b70465SJoe Wallwork Level: beginner 172*31b70465SJoe Wallwork 173*31b70465SJoe Wallwork .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude() 174*31b70465SJoe Wallwork */ 175*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min) 176*31b70465SJoe Wallwork { 177*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 178*31b70465SJoe Wallwork PetscErrorCode ierr; 179*31b70465SJoe Wallwork 180*31b70465SJoe Wallwork PetscFunctionBegin; 181*31b70465SJoe Wallwork if (!plex->metricCtx) { 182*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 183*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 184*31b70465SJoe Wallwork } 185*31b70465SJoe Wallwork *h_min = plex->metricCtx->h_min; 186*31b70465SJoe Wallwork PetscFunctionReturn(0); 187*31b70465SJoe Wallwork } 188*31b70465SJoe Wallwork 189*31b70465SJoe Wallwork /* 190*31b70465SJoe Wallwork DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude 191*31b70465SJoe Wallwork 192*31b70465SJoe Wallwork Input parameters: 193*31b70465SJoe Wallwork + dm - The DM 194*31b70465SJoe Wallwork - h_max - The maximum tolerated metric magnitude 195*31b70465SJoe Wallwork 196*31b70465SJoe Wallwork Level: beginner 197*31b70465SJoe Wallwork 198*31b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude() 199*31b70465SJoe Wallwork */ 200*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max) 201*31b70465SJoe Wallwork { 202*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 203*31b70465SJoe Wallwork PetscErrorCode ierr; 204*31b70465SJoe Wallwork 205*31b70465SJoe Wallwork PetscFunctionBegin; 206*31b70465SJoe Wallwork if (!plex->metricCtx) { 207*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 208*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 209*31b70465SJoe Wallwork } 210*31b70465SJoe Wallwork if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max); 211*31b70465SJoe Wallwork plex->metricCtx->h_max = h_max; 212*31b70465SJoe Wallwork PetscFunctionReturn(0); 213*31b70465SJoe Wallwork } 214*31b70465SJoe Wallwork 215*31b70465SJoe Wallwork /* 216*31b70465SJoe Wallwork DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude 217*31b70465SJoe Wallwork 218*31b70465SJoe Wallwork Input parameters: 219*31b70465SJoe Wallwork . dm - The DM 220*31b70465SJoe Wallwork 221*31b70465SJoe Wallwork Input parameters: 222*31b70465SJoe Wallwork . h_max - The maximum tolerated metric magnitude 223*31b70465SJoe Wallwork 224*31b70465SJoe Wallwork Level: beginner 225*31b70465SJoe Wallwork 226*31b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude() 227*31b70465SJoe Wallwork */ 228*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max) 229*31b70465SJoe Wallwork { 230*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 231*31b70465SJoe Wallwork PetscErrorCode ierr; 232*31b70465SJoe Wallwork 233*31b70465SJoe Wallwork PetscFunctionBegin; 234*31b70465SJoe Wallwork if (!plex->metricCtx) { 235*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 236*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 237*31b70465SJoe Wallwork } 238*31b70465SJoe Wallwork *h_max = plex->metricCtx->h_max; 239*31b70465SJoe Wallwork PetscFunctionReturn(0); 240*31b70465SJoe Wallwork } 241*31b70465SJoe Wallwork 242*31b70465SJoe Wallwork /* 243*31b70465SJoe Wallwork DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy 244*31b70465SJoe Wallwork 245*31b70465SJoe Wallwork Input parameters: 246*31b70465SJoe Wallwork + dm - The DM 247*31b70465SJoe Wallwork - a_max - The maximum tolerated metric anisotropy 248*31b70465SJoe Wallwork 249*31b70465SJoe Wallwork Level: beginner 250*31b70465SJoe Wallwork 251*31b70465SJoe Wallwork Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one. 252*31b70465SJoe Wallwork 253*31b70465SJoe Wallwork .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude() 254*31b70465SJoe Wallwork */ 255*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max) 256*31b70465SJoe Wallwork { 257*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 258*31b70465SJoe Wallwork PetscErrorCode ierr; 259*31b70465SJoe Wallwork 260*31b70465SJoe Wallwork PetscFunctionBegin; 261*31b70465SJoe Wallwork if (!plex->metricCtx) { 262*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 263*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 264*31b70465SJoe Wallwork } 265*31b70465SJoe Wallwork if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max); 266*31b70465SJoe Wallwork plex->metricCtx->a_max = a_max; 267*31b70465SJoe Wallwork PetscFunctionReturn(0); 268*31b70465SJoe Wallwork } 269*31b70465SJoe Wallwork 270*31b70465SJoe Wallwork /* 271*31b70465SJoe Wallwork DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy 272*31b70465SJoe Wallwork 273*31b70465SJoe Wallwork Input parameters: 274*31b70465SJoe Wallwork . dm - The DM 275*31b70465SJoe Wallwork 276*31b70465SJoe Wallwork Input parameters: 277*31b70465SJoe Wallwork . a_max - The maximum tolerated metric anisotropy 278*31b70465SJoe Wallwork 279*31b70465SJoe Wallwork Level: beginner 280*31b70465SJoe Wallwork 281*31b70465SJoe Wallwork .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude() 282*31b70465SJoe Wallwork */ 283*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max) 284*31b70465SJoe Wallwork { 285*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 286*31b70465SJoe Wallwork PetscErrorCode ierr; 287*31b70465SJoe Wallwork 288*31b70465SJoe Wallwork PetscFunctionBegin; 289*31b70465SJoe Wallwork if (!plex->metricCtx) { 290*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 291*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 292*31b70465SJoe Wallwork } 293*31b70465SJoe Wallwork *a_max = plex->metricCtx->a_max; 294*31b70465SJoe Wallwork PetscFunctionReturn(0); 295*31b70465SJoe Wallwork } 296*31b70465SJoe Wallwork 297*31b70465SJoe Wallwork /* 298*31b70465SJoe Wallwork DMPlexMetricSetTargetComplexity - Set the target metric complexity 299*31b70465SJoe Wallwork 300*31b70465SJoe Wallwork Input parameters: 301*31b70465SJoe Wallwork + dm - The DM 302*31b70465SJoe Wallwork - targetComplexity - The target metric complexity 303*31b70465SJoe Wallwork 304*31b70465SJoe Wallwork Level: beginner 305*31b70465SJoe Wallwork 306*31b70465SJoe Wallwork .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder() 307*31b70465SJoe Wallwork */ 308*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity) 309*31b70465SJoe Wallwork { 310*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 311*31b70465SJoe Wallwork PetscErrorCode ierr; 312*31b70465SJoe Wallwork 313*31b70465SJoe Wallwork PetscFunctionBegin; 314*31b70465SJoe Wallwork if (!plex->metricCtx) { 315*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 316*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 317*31b70465SJoe Wallwork } 318*31b70465SJoe Wallwork if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive"); 319*31b70465SJoe Wallwork plex->metricCtx->targetComplexity = targetComplexity; 320*31b70465SJoe Wallwork PetscFunctionReturn(0); 321*31b70465SJoe Wallwork } 322*31b70465SJoe Wallwork 323*31b70465SJoe Wallwork /* 324*31b70465SJoe Wallwork DMPlexMetricGetTargetComplexity - Get the target metric complexity 325*31b70465SJoe Wallwork 326*31b70465SJoe Wallwork Input parameters: 327*31b70465SJoe Wallwork . dm - The DM 328*31b70465SJoe Wallwork 329*31b70465SJoe Wallwork Input parameters: 330*31b70465SJoe Wallwork . targetComplexity - The target metric complexity 331*31b70465SJoe Wallwork 332*31b70465SJoe Wallwork Level: beginner 333*31b70465SJoe Wallwork 334*31b70465SJoe Wallwork .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder() 335*31b70465SJoe Wallwork */ 336*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity) 337*31b70465SJoe Wallwork { 338*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 339*31b70465SJoe Wallwork PetscErrorCode ierr; 340*31b70465SJoe Wallwork 341*31b70465SJoe Wallwork PetscFunctionBegin; 342*31b70465SJoe Wallwork if (!plex->metricCtx) { 343*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 344*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 345*31b70465SJoe Wallwork } 346*31b70465SJoe Wallwork *targetComplexity = plex->metricCtx->targetComplexity; 347*31b70465SJoe Wallwork PetscFunctionReturn(0); 348*31b70465SJoe Wallwork } 349*31b70465SJoe Wallwork 350*31b70465SJoe Wallwork /* 351*31b70465SJoe Wallwork DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization 352*31b70465SJoe Wallwork 353*31b70465SJoe Wallwork Input parameters: 354*31b70465SJoe Wallwork + dm - The DM 355*31b70465SJoe Wallwork - p - The normalization order 356*31b70465SJoe Wallwork 357*31b70465SJoe Wallwork Level: beginner 358*31b70465SJoe Wallwork 359*31b70465SJoe Wallwork .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity() 360*31b70465SJoe Wallwork */ 361*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p) 362*31b70465SJoe Wallwork { 363*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 364*31b70465SJoe Wallwork PetscErrorCode ierr; 365*31b70465SJoe Wallwork 366*31b70465SJoe Wallwork PetscFunctionBegin; 367*31b70465SJoe Wallwork if (!plex->metricCtx) { 368*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 369*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 370*31b70465SJoe Wallwork } 371*31b70465SJoe Wallwork if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater"); 372*31b70465SJoe Wallwork plex->metricCtx->p = p; 373*31b70465SJoe Wallwork PetscFunctionReturn(0); 374*31b70465SJoe Wallwork } 375*31b70465SJoe Wallwork 376*31b70465SJoe Wallwork /* 377*31b70465SJoe Wallwork DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization 378*31b70465SJoe Wallwork 379*31b70465SJoe Wallwork Input parameters: 380*31b70465SJoe Wallwork . dm - The DM 381*31b70465SJoe Wallwork 382*31b70465SJoe Wallwork Input parameters: 383*31b70465SJoe Wallwork . p - The normalization order 384*31b70465SJoe Wallwork 385*31b70465SJoe Wallwork Level: beginner 386*31b70465SJoe Wallwork 387*31b70465SJoe Wallwork .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity() 388*31b70465SJoe Wallwork */ 389*31b70465SJoe Wallwork PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p) 390*31b70465SJoe Wallwork { 391*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 392*31b70465SJoe Wallwork PetscErrorCode ierr; 393*31b70465SJoe Wallwork 394*31b70465SJoe Wallwork PetscFunctionBegin; 395*31b70465SJoe Wallwork if (!plex->metricCtx) { 396*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 397*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 398*31b70465SJoe Wallwork } 399*31b70465SJoe Wallwork *p = plex->metricCtx->p; 400*31b70465SJoe Wallwork PetscFunctionReturn(0); 401*31b70465SJoe Wallwork } 402*31b70465SJoe Wallwork 40320282da2SJoe Wallwork PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric) 40420282da2SJoe Wallwork { 40520282da2SJoe Wallwork MPI_Comm comm; 40620282da2SJoe Wallwork PetscErrorCode ierr; 40720282da2SJoe Wallwork PetscFE fe; 40820282da2SJoe Wallwork PetscInt dim; 40920282da2SJoe Wallwork 41020282da2SJoe Wallwork PetscFunctionBegin; 41120282da2SJoe Wallwork 41220282da2SJoe Wallwork /* Extract metadata from dm */ 41320282da2SJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 41420282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 41520282da2SJoe Wallwork 41620282da2SJoe Wallwork /* Create a P1 field of the requested size */ 41720282da2SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 41820282da2SJoe Wallwork ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr); 41920282da2SJoe Wallwork ierr = DMCreateDS(dm);CHKERRQ(ierr); 42020282da2SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 42120282da2SJoe Wallwork ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr); 42220282da2SJoe Wallwork 42320282da2SJoe Wallwork PetscFunctionReturn(0); 42420282da2SJoe Wallwork } 42520282da2SJoe Wallwork 42620282da2SJoe Wallwork /* 42720282da2SJoe Wallwork DMPlexMetricCreate - Create a Riemannian metric field 42820282da2SJoe Wallwork 42920282da2SJoe Wallwork Input parameters: 43020282da2SJoe Wallwork + dm - The DM 43120282da2SJoe Wallwork - f - The field number to use 43220282da2SJoe Wallwork 43320282da2SJoe Wallwork Output parameter: 43420282da2SJoe Wallwork . metric - The metric 43520282da2SJoe Wallwork 43620282da2SJoe Wallwork Level: beginner 43720282da2SJoe Wallwork 438*31b70465SJoe Wallwork Notes: 439*31b70465SJoe Wallwork 440*31b70465SJoe Wallwork It is assumed that the DM is comprised of simplices. 441*31b70465SJoe Wallwork 442*31b70465SJoe Wallwork Command line options for Riemannian metrics: 443*31b70465SJoe Wallwork 444*31b70465SJoe Wallwork -dm_plex_metric_isotropic - Is the metric isotropic? 445*31b70465SJoe Wallwork -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 446*31b70465SJoe Wallwork -dm_plex_metric_h_min - Minimum tolerated metric magnitude 447*31b70465SJoe Wallwork -dm_plex_metric_h_max - Maximum tolerated metric magnitude 448*31b70465SJoe Wallwork -dm_plex_metric_a_max - Maximum tolerated anisotropy 449*31b70465SJoe Wallwork -dm_plex_metric_p - L-p normalization order 450*31b70465SJoe Wallwork -dm_plex_metric_target_complexity - Target metric complexity 45120282da2SJoe Wallwork 45220282da2SJoe Wallwork .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic() 45320282da2SJoe Wallwork */ 45420282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric) 45520282da2SJoe Wallwork { 456*31b70465SJoe Wallwork DM_Plex *plex = (DM_Plex *) dm->data; 45720282da2SJoe Wallwork PetscErrorCode ierr; 45820282da2SJoe Wallwork PetscInt coordDim, Nd; 45920282da2SJoe Wallwork 46020282da2SJoe Wallwork PetscFunctionBegin; 461*31b70465SJoe Wallwork if (!plex->metricCtx) { 462*31b70465SJoe Wallwork ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 463*31b70465SJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 464*31b70465SJoe Wallwork } 46520282da2SJoe Wallwork ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 46620282da2SJoe Wallwork Nd = coordDim*coordDim; 46720282da2SJoe Wallwork ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr); 46820282da2SJoe Wallwork PetscFunctionReturn(0); 46920282da2SJoe Wallwork } 47020282da2SJoe Wallwork 47120282da2SJoe Wallwork typedef struct { 47220282da2SJoe Wallwork PetscReal scaling; /* Scaling for uniform metric diagonal */ 47320282da2SJoe Wallwork } DMPlexMetricUniformCtx; 47420282da2SJoe Wallwork 47520282da2SJoe Wallwork static PetscErrorCode diagonal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 47620282da2SJoe Wallwork { 47720282da2SJoe Wallwork DMPlexMetricUniformCtx *user = (DMPlexMetricUniformCtx*)ctx; 47820282da2SJoe Wallwork PetscInt i, j; 47920282da2SJoe Wallwork 48020282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 48120282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 48220282da2SJoe Wallwork if (i == j) u[i+dim*j] = user->scaling; 48320282da2SJoe Wallwork else u[i+dim*j] = 0.0; 48420282da2SJoe Wallwork } 48520282da2SJoe Wallwork } 48620282da2SJoe Wallwork return 0; 48720282da2SJoe Wallwork } 48820282da2SJoe Wallwork 48920282da2SJoe Wallwork /* 49020282da2SJoe Wallwork DMPlexMetricCreateUniform - Construct a uniform isotropic metric 49120282da2SJoe Wallwork 49220282da2SJoe Wallwork Input parameters: 49320282da2SJoe Wallwork + dm - The DM 49420282da2SJoe Wallwork . f - The field number to use 49520282da2SJoe Wallwork - alpha - Scaling parameter for the diagonal 49620282da2SJoe Wallwork 49720282da2SJoe Wallwork Output parameter: 49820282da2SJoe Wallwork . metric - The uniform metric 49920282da2SJoe Wallwork 50020282da2SJoe Wallwork Level: beginner 50120282da2SJoe Wallwork 50220282da2SJoe Wallwork Note: It is assumed that the DM is comprised of simplices. 50320282da2SJoe Wallwork 50420282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic() 50520282da2SJoe Wallwork */ 50620282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric) 50720282da2SJoe Wallwork { 50820282da2SJoe Wallwork DMPlexMetricUniformCtx user; 50920282da2SJoe Wallwork PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 51020282da2SJoe Wallwork PetscErrorCode ierr; 51120282da2SJoe Wallwork void *ctxs[1]; 51220282da2SJoe Wallwork 51320282da2SJoe Wallwork PetscFunctionBegin; 51420282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr); 51520282da2SJoe Wallwork if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined"); 51620282da2SJoe Wallwork if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha); 51720282da2SJoe Wallwork else user.scaling = alpha; 51820282da2SJoe Wallwork funcs[0] = diagonal; 51920282da2SJoe Wallwork ctxs[0] = &user; 52020282da2SJoe Wallwork ierr = DMProjectFunctionLocal(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, *metric);CHKERRQ(ierr); 52120282da2SJoe Wallwork PetscFunctionReturn(0); 52220282da2SJoe Wallwork } 52320282da2SJoe Wallwork 52420282da2SJoe Wallwork /* 52520282da2SJoe Wallwork DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator 52620282da2SJoe Wallwork 52720282da2SJoe Wallwork Input parameters: 52820282da2SJoe Wallwork + dm - The DM 52920282da2SJoe Wallwork . f - The field number to use 53020282da2SJoe Wallwork - indicator - The error indicator 53120282da2SJoe Wallwork 53220282da2SJoe Wallwork Output parameter: 53320282da2SJoe Wallwork . metric - The isotropic metric 53420282da2SJoe Wallwork 53520282da2SJoe Wallwork Level: beginner 53620282da2SJoe Wallwork 53720282da2SJoe Wallwork Notes: 53820282da2SJoe Wallwork 53920282da2SJoe Wallwork It is assumed that the DM is comprised of simplices. 54020282da2SJoe Wallwork 54120282da2SJoe Wallwork The indicator needs to be a scalar field defined at *vertices*. 54220282da2SJoe Wallwork 54320282da2SJoe Wallwork .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform() 54420282da2SJoe Wallwork */ 54520282da2SJoe Wallwork PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric) 54620282da2SJoe Wallwork { 54720282da2SJoe Wallwork DM dmIndi; 54820282da2SJoe Wallwork PetscErrorCode ierr; 54920282da2SJoe Wallwork PetscInt dim, d, vStart, vEnd, v; 55020282da2SJoe Wallwork const PetscScalar *indi; 55120282da2SJoe Wallwork PetscScalar *met; 55220282da2SJoe Wallwork 55320282da2SJoe Wallwork PetscFunctionBegin; 55420282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 55520282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr); 55620282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 55720282da2SJoe Wallwork ierr = VecGetArrayRead(indicator, &indi);CHKERRQ(ierr); 55820282da2SJoe Wallwork ierr = VecGetArrayWrite(*metric, &met);CHKERRQ(ierr); 55920282da2SJoe Wallwork ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr); 56020282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 56120282da2SJoe Wallwork PetscScalar *vindi, *vmet; 56220282da2SJoe Wallwork ierr = DMPlexPointLocalRead(dmIndi, v, indi, &vindi);CHKERRQ(ierr); 56320282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr); 56420282da2SJoe Wallwork for (d = 0; d < dim; ++d) vmet[d*(dim+1)] = vindi[0]; 56520282da2SJoe Wallwork } 56620282da2SJoe Wallwork ierr = VecRestoreArrayWrite(*metric, &met);CHKERRQ(ierr); 56720282da2SJoe Wallwork ierr = VecRestoreArrayRead(indicator, &indi);CHKERRQ(ierr); 56820282da2SJoe Wallwork PetscFunctionReturn(0); 56920282da2SJoe Wallwork } 57020282da2SJoe Wallwork 57120282da2SJoe Wallwork static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[]) 57220282da2SJoe Wallwork { 57320282da2SJoe Wallwork PetscErrorCode ierr; 57420282da2SJoe Wallwork PetscInt i, j, k; 57520282da2SJoe 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); 57620282da2SJoe Wallwork PetscScalar *Mpos; 57720282da2SJoe Wallwork 57820282da2SJoe Wallwork PetscFunctionBegin; 57920282da2SJoe Wallwork ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr); 58020282da2SJoe Wallwork 58120282da2SJoe Wallwork /* Symmetrize */ 58220282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 58320282da2SJoe Wallwork Mpos[i*dim+i] = Mp[i*dim+i]; 58420282da2SJoe Wallwork for (j = i+1; j < dim; ++j) { 58520282da2SJoe Wallwork Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 58620282da2SJoe Wallwork Mpos[j*dim+i] = Mpos[i*dim+j]; 58720282da2SJoe Wallwork } 58820282da2SJoe Wallwork } 58920282da2SJoe Wallwork 59020282da2SJoe Wallwork /* Compute eigendecomposition */ 59120282da2SJoe Wallwork { 59220282da2SJoe Wallwork PetscScalar *work; 59320282da2SJoe Wallwork PetscBLASInt lwork; 59420282da2SJoe Wallwork 59520282da2SJoe Wallwork lwork = 5*dim; 59620282da2SJoe Wallwork ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr); 59720282da2SJoe Wallwork { 59820282da2SJoe Wallwork PetscBLASInt lierr; 59920282da2SJoe Wallwork PetscBLASInt nb; 60020282da2SJoe Wallwork 60120282da2SJoe Wallwork ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr); 60220282da2SJoe Wallwork ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 60320282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 60420282da2SJoe Wallwork { 60520282da2SJoe Wallwork PetscReal *rwork; 60620282da2SJoe Wallwork ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 60720282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr)); 60820282da2SJoe Wallwork ierr = PetscFree(rwork);CHKERRQ(ierr); 60920282da2SJoe Wallwork } 61020282da2SJoe Wallwork #else 61120282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr)); 61220282da2SJoe Wallwork #endif 61320282da2SJoe Wallwork if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 61420282da2SJoe Wallwork ierr = PetscFPTrapPop();CHKERRQ(ierr); 61520282da2SJoe Wallwork } 61620282da2SJoe Wallwork ierr = PetscFree(work);CHKERRQ(ierr); 61720282da2SJoe Wallwork } 61820282da2SJoe Wallwork 61920282da2SJoe Wallwork /* Reflect to positive orthant and enforce maximum and minimum size */ 62020282da2SJoe Wallwork max_eig = 0.0; 62120282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 62220282da2SJoe Wallwork eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i]))); 62320282da2SJoe Wallwork max_eig = PetscMax(eigs[i], max_eig); 62420282da2SJoe Wallwork } 62520282da2SJoe Wallwork 62620282da2SJoe Wallwork /* Enforce maximum anisotropy */ 62720282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 62820282da2SJoe Wallwork if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min); 62920282da2SJoe Wallwork } 63020282da2SJoe Wallwork 63120282da2SJoe Wallwork /* Reconstruct Hessian */ 63220282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 63320282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 63420282da2SJoe Wallwork Mp[i*dim+j] = 0.0; 63520282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 63620282da2SJoe Wallwork Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j]; 63720282da2SJoe Wallwork } 63820282da2SJoe Wallwork } 63920282da2SJoe Wallwork } 64020282da2SJoe Wallwork ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr); 64120282da2SJoe Wallwork 64220282da2SJoe Wallwork PetscFunctionReturn(0); 64320282da2SJoe Wallwork } 64420282da2SJoe Wallwork 64520282da2SJoe Wallwork /* 64620282da2SJoe Wallwork DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric 64720282da2SJoe Wallwork 64820282da2SJoe Wallwork Input parameters: 64920282da2SJoe Wallwork + dm - The DM 65020282da2SJoe Wallwork . restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced? 65120282da2SJoe Wallwork - metric - The metric 65220282da2SJoe Wallwork 65320282da2SJoe Wallwork Output parameter: 65420282da2SJoe Wallwork . metric - The metric 65520282da2SJoe Wallwork 65620282da2SJoe Wallwork Level: beginner 65720282da2SJoe Wallwork 65820282da2SJoe Wallwork .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection() 65920282da2SJoe Wallwork */ 66020282da2SJoe Wallwork PetscErrorCode DMPlexMetricEnforceSPD(DM dm, PetscBool restrictSizes, Vec metric) 66120282da2SJoe Wallwork { 66220282da2SJoe Wallwork DMPlexMetricCtx *user; 66320282da2SJoe Wallwork PetscErrorCode ierr; 66420282da2SJoe Wallwork PetscInt dim, vStart, vEnd, v; 66520282da2SJoe Wallwork PetscScalar *met; 66620282da2SJoe Wallwork PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0; 66720282da2SJoe Wallwork 66820282da2SJoe Wallwork PetscFunctionBegin; 66920282da2SJoe Wallwork 67020282da2SJoe Wallwork /* Extract metadata from dm */ 67120282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 67220282da2SJoe Wallwork ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr); 67320282da2SJoe Wallwork if (restrictSizes) { 67420282da2SJoe Wallwork if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max); 67520282da2SJoe Wallwork if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min); 67620282da2SJoe Wallwork if (user->a_max > 1.0) a_max = user->a_max; 67720282da2SJoe Wallwork } 67820282da2SJoe Wallwork 67920282da2SJoe Wallwork /* Enforce SPD */ 68020282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 68120282da2SJoe Wallwork ierr = VecGetArray(metric, &met);CHKERRQ(ierr); 68220282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 68320282da2SJoe Wallwork PetscScalar *vmet; 68420282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr); 68520282da2SJoe Wallwork ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, vmet);CHKERRQ(ierr); 68620282da2SJoe Wallwork } 68720282da2SJoe Wallwork ierr = VecRestoreArray(metric, &met);CHKERRQ(ierr); 68820282da2SJoe Wallwork 68920282da2SJoe Wallwork PetscFunctionReturn(0); 69020282da2SJoe Wallwork } 69120282da2SJoe Wallwork 69220282da2SJoe Wallwork static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 69320282da2SJoe Wallwork const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 69420282da2SJoe Wallwork const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 69520282da2SJoe Wallwork PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 69620282da2SJoe Wallwork { 69720282da2SJoe Wallwork const PetscScalar p = constants[0]; 69820282da2SJoe Wallwork PetscReal detH = 0.0; 69920282da2SJoe Wallwork 70020282da2SJoe Wallwork if (dim == 2) DMPlex_Det2D_Scalar_Internal(&detH, u); 70120282da2SJoe Wallwork else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detH, u); 70220282da2SJoe Wallwork f0[0] = PetscPowReal(detH, p/(2.0*p + dim)); 70320282da2SJoe Wallwork } 70420282da2SJoe Wallwork 70520282da2SJoe Wallwork /* 70620282da2SJoe Wallwork DMPlexMetricNormalize - Apply L-p normalization to a metric 70720282da2SJoe Wallwork 70820282da2SJoe Wallwork Input parameters: 70920282da2SJoe Wallwork + dm - The DM 71020282da2SJoe Wallwork . metricIn - The unnormalized metric 71120282da2SJoe Wallwork - restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced? 71220282da2SJoe Wallwork 71320282da2SJoe Wallwork Output parameter: 71420282da2SJoe Wallwork . metricOut - The normalized metric 71520282da2SJoe Wallwork 71620282da2SJoe Wallwork Level: beginner 71720282da2SJoe Wallwork 71820282da2SJoe Wallwork .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection() 71920282da2SJoe Wallwork */ 72020282da2SJoe Wallwork PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, Vec *metricOut) 72120282da2SJoe Wallwork { 72220282da2SJoe Wallwork DMPlexMetricCtx *user; 72320282da2SJoe Wallwork MPI_Comm comm; 72420282da2SJoe Wallwork PetscDS ds; 72520282da2SJoe Wallwork PetscErrorCode ierr; 72620282da2SJoe Wallwork PetscInt dim, Nd, vStart, vEnd, v, i; 72720282da2SJoe Wallwork PetscScalar *met, integral, constants[1]; 72820282da2SJoe Wallwork PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, target; 72920282da2SJoe Wallwork 73020282da2SJoe Wallwork PetscFunctionBegin; 73120282da2SJoe Wallwork 73220282da2SJoe Wallwork /* Extract metadata from dm */ 73320282da2SJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 73420282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 73520282da2SJoe Wallwork Nd = dim*dim; 73620282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 73720282da2SJoe Wallwork ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr); 73820282da2SJoe Wallwork if (restrictSizes && user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max; 73920282da2SJoe Wallwork if (PetscAbsReal(user->p) >= 1.0) p = user->p; 74020282da2SJoe Wallwork else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric normalization order %f should be greater than one.", user->p); 74120282da2SJoe Wallwork constants[0] = p; 74220282da2SJoe Wallwork if (user->targetComplexity > 0.0) target = user->targetComplexity; 74320282da2SJoe Wallwork else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Target metric complexity %f should be positive.", user->targetComplexity); 74420282da2SJoe Wallwork 74520282da2SJoe Wallwork /* Set up metric and ensure it is SPD */ 74620282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr); 74720282da2SJoe Wallwork ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr); 74820282da2SJoe Wallwork ierr = DMPlexMetricEnforceSPD(dm, PETSC_FALSE, *metricOut);CHKERRQ(ierr); 74920282da2SJoe Wallwork 75020282da2SJoe Wallwork /* Compute global normalization factor */ 75120282da2SJoe Wallwork ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 75220282da2SJoe Wallwork ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr); 75320282da2SJoe Wallwork ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr); 75420282da2SJoe Wallwork ierr = DMPlexComputeIntegralFEM(dm, *metricOut, &integral, NULL);CHKERRQ(ierr); 75520282da2SJoe Wallwork factGlob = PetscPowReal(target/PetscRealPart(integral), 2.0/dim); 75620282da2SJoe Wallwork 75720282da2SJoe Wallwork /* Apply local scaling */ 75820282da2SJoe Wallwork a_max = 0.0; 75920282da2SJoe Wallwork if (restrictSizes) { 76020282da2SJoe Wallwork if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max); 76120282da2SJoe Wallwork if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min); 76220282da2SJoe Wallwork if (!user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max; 76320282da2SJoe Wallwork } 76420282da2SJoe Wallwork ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr); 76520282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 76620282da2SJoe Wallwork PetscScalar *Mp; 76720282da2SJoe Wallwork PetscReal detM, fact; 76820282da2SJoe Wallwork 76920282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr); 77020282da2SJoe Wallwork if (dim == 2) DMPlex_Det2D_Scalar_Internal(&detM, Mp); 77120282da2SJoe Wallwork else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detM, Mp); 77220282da2SJoe Wallwork else SETERRQ1(comm, PETSC_ERR_SUP, "Dimension %d not supported", dim); 77320282da2SJoe Wallwork fact = factGlob * PetscPowReal(PetscAbsReal(detM), -1.0/(2*p+dim)); 77420282da2SJoe Wallwork for (i = 0; i < Nd; ++i) Mp[i] *= fact; 77520282da2SJoe Wallwork if (restrictSizes) { ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, Mp);CHKERRQ(ierr); } 77620282da2SJoe Wallwork } 77720282da2SJoe Wallwork ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr); 77820282da2SJoe Wallwork 77920282da2SJoe Wallwork PetscFunctionReturn(0); 78020282da2SJoe Wallwork } 78120282da2SJoe Wallwork 78220282da2SJoe Wallwork /* 78320282da2SJoe Wallwork DMPlexMetricAverage - Compute the average of a list of metrics 78420282da2SJoe Wallwork 78520282da2SJoe Wallwork Input Parameter: 78620282da2SJoe Wallwork + dm - The DM 78720282da2SJoe Wallwork . numMetrics - The number of metrics to be averaged 78820282da2SJoe Wallwork . weights - Weights for the average 78920282da2SJoe Wallwork - metrics - The metrics to be averaged 79020282da2SJoe Wallwork 79120282da2SJoe Wallwork Output Parameter: 79220282da2SJoe Wallwork . metricAvg - The averaged metric 79320282da2SJoe Wallwork 79420282da2SJoe Wallwork Level: beginner 79520282da2SJoe Wallwork 79620282da2SJoe Wallwork Notes: 79720282da2SJoe Wallwork The weights should sum to unity. 79820282da2SJoe Wallwork 79920282da2SJoe Wallwork If weights are not provided then an unweighted average is used. 80020282da2SJoe Wallwork 80120282da2SJoe Wallwork .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection() 80220282da2SJoe Wallwork */ 80320282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg) 80420282da2SJoe Wallwork { 80520282da2SJoe Wallwork PetscBool haveWeights = PETSC_TRUE; 80620282da2SJoe Wallwork PetscErrorCode ierr; 80720282da2SJoe Wallwork PetscInt i; 80820282da2SJoe Wallwork PetscReal sum = 0.0, tol = 1.0e-10; 80920282da2SJoe Wallwork 81020282da2SJoe Wallwork PetscFunctionBegin; 81120282da2SJoe Wallwork if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); } 81220282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr); 81320282da2SJoe Wallwork ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr); 81420282da2SJoe Wallwork 81520282da2SJoe Wallwork /* Default to the unweighted case */ 81620282da2SJoe Wallwork if (!weights) { 81720282da2SJoe Wallwork ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr); 81820282da2SJoe Wallwork haveWeights = PETSC_FALSE; 81920282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; } 82020282da2SJoe Wallwork } 82120282da2SJoe Wallwork 82220282da2SJoe Wallwork /* Check weights sum to unity */ 82320282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) { sum += weights[i]; } 82420282da2SJoe Wallwork if (PetscAbsReal(sum - 1) > tol) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); } 82520282da2SJoe Wallwork 82620282da2SJoe Wallwork /* Compute metric average */ 82720282da2SJoe Wallwork for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); } 82820282da2SJoe Wallwork if (!haveWeights) {ierr = PetscFree(weights); } 82920282da2SJoe Wallwork PetscFunctionReturn(0); 83020282da2SJoe Wallwork } 83120282da2SJoe Wallwork 83220282da2SJoe Wallwork /* 83320282da2SJoe Wallwork DMPlexMetricAverage2 - Compute the unweighted average of two metrics 83420282da2SJoe Wallwork 83520282da2SJoe Wallwork Input Parameter: 83620282da2SJoe Wallwork + dm - The DM 83720282da2SJoe Wallwork . metric1 - The first metric to be averaged 83820282da2SJoe Wallwork - metric2 - The second metric to be averaged 83920282da2SJoe Wallwork 84020282da2SJoe Wallwork Output Parameter: 84120282da2SJoe Wallwork . metricAvg - The averaged metric 84220282da2SJoe Wallwork 84320282da2SJoe Wallwork Level: beginner 84420282da2SJoe Wallwork 84520282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3() 84620282da2SJoe Wallwork */ 84720282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg) 84820282da2SJoe Wallwork { 84920282da2SJoe Wallwork PetscErrorCode ierr; 85020282da2SJoe Wallwork PetscReal weights[2] = {0.5, 0.5}; 85120282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 85220282da2SJoe Wallwork 85320282da2SJoe Wallwork PetscFunctionBegin; 85420282da2SJoe Wallwork ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr); 85520282da2SJoe Wallwork PetscFunctionReturn(0); 85620282da2SJoe Wallwork } 85720282da2SJoe Wallwork 85820282da2SJoe Wallwork /* 85920282da2SJoe Wallwork DMPlexMetricAverage3 - Compute the unweighted average of three metrics 86020282da2SJoe Wallwork 86120282da2SJoe Wallwork Input Parameter: 86220282da2SJoe Wallwork + dm - The DM 86320282da2SJoe Wallwork . metric1 - The first metric to be averaged 86420282da2SJoe Wallwork . metric2 - The second metric to be averaged 86520282da2SJoe Wallwork - metric3 - The third metric to be averaged 86620282da2SJoe Wallwork 86720282da2SJoe Wallwork Output Parameter: 86820282da2SJoe Wallwork . metricAvg - The averaged metric 86920282da2SJoe Wallwork 87020282da2SJoe Wallwork Level: beginner 87120282da2SJoe Wallwork 87220282da2SJoe Wallwork .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2() 87320282da2SJoe Wallwork */ 87420282da2SJoe Wallwork PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg) 87520282da2SJoe Wallwork { 87620282da2SJoe Wallwork PetscErrorCode ierr; 87720282da2SJoe Wallwork PetscReal weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0}; 87820282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 87920282da2SJoe Wallwork 88020282da2SJoe Wallwork PetscFunctionBegin; 88120282da2SJoe Wallwork ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr); 88220282da2SJoe Wallwork PetscFunctionReturn(0); 88320282da2SJoe Wallwork } 88420282da2SJoe Wallwork 88520282da2SJoe Wallwork static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[]) 88620282da2SJoe Wallwork { 88720282da2SJoe Wallwork PetscErrorCode ierr; 88820282da2SJoe Wallwork PetscInt i, j, k, l, m; 88920282da2SJoe Wallwork PetscReal *evals, *evals1; 89020282da2SJoe Wallwork PetscScalar *evecs, *sqrtM1, *isqrtM1; 89120282da2SJoe Wallwork 89220282da2SJoe Wallwork PetscFunctionBegin; 89320282da2SJoe Wallwork ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr); 89420282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 89520282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 89620282da2SJoe Wallwork evecs[i*dim+j] = M1[i*dim+j]; 89720282da2SJoe Wallwork } 89820282da2SJoe Wallwork } 89920282da2SJoe Wallwork { 90020282da2SJoe Wallwork PetscScalar *work; 90120282da2SJoe Wallwork PetscBLASInt lwork; 90220282da2SJoe Wallwork 90320282da2SJoe Wallwork lwork = 5*dim; 90420282da2SJoe Wallwork ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr); 90520282da2SJoe Wallwork { 90620282da2SJoe Wallwork PetscBLASInt lierr, nb; 90720282da2SJoe Wallwork PetscReal sqrtk; 90820282da2SJoe Wallwork 90920282da2SJoe Wallwork /* Compute eigendecomposition of M1 */ 91020282da2SJoe Wallwork ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr); 91120282da2SJoe Wallwork ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 91220282da2SJoe Wallwork #if defined(PETSC_USE_COMPLEX) 91320282da2SJoe Wallwork { 91420282da2SJoe Wallwork PetscReal *rwork; 91520282da2SJoe Wallwork ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 91620282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr)); 91720282da2SJoe Wallwork ierr = PetscFree(rwork);CHKERRQ(ierr); 91820282da2SJoe Wallwork } 91920282da2SJoe Wallwork #else 92020282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr)); 92120282da2SJoe Wallwork #endif 92220282da2SJoe Wallwork if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 92320282da2SJoe Wallwork ierr = PetscFPTrapPop(); 92420282da2SJoe Wallwork 92520282da2SJoe Wallwork /* Compute square root and reciprocal */ 92620282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 92720282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 92820282da2SJoe Wallwork sqrtM1[i*dim+j] = 0.0; 92920282da2SJoe Wallwork isqrtM1[i*dim+j] = 0.0; 93020282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 93120282da2SJoe Wallwork sqrtk = PetscSqrtReal(evals1[k]); 93220282da2SJoe Wallwork sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j]; 93320282da2SJoe Wallwork isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j]; 93420282da2SJoe Wallwork } 93520282da2SJoe Wallwork } 93620282da2SJoe Wallwork } 93720282da2SJoe Wallwork 93820282da2SJoe Wallwork /* Map into the space spanned by the eigenvectors of M1 */ 93920282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 94020282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 94120282da2SJoe Wallwork evecs[i*dim+j] = 0.0; 94220282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 94320282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 94420282da2SJoe Wallwork evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 94520282da2SJoe Wallwork } 94620282da2SJoe Wallwork } 94720282da2SJoe Wallwork } 94820282da2SJoe Wallwork } 94920282da2SJoe Wallwork 95020282da2SJoe Wallwork /* Compute eigendecomposition */ 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, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 95720282da2SJoe Wallwork ierr = PetscFree(rwork);CHKERRQ(ierr); 95820282da2SJoe Wallwork } 95920282da2SJoe Wallwork #else 96020282da2SJoe Wallwork PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 96120282da2SJoe Wallwork #endif 96220282da2SJoe Wallwork if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 96320282da2SJoe Wallwork ierr = PetscFPTrapPop(); 96420282da2SJoe Wallwork 96520282da2SJoe Wallwork /* Modify eigenvalues */ 96620282da2SJoe Wallwork for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]); 96720282da2SJoe Wallwork 96820282da2SJoe Wallwork /* Map back to get the intersection */ 96920282da2SJoe Wallwork for (i = 0; i < dim; ++i) { 97020282da2SJoe Wallwork for (j = 0; j < dim; ++j) { 97120282da2SJoe Wallwork M2[i*dim+j] = 0.0; 97220282da2SJoe Wallwork for (k = 0; k < dim; ++k) { 97320282da2SJoe Wallwork for (l = 0; l < dim; ++l) { 97420282da2SJoe Wallwork for (m = 0; m < dim; ++m) { 97520282da2SJoe Wallwork M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m]; 97620282da2SJoe Wallwork } 97720282da2SJoe Wallwork } 97820282da2SJoe Wallwork } 97920282da2SJoe Wallwork } 98020282da2SJoe Wallwork } 98120282da2SJoe Wallwork } 98220282da2SJoe Wallwork ierr = PetscFree(work);CHKERRQ(ierr); 98320282da2SJoe Wallwork } 98420282da2SJoe Wallwork ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr); 98520282da2SJoe Wallwork PetscFunctionReturn(0); 98620282da2SJoe Wallwork } 98720282da2SJoe Wallwork 98820282da2SJoe Wallwork /* 98920282da2SJoe Wallwork DMPlexMetricIntersection - Compute the intersection of a list of metrics 99020282da2SJoe Wallwork 99120282da2SJoe Wallwork Input Parameter: 99220282da2SJoe Wallwork + dm - The DM 99320282da2SJoe Wallwork . numMetrics - The number of metrics to be intersected 99420282da2SJoe Wallwork - metrics - The metrics to be intersected 99520282da2SJoe Wallwork 99620282da2SJoe Wallwork Output Parameter: 99720282da2SJoe Wallwork . metricInt - The intersected metric 99820282da2SJoe Wallwork 99920282da2SJoe Wallwork Level: beginner 100020282da2SJoe Wallwork 100120282da2SJoe Wallwork Notes: 100220282da2SJoe Wallwork 100320282da2SJoe Wallwork The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics. 100420282da2SJoe Wallwork 100520282da2SJoe Wallwork The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2. 100620282da2SJoe Wallwork 100720282da2SJoe Wallwork .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage() 100820282da2SJoe Wallwork */ 100920282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt) 101020282da2SJoe Wallwork { 101120282da2SJoe Wallwork PetscErrorCode ierr; 101220282da2SJoe Wallwork PetscInt dim, vStart, vEnd, v, i; 101320282da2SJoe Wallwork PetscScalar *met, *meti, *M, *Mi; 101420282da2SJoe Wallwork 101520282da2SJoe Wallwork PetscFunctionBegin; 101620282da2SJoe Wallwork if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); } 101720282da2SJoe Wallwork 101820282da2SJoe Wallwork /* Extract metadata from dm */ 101920282da2SJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 102020282da2SJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr); 102120282da2SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 102220282da2SJoe Wallwork 102320282da2SJoe Wallwork /* Copy over the first metric */ 102420282da2SJoe Wallwork ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr); 102520282da2SJoe Wallwork 102620282da2SJoe Wallwork /* Intersect subsequent metrics in turn */ 102720282da2SJoe Wallwork if (numMetrics > 1) { 102820282da2SJoe Wallwork ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr); 102920282da2SJoe Wallwork for (i = 1; i < numMetrics; ++i) { 103020282da2SJoe Wallwork ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr); 103120282da2SJoe Wallwork for (v = vStart; v < vEnd; ++v) { 103220282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr); 103320282da2SJoe Wallwork ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr); 103420282da2SJoe Wallwork ierr = DMPlexMetricIntersection_Private(dim, Mi, M);CHKERRQ(ierr); 103520282da2SJoe Wallwork } 103620282da2SJoe Wallwork ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr); 103720282da2SJoe Wallwork } 103820282da2SJoe Wallwork ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr); 103920282da2SJoe Wallwork } 104020282da2SJoe Wallwork 104120282da2SJoe Wallwork PetscFunctionReturn(0); 104220282da2SJoe Wallwork } 104320282da2SJoe Wallwork 104420282da2SJoe Wallwork /* 104520282da2SJoe Wallwork DMPlexMetricIntersection2 - Compute the intersection of two metrics 104620282da2SJoe Wallwork 104720282da2SJoe Wallwork Input Parameter: 104820282da2SJoe Wallwork + dm - The DM 104920282da2SJoe Wallwork . metric1 - The first metric to be intersected 105020282da2SJoe Wallwork - metric2 - The second metric to be intersected 105120282da2SJoe Wallwork 105220282da2SJoe Wallwork Output Parameter: 105320282da2SJoe Wallwork . metricInt - The intersected metric 105420282da2SJoe Wallwork 105520282da2SJoe Wallwork Level: beginner 105620282da2SJoe Wallwork 105720282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3() 105820282da2SJoe Wallwork */ 105920282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt) 106020282da2SJoe Wallwork { 106120282da2SJoe Wallwork PetscErrorCode ierr; 106220282da2SJoe Wallwork Vec metrics[2] = {metric1, metric2}; 106320282da2SJoe Wallwork 106420282da2SJoe Wallwork PetscFunctionBegin; 106520282da2SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr); 106620282da2SJoe Wallwork PetscFunctionReturn(0); 106720282da2SJoe Wallwork } 106820282da2SJoe Wallwork 106920282da2SJoe Wallwork /* 107020282da2SJoe Wallwork DMPlexMetricIntersection3 - Compute the intersection of three metrics 107120282da2SJoe Wallwork 107220282da2SJoe Wallwork Input Parameter: 107320282da2SJoe Wallwork + dm - The DM 107420282da2SJoe Wallwork . metric1 - The first metric to be intersected 107520282da2SJoe Wallwork . metric2 - The second metric to be intersected 107620282da2SJoe Wallwork - metric3 - The third metric to be intersected 107720282da2SJoe Wallwork 107820282da2SJoe Wallwork Output Parameter: 107920282da2SJoe Wallwork . metricInt - The intersected metric 108020282da2SJoe Wallwork 108120282da2SJoe Wallwork Level: beginner 108220282da2SJoe Wallwork 108320282da2SJoe Wallwork .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2() 108420282da2SJoe Wallwork */ 108520282da2SJoe Wallwork PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt) 108620282da2SJoe Wallwork { 108720282da2SJoe Wallwork PetscErrorCode ierr; 108820282da2SJoe Wallwork Vec metrics[3] = {metric1, metric2, metric3}; 108920282da2SJoe Wallwork 109020282da2SJoe Wallwork PetscFunctionBegin; 109120282da2SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr); 109220282da2SJoe Wallwork PetscFunctionReturn(0); 109320282da2SJoe Wallwork } 1094