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