1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscblaslapack.h> 3 4 PetscErrorCode DMPlexMetricSetFromOptions(DM dm) 5 { 6 MPI_Comm comm; 7 PetscBool isotropic = PETSC_FALSE, restrictAnisotropyFirst = PETSC_FALSE; 8 PetscErrorCode ierr; 9 PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 1.0e+05, p = 1.0, target = 1000.0; 10 11 PetscFunctionBegin; 12 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 13 ierr = PetscOptionsBegin(comm, "", "Riemannian metric options", "DMPlexMetric");CHKERRQ(ierr); 14 ierr = PetscOptionsBool("-dm_plex_metric_isotropic", "Is the metric isotropic?", "DMPlexMetricCreateIsotropic", isotropic, &isotropic, NULL);CHKERRQ(ierr); 15 ierr = DMPlexMetricSetIsotropic(dm, isotropic); 16 ierr = PetscOptionsBool("-dm_plex_metric_restrict_anisotropy_first", "Should anisotropy be restricted before normalization?", "DMPlexNormalize", restrictAnisotropyFirst, &restrictAnisotropyFirst, NULL);CHKERRQ(ierr); 17 ierr = DMPlexMetricSetRestrictAnisotropyFirst(dm, restrictAnisotropyFirst); 18 ierr = PetscOptionsReal("-dm_plex_metric_h_min", "Minimum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_min, &h_min, NULL);CHKERRQ(ierr); 19 ierr = DMPlexMetricSetMinimumMagnitude(dm, h_min);CHKERRQ(ierr); 20 ierr = PetscOptionsReal("-dm_plex_metric_h_max", "Maximum tolerated metric magnitude", "DMPlexMetricEnforceSPD", h_max, &h_max, NULL);CHKERRQ(ierr); 21 ierr = DMPlexMetricSetMaximumMagnitude(dm, h_max);CHKERRQ(ierr); 22 ierr = PetscOptionsReal("-dm_plex_metric_a_max", "Maximum tolerated anisotropy", "DMPlexMetricEnforceSPD", a_max, &a_max, NULL);CHKERRQ(ierr); 23 ierr = DMPlexMetricSetMaximumAnisotropy(dm, a_max);CHKERRQ(ierr); 24 ierr = PetscOptionsReal("-dm_plex_metric_p", "L-p normalization order", "DMPlexMetricNormalize", p, &p, NULL);CHKERRQ(ierr); 25 ierr = DMPlexMetricSetNormalizationOrder(dm, p);CHKERRQ(ierr); 26 ierr = PetscOptionsReal("-dm_plex_metric_target_complexity", "Target metric complexity", "DMPlexMetricNormalize", target, &target, NULL);CHKERRQ(ierr); 27 ierr = DMPlexMetricSetTargetComplexity(dm, target);CHKERRQ(ierr); 28 ierr = PetscOptionsEnd();CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 /* 33 DMPlexMetricSetIsotropic - Record whether a metric is isotropic 34 35 Input parameters: 36 + dm - The DM 37 - isotropic - Is the metric isotropic? 38 39 Level: beginner 40 41 .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 42 */ 43 PetscErrorCode DMPlexMetricSetIsotropic(DM dm, PetscBool isotropic) 44 { 45 DM_Plex *plex = (DM_Plex *) dm->data; 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 if (!plex->metricCtx) { 50 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 51 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 52 } 53 plex->metricCtx->isotropic = isotropic; 54 PetscFunctionReturn(0); 55 } 56 57 /* 58 DMPlexMetricIsIsotropic - Is a metric is isotropic? 59 60 Input parameters: 61 . dm - The DM 62 63 Output parameters: 64 . isotropic - Is the metric isotropic? 65 66 Level: beginner 67 68 .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 69 */ 70 PetscErrorCode DMPlexMetricIsIsotropic(DM dm, PetscBool *isotropic) 71 { 72 DM_Plex *plex = (DM_Plex *) dm->data; 73 PetscErrorCode ierr; 74 75 PetscFunctionBegin; 76 if (!plex->metricCtx) { 77 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 78 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 79 } 80 *isotropic = plex->metricCtx->isotropic; 81 PetscFunctionReturn(0); 82 } 83 84 /* 85 DMPlexMetricSetRestrictAnisotropyFirst - Record whether anisotropy should be restricted before normalization 86 87 Input parameters: 88 + dm - The DM 89 - restrictAnisotropyFirst - Should anisotropy be normalized first? 90 91 Level: beginner 92 93 .seealso: DMPlexMetricSetIsotropic(), DMPlexMetricRestrictAnisotropyFirst() 94 */ 95 PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM dm, PetscBool restrictAnisotropyFirst) 96 { 97 DM_Plex *plex = (DM_Plex *) dm->data; 98 PetscErrorCode ierr; 99 100 PetscFunctionBegin; 101 if (!plex->metricCtx) { 102 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 103 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 104 } 105 plex->metricCtx->restrictAnisotropyFirst = restrictAnisotropyFirst; 106 PetscFunctionReturn(0); 107 } 108 109 /* 110 DMPlexMetricRestrictAnisotropyFirst - Is anisotropy restricted before normalization or after? 111 112 Input parameters: 113 . dm - The DM 114 115 Output parameters: 116 . restrictAnisotropyFirst - Is anisotropy be normalized first? 117 118 Level: beginner 119 120 .seealso: DMPlexMetricIsIsotropic(), DMPlexMetricSetRestrictAnisotropyFirst() 121 */ 122 PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM dm, PetscBool *restrictAnisotropyFirst) 123 { 124 DM_Plex *plex = (DM_Plex *) dm->data; 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 if (!plex->metricCtx) { 129 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 130 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 131 } 132 *restrictAnisotropyFirst = plex->metricCtx->restrictAnisotropyFirst; 133 PetscFunctionReturn(0); 134 } 135 136 /* 137 DMPlexMetricSetMinimumMagnitude - Set the minimum tolerated metric magnitude 138 139 Input parameters: 140 + dm - The DM 141 - h_min - The minimum tolerated metric magnitude 142 143 Level: beginner 144 145 .seealso: DMPlexMetricGetMinimumMagnitude(), DMPlexMetricSetMaximumMagnitude() 146 */ 147 PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min) 148 { 149 DM_Plex *plex = (DM_Plex *) dm->data; 150 PetscErrorCode ierr; 151 152 PetscFunctionBegin; 153 if (!plex->metricCtx) { 154 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 155 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 156 } 157 if (h_min <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_min); 158 plex->metricCtx->h_min = h_min; 159 PetscFunctionReturn(0); 160 } 161 162 /* 163 DMPlexMetricGetMinimumMagnitude - Get the minimum tolerated metric magnitude 164 165 Input parameters: 166 . dm - The DM 167 168 Input parameters: 169 . h_min - The minimum tolerated metric magnitude 170 171 Level: beginner 172 173 .seealso: DMPlexMetricSetMinimumMagnitude(), DMPlexMetricGetMaximumMagnitude() 174 */ 175 PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM dm, PetscReal *h_min) 176 { 177 DM_Plex *plex = (DM_Plex *) dm->data; 178 PetscErrorCode ierr; 179 180 PetscFunctionBegin; 181 if (!plex->metricCtx) { 182 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 183 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 184 } 185 *h_min = plex->metricCtx->h_min; 186 PetscFunctionReturn(0); 187 } 188 189 /* 190 DMPlexMetricSetMaximumMagnitude - Set the maximum tolerated metric magnitude 191 192 Input parameters: 193 + dm - The DM 194 - h_max - The maximum tolerated metric magnitude 195 196 Level: beginner 197 198 .seealso: DMPlexMetricGetMaximumMagnitude(), DMPlexMetricSetMinimumMagnitude() 199 */ 200 PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max) 201 { 202 DM_Plex *plex = (DM_Plex *) dm->data; 203 PetscErrorCode ierr; 204 205 PetscFunctionBegin; 206 if (!plex->metricCtx) { 207 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 208 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 209 } 210 if (h_max <= 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric magnitudes must be positive, not %.4e", h_max); 211 plex->metricCtx->h_max = h_max; 212 PetscFunctionReturn(0); 213 } 214 215 /* 216 DMPlexMetricGetMaximumMagnitude - Get the maximum tolerated metric magnitude 217 218 Input parameters: 219 . dm - The DM 220 221 Input parameters: 222 . h_max - The maximum tolerated metric magnitude 223 224 Level: beginner 225 226 .seealso: DMPlexMetricSetMaximumMagnitude(), DMPlexMetricGetMinimumMagnitude() 227 */ 228 PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM dm, PetscReal *h_max) 229 { 230 DM_Plex *plex = (DM_Plex *) dm->data; 231 PetscErrorCode ierr; 232 233 PetscFunctionBegin; 234 if (!plex->metricCtx) { 235 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 236 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 237 } 238 *h_max = plex->metricCtx->h_max; 239 PetscFunctionReturn(0); 240 } 241 242 /* 243 DMPlexMetricSetMaximumAnisotropy - Set the maximum tolerated metric anisotropy 244 245 Input parameters: 246 + dm - The DM 247 - a_max - The maximum tolerated metric anisotropy 248 249 Level: beginner 250 251 Note: If the value zero is given then anisotropy will not be restricted. Otherwise, it should be at least one. 252 253 .seealso: DMPlexMetricGetMaximumAnisotropy(), DMPlexMetricSetMaximumMagnitude() 254 */ 255 PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max) 256 { 257 DM_Plex *plex = (DM_Plex *) dm->data; 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 if (!plex->metricCtx) { 262 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 263 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 264 } 265 if (a_max < 1.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Anisotropy must be at least one, not %.4e", a_max); 266 plex->metricCtx->a_max = a_max; 267 PetscFunctionReturn(0); 268 } 269 270 /* 271 DMPlexMetricGetMaximumAnisotropy - Get the maximum tolerated metric anisotropy 272 273 Input parameters: 274 . dm - The DM 275 276 Input parameters: 277 . a_max - The maximum tolerated metric anisotropy 278 279 Level: beginner 280 281 .seealso: DMPlexMetricSetMaximumAnisotropy(), DMPlexMetricGetMaximumMagnitude() 282 */ 283 PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM dm, PetscReal *a_max) 284 { 285 DM_Plex *plex = (DM_Plex *) dm->data; 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 if (!plex->metricCtx) { 290 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 291 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 292 } 293 *a_max = plex->metricCtx->a_max; 294 PetscFunctionReturn(0); 295 } 296 297 /* 298 DMPlexMetricSetTargetComplexity - Set the target metric complexity 299 300 Input parameters: 301 + dm - The DM 302 - targetComplexity - The target metric complexity 303 304 Level: beginner 305 306 .seealso: DMPlexMetricGetTargetComplexity(), DMPlexMetricSetNormalizationOrder() 307 */ 308 PetscErrorCode DMPlexMetricSetTargetComplexity(DM dm, PetscReal targetComplexity) 309 { 310 DM_Plex *plex = (DM_Plex *) dm->data; 311 PetscErrorCode ierr; 312 313 PetscFunctionBegin; 314 if (!plex->metricCtx) { 315 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 316 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 317 } 318 if (targetComplexity <= 0.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric complexity must be positive"); 319 plex->metricCtx->targetComplexity = targetComplexity; 320 PetscFunctionReturn(0); 321 } 322 323 /* 324 DMPlexMetricGetTargetComplexity - Get the target metric complexity 325 326 Input parameters: 327 . dm - The DM 328 329 Input parameters: 330 . targetComplexity - The target metric complexity 331 332 Level: beginner 333 334 .seealso: DMPlexMetricSetTargetComplexity(), DMPlexMetricGetNormalizationOrder() 335 */ 336 PetscErrorCode DMPlexMetricGetTargetComplexity(DM dm, PetscReal *targetComplexity) 337 { 338 DM_Plex *plex = (DM_Plex *) dm->data; 339 PetscErrorCode ierr; 340 341 PetscFunctionBegin; 342 if (!plex->metricCtx) { 343 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 344 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 345 } 346 *targetComplexity = plex->metricCtx->targetComplexity; 347 PetscFunctionReturn(0); 348 } 349 350 /* 351 DMPlexMetricSetNormalizationOrder - Set the order p for L-p normalization 352 353 Input parameters: 354 + dm - The DM 355 - p - The normalization order 356 357 Level: beginner 358 359 .seealso: DMPlexMetricGetNormalizationOrder(), DMPlexMetricSetTargetComplexity() 360 */ 361 PetscErrorCode DMPlexMetricSetNormalizationOrder(DM dm, PetscReal p) 362 { 363 DM_Plex *plex = (DM_Plex *) dm->data; 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 if (!plex->metricCtx) { 368 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 369 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 370 } 371 if (p < 1.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normalization order must be one or greater"); 372 plex->metricCtx->p = p; 373 PetscFunctionReturn(0); 374 } 375 376 /* 377 DMPlexMetricGetNormalizationOrder - Get the order p for L-p normalization 378 379 Input parameters: 380 . dm - The DM 381 382 Input parameters: 383 . p - The normalization order 384 385 Level: beginner 386 387 .seealso: DMPlexMetricSetNormalizationOrder(), DMPlexMetricGetTargetComplexity() 388 */ 389 PetscErrorCode DMPlexMetricGetNormalizationOrder(DM dm, PetscReal *p) 390 { 391 DM_Plex *plex = (DM_Plex *) dm->data; 392 PetscErrorCode ierr; 393 394 PetscFunctionBegin; 395 if (!plex->metricCtx) { 396 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 397 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 398 } 399 *p = plex->metricCtx->p; 400 PetscFunctionReturn(0); 401 } 402 403 PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric) 404 { 405 MPI_Comm comm; 406 PetscErrorCode ierr; 407 PetscFE fe; 408 PetscInt dim; 409 410 PetscFunctionBegin; 411 412 /* Extract metadata from dm */ 413 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 414 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 415 416 /* Create a P1 field of the requested size */ 417 ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 418 ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr); 419 ierr = DMCreateDS(dm);CHKERRQ(ierr); 420 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 421 ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr); 422 423 PetscFunctionReturn(0); 424 } 425 426 /* 427 DMPlexMetricCreate - Create a Riemannian metric field 428 429 Input parameters: 430 + dm - The DM 431 - f - The field number to use 432 433 Output parameter: 434 . metric - The metric 435 436 Level: beginner 437 438 Notes: 439 440 It is assumed that the DM is comprised of simplices. 441 442 Command line options for Riemannian metrics: 443 444 -dm_plex_metric_isotropic - Is the metric isotropic? 445 -dm_plex_metric_restrict_anisotropy_first - Should anisotropy be restricted before normalization? 446 -dm_plex_metric_h_min - Minimum tolerated metric magnitude 447 -dm_plex_metric_h_max - Maximum tolerated metric magnitude 448 -dm_plex_metric_a_max - Maximum tolerated anisotropy 449 -dm_plex_metric_p - L-p normalization order 450 -dm_plex_metric_target_complexity - Target metric complexity 451 452 .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic() 453 */ 454 PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric) 455 { 456 DM_Plex *plex = (DM_Plex *) dm->data; 457 PetscErrorCode ierr; 458 PetscInt coordDim, Nd; 459 460 PetscFunctionBegin; 461 if (!plex->metricCtx) { 462 ierr = PetscNew(&plex->metricCtx);CHKERRQ(ierr); 463 ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 464 } 465 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 466 Nd = coordDim*coordDim; 467 ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr); 468 PetscFunctionReturn(0); 469 } 470 471 typedef struct { 472 PetscReal scaling; /* Scaling for uniform metric diagonal */ 473 } DMPlexMetricUniformCtx; 474 475 static PetscErrorCode diagonal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 476 { 477 DMPlexMetricUniformCtx *user = (DMPlexMetricUniformCtx*)ctx; 478 PetscInt i, j; 479 480 for (i = 0; i < dim; ++i) { 481 for (j = 0; j < dim; ++j) { 482 if (i == j) u[i+dim*j] = user->scaling; 483 else u[i+dim*j] = 0.0; 484 } 485 } 486 return 0; 487 } 488 489 /* 490 DMPlexMetricCreateUniform - Construct a uniform isotropic metric 491 492 Input parameters: 493 + dm - The DM 494 . f - The field number to use 495 - alpha - Scaling parameter for the diagonal 496 497 Output parameter: 498 . metric - The uniform metric 499 500 Level: beginner 501 502 Note: It is assumed that the DM is comprised of simplices. 503 504 .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic() 505 */ 506 PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric) 507 { 508 DMPlexMetricUniformCtx user; 509 PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 510 PetscErrorCode ierr; 511 void *ctxs[1]; 512 513 PetscFunctionBegin; 514 ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr); 515 if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined"); 516 if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha); 517 else user.scaling = alpha; 518 funcs[0] = diagonal; 519 ctxs[0] = &user; 520 ierr = DMProjectFunctionLocal(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, *metric);CHKERRQ(ierr); 521 PetscFunctionReturn(0); 522 } 523 524 /* 525 DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator 526 527 Input parameters: 528 + dm - The DM 529 . f - The field number to use 530 - indicator - The error indicator 531 532 Output parameter: 533 . metric - The isotropic metric 534 535 Level: beginner 536 537 Notes: 538 539 It is assumed that the DM is comprised of simplices. 540 541 The indicator needs to be a scalar field defined at *vertices*. 542 543 .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform() 544 */ 545 PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric) 546 { 547 DM dmIndi; 548 PetscErrorCode ierr; 549 PetscInt dim, d, vStart, vEnd, v; 550 const PetscScalar *indi; 551 PetscScalar *met; 552 553 PetscFunctionBegin; 554 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 555 ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr); 556 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 557 ierr = VecGetArrayRead(indicator, &indi);CHKERRQ(ierr); 558 ierr = VecGetArrayWrite(*metric, &met);CHKERRQ(ierr); 559 ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr); 560 for (v = vStart; v < vEnd; ++v) { 561 PetscScalar *vindi, *vmet; 562 ierr = DMPlexPointLocalRead(dmIndi, v, indi, &vindi);CHKERRQ(ierr); 563 ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr); 564 for (d = 0; d < dim; ++d) vmet[d*(dim+1)] = vindi[0]; 565 } 566 ierr = VecRestoreArrayWrite(*metric, &met);CHKERRQ(ierr); 567 ierr = VecRestoreArrayRead(indicator, &indi);CHKERRQ(ierr); 568 PetscFunctionReturn(0); 569 } 570 571 static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[]) 572 { 573 PetscErrorCode ierr; 574 PetscInt i, j, k; 575 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); 576 PetscScalar *Mpos; 577 578 PetscFunctionBegin; 579 ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr); 580 581 /* Symmetrize */ 582 for (i = 0; i < dim; ++i) { 583 Mpos[i*dim+i] = Mp[i*dim+i]; 584 for (j = i+1; j < dim; ++j) { 585 Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]); 586 Mpos[j*dim+i] = Mpos[i*dim+j]; 587 } 588 } 589 590 /* Compute eigendecomposition */ 591 { 592 PetscScalar *work; 593 PetscBLASInt lwork; 594 595 lwork = 5*dim; 596 ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr); 597 { 598 PetscBLASInt lierr; 599 PetscBLASInt nb; 600 601 ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr); 602 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 603 #if defined(PETSC_USE_COMPLEX) 604 { 605 PetscReal *rwork; 606 ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 607 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr)); 608 ierr = PetscFree(rwork);CHKERRQ(ierr); 609 } 610 #else 611 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr)); 612 #endif 613 if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 614 ierr = PetscFPTrapPop();CHKERRQ(ierr); 615 } 616 ierr = PetscFree(work);CHKERRQ(ierr); 617 } 618 619 /* Reflect to positive orthant and enforce maximum and minimum size */ 620 max_eig = 0.0; 621 for (i = 0; i < dim; ++i) { 622 eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i]))); 623 max_eig = PetscMax(eigs[i], max_eig); 624 } 625 626 /* Enforce maximum anisotropy */ 627 for (i = 0; i < dim; ++i) { 628 if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min); 629 } 630 631 /* Reconstruct Hessian */ 632 for (i = 0; i < dim; ++i) { 633 for (j = 0; j < dim; ++j) { 634 Mp[i*dim+j] = 0.0; 635 for (k = 0; k < dim; ++k) { 636 Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j]; 637 } 638 } 639 } 640 ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr); 641 642 PetscFunctionReturn(0); 643 } 644 645 /* 646 DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric 647 648 Input parameters: 649 + dm - The DM 650 . restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced? 651 - metric - The metric 652 653 Output parameter: 654 . metric - The metric 655 656 Level: beginner 657 658 .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection() 659 */ 660 PetscErrorCode DMPlexMetricEnforceSPD(DM dm, PetscBool restrictSizes, Vec metric) 661 { 662 DMPlexMetricCtx *user; 663 PetscErrorCode ierr; 664 PetscInt dim, vStart, vEnd, v; 665 PetscScalar *met; 666 PetscReal h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0; 667 668 PetscFunctionBegin; 669 670 /* Extract metadata from dm */ 671 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 672 ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr); 673 if (restrictSizes) { 674 if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max); 675 if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min); 676 if (user->a_max > 1.0) a_max = user->a_max; 677 } 678 679 /* Enforce SPD */ 680 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 681 ierr = VecGetArray(metric, &met);CHKERRQ(ierr); 682 for (v = vStart; v < vEnd; ++v) { 683 PetscScalar *vmet; 684 ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr); 685 ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, vmet);CHKERRQ(ierr); 686 } 687 ierr = VecRestoreArray(metric, &met);CHKERRQ(ierr); 688 689 PetscFunctionReturn(0); 690 } 691 692 static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 693 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 694 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 695 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 696 { 697 const PetscScalar p = constants[0]; 698 PetscReal detH = 0.0; 699 700 if (dim == 2) DMPlex_Det2D_Scalar_Internal(&detH, u); 701 else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detH, u); 702 f0[0] = PetscPowReal(detH, p/(2.0*p + dim)); 703 } 704 705 /* 706 DMPlexMetricNormalize - Apply L-p normalization to a metric 707 708 Input parameters: 709 + dm - The DM 710 . metricIn - The unnormalized metric 711 - restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced? 712 713 Output parameter: 714 . metricOut - The normalized metric 715 716 Level: beginner 717 718 .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection() 719 */ 720 PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, Vec *metricOut) 721 { 722 DMPlexMetricCtx *user; 723 MPI_Comm comm; 724 PetscDS ds; 725 PetscErrorCode ierr; 726 PetscInt dim, Nd, vStart, vEnd, v, i; 727 PetscScalar *met, integral, constants[1]; 728 PetscReal p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, target; 729 730 PetscFunctionBegin; 731 732 /* Extract metadata from dm */ 733 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 734 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 735 Nd = dim*dim; 736 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 737 ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr); 738 if (restrictSizes && user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max; 739 if (PetscAbsReal(user->p) >= 1.0) p = user->p; 740 else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric normalization order %f should be greater than one.", user->p); 741 constants[0] = p; 742 if (user->targetComplexity > 0.0) target = user->targetComplexity; 743 else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Target metric complexity %f should be positive.", user->targetComplexity); 744 745 /* Set up metric and ensure it is SPD */ 746 ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr); 747 ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr); 748 ierr = DMPlexMetricEnforceSPD(dm, PETSC_FALSE, *metricOut);CHKERRQ(ierr); 749 750 /* Compute global normalization factor */ 751 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 752 ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr); 753 ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr); 754 ierr = DMPlexComputeIntegralFEM(dm, *metricOut, &integral, NULL);CHKERRQ(ierr); 755 factGlob = PetscPowReal(target/PetscRealPart(integral), 2.0/dim); 756 757 /* Apply local scaling */ 758 a_max = 0.0; 759 if (restrictSizes) { 760 if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max); 761 if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min); 762 if (!user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max; 763 } 764 ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr); 765 for (v = vStart; v < vEnd; ++v) { 766 PetscScalar *Mp; 767 PetscReal detM, fact; 768 769 ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr); 770 if (dim == 2) DMPlex_Det2D_Scalar_Internal(&detM, Mp); 771 else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detM, Mp); 772 else SETERRQ1(comm, PETSC_ERR_SUP, "Dimension %d not supported", dim); 773 fact = factGlob * PetscPowReal(PetscAbsReal(detM), -1.0/(2*p+dim)); 774 for (i = 0; i < Nd; ++i) Mp[i] *= fact; 775 if (restrictSizes) { ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, Mp);CHKERRQ(ierr); } 776 } 777 ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr); 778 779 PetscFunctionReturn(0); 780 } 781 782 /* 783 DMPlexMetricAverage - Compute the average of a list of metrics 784 785 Input Parameter: 786 + dm - The DM 787 . numMetrics - The number of metrics to be averaged 788 . weights - Weights for the average 789 - metrics - The metrics to be averaged 790 791 Output Parameter: 792 . metricAvg - The averaged metric 793 794 Level: beginner 795 796 Notes: 797 The weights should sum to unity. 798 799 If weights are not provided then an unweighted average is used. 800 801 .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection() 802 */ 803 PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg) 804 { 805 PetscBool haveWeights = PETSC_TRUE; 806 PetscErrorCode ierr; 807 PetscInt i; 808 PetscReal sum = 0.0, tol = 1.0e-10; 809 810 PetscFunctionBegin; 811 if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); } 812 ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr); 813 ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr); 814 815 /* Default to the unweighted case */ 816 if (!weights) { 817 ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr); 818 haveWeights = PETSC_FALSE; 819 for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; } 820 } 821 822 /* Check weights sum to unity */ 823 for (i = 0; i < numMetrics; ++i) { sum += weights[i]; } 824 if (PetscAbsReal(sum - 1) > tol) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); } 825 826 /* Compute metric average */ 827 for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); } 828 if (!haveWeights) {ierr = PetscFree(weights); } 829 PetscFunctionReturn(0); 830 } 831 832 /* 833 DMPlexMetricAverage2 - Compute the unweighted average of two metrics 834 835 Input Parameter: 836 + dm - The DM 837 . metric1 - The first metric to be averaged 838 - metric2 - The second metric to be averaged 839 840 Output Parameter: 841 . metricAvg - The averaged metric 842 843 Level: beginner 844 845 .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3() 846 */ 847 PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg) 848 { 849 PetscErrorCode ierr; 850 PetscReal weights[2] = {0.5, 0.5}; 851 Vec metrics[2] = {metric1, metric2}; 852 853 PetscFunctionBegin; 854 ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr); 855 PetscFunctionReturn(0); 856 } 857 858 /* 859 DMPlexMetricAverage3 - Compute the unweighted average of three metrics 860 861 Input Parameter: 862 + dm - The DM 863 . metric1 - The first metric to be averaged 864 . metric2 - The second metric to be averaged 865 - metric3 - The third metric to be averaged 866 867 Output Parameter: 868 . metricAvg - The averaged metric 869 870 Level: beginner 871 872 .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2() 873 */ 874 PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg) 875 { 876 PetscErrorCode ierr; 877 PetscReal weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0}; 878 Vec metrics[3] = {metric1, metric2, metric3}; 879 880 PetscFunctionBegin; 881 ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr); 882 PetscFunctionReturn(0); 883 } 884 885 static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[]) 886 { 887 PetscErrorCode ierr; 888 PetscInt i, j, k, l, m; 889 PetscReal *evals, *evals1; 890 PetscScalar *evecs, *sqrtM1, *isqrtM1; 891 892 PetscFunctionBegin; 893 ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr); 894 for (i = 0; i < dim; ++i) { 895 for (j = 0; j < dim; ++j) { 896 evecs[i*dim+j] = M1[i*dim+j]; 897 } 898 } 899 { 900 PetscScalar *work; 901 PetscBLASInt lwork; 902 903 lwork = 5*dim; 904 ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr); 905 { 906 PetscBLASInt lierr, nb; 907 PetscReal sqrtk; 908 909 /* Compute eigendecomposition of M1 */ 910 ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr); 911 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 912 #if defined(PETSC_USE_COMPLEX) 913 { 914 PetscReal *rwork; 915 ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 916 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr)); 917 ierr = PetscFree(rwork);CHKERRQ(ierr); 918 } 919 #else 920 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr)); 921 #endif 922 if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 923 ierr = PetscFPTrapPop(); 924 925 /* Compute square root and reciprocal */ 926 for (i = 0; i < dim; ++i) { 927 for (j = 0; j < dim; ++j) { 928 sqrtM1[i*dim+j] = 0.0; 929 isqrtM1[i*dim+j] = 0.0; 930 for (k = 0; k < dim; ++k) { 931 sqrtk = PetscSqrtReal(evals1[k]); 932 sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j]; 933 isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j]; 934 } 935 } 936 } 937 938 /* Map into the space spanned by the eigenvectors of M1 */ 939 for (i = 0; i < dim; ++i) { 940 for (j = 0; j < dim; ++j) { 941 evecs[i*dim+j] = 0.0; 942 for (k = 0; k < dim; ++k) { 943 for (l = 0; l < dim; ++l) { 944 evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l]; 945 } 946 } 947 } 948 } 949 950 /* Compute eigendecomposition */ 951 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 952 #if defined(PETSC_USE_COMPLEX) 953 { 954 PetscReal *rwork; 955 ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr); 956 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr)); 957 ierr = PetscFree(rwork);CHKERRQ(ierr); 958 } 959 #else 960 PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr)); 961 #endif 962 if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr); 963 ierr = PetscFPTrapPop(); 964 965 /* Modify eigenvalues */ 966 for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]); 967 968 /* Map back to get the intersection */ 969 for (i = 0; i < dim; ++i) { 970 for (j = 0; j < dim; ++j) { 971 M2[i*dim+j] = 0.0; 972 for (k = 0; k < dim; ++k) { 973 for (l = 0; l < dim; ++l) { 974 for (m = 0; m < dim; ++m) { 975 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m]; 976 } 977 } 978 } 979 } 980 } 981 } 982 ierr = PetscFree(work);CHKERRQ(ierr); 983 } 984 ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr); 985 PetscFunctionReturn(0); 986 } 987 988 /* 989 DMPlexMetricIntersection - Compute the intersection of a list of metrics 990 991 Input Parameter: 992 + dm - The DM 993 . numMetrics - The number of metrics to be intersected 994 - metrics - The metrics to be intersected 995 996 Output Parameter: 997 . metricInt - The intersected metric 998 999 Level: beginner 1000 1001 Notes: 1002 1003 The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics. 1004 1005 The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2. 1006 1007 .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage() 1008 */ 1009 PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt) 1010 { 1011 PetscErrorCode ierr; 1012 PetscInt dim, vStart, vEnd, v, i; 1013 PetscScalar *met, *meti, *M, *Mi; 1014 1015 PetscFunctionBegin; 1016 if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); } 1017 1018 /* Extract metadata from dm */ 1019 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1020 ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr); 1021 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1022 1023 /* Copy over the first metric */ 1024 ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr); 1025 1026 /* Intersect subsequent metrics in turn */ 1027 if (numMetrics > 1) { 1028 ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr); 1029 for (i = 1; i < numMetrics; ++i) { 1030 ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr); 1031 for (v = vStart; v < vEnd; ++v) { 1032 ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr); 1033 ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr); 1034 ierr = DMPlexMetricIntersection_Private(dim, Mi, M);CHKERRQ(ierr); 1035 } 1036 ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr); 1037 } 1038 ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr); 1039 } 1040 1041 PetscFunctionReturn(0); 1042 } 1043 1044 /* 1045 DMPlexMetricIntersection2 - Compute the intersection of two metrics 1046 1047 Input Parameter: 1048 + dm - The DM 1049 . metric1 - The first metric to be intersected 1050 - metric2 - The second metric to be intersected 1051 1052 Output Parameter: 1053 . metricInt - The intersected metric 1054 1055 Level: beginner 1056 1057 .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3() 1058 */ 1059 PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt) 1060 { 1061 PetscErrorCode ierr; 1062 Vec metrics[2] = {metric1, metric2}; 1063 1064 PetscFunctionBegin; 1065 ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr); 1066 PetscFunctionReturn(0); 1067 } 1068 1069 /* 1070 DMPlexMetricIntersection3 - Compute the intersection of three metrics 1071 1072 Input Parameter: 1073 + dm - The DM 1074 . metric1 - The first metric to be intersected 1075 . metric2 - The second metric to be intersected 1076 - metric3 - The third metric to be intersected 1077 1078 Output Parameter: 1079 . metricInt - The intersected metric 1080 1081 Level: beginner 1082 1083 .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2() 1084 */ 1085 PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt) 1086 { 1087 PetscErrorCode ierr; 1088 Vec metrics[3] = {metric1, metric2, metric3}; 1089 1090 PetscFunctionBegin; 1091 ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr); 1092 PetscFunctionReturn(0); 1093 } 1094