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