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