1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscblaslapack.h> 3 4 PetscLogEvent DMPLEX_MetricEnforceSPD, DMPLEX_MetricNormalize, DMPLEX_MetricAverage, DMPLEX_MetricIntersection; 5 6 PetscErrorCode DMPlexMetricSetFromOptions(DM dm) 7 { 8 MPI_Comm comm; 9 PetscBool isotropic = PETSC_FALSE, uniform = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE; 10 PetscBool noInsert = PETSC_FALSE, noSwap = PETSC_FALSE, noMove = PETSC_FALSE, noSurf = PETSC_FALSE; 11 PetscErrorCode ierr; 12 PetscInt verbosity = -1, numIter = 3; 13 PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0, beta = 1.3, hausd = 0.01; 14 15 PetscFunctionBegin; 16 CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 17 ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr); 18 CHKERRQ(PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL)); 19 CHKERRQ(DMPlexMetricSetIsotropic(dm, isotropic)); 20 CHKERRQ(PetscOptionsBool("-dm_plex_metric_uniform", "Is the metric uniform?", "DMPlexMetricCreateUniform", uniform, &uniform, NULL)); 21 CHKERRQ(DMPlexMetricSetUniform(dm, uniform)); 22 CHKERRQ(PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL)); 23 CHKERRQ(DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst)); 24 CHKERRQ(PetscOptionsBool("-dm_plex_metric_no_insert", "Turn off node insertion and deletion", "DMAdaptMetric", noInsert, &noInsert, NULL)); 25 CHKERRQ(DMPlexMetricSetNoInsertion(dm, noInsert)); 26 CHKERRQ(PetscOptionsBool("-dm_plex_metric_no_swap", "Turn off facet swapping", "DMAdaptMetric", noSwap, &noSwap, NULL)); 27 CHKERRQ(DMPlexMetricSetNoSwapping(dm, noSwap)); 28 CHKERRQ(PetscOptionsBool("-dm_plex_metric_no_move", "Turn off facet node movement", "DMAdaptMetric", noMove, &noMove, NULL)); 29 CHKERRQ(DMPlexMetricSetNoMovement(dm, noMove)); 30 CHKERRQ(PetscOptionsBool("-dm_plex_metric_no_surf", "Turn off surface modification", "DMAdaptMetric", noSurf, &noSurf, NULL)); 31 CHKERRQ(DMPlexMetricSetNoSurf(dm, noSurf)); 32 CHKERRQ(PetscOptionsBoundedInt("-dm_plex_metric_num_iterations", "Number of ParMmg adaptation iterations", "DMAdaptMetric", numIter, &numIter, NULL, 0)); 33 CHKERRQ(DMPlexMetricSetNumIterations(dm, numIter)); 34 CHKERRQ(PetscOptionsRangeInt("-dm_plex_metric_verbosity", "Verbosity of metric-based mesh adaptation package (-1 = silent, 10 = maximum)", "DMAdaptMetric", verbosity, &verbosity, NULL, -1, 10)); 35 CHKERRQ(DMPlexMetricSetVerbosity(dm, verbosity)); 36 CHKERRQ(PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL)); 37 CHKERRQ(DMPlexMetricSetMinimumMagnitude(dm, h_min)); 38 CHKERRQ(PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL)); 39 CHKERRQ(DMPlexMetricSetMaximumMagnitude(dm, h_max)); 40 CHKERRQ(PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL)); 41 CHKERRQ(DMPlexMetricSetMaximumAnisotropy(dm, a_max)); 42 CHKERRQ(PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL)); 43 CHKERRQ(DMPlexMetricSetNormalizationOrder(dm, p)); 44 CHKERRQ(PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL)); 45 CHKERRQ(DMPlexMetricSetTargetComplexity(dm, target)); 46 CHKERRQ(PetscOptionsReal("-dm_plex_metric_gradation_factor", "Metric gradation factor", "DMAdaptMetric", beta, &beta, NULL)); 47 CHKERRQ(DMPlexMetricSetGradationFactor(dm, beta)); 48 CHKERRQ(PetscOptionsReal("-dm_plex_metric_hausdorff_number", "Metric Hausdorff number", "DMAdaptMetric", hausd, &hausd, NULL)); 49 CHKERRQ(DMPlexMetricSetHausdorffNumber(dm, hausd)); 50 ierr = PetscOptionsEnd();CHKERRQ(ierr); 51 PetscFunctionReturn(0); 52 } 53 54 /*@ 55 DMPlexMetricSetIsotropic - Record whether a metric is isotropic 56 57 Input parameters: 58 + dm - The DM 59 - isotropic - Is the metric isotropic? 60 61 Level: beginner 62 63 .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetUniform(), DMPlexMetricSetRestrictAnisotropyFirst() 64 @*/ 65 PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic) 66 { 67 DM_Plex *plex = (DM_Plex *) dm->data; 68 69 PetscFunctionBegin; 70 if (!plex->metricCtx) { 71 CHKERRQ(PetscNew(&plex->metricCtx)); 72 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 73 } 74 plex->metricCtx->isotropic = isotropic; 75 PetscFunctionReturn(0); 76 } 77 78 /*@ 79 DMPlexMetricIsIsotropic - Is a metric isotropic? 80 81 Input parameters: 82 . dm - The DM 83 84 Output parameters: 85 . isotropic - Is the metric isotropic? 86 87 Level: beginner 88 89 .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricIsUniform(), DMPlexMetricRestrictAnisotropyFirst() 90 @*/ 91 PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic) 92 { 93 DM_Plex *plex = (DM_Plex *) dm->data; 94 95 PetscFunctionBegin; 96 if (!plex->metricCtx) { 97 CHKERRQ(PetscNew(&plex->metricCtx)); 98 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 99 } 100 *isotropic = plex->metricCtx->isotropic; 101 PetscFunctionReturn(0); 102 } 103 104 /*@ 105 DMPlexMetricSetUniform - Record whether a metric is uniform 106 107 Input parameters: 108 + dm - The DM 109 - uniform - Is the metric uniform? 110 111 Level: beginner 112 113 Notes: 114 115 If the metric is specified as uniform then it is assumed to be isotropic, too. 116 117 .seealso: DMPlexMetricIsUniform(), DMPlexMetricSetIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 118 @*/ 119 PetscErrorCode DMPlexMetricSetUniform(DM dm, PetscBool uniform) 120 { 121 DM_Plex *plex = (DM_Plex *) dm->data; 122 123 PetscFunctionBegin; 124 if (!plex->metricCtx) { 125 CHKERRQ(PetscNew(&plex->metricCtx)); 126 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 127 } 128 plex->metricCtx->uniform = uniform; 129 if (uniform) plex->metricCtx->isotropic = uniform; 130 PetscFunctionReturn(0); 131 } 132 133 /*@ 134 DMPlexMetricIsUniform - Is a metric uniform? 135 136 Input parameters: 137 . dm - The DM 138 139 Output parameters: 140 . uniform - Is the metric uniform? 141 142 Level: beginner 143 144 .seealso: DMPlexMetricSetUniform(), DMPlexMetricIsIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 145 @*/ 146 PetscErrorCode DMPlexMetricIsUniform(DM dm, PetscBool *uniform) 147 { 148 DM_Plex *plex = (DM_Plex *) dm->data; 149 150 PetscFunctionBegin; 151 if (!plex->metricCtx) { 152 CHKERRQ(PetscNew(&plex->metricCtx)); 153 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 154 } 155 *uniform = plex->metricCtx->uniform; 156 PetscFunctionReturn(0); 157 } 158 159 /*@ 160 DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization 161 162 Input parameters: 163 + dm - The DM 164 - restrictAnisotropyFirst - Should anisotropy be normalized first? 165 166 Level: beginner 167 168 .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 169 @*/ 170 PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst) 171 { 172 DM_Plex *plex = (DM_Plex *) dm->data; 173 174 PetscFunctionBegin; 175 if (!plex->metricCtx) { 176 CHKERRQ(PetscNew(&plex->metricCtx)); 177 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 178 } 179 plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst; 180 PetscFunctionReturn(0); 181 } 182 183 /*@ 184 DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after? 185 186 Input parameters: 187 . dm - The DM 188 189 Output parameters: 190 . restrictAnisotropyFirst - Is anisotropy be normalized first? 191 192 Level: beginner 193 194 .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 195 @*/ 196 PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst) 197 { 198 DM_Plex *plex = (DM_Plex *) dm->data; 199 200 PetscFunctionBegin; 201 if (!plex->metricCtx) { 202 CHKERRQ(PetscNew(&plex->metricCtx)); 203 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 204 } 205 *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst; 206 PetscFunctionReturn(0); 207 } 208 209 /*@ 210 DMPlexMetricSetNoInsertion - Should node insertion and deletion be turned off? 211 212 Input parameters: 213 + dm - The DM 214 - noInsert - Should node insertion and deletion be turned off? 215 216 Level: beginner 217 218 Notes: 219 This is only used by Mmg and ParMmg (not Pragmatic). 220 221 .seealso: DMPlexMetricNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf() 222 @*/ 223 PetscErrorCode DMPlexMetricSetNoInsertion(DM dm, PetscBool noInsert) 224 { 225 DM_Plex *plex = (DM_Plex *) dm->data; 226 227 PetscFunctionBegin; 228 if (!plex->metricCtx) { 229 CHKERRQ(PetscNew(&plex->metricCtx)); 230 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 231 } 232 plex->metricCtx->noInsert = noInsert; 233 PetscFunctionReturn(0); 234 } 235 236 /*@ 237 DMPlexMetricNoInsertion - Are node insertion and deletion turned off? 238 239 Input parameters: 240 . dm - The DM 241 242 Output parameters: 243 . noInsert - Are node insertion and deletion turned off? 244 245 Level: beginner 246 247 Notes: 248 This is only used by Mmg and ParMmg (not Pragmatic). 249 250 .seealso: DMPlexMetricSetNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf() 251 @*/ 252 PetscErrorCode DMPlexMetricNoInsertion(DM dm, PetscBool *noInsert) 253 { 254 DM_Plex *plex = (DM_Plex *) dm->data; 255 256 PetscFunctionBegin; 257 if (!plex->metricCtx) { 258 CHKERRQ(PetscNew(&plex->metricCtx)); 259 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 260 } 261 *noInsert = plex->metricCtx->noInsert; 262 PetscFunctionReturn(0); 263 } 264 265 /*@ 266 DMPlexMetricSetNoSwapping - Should facet swapping be turned off? 267 268 Input parameters: 269 + dm - The DM 270 - noSwap - Should facet swapping be turned off? 271 272 Level: beginner 273 274 Notes: 275 This is only used by Mmg and ParMmg (not Pragmatic). 276 277 .seealso: DMPlexMetricNoSwapping(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoSurf() 278 @*/ 279 PetscErrorCode DMPlexMetricSetNoSwapping(DM dm, PetscBool noSwap) 280 { 281 DM_Plex *plex = (DM_Plex *) dm->data; 282 283 PetscFunctionBegin; 284 if (!plex->metricCtx) { 285 CHKERRQ(PetscNew(&plex->metricCtx)); 286 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 287 } 288 plex->metricCtx->noSwap = noSwap; 289 PetscFunctionReturn(0); 290 } 291 292 /*@ 293 DMPlexMetricNoSwapping - Is facet swapping turned off? 294 295 Input parameters: 296 . dm - The DM 297 298 Output parameters: 299 . noSwap - Is facet swapping turned off? 300 301 Level: beginner 302 303 Notes: 304 This is only used by Mmg and ParMmg (not Pragmatic). 305 306 .seealso: DMPlexMetricSetNoSwapping(), DMPlexMetricNoInsertion(), DMPlexMetricNoMovement(), DMPlexMetricNoSurf() 307 @*/ 308 PetscErrorCode DMPlexMetricNoSwapping(DM dm, PetscBool *noSwap) 309 { 310 DM_Plex *plex = (DM_Plex *) dm->data; 311 312 PetscFunctionBegin; 313 if (!plex->metricCtx) { 314 CHKERRQ(PetscNew(&plex->metricCtx)); 315 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 316 } 317 *noSwap = plex->metricCtx->noSwap; 318 PetscFunctionReturn(0); 319 } 320 321 /*@ 322 DMPlexMetricSetNoMovement - Should node movement be turned off? 323 324 Input parameters: 325 + dm - The DM 326 - noMove - Should node movement be turned off? 327 328 Level: beginner 329 330 Notes: 331 This is only used by Mmg and ParMmg (not Pragmatic). 332 333 .seealso: DMPlexMetricNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping(), DMPlexMetricSetNoSurf() 334 @*/ 335 PetscErrorCode DMPlexMetricSetNoMovement(DM dm, PetscBool noMove) 336 { 337 DM_Plex *plex = (DM_Plex *) dm->data; 338 339 PetscFunctionBegin; 340 if (!plex->metricCtx) { 341 CHKERRQ(PetscNew(&plex->metricCtx)); 342 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 343 } 344 plex->metricCtx->noMove = noMove; 345 PetscFunctionReturn(0); 346 } 347 348 /*@ 349 DMPlexMetricNoMovement - Is node movement turned off? 350 351 Input parameters: 352 . dm - The DM 353 354 Output parameters: 355 . noMove - Is node movement turned off? 356 357 Level: beginner 358 359 Notes: 360 This is only used by Mmg and ParMmg (not Pragmatic). 361 362 .seealso: DMPlexMetricSetNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping(), DMPlexMetricNoSurf() 363 @*/ 364 PetscErrorCode DMPlexMetricNoMovement(DM dm, PetscBool *noMove) 365 { 366 DM_Plex *plex = (DM_Plex *) dm->data; 367 368 PetscFunctionBegin; 369 if (!plex->metricCtx) { 370 CHKERRQ(PetscNew(&plex->metricCtx)); 371 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 372 } 373 *noMove = plex->metricCtx->noMove; 374 PetscFunctionReturn(0); 375 } 376 377 /*@ 378 DMPlexMetricSetNoSurf - Should surface modification be turned off? 379 380 Input parameters: 381 + dm - The DM 382 - noSurf - Should surface modification be turned off? 383 384 Level: beginner 385 386 Notes: 387 This is only used by Mmg and ParMmg (not Pragmatic). 388 389 .seealso: DMPlexMetricNoSurf(), DMPlexMetricSetNoMovement(), DMPlexMetricSetNoInsertion(), DMPlexMetricSetNoSwapping() 390 @*/ 391 PetscErrorCode DMPlexMetricSetNoSurf(DM dm, PetscBool noSurf) 392 { 393 DM_Plex *plex = (DM_Plex *) dm->data; 394 395 PetscFunctionBegin; 396 if (!plex->metricCtx) { 397 CHKERRQ(PetscNew(&plex->metricCtx)); 398 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 399 } 400 plex->metricCtx->noSurf = noSurf; 401 PetscFunctionReturn(0); 402 } 403 404 /*@ 405 DMPlexMetricNoSurf - Is surface modification turned off? 406 407 Input parameters: 408 . dm - The DM 409 410 Output parameters: 411 . noSurf - Is surface modification turned off? 412 413 Level: beginner 414 415 Notes: 416 This is only used by Mmg and ParMmg (not Pragmatic). 417 418 .seealso: DMPlexMetricSetNoSurf(), DMPlexMetricNoMovement(), DMPlexMetricNoInsertion(), DMPlexMetricNoSwapping() 419 @*/ 420 PetscErrorCode DMPlexMetricNoSurf(DM dm, PetscBool *noSurf) 421 { 422 DM_Plex *plex = (DM_Plex *) dm->data; 423 424 PetscFunctionBegin; 425 if (!plex->metricCtx) { 426 CHKERRQ(PetscNew(&plex->metricCtx)); 427 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 428 } 429 *noSurf = plex->metricCtx->noSurf; 430 PetscFunctionReturn(0); 431 } 432 433 /*@ 434 DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude 435 436 Input parameters: 437 + dm - The DM 438 - h_min - The minimum tolerated metric magnitude 439 440 Level: beginner 441 442 .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude() 443 @*/ 444 PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min) 445 { 446 DM_Plex *plex = (DM_Plex *) dm->data; 447 448 PetscFunctionBegin; 449 if (!plex->metricCtx) { 450 CHKERRQ(PetscNew(&plex->metricCtx)); 451 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 452 } 453 PetscCheckFalse(h_min <= 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min); 454 plex->metricCtx->h_min = h_min; 455 PetscFunctionReturn(0); 456 } 457 458 /*@ 459 DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude 460 461 Input parameters: 462 . dm - The DM 463 464 Output parameters: 465 . h_min - The minimum tolerated metric magnitude 466 467 Level: beginner 468 469 .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude() 470 @*/ 471 PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min) 472 { 473 DM_Plex *plex = (DM_Plex *) dm->data; 474 475 PetscFunctionBegin; 476 if (!plex->metricCtx) { 477 CHKERRQ(PetscNew(&plex->metricCtx)); 478 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 479 } 480 *h_min = plex->metricCtx->h_min; 481 PetscFunctionReturn(0); 482 } 483 484 /*@ 485 DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude 486 487 Input parameters: 488 + dm - The DM 489 - h_max - The maximum tolerated metric magnitude 490 491 Level: beginner 492 493 .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude() 494 @*/ 495 PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max) 496 { 497 DM_Plex *plex = (DM_Plex *) dm->data; 498 499 PetscFunctionBegin; 500 if (!plex->metricCtx) { 501 CHKERRQ(PetscNew(&plex->metricCtx)); 502 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 503 } 504 PetscCheckFalse(h_max <= 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max); 505 plex->metricCtx->h_max = h_max; 506 PetscFunctionReturn(0); 507 } 508 509 /*@ 510 DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude 511 512 Input parameters: 513 . dm - The DM 514 515 Output parameters: 516 . h_max - The maximum tolerated metric magnitude 517 518 Level: beginner 519 520 .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude() 521 @*/ 522 PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max) 523 { 524 DM_Plex *plex = (DM_Plex *) dm->data; 525 526 PetscFunctionBegin; 527 if (!plex->metricCtx) { 528 CHKERRQ(PetscNew(&plex->metricCtx)); 529 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 530 } 531 *h_max = plex->metricCtx->h_max; 532 PetscFunctionReturn(0); 533 } 534 535 /*@ 536 DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy 537 538 Input parameters: 539 + dm - The DM 540 - a_max - The maximum tolerated metric anisotropy 541 542 Level: beginner 543 544 Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one. 545 546 .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude() 547 @*/ 548 PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max) 549 { 550 DM_Plex *plex = (DM_Plex *) dm->data; 551 552 PetscFunctionBegin; 553 if (!plex->metricCtx) { 554 CHKERRQ(PetscNew(&plex->metricCtx)); 555 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 556 } 557 PetscCheckFalse(a_max < 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max); 558 plex->metricCtx->a_max = a_max; 559 PetscFunctionReturn(0); 560 } 561 562 /*@ 563 DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy 564 565 Input parameters: 566 . dm - The DM 567 568 Output parameters: 569 . a_max - The maximum tolerated metric anisotropy 570 571 Level: beginner 572 573 .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude() 574 @*/ 575 PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max) 576 { 577 DM_Plex *plex = (DM_Plex *) dm->data; 578 579 PetscFunctionBegin; 580 if (!plex->metricCtx) { 581 CHKERRQ(PetscNew(&plex->metricCtx)); 582 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 583 } 584 *a_max = plex->metricCtx->a_max; 585 PetscFunctionReturn(0); 586 } 587 588 /*@ 589 DMPlexMetricSetTargetComplexity - Set the target metric complexity 590 591 Input parameters: 592 + dm - The DM 593 - targetComplexity - The target metric complexity 594 595 Level: beginner 596 597 .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder() 598 @*/ 599 PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity) 600 { 601 DM_Plex *plex = (DM_Plex *) dm->data; 602 603 PetscFunctionBegin; 604 if (!plex->metricCtx) { 605 CHKERRQ(PetscNew(&plex->metricCtx)); 606 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 607 } 608 PetscCheckFalse(targetComplexity <= 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive"); 609 plex->metricCtx->targetComplexity = targetComplexity; 610 PetscFunctionReturn(0); 611 } 612 613 /*@ 614 DMPlexMetricGetTargetComplexity - Get the target metric complexity 615 616 Input parameters: 617 . dm - The DM 618 619 Output parameters: 620 . targetComplexity - The target metric complexity 621 622 Level: beginner 623 624 .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder() 625 @*/ 626 PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity) 627 { 628 DM_Plex *plex = (DM_Plex *) dm->data; 629 630 PetscFunctionBegin; 631 if (!plex->metricCtx) { 632 CHKERRQ(PetscNew(&plex->metricCtx)); 633 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 634 } 635 *targetComplexity = plex->metricCtx->targetComplexity; 636 PetscFunctionReturn(0); 637 } 638 639 /*@ 640 DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization 641 642 Input parameters: 643 + dm - The DM 644 - p - The normalization order 645 646 Level: beginner 647 648 .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity() 649 @*/ 650 PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p) 651 { 652 DM_Plex *plex = (DM_Plex *) dm->data; 653 654 PetscFunctionBegin; 655 if (!plex->metricCtx) { 656 CHKERRQ(PetscNew(&plex->metricCtx)); 657 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 658 } 659 PetscCheckFalse(p < 1.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater"); 660 plex->metricCtx->p = p; 661 PetscFunctionReturn(0); 662 } 663 664 /*@ 665 DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization 666 667 Input parameters: 668 . dm - The DM 669 670 Output parameters: 671 . p - The normalization order 672 673 Level: beginner 674 675 .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity() 676 @*/ 677 PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p) 678 { 679 DM_Plex *plex = (DM_Plex *) dm->data; 680 681 PetscFunctionBegin; 682 if (!plex->metricCtx) { 683 CHKERRQ(PetscNew(&plex->metricCtx)); 684 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 685 } 686 *p = plex->metricCtx->p; 687 PetscFunctionReturn(0); 688 } 689 690 /*@ 691 DMPlexMetricSetGradationFactor - Set the metric gradation factor 692 693 Input parameters: 694 + dm - The DM 695 - beta - The metric gradation factor 696 697 Level: beginner 698 699 Notes: 700 701 The gradation factor is the maximum tolerated length ratio between adjacent edges. 702 703 Turn off gradation by passing the value -1. Otherwise, pass a positive value. 704 705 This is only used by Mmg and ParMmg (not Pragmatic). 706 707 .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber() 708 @*/ 709 PetscErrorCode DMPlexMetricSetGradationFactor(DM dm, PetscReal beta) 710 { 711 DM_Plex *plex = (DM_Plex *) dm->data; 712 713 PetscFunctionBegin; 714 if (!plex->metricCtx) { 715 CHKERRQ(PetscNew(&plex->metricCtx)); 716 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 717 } 718 plex->metricCtx->gradationFactor = beta; 719 PetscFunctionReturn(0); 720 } 721 722 /*@ 723 DMPlexMetricGetGradationFactor - Get the metric gradation factor 724 725 Input parameters: 726 . dm - The DM 727 728 Output parameters: 729 . beta - The metric gradation factor 730 731 Level: beginner 732 733 Notes: 734 735 The gradation factor is the maximum tolerated length ratio between adjacent edges. 736 737 The value -1 implies that gradation is turned off. 738 739 This is only used by Mmg and ParMmg (not Pragmatic). 740 741 .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber() 742 @*/ 743 PetscErrorCode DMPlexMetricGetGradationFactor(DM dm, PetscReal *beta) 744 { 745 DM_Plex *plex = (DM_Plex *) dm->data; 746 747 PetscFunctionBegin; 748 if (!plex->metricCtx) { 749 CHKERRQ(PetscNew(&plex->metricCtx)); 750 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 751 } 752 *beta = plex->metricCtx->gradationFactor; 753 PetscFunctionReturn(0); 754 } 755 756 /*@ 757 DMPlexMetricSetHausdorffNumber - Set the metric Hausdorff number 758 759 Input parameters: 760 + dm - The DM 761 - hausd - The metric Hausdorff number 762 763 Level: beginner 764 765 Notes: 766 767 The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the 768 boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the 769 high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for 770 an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease 771 (resp. increase) the Hausdorff parameter. (Taken from 772 https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd). 773 774 This is only used by Mmg and ParMmg (not Pragmatic). 775 776 .seealso: DMPlexMetricSetGradationFactor(), DMPlexMetricGetHausdorffNumber() 777 @*/ 778 PetscErrorCode DMPlexMetricSetHausdorffNumber(DM dm, PetscReal hausd) 779 { 780 DM_Plex *plex = (DM_Plex *) dm->data; 781 782 PetscFunctionBegin; 783 if (!plex->metricCtx) { 784 CHKERRQ(PetscNew(&plex->metricCtx)); 785 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 786 } 787 plex->metricCtx->hausdorffNumber = hausd; 788 PetscFunctionReturn(0); 789 } 790 791 /*@ 792 DMPlexMetricGetHausdorffNumber - Get the metric Hausdorff number 793 794 Input parameters: 795 . dm - The DM 796 797 Output parameters: 798 . hausd - The metric Hausdorff number 799 800 Level: beginner 801 802 Notes: 803 804 The Hausdorff number imposes the maximal distance between the piecewise linear approximation of the 805 boundary and the reconstructed ideal boundary. Thus, a low Hausdorff parameter leads to refine the 806 high curvature areas. By default, the Hausdorff value is set to 0.01, which is a suitable value for 807 an object of size 1 in each direction. For smaller (resp. larger) objects, you may need to decrease 808 (resp. increase) the Hausdorff parameter. (Taken from 809 https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd). 810 811 This is only used by Mmg and ParMmg (not Pragmatic). 812 813 .seealso: DMPlexMetricGetGradationFactor(), DMPlexMetricSetHausdorffNumber() 814 @*/ 815 PetscErrorCode DMPlexMetricGetHausdorffNumber(DM dm, PetscReal *hausd) 816 { 817 DM_Plex *plex = (DM_Plex *) dm->data; 818 819 PetscFunctionBegin; 820 if (!plex->metricCtx) { 821 CHKERRQ(PetscNew(&plex->metricCtx)); 822 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 823 } 824 *hausd = plex->metricCtx->hausdorffNumber; 825 PetscFunctionReturn(0); 826 } 827 828 /*@ 829 DMPlexMetricSetVerbosity - Set the verbosity of the mesh adaptation package 830 831 Input parameters: 832 + dm - The DM 833 - verbosity - The verbosity, where -1 is silent and 10 is maximum 834 835 Level: beginner 836 837 Notes: 838 This is only used by Mmg and ParMmg (not Pragmatic). 839 840 .seealso: DMPlexMetricGetVerbosity(), DMPlexMetricSetNumIterations() 841 @*/ 842 PetscErrorCode DMPlexMetricSetVerbosity(DM dm, PetscInt verbosity) 843 { 844 DM_Plex *plex = (DM_Plex *) dm->data; 845 846 PetscFunctionBegin; 847 if (!plex->metricCtx) { 848 CHKERRQ(PetscNew(&plex->metricCtx)); 849 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 850 } 851 plex->metricCtx->verbosity = verbosity; 852 PetscFunctionReturn(0); 853 } 854 855 /*@ 856 DMPlexMetricGetVerbosity - Get the verbosity of the mesh adaptation package 857 858 Input parameters: 859 . dm - The DM 860 861 Output parameters: 862 . verbosity - The verbosity, where -1 is silent and 10 is maximum 863 864 Level: beginner 865 866 Notes: 867 This is only used by Mmg and ParMmg (not Pragmatic). 868 869 .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations() 870 @*/ 871 PetscErrorCode DMPlexMetricGetVerbosity(DM dm, PetscInt *verbosity) 872 { 873 DM_Plex *plex = (DM_Plex *) dm->data; 874 875 PetscFunctionBegin; 876 if (!plex->metricCtx) { 877 CHKERRQ(PetscNew(&plex->metricCtx)); 878 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 879 } 880 *verbosity = plex->metricCtx->verbosity; 881 PetscFunctionReturn(0); 882 } 883 884 /*@ 885 DMPlexMetricSetNumIterations - Set the number of parallel adaptation iterations 886 887 Input parameters: 888 + dm - The DM 889 - numIter - the number of parallel adaptation iterations 890 891 Level: beginner 892 893 Notes: 894 This is only used by ParMmg (not Pragmatic or Mmg). 895 896 .seealso: DMPlexMetricSetVerbosity(), DMPlexMetricGetNumIterations() 897 @*/ 898 PetscErrorCode DMPlexMetricSetNumIterations(DM dm, PetscInt numIter) 899 { 900 DM_Plex *plex = (DM_Plex *) dm->data; 901 902 PetscFunctionBegin; 903 if (!plex->metricCtx) { 904 CHKERRQ(PetscNew(&plex->metricCtx)); 905 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 906 } 907 plex->metricCtx->numIter = numIter; 908 PetscFunctionReturn(0); 909 } 910 911 /*@ 912 DMPlexMetricGetNumIterations - Get the number of parallel adaptation iterations 913 914 Input parameters: 915 . dm - The DM 916 917 Output parameters: 918 . numIter - the number of parallel adaptation iterations 919 920 Level: beginner 921 922 Notes: 923 This is only used by Mmg and ParMmg (not Pragmatic or Mmg). 924 925 .seealso: DMPlexMetricSetNumIterations(), DMPlexMetricGetVerbosity() 926 @*/ 927 PetscErrorCode DMPlexMetricGetNumIterations(DM dm, PetscInt *numIter) 928 { 929 DM_Plex *plex = (DM_Plex *) dm->data; 930 931 PetscFunctionBegin; 932 if (!plex->metricCtx) { 933 CHKERRQ(PetscNew(&plex->metricCtx)); 934 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 935 } 936 *numIter = plex->metricCtx->numIter; 937 PetscFunctionReturn(0); 938 } 939 940 PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric) 941 { 942 MPI_Comm comm; 943 PetscFE fe; 944 PetscInt dim; 945 946 PetscFunctionBegin; 947 948 /* Extract metadata from dm */ 949 CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 950 CHKERRQ(DMGetDimension(dm, &dim)); 951 952 /* Create a P1 field of the requested size */ 953 CHKERRQ(PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 954 CHKERRQ(DMSetField(dm, f, NULL, (PetscObject)fe)); 955 CHKERRQ(DMCreateDS(dm)); 956 CHKERRQ(PetscFEDestroy(&fe)); 957 CHKERRQ(DMCreateLocalVector(dm, metric)); 958 959 PetscFunctionReturn(0); 960 } 961 962 /*@ 963 DMPlexMetricCreate - Create a Riemannian metric field 964 965 Input parameters: 966 + dm - The DM 967 - f - The field number to use 968 969 Output parameter: 970 . metric - The metric 971 972 Level: beginner 973 974 Notes: 975 976 It is assumed that the DM is comprised of simplices. 977 978 Command line options for Riemannian metrics: 979 980 + -dm_plex_metric_isotropic - Is the metric isotropic? 981 . -dm_plex_metric_uniform - Is the metric uniform? 982 . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 983 . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 984 . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 985 . -dm_plex_metric_a_max - Maximum tolerated anisotropy 986 . -dm_plex_metric_p - L-p normalization order 987 - -dm_plex_metric_target_complexity - Target metric complexity 988 989 Switching between remeshers can be achieved using 990 991 . -dm_adaptor <pragmatic/mmg/parmmg> - specify dm adaptor to use 992 993 Further options that are only relevant to Mmg and ParMmg: 994 995 + -dm_plex_metric_gradation_factor - Maximum ratio by which edge lengths may grow during gradation 996 . -dm_plex_metric_num_iterations - Number of parallel mesh adaptation iterations for ParMmg 997 . -dm_plex_metric_no_insert - Should node insertion/deletion be turned off? 998 . -dm_plex_metric_no_swap - Should facet swapping be turned off? 999 . -dm_plex_metric_no_move - Should node movement be turned off? 1000 - -dm_plex_metric_verbosity - Choose a verbosity level from -1 (silent) to 10 (maximum). 1001 1002 .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic() 1003 @*/ 1004 PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric) 1005 { 1006 DM_Plex *plex = (DM_Plex *) dm->data; 1007 PetscBool isotropic, uniform; 1008 PetscInt coordDim, Nd; 1009 1010 PetscFunctionBegin; 1011 if (!plex->metricCtx) { 1012 CHKERRQ(PetscNew(&plex->metricCtx)); 1013 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 1014 } 1015 CHKERRQ(DMGetCoordinateDim(dm, &coordDim)); 1016 Nd = coordDim*coordDim; 1017 CHKERRQ(DMPlexMetricIsUniform(dm, &uniform)); 1018 CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic)); 1019 if (uniform) { 1020 MPI_Comm comm; 1021 1022 CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 1023 PetscCheck(isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported"); 1024 CHKERRQ(VecCreate(comm, metric)); 1025 CHKERRQ(VecSetSizes(*metric, 1, PETSC_DECIDE)); 1026 CHKERRQ(VecSetFromOptions(*metric)); 1027 } else if (isotropic) { 1028 CHKERRQ(DMPlexP1FieldCreate_Private(dm, f, 1, metric)); 1029 } else { 1030 CHKERRQ(DMPlexP1FieldCreate_Private(dm, f, Nd, metric)); 1031 } 1032 PetscFunctionReturn(0); 1033 } 1034 1035 /*@ 1036 DMPlexMetricCreateUniform - Construct a uniform isotropic metric 1037 1038 Input parameters: 1039 + dm - The DM 1040 . f - The field number to use 1041 - alpha - Scaling parameter for the diagonal 1042 1043 Output parameter: 1044 . metric - The uniform metric 1045 1046 Level: beginner 1047 1048 Note: In this case, the metric is constant in space and so is only specified for a single vertex. 1049 1050 .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic() 1051 @*/ 1052 PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric) 1053 { 1054 PetscFunctionBegin; 1055 CHKERRQ(DMPlexMetricSetUniform(dm, PETSC_TRUE)); 1056 CHKERRQ(DMPlexMetricCreate(dm, f, metric)); 1057 PetscCheck(alpha,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined"); 1058 PetscCheck(alpha >= 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha); 1059 CHKERRQ(VecSet(*metric, alpha)); 1060 CHKERRQ(VecAssemblyBegin(*metric)); 1061 CHKERRQ(VecAssemblyEnd(*metric)); 1062 PetscFunctionReturn(0); 1063 } 1064 1065 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1066 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1067 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1068 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 1069 { 1070 f0[0] = u[0]; 1071 } 1072 1073 /*@ 1074 DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator 1075 1076 Input parameters: 1077 + dm - The DM 1078 . f - The field number to use 1079 - indicator - The error indicator 1080 1081 Output parameter: 1082 . metric - The isotropic metric 1083 1084 Level: beginner 1085 1086 Notes: 1087 1088 It is assumed that the DM is comprised of simplices. 1089 1090 The indicator needs to be a scalar field. If it is not defined vertex-wise, then it is projected appropriately. 1091 1092 .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform() 1093 @*/ 1094 PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric) 1095 { 1096 PetscInt m, n; 1097 1098 PetscFunctionBegin; 1099 CHKERRQ(DMPlexMetricSetIsotropic(dm, PETSC_TRUE)); 1100 CHKERRQ(DMPlexMetricCreate(dm, f, metric)); 1101 CHKERRQ(VecGetSize(indicator, &m)); 1102 CHKERRQ(VecGetSize(*metric, &n)); 1103 if (m == n) CHKERRQ(VecCopy(indicator, *metric)); 1104 else { 1105 DM dmIndi; 1106 void (*funcs[1])(PetscInt, PetscInt, PetscInt, 1107 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1108 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1109 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 1110 1111 CHKERRQ(VecGetDM(indicator, &dmIndi)); 1112 funcs[0] = identity; 1113 CHKERRQ(DMProjectFieldLocal(dmIndi, 0.0, indicator, funcs, INSERT_VALUES, *metric)); 1114 } 1115 PetscFunctionReturn(0); 1116 } 1117 1118 static PetscErrorCode LAPACKsyevFail(PetscInt dim, PetscScalar Mpos[]) 1119 { 1120 PetscInt i, j; 1121 1122 PetscFunctionBegin; 1123 PetscPrintf(PETSC_COMM_SELF, "Failed to apply LAPACKsyev to the matrix\n"); 1124 for (i = 0; i < dim; ++i) { 1125 if (i == 0) PetscPrintf(PETSC_COMM_SELF, " [["); 1126 else PetscPrintf(PETSC_COMM_SELF, " ["); 1127 for (j = 0; j < dim; ++j) { 1128 if (j < dim-1) PetscPrintf(PETSC_COMM_SELF, "%15.8e, ", Mpos[i*dim+j]); 1129 else PetscPrintf(PETSC_COMM_SELF, "%15.8e", Mpos[i*dim+j]); 1130 } 1131 if (i < dim-1) PetscPrintf(PETSC_COMM_SELF, "]\n"); 1132 else PetscPrintf(PETSC_COMM_SELF, "]]\n"); 1133 } 1134 PetscFunctionReturn(0); 1135 } 1136 1137 static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[], PetscScalar *detMp) 1138 { 1139 PetscInt i, j, k; 1140 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); 1141 PetscScalar *Mpos; 1142 1143 PetscFunctionBegin; 1144 CHKERRQ(PetscMalloc2(dim*dim, &Mpos, dim, &eigs)); 1145 1146 /* Symmetrize */ 1147 for (i = 0; i < dim; ++i) { 1148 Mpos[i*dim+i] = Mp[i*dim+i]; 1149 for (j = i+1; j < dim; ++j) { 1150 Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 1151 Mpos[j*dim+i] = Mpos[i*dim+j]; 1152 } 1153 } 1154 1155 /* Compute eigendecomposition */ 1156 if (dim == 1) { 1157 1158 /* Isotropic case */ 1159 eigs[0] = PetscRealPart(Mpos[0]); 1160 Mpos[0] = 1.0; 1161 } else { 1162 1163 /* Anisotropic case */ 1164 PetscScalar *work; 1165 PetscBLASInt lwork; 1166 1167 lwork = 5*dim; 1168 CHKERRQ(PetscMalloc1(5*dim, &work)); 1169 { 1170 PetscBLASInt lierr; 1171 PetscBLASInt nb; 1172 1173 CHKERRQ(PetscBLASIntCast(dim, &nb)); 1174 CHKERRQ(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1175 #if defined(PETSC_USE_COMPLEX) 1176 { 1177 PetscReal *rwork; 1178 CHKERRQ(PetscMalloc1(3*dim, &rwork)); 1179 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr)); 1180 CHKERRQ(PetscFree(rwork)); 1181 } 1182 #else 1183 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr)); 1184 #endif 1185 if (lierr) { 1186 for (i = 0; i < dim; ++i) { 1187 Mpos[i*dim+i] = Mp[i*dim+i]; 1188 for (j = i+1; j < dim; ++j) { 1189 Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 1190 Mpos[j*dim+i] = Mpos[i*dim+j]; 1191 } 1192 } 1193 LAPACKsyevFail(dim, Mpos); 1194 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 1195 } 1196 CHKERRQ(PetscFPTrapPop()); 1197 } 1198 CHKERRQ(PetscFree(work)); 1199 } 1200 1201 /* Reflect to positive orthant and enforce maximum and minimum size */ 1202 max_eig = 0.0; 1203 for (i = 0; i < dim; ++i) { 1204 eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i]))); 1205 max_eig = PetscMax(eigs[i], max_eig); 1206 } 1207 1208 /* Enforce maximum anisotropy and compute determinant */ 1209 *detMp = 1.0; 1210 for (i = 0; i < dim; ++i) { 1211 if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min); 1212 *detMp *= eigs[i]; 1213 } 1214 1215 /* Reconstruct Hessian */ 1216 for (i = 0; i < dim; ++i) { 1217 for (j = 0; j < dim; ++j) { 1218 Mp[i*dim+j] = 0.0; 1219 for (k = 0; k < dim; ++k) { 1220 Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j]; 1221 } 1222 } 1223 } 1224 CHKERRQ(PetscFree2(Mpos, eigs)); 1225 1226 PetscFunctionReturn(0); 1227 } 1228 1229 /*@ 1230 DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric 1231 1232 Input parameters: 1233 + dm - The DM 1234 . metricIn - The metric 1235 . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 1236 - restrictAnisotropy - Should maximum anisotropy be enforced? 1237 1238 Output parameter: 1239 + metricOut - The metric 1240 - determinant - Its determinant 1241 1242 Level: beginner 1243 1244 Notes: 1245 1246 Relevant command line options: 1247 1248 + -dm_plex_metric_isotropic - Is the metric isotropic? 1249 . -dm_plex_metric_uniform - Is the metric uniform? 1250 . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 1251 . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 1252 - -dm_plex_metric_a_max - Maximum tolerated anisotropy 1253 1254 .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection() 1255 @*/ 1256 PetscErrorCode DMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut, Vec *determinant) 1257 { 1258 DM dmDet; 1259 PetscBool isotropic, uniform; 1260 PetscInt dim, vStart, vEnd, v; 1261 PetscScalar *met, *det; 1262 PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0; 1263 1264 PetscFunctionBegin; 1265 CHKERRQ(PetscLogEventBegin(DMPLEX_MetricEnforceSPD,0,0,0,0)); 1266 1267 /* Extract metadata from dm */ 1268 CHKERRQ(DMGetDimension(dm, &dim)); 1269 if (restrictSizes) { 1270 CHKERRQ(DMPlexMetricGetMinimumMagnitude(dm, &h_min)); 1271 CHKERRQ(DMPlexMetricGetMaximumMagnitude(dm, &h_max)); 1272 h_min = PetscMax(h_min, 1.0e-30); 1273 h_max = PetscMin(h_max, 1.0e+30); 1274 PetscCheck(h_min < h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max); 1275 } 1276 if (restrictAnisotropy) { 1277 CHKERRQ(DMPlexMetricGetMaximumAnisotropy(dm, &a_max)); 1278 a_max = PetscMin(a_max, 1.0e+30); 1279 } 1280 1281 /* Setup output metric */ 1282 CHKERRQ(DMPlexMetricCreate(dm, 0, metricOut)); 1283 CHKERRQ(VecCopy(metricIn, *metricOut)); 1284 1285 /* Enforce SPD and extract determinant */ 1286 CHKERRQ(VecGetArray(*metricOut, &met)); 1287 CHKERRQ(DMPlexMetricIsUniform(dm, &uniform)); 1288 CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic)); 1289 if (uniform) { 1290 PetscCheck(isotropic, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics cannot exist"); 1291 1292 /* Uniform case */ 1293 CHKERRQ(VecDuplicate(metricIn, determinant)); 1294 CHKERRQ(VecGetArray(*determinant, &det)); 1295 CHKERRQ(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det)); 1296 CHKERRQ(VecRestoreArray(*determinant, &det)); 1297 } else { 1298 1299 /* Spatially varying case */ 1300 PetscInt nrow; 1301 1302 if (isotropic) nrow = 1; 1303 else nrow = dim; 1304 CHKERRQ(DMClone(dm, &dmDet)); 1305 CHKERRQ(DMPlexP1FieldCreate_Private(dmDet, 0, 1, determinant)); 1306 CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1307 CHKERRQ(VecGetArray(*determinant, &det)); 1308 for (v = vStart; v < vEnd; ++v) { 1309 PetscScalar *vmet, *vdet; 1310 CHKERRQ(DMPlexPointLocalRef(dm, v, met, &vmet)); 1311 CHKERRQ(DMPlexPointLocalRef(dmDet, v, det, &vdet)); 1312 CHKERRQ(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, vmet, vdet)); 1313 } 1314 CHKERRQ(VecRestoreArray(*determinant, &det)); 1315 } 1316 CHKERRQ(VecRestoreArray(*metricOut, &met)); 1317 1318 CHKERRQ(PetscLogEventEnd(DMPLEX_MetricEnforceSPD,0,0,0,0)); 1319 PetscFunctionReturn(0); 1320 } 1321 1322 static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1323 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1324 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1325 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 1326 { 1327 const PetscScalar p = constants[0]; 1328 1329 f0[0] = PetscPowReal(u[0], p/(2.0*p + dim)); 1330 } 1331 1332 /*@ 1333 DMPlexMetricNormalize - Apply L-p normalization to a metric 1334 1335 Input parameters: 1336 + dm - The DM 1337 . metricIn - The unnormalized metric 1338 . restrictSizes - Should maximum/minimum metric magnitudes be enforced? 1339 - restrictAnisotropy - Should maximum metric anisotropy be enforced? 1340 1341 Output parameter: 1342 . metricOut - The normalized metric 1343 1344 Level: beginner 1345 1346 Notes: 1347 1348 Relevant command line options: 1349 1350 + -dm_plex_metric_isotropic - Is the metric isotropic? 1351 . -dm_plex_metric_uniform - Is the metric uniform? 1352 . -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 1353 . -dm_plex_metric_h_min - Minimum tolerated metric magnitude 1354 . -dm_plex_metric_h_max - Maximum tolerated metric magnitude 1355 . -dm_plex_metric_a_max - Maximum tolerated anisotropy 1356 . -dm_plex_metric_p - L-p normalization order 1357 - -dm_plex_metric_target_complexity - Target metric complexity 1358 1359 .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection() 1360 @*/ 1361 PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec *metricOut) 1362 { 1363 DM dmDet; 1364 MPI_Comm comm; 1365 PetscBool restrictAnisotropyFirst, isotropic, uniform; 1366 PetscDS ds; 1367 PetscInt dim, Nd, vStart, vEnd, v, i; 1368 PetscScalar *met, *det, integral, constants[1]; 1369 PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, fact, target, realIntegral; 1370 Vec determinant; 1371 1372 PetscFunctionBegin; 1373 CHKERRQ(PetscLogEventBegin(DMPLEX_MetricNormalize,0,0,0,0)); 1374 1375 /* Extract metadata from dm */ 1376 CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 1377 CHKERRQ(DMGetDimension(dm, &dim)); 1378 CHKERRQ(DMPlexMetricIsUniform(dm, &uniform)); 1379 CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic)); 1380 if (isotropic) Nd = 1; 1381 else Nd = dim*dim; 1382 1383 /* Set up metric and ensure it is SPD */ 1384 CHKERRQ(DMPlexMetricRestrictAnisotropyFirst(dm, &restrictAnisotropyFirst)); 1385 CHKERRQ(DMPlexMetricEnforceSPD(dm, metricIn, PETSC_FALSE, (PetscBool)(restrictAnisotropy && restrictAnisotropyFirst), metricOut, &determinant)); 1386 1387 /* Compute global normalization factor */ 1388 CHKERRQ(DMPlexMetricGetTargetComplexity(dm, &target)); 1389 CHKERRQ(DMPlexMetricGetNormalizationOrder(dm, &p)); 1390 constants[0] = p; 1391 if (uniform) { 1392 PetscCheck(isotropic,comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported"); 1393 DM dmTmp; 1394 Vec tmp; 1395 1396 CHKERRQ(DMClone(dm, &dmTmp)); 1397 CHKERRQ(DMPlexP1FieldCreate_Private(dmTmp, 0, 1, &tmp)); 1398 CHKERRQ(VecGetArray(determinant, &det)); 1399 CHKERRQ(VecSet(tmp, det[0])); 1400 CHKERRQ(VecRestoreArray(determinant, &det)); 1401 CHKERRQ(DMGetDS(dmTmp, &ds)); 1402 CHKERRQ(PetscDSSetConstants(ds, 1, constants)); 1403 CHKERRQ(PetscDSSetObjective(ds, 0, detMFunc)); 1404 CHKERRQ(DMPlexComputeIntegralFEM(dmTmp, tmp, &integral, NULL)); 1405 CHKERRQ(VecDestroy(&tmp)); 1406 CHKERRQ(DMDestroy(&dmTmp)); 1407 } else { 1408 CHKERRQ(VecGetDM(determinant, &dmDet)); 1409 CHKERRQ(DMGetDS(dmDet, &ds)); 1410 CHKERRQ(PetscDSSetConstants(ds, 1, constants)); 1411 CHKERRQ(PetscDSSetObjective(ds, 0, detMFunc)); 1412 CHKERRQ(DMPlexComputeIntegralFEM(dmDet, determinant, &integral, NULL)); 1413 } 1414 realIntegral = PetscRealPart(integral); 1415 PetscCheckFalse(realIntegral < 1.0e-30,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global metric normalization factor should be strictly positive, not %.4e Is the input metric positive-definite?", realIntegral); 1416 factGlob = PetscPowReal(target/realIntegral, 2.0/dim); 1417 1418 /* Apply local scaling */ 1419 if (restrictSizes) { 1420 CHKERRQ(DMPlexMetricGetMinimumMagnitude(dm, &h_min)); 1421 CHKERRQ(DMPlexMetricGetMaximumMagnitude(dm, &h_max)); 1422 h_min = PetscMax(h_min, 1.0e-30); 1423 h_max = PetscMin(h_max, 1.0e+30); 1424 PetscCheckFalse(h_min >= h_max,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Incompatible min/max metric magnitudes (%.4e not smaller than %.4e)", h_min, h_max); 1425 } 1426 if (restrictAnisotropy && !restrictAnisotropyFirst) { 1427 CHKERRQ(DMPlexMetricGetMaximumAnisotropy(dm, &a_max)); 1428 a_max = PetscMin(a_max, 1.0e+30); 1429 } 1430 CHKERRQ(VecGetArray(*metricOut, &met)); 1431 CHKERRQ(VecGetArray(determinant, &det)); 1432 if (uniform) { 1433 1434 /* Uniform case */ 1435 met[0] *= factGlob * PetscPowReal(PetscAbsScalar(det[0]), -1.0/(2*p+dim)); 1436 if (restrictSizes) CHKERRQ(DMPlexMetricModify_Private(1, h_min, h_max, a_max, met, det)); 1437 } else { 1438 1439 /* Spatially varying case */ 1440 PetscInt nrow; 1441 1442 if (isotropic) nrow = 1; 1443 else nrow = dim; 1444 CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1445 for (v = vStart; v < vEnd; ++v) { 1446 PetscScalar *Mp, *detM; 1447 1448 CHKERRQ(DMPlexPointLocalRef(dm, v, met, &Mp)); 1449 CHKERRQ(DMPlexPointLocalRef(dmDet, v, det, &detM)); 1450 fact = factGlob * PetscPowReal(PetscAbsScalar(detM[0]), -1.0/(2*p+dim)); 1451 for (i = 0; i < Nd; ++i) Mp[i] *= fact; 1452 if (restrictSizes) CHKERRQ(DMPlexMetricModify_Private(nrow, h_min, h_max, a_max, Mp, detM)); 1453 } 1454 } 1455 CHKERRQ(VecRestoreArray(determinant, &det)); 1456 CHKERRQ(VecRestoreArray(*metricOut, &met)); 1457 CHKERRQ(VecDestroy(&determinant)); 1458 if (!uniform) CHKERRQ(DMDestroy(&dmDet)); 1459 1460 CHKERRQ(PetscLogEventEnd(DMPLEX_MetricNormalize,0,0,0,0)); 1461 PetscFunctionReturn(0); 1462 } 1463 1464 /*@ 1465 DMPlexMetricAverage - Compute the average of a list of metrics 1466 1467 Input Parameters: 1468 + dm - The DM 1469 . numMetrics - The number of metrics to be averaged 1470 . weights - Weights for the average 1471 - metrics - The metrics to be averaged 1472 1473 Output Parameter: 1474 . metricAvg - The averaged metric 1475 1476 Level: beginner 1477 1478 Notes: 1479 The weights should sum to unity. 1480 1481 If weights are not provided then an unweighted average is used. 1482 1483 .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection() 1484 @*/ 1485 PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg) 1486 { 1487 PetscBool haveWeights = PETSC_TRUE; 1488 PetscInt i, m, n; 1489 PetscReal sum = 0.0, tol = 1.0e-10; 1490 1491 PetscFunctionBegin; 1492 CHKERRQ(PetscLogEventBegin(DMPLEX_MetricAverage,0,0,0,0)); 1493 PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); 1494 CHKERRQ(DMPlexMetricCreate(dm, 0, metricAvg)); 1495 CHKERRQ(VecSet(*metricAvg, 0.0)); 1496 CHKERRQ(VecGetSize(*metricAvg, &m)); 1497 for (i = 0; i < numMetrics; ++i) { 1498 CHKERRQ(VecGetSize(metrics[i], &n)); 1499 PetscCheck(m == n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Averaging different metric types not implemented"); 1500 } 1501 1502 /* Default to the unweighted case */ 1503 if (!weights) { 1504 CHKERRQ(PetscMalloc1(numMetrics, &weights)); 1505 haveWeights = PETSC_FALSE; 1506 for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; } 1507 } 1508 1509 /* Check weights sum to unity */ 1510 for (i = 0; i < numMetrics; ++i) sum += weights[i]; 1511 PetscCheck(PetscAbsReal(sum - 1) <= tol,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); 1512 1513 /* Compute metric average */ 1514 for (i = 0; i < numMetrics; ++i) CHKERRQ(VecAXPY(*metricAvg, weights[i], metrics[i])); 1515 if (!haveWeights) CHKERRQ(PetscFree(weights)); 1516 1517 CHKERRQ(PetscLogEventEnd(DMPLEX_MetricAverage,0,0,0,0)); 1518 PetscFunctionReturn(0); 1519 } 1520 1521 /*@ 1522 DMPlexMetricAverage2 - Compute the unweighted average of two metrics 1523 1524 Input Parameters: 1525 + dm - The DM 1526 . metric1 - The first metric to be averaged 1527 - metric2 - The second metric to be averaged 1528 1529 Output Parameter: 1530 . metricAvg - The averaged metric 1531 1532 Level: beginner 1533 1534 .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3() 1535 @*/ 1536 PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg) 1537 { 1538 PetscReal weights[2] = {0.5, 0.5}; 1539 Vec metrics[2] = {metric1, metric2}; 1540 1541 PetscFunctionBegin; 1542 CHKERRQ(DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg)); 1543 PetscFunctionReturn(0); 1544 } 1545 1546 /*@ 1547 DMPlexMetricAverage3 - Compute the unweighted average of three metrics 1548 1549 Input Parameters: 1550 + dm - The DM 1551 . metric1 - The first metric to be averaged 1552 . metric2 - The second metric to be averaged 1553 - metric3 - The third metric to be averaged 1554 1555 Output Parameter: 1556 . metricAvg - The averaged metric 1557 1558 Level: beginner 1559 1560 .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2() 1561 @*/ 1562 PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg) 1563 { 1564 PetscReal weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0}; 1565 Vec metrics[3] = {metric1, metric2, metric3}; 1566 1567 PetscFunctionBegin; 1568 CHKERRQ(DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg)); 1569 PetscFunctionReturn(0); 1570 } 1571 1572 static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[]) 1573 { 1574 PetscInt i, j, k, l, m; 1575 PetscReal *evals, *evals1; 1576 PetscScalar *evecs, *sqrtM1, *isqrtM1; 1577 1578 PetscFunctionBegin; 1579 1580 /* Isotropic case */ 1581 if (dim == 1) { 1582 M2[0] = (PetscScalar)PetscMin(PetscRealPart(M1[0]), PetscRealPart(M2[0])); 1583 PetscFunctionReturn(0); 1584 } 1585 1586 /* Anisotropic case */ 1587 CHKERRQ(PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1)); 1588 for (i = 0; i < dim; ++i) { 1589 for (j = 0; j < dim; ++j) { 1590 evecs[i*dim+j] = M1[i*dim+j]; 1591 } 1592 } 1593 { 1594 PetscScalar *work; 1595 PetscBLASInt lwork; 1596 1597 lwork = 5*dim; 1598 CHKERRQ(PetscMalloc1(5*dim, &work)); 1599 { 1600 PetscBLASInt lierr, nb; 1601 PetscReal sqrtk; 1602 1603 /* Compute eigendecomposition of M1 */ 1604 CHKERRQ(PetscBLASIntCast(dim, &nb)); 1605 CHKERRQ(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1606 #if defined(PETSC_USE_COMPLEX) 1607 { 1608 PetscReal *rwork; 1609 CHKERRQ(PetscMalloc1(3*dim, &rwork)); 1610 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr)); 1611 CHKERRQ(PetscFree(rwork)); 1612 } 1613 #else 1614 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr)); 1615 #endif 1616 if (lierr) { 1617 LAPACKsyevFail(dim, M1); 1618 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 1619 } 1620 CHKERRQ(PetscFPTrapPop()); 1621 1622 /* Compute square root and reciprocal */ 1623 for (i = 0; i < dim; ++i) { 1624 for (j = 0; j < dim; ++j) { 1625 sqrtM1[i*dim+j] = 0.0; 1626 isqrtM1[i*dim+j] = 0.0; 1627 for (k = 0; k < dim; ++k) { 1628 sqrtk = PetscSqrtReal(evals1[k]); 1629 sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j]; 1630 isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j]; 1631 } 1632 } 1633 } 1634 1635 /* Map into the space spanned by the eigenvectors of M1 */ 1636 for (i = 0; i < dim; ++i) { 1637 for (j = 0; j < dim; ++j) { 1638 evecs[i*dim+j] = 0.0; 1639 for (k = 0; k < dim; ++k) { 1640 for (l = 0; l < dim; ++l) { 1641 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 1642 } 1643 } 1644 } 1645 } 1646 1647 /* Compute eigendecomposition */ 1648 CHKERRQ(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1649 #if defined(PETSC_USE_COMPLEX) 1650 { 1651 PetscReal *rwork; 1652 CHKERRQ(PetscMalloc1(3*dim, &rwork)); 1653 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 1654 CHKERRQ(PetscFree(rwork)); 1655 } 1656 #else 1657 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 1658 #endif 1659 if (lierr) { 1660 for (i = 0; i < dim; ++i) { 1661 for (j = 0; j < dim; ++j) { 1662 evecs[i*dim+j] = 0.0; 1663 for (k = 0; k < dim; ++k) { 1664 for (l = 0; l < dim; ++l) { 1665 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 1666 } 1667 } 1668 } 1669 } 1670 LAPACKsyevFail(dim, evecs); 1671 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 1672 } 1673 CHKERRQ(PetscFPTrapPop()); 1674 1675 /* Modify eigenvalues */ 1676 for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]); 1677 1678 /* Map back to get the intersection */ 1679 for (i = 0; i < dim; ++i) { 1680 for (j = 0; j < dim; ++j) { 1681 M2[i*dim+j] = 0.0; 1682 for (k = 0; k < dim; ++k) { 1683 for (l = 0; l < dim; ++l) { 1684 for (m = 0; m < dim; ++m) { 1685 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m]; 1686 } 1687 } 1688 } 1689 } 1690 } 1691 } 1692 CHKERRQ(PetscFree(work)); 1693 } 1694 CHKERRQ(PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1)); 1695 PetscFunctionReturn(0); 1696 } 1697 1698 /*@ 1699 DMPlexMetricIntersection - Compute the intersection of a list of metrics 1700 1701 Input Parameters: 1702 + dm - The DM 1703 . numMetrics - The number of metrics to be intersected 1704 - metrics - The metrics to be intersected 1705 1706 Output Parameter: 1707 . metricInt - The intersected metric 1708 1709 Level: beginner 1710 1711 Notes: 1712 1713 The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics. 1714 1715 The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2. 1716 1717 .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage() 1718 @*/ 1719 PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt) 1720 { 1721 PetscBool isotropic, uniform; 1722 PetscInt v, i, m, n; 1723 PetscScalar *met, *meti; 1724 1725 PetscFunctionBegin; 1726 CHKERRQ(PetscLogEventBegin(DMPLEX_MetricIntersection,0,0,0,0)); 1727 PetscCheck(numMetrics >= 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); 1728 1729 /* Copy over the first metric */ 1730 CHKERRQ(DMPlexMetricCreate(dm, 0, metricInt)); 1731 CHKERRQ(VecCopy(metrics[0], *metricInt)); 1732 if (numMetrics == 1) PetscFunctionReturn(0); 1733 CHKERRQ(VecGetSize(*metricInt, &m)); 1734 for (i = 0; i < numMetrics; ++i) { 1735 CHKERRQ(VecGetSize(metrics[i], &n)); 1736 PetscCheckFalse(m != n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Intersecting different metric types not implemented"); 1737 } 1738 1739 /* Intersect subsequent metrics in turn */ 1740 CHKERRQ(DMPlexMetricIsUniform(dm, &uniform)); 1741 CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic)); 1742 if (uniform) { 1743 1744 /* Uniform case */ 1745 CHKERRQ(VecGetArray(*metricInt, &met)); 1746 for (i = 1; i < numMetrics; ++i) { 1747 CHKERRQ(VecGetArray(metrics[i], &meti)); 1748 CHKERRQ(DMPlexMetricIntersection_Private(1, meti, met)); 1749 CHKERRQ(VecRestoreArray(metrics[i], &meti)); 1750 } 1751 CHKERRQ(VecRestoreArray(*metricInt, &met)); 1752 } else { 1753 1754 /* Spatially varying case */ 1755 PetscInt dim, vStart, vEnd, nrow; 1756 PetscScalar *M, *Mi; 1757 1758 CHKERRQ(DMGetDimension(dm, &dim)); 1759 if (isotropic) nrow = 1; 1760 else nrow = dim; 1761 CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1762 CHKERRQ(VecGetArray(*metricInt, &met)); 1763 for (i = 1; i < numMetrics; ++i) { 1764 CHKERRQ(VecGetArray(metrics[i], &meti)); 1765 for (v = vStart; v < vEnd; ++v) { 1766 CHKERRQ(DMPlexPointLocalRef(dm, v, met, &M)); 1767 CHKERRQ(DMPlexPointLocalRef(dm, v, meti, &Mi)); 1768 CHKERRQ(DMPlexMetricIntersection_Private(nrow, Mi, M)); 1769 } 1770 CHKERRQ(VecRestoreArray(metrics[i], &meti)); 1771 } 1772 CHKERRQ(VecRestoreArray(*metricInt, &met)); 1773 } 1774 1775 CHKERRQ(PetscLogEventEnd(DMPLEX_MetricIntersection,0,0,0,0)); 1776 PetscFunctionReturn(0); 1777 } 1778 1779 /*@ 1780 DMPlexMetricIntersection2 - Compute the intersection of two metrics 1781 1782 Input Parameters: 1783 + dm - The DM 1784 . metric1 - The first metric to be intersected 1785 - metric2 - The second metric to be intersected 1786 1787 Output Parameter: 1788 . metricInt - The intersected metric 1789 1790 Level: beginner 1791 1792 .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3() 1793 @*/ 1794 PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt) 1795 { 1796 Vec metrics[2] = {metric1, metric2}; 1797 1798 PetscFunctionBegin; 1799 CHKERRQ(DMPlexMetricIntersection(dm, 2, metrics, metricInt)); 1800 PetscFunctionReturn(0); 1801 } 1802 1803 /*@ 1804 DMPlexMetricIntersection3 - Compute the intersection of three metrics 1805 1806 Input Parameters: 1807 + dm - The DM 1808 . metric1 - The first metric to be intersected 1809 . metric2 - The second metric to be intersected 1810 - metric3 - The third metric to be intersected 1811 1812 Output Parameter: 1813 . metricInt - The intersected metric 1814 1815 Level: beginner 1816 1817 .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2() 1818 @*/ 1819 PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt) 1820 { 1821 Vec metrics[3] = {metric1, metric2, metric3}; 1822 1823 PetscFunctionBegin; 1824 CHKERRQ(DMPlexMetricIntersection(dm, 3, metrics, metricInt)); 1825 PetscFunctionReturn(0); 1826 } 1827