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