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