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