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