1 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 2 3 #include <petscdmplex.h> /* For DMProjectCoordinates() */ 4 #include <petscsf.h> /* For DMLocatePoints() */ 5 6 PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, void *ctx) 7 { 8 DM dm_coord, dmc_coord; 9 Vec coords, ccoords; 10 Mat inject; 11 PetscFunctionBegin; 12 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 13 PetscCall(DMGetCoordinateDM(dmc, &dmc_coord)); 14 PetscCall(DMGetCoordinates(dm, &coords)); 15 PetscCall(DMGetCoordinates(dmc, &ccoords)); 16 if (coords && !ccoords) { 17 PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords)); 18 PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates")); 19 PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject)); 20 PetscCall(MatRestrict(inject, coords, ccoords)); 21 PetscCall(MatDestroy(&inject)); 22 PetscCall(DMSetCoordinates(dmc, ccoords)); 23 PetscCall(VecDestroy(&ccoords)); 24 } 25 PetscFunctionReturn(0); 26 } 27 28 static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx) 29 { 30 DM dm_coord, subdm_coord; 31 Vec coords, ccoords, clcoords; 32 VecScatter *scat_i, *scat_g; 33 PetscFunctionBegin; 34 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 35 PetscCall(DMGetCoordinateDM(subdm, &subdm_coord)); 36 PetscCall(DMGetCoordinates(dm, &coords)); 37 PetscCall(DMGetCoordinates(subdm, &ccoords)); 38 if (coords && !ccoords) { 39 PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords)); 40 PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates")); 41 PetscCall(DMCreateLocalVector(subdm_coord, &clcoords)); 42 PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates")); 43 PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g)); 44 PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD)); 45 PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD)); 46 PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD)); 47 PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD)); 48 PetscCall(DMSetCoordinates(subdm, ccoords)); 49 PetscCall(DMSetCoordinatesLocal(subdm, clcoords)); 50 PetscCall(VecScatterDestroy(&scat_i[0])); 51 PetscCall(VecScatterDestroy(&scat_g[0])); 52 PetscCall(VecDestroy(&ccoords)); 53 PetscCall(VecDestroy(&clcoords)); 54 PetscCall(PetscFree(scat_i)); 55 PetscCall(PetscFree(scat_g)); 56 } 57 PetscFunctionReturn(0); 58 } 59 60 /*@ 61 DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates 62 63 Collective on dm 64 65 Input Parameter: 66 . dm - the DM 67 68 Output Parameter: 69 . cdm - coordinate DM 70 71 Level: intermediate 72 73 .seealso: `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 74 @*/ 75 PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm) 76 { 77 PetscFunctionBegin; 78 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 79 PetscValidPointer(cdm, 2); 80 if (!dm->coordinates[0].dm) { 81 DM cdm; 82 83 PetscUseTypeMethod(dm, createcoordinatedm, &cdm); 84 PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM")); 85 /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup 86 * until the call to CreateCoordinateDM) */ 87 PetscCall(DMDestroy(&dm->coordinates[0].dm)); 88 dm->coordinates[0].dm = cdm; 89 } 90 *cdm = dm->coordinates[0].dm; 91 PetscFunctionReturn(0); 92 } 93 94 /*@ 95 DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates 96 97 Logically Collective on dm 98 99 Input Parameters: 100 + dm - the DM 101 - cdm - coordinate DM 102 103 Level: intermediate 104 105 .seealso: `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 106 @*/ 107 PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm) 108 { 109 PetscFunctionBegin; 110 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 111 PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 112 PetscCall(PetscObjectReference((PetscObject)cdm)); 113 PetscCall(DMDestroy(&dm->coordinates[0].dm)); 114 dm->coordinates[0].dm = cdm; 115 PetscFunctionReturn(0); 116 } 117 118 /*@ 119 DMGetCellCoordinateDM - Gets the DM that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates 120 121 Collective on dm 122 123 Input Parameter: 124 . dm - the DM 125 126 Output Parameter: 127 . cdm - cellwise coordinate DM, or NULL if they are not defined 128 129 Note: 130 Call DMLocalizeCoordinates() to automatically create cellwise coordinates for periodic geometries. 131 132 Level: intermediate 133 134 .seealso: `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMLocalizeCoordinates()` 135 @*/ 136 PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm) 137 { 138 PetscFunctionBegin; 139 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 140 PetscValidPointer(cdm, 2); 141 *cdm = dm->coordinates[1].dm; 142 PetscFunctionReturn(0); 143 } 144 145 /*@ 146 DMSetCellCoordinateDM - Sets the DM that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates 147 148 Logically Collective on dm 149 150 Input Parameters: 151 + dm - the DM 152 - cdm - cellwise coordinate DM 153 154 Level: intermediate 155 156 .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()` 157 @*/ 158 PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm) 159 { 160 PetscInt dim; 161 162 PetscFunctionBegin; 163 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 164 if (cdm) { 165 PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 166 PetscCall(DMGetCoordinateDim(dm, &dim)); 167 dm->coordinates[1].dim = dim; 168 } 169 PetscCall(PetscObjectReference((PetscObject)cdm)); 170 PetscCall(DMDestroy(&dm->coordinates[1].dm)); 171 dm->coordinates[1].dm = cdm; 172 PetscFunctionReturn(0); 173 } 174 175 /*@ 176 DMGetCoordinateDim - Retrieve the dimension of embedding space for coordinate values. 177 178 Not Collective 179 180 Input Parameter: 181 . dm - The DM object 182 183 Output Parameter: 184 . dim - The embedding dimension 185 186 Level: intermediate 187 188 .seealso: `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 189 @*/ 190 PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim) 191 { 192 PetscFunctionBegin; 193 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 194 PetscValidIntPointer(dim, 2); 195 if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim; 196 *dim = dm->coordinates[0].dim; 197 PetscFunctionReturn(0); 198 } 199 200 /*@ 201 DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values. 202 203 Not Collective 204 205 Input Parameters: 206 + dm - The DM object 207 - dim - The embedding dimension 208 209 Level: intermediate 210 211 .seealso: `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()` 212 @*/ 213 PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim) 214 { 215 PetscDS ds; 216 PetscInt Nds, n; 217 218 PetscFunctionBegin; 219 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 220 dm->coordinates[0].dim = dim; 221 if (dm->dim >= 0) { 222 PetscCall(DMGetNumDS(dm, &Nds)); 223 for (n = 0; n < Nds; ++n) { 224 PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds)); 225 PetscCall(PetscDSSetCoordinateDimension(ds, dim)); 226 } 227 } 228 PetscFunctionReturn(0); 229 } 230 231 /*@ 232 DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh. 233 234 Collective on dm 235 236 Input Parameter: 237 . dm - The DM object 238 239 Output Parameter: 240 . section - The PetscSection object 241 242 Level: intermediate 243 244 .seealso: `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 245 @*/ 246 PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section) 247 { 248 DM cdm; 249 250 PetscFunctionBegin; 251 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 252 PetscValidPointer(section, 2); 253 PetscCall(DMGetCoordinateDM(dm, &cdm)); 254 PetscCall(DMGetLocalSection(cdm, section)); 255 PetscFunctionReturn(0); 256 } 257 258 /*@ 259 DMSetCoordinateSection - Set the layout of coordinate values over the mesh. 260 261 Not Collective 262 263 Input Parameters: 264 + dm - The DM object 265 . dim - The embedding dimension, or PETSC_DETERMINE 266 - section - The PetscSection object 267 268 Level: intermediate 269 270 .seealso: `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()` 271 @*/ 272 PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section) 273 { 274 DM cdm; 275 276 PetscFunctionBegin; 277 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 278 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 279 PetscCall(DMGetCoordinateDM(dm, &cdm)); 280 PetscCall(DMSetLocalSection(cdm, section)); 281 if (dim == PETSC_DETERMINE) { 282 PetscInt d = PETSC_DEFAULT; 283 PetscInt pStart, pEnd, vStart, vEnd, v, dd; 284 285 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 286 PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd)); 287 pStart = PetscMax(vStart, pStart); 288 pEnd = PetscMin(vEnd, pEnd); 289 for (v = pStart; v < pEnd; ++v) { 290 PetscCall(PetscSectionGetDof(section, v, &dd)); 291 if (dd) { 292 d = dd; 293 break; 294 } 295 } 296 if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d)); 297 } 298 PetscFunctionReturn(0); 299 } 300 301 /*@ 302 DMGetCellCoordinateSection - Retrieve the layout of cellwise coordinate values over the mesh. 303 304 Collective on dm 305 306 Input Parameter: 307 . dm - The DM object 308 309 Output Parameter: 310 . section - The PetscSection object, or NULL if no cellwise coordinates are defined 311 312 Level: intermediate 313 314 .seealso: `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 315 @*/ 316 PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section) 317 { 318 DM cdm; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 322 PetscValidPointer(section, 2); 323 *section = NULL; 324 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 325 if (cdm) PetscCall(DMGetLocalSection(cdm, section)); 326 PetscFunctionReturn(0); 327 } 328 329 /*@ 330 DMSetCellCoordinateSection - Set the layout of cellwise coordinate values over the mesh. 331 332 Not Collective 333 334 Input Parameters: 335 + dm - The DM object 336 . dim - The embedding dimension, or PETSC_DETERMINE 337 - section - The PetscSection object for a cellwise layout 338 339 Level: intermediate 340 341 .seealso: `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 342 @*/ 343 PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section) 344 { 345 DM cdm; 346 347 PetscFunctionBegin; 348 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 349 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 350 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 351 PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates"); 352 PetscCall(DMSetLocalSection(cdm, section)); 353 if (dim == PETSC_DETERMINE) { 354 PetscInt d = PETSC_DEFAULT; 355 PetscInt pStart, pEnd, vStart, vEnd, v, dd; 356 357 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 358 PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd)); 359 pStart = PetscMax(vStart, pStart); 360 pEnd = PetscMin(vEnd, pEnd); 361 for (v = pStart; v < pEnd; ++v) { 362 PetscCall(PetscSectionGetDof(section, v, &dd)); 363 if (dd) { 364 d = dd; 365 break; 366 } 367 } 368 if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d)); 369 } 370 PetscFunctionReturn(0); 371 } 372 373 /*@ 374 DMGetCoordinates - Gets a global vector with the coordinates associated with the DM. 375 376 Collective on dm 377 378 Input Parameter: 379 . dm - the DM 380 381 Output Parameter: 382 . c - global coordinate vector 383 384 Note: 385 This is a borrowed reference, so the user should NOT destroy this vector. When the DM is 386 destroyed the array will no longer be valid. 387 388 Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates). 389 390 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 391 and (x_0,y_0,z_0,x_1,y_1,z_1...) 392 393 Level: intermediate 394 395 .seealso: `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()` 396 @*/ 397 PetscErrorCode DMGetCoordinates(DM dm, Vec *c) 398 { 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 401 PetscValidPointer(c, 2); 402 if (!dm->coordinates[0].x && dm->coordinates[0].xl) { 403 DM cdm = NULL; 404 405 PetscCall(DMGetCoordinateDM(dm, &cdm)); 406 PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x)); 407 PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates")); 408 PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x)); 409 PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x)); 410 } 411 *c = dm->coordinates[0].x; 412 PetscFunctionReturn(0); 413 } 414 415 /*@ 416 DMSetCoordinates - Sets into the DM a global vector that holds the coordinates 417 418 Collective on dm 419 420 Input Parameters: 421 + dm - the DM 422 - c - coordinate vector 423 424 Notes: 425 The coordinates do include those for ghost points, which are in the local vector. 426 427 The vector c should be destroyed by the caller. 428 429 Level: intermediate 430 431 .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()` 432 @*/ 433 PetscErrorCode DMSetCoordinates(DM dm, Vec c) 434 { 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 437 if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 438 PetscCall(PetscObjectReference((PetscObject)c)); 439 PetscCall(VecDestroy(&dm->coordinates[0].x)); 440 dm->coordinates[0].x = c; 441 PetscCall(VecDestroy(&dm->coordinates[0].xl)); 442 PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL)); 443 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL)); 444 PetscFunctionReturn(0); 445 } 446 447 /*@ 448 DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the DM. 449 450 Collective on dm 451 452 Input Parameter: 453 . dm - the DM 454 455 Output Parameter: 456 . c - global coordinate vector 457 458 Note: 459 This is a borrowed reference, so the user should NOT destroy this vector. When the DM is 460 destroyed the array will no longer be valid. 461 462 Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates). 463 464 Level: intermediate 465 466 .seealso: `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()` 467 @*/ 468 PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c) 469 { 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 472 PetscValidPointer(c, 2); 473 if (!dm->coordinates[1].x && dm->coordinates[1].xl) { 474 DM cdm = NULL; 475 476 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 477 PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x)); 478 PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates")); 479 PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x)); 480 PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x)); 481 } 482 *c = dm->coordinates[1].x; 483 PetscFunctionReturn(0); 484 } 485 486 /*@ 487 DMSetCellCoordinates - Sets into the DM a global vector that holds the cellwise coordinates 488 489 Collective on dm 490 491 Input Parameters: 492 + dm - the DM 493 - c - cellwise coordinate vector 494 495 Notes: 496 The coordinates do include those for ghost points, which are in the local vector. 497 498 The vector c should be destroyed by the caller. 499 500 Level: intermediate 501 502 .seealso: `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()` 503 @*/ 504 PetscErrorCode DMSetCellCoordinates(DM dm, Vec c) 505 { 506 PetscFunctionBegin; 507 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 508 if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 509 PetscCall(PetscObjectReference((PetscObject)c)); 510 PetscCall(VecDestroy(&dm->coordinates[1].x)); 511 dm->coordinates[1].x = c; 512 PetscCall(VecDestroy(&dm->coordinates[1].xl)); 513 PetscFunctionReturn(0); 514 } 515 516 /*@ 517 DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that DMGetCoordinatesLocalNoncollective() can be used as non-collective afterwards. 518 519 Collective on dm 520 521 Input Parameter: 522 . dm - the DM 523 524 Level: advanced 525 526 .seealso: `DMGetCoordinatesLocalNoncollective()` 527 @*/ 528 PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm) 529 { 530 PetscFunctionBegin; 531 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 532 if (!dm->coordinates[0].xl && dm->coordinates[0].x) { 533 DM cdm = NULL; 534 535 PetscCall(DMGetCoordinateDM(dm, &cdm)); 536 PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl)); 537 PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates")); 538 PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl)); 539 PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl)); 540 } 541 PetscFunctionReturn(0); 542 } 543 544 /*@ 545 DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM. 546 547 Collective on dm 548 549 Input Parameter: 550 . dm - the DM 551 552 Output Parameter: 553 . c - coordinate vector 554 555 Note: 556 This is a borrowed reference, so the user should NOT destroy this vector 557 558 Each process has the local and ghost coordinates 559 560 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 561 and (x_0,y_0,z_0,x_1,y_1,z_1...) 562 563 Level: intermediate 564 565 .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()` 566 @*/ 567 PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c) 568 { 569 PetscFunctionBegin; 570 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 571 PetscValidPointer(c, 2); 572 PetscCall(DMGetCoordinatesLocalSetUp(dm)); 573 *c = dm->coordinates[0].xl; 574 PetscFunctionReturn(0); 575 } 576 577 /*@ 578 DMGetCoordinatesLocalNoncollective - Non-collective version of DMGetCoordinatesLocal(). Fails if global coordinates have been set and DMGetCoordinatesLocalSetUp() not called. 579 580 Not collective 581 582 Input Parameter: 583 . dm - the DM 584 585 Output Parameter: 586 . c - coordinate vector 587 588 Level: advanced 589 590 .seealso: `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()` 591 @*/ 592 PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c) 593 { 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 596 PetscValidPointer(c, 2); 597 PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called"); 598 *c = dm->coordinates[0].xl; 599 PetscFunctionReturn(0); 600 } 601 602 /*@ 603 DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and section describing its layout. 604 605 Not collective 606 607 Input Parameters: 608 + dm - the DM 609 - p - the IS of points whose coordinates will be returned 610 611 Output Parameters: 612 + pCoordSection - the PetscSection describing the layout of pCoord, i.e. each point corresponds to one point in p, and DOFs correspond to coordinates 613 - pCoord - the Vec with coordinates of points in p 614 615 Note: 616 DMGetCoordinatesLocalSetUp() must be called first. This function employs DMGetCoordinatesLocalNoncollective() so it is not collective. 617 618 This creates a new vector, so the user SHOULD destroy this vector 619 620 Each process has the local and ghost coordinates 621 622 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 623 and (x_0,y_0,z_0,x_1,y_1,z_1...) 624 625 Level: advanced 626 627 .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()` 628 @*/ 629 PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord) 630 { 631 DM cdm; 632 PetscSection cs, newcs; 633 Vec coords; 634 const PetscScalar *arr; 635 PetscScalar *newarr = NULL; 636 PetscInt n; 637 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 640 PetscValidHeaderSpecific(p, IS_CLASSID, 2); 641 if (pCoordSection) PetscValidPointer(pCoordSection, 3); 642 if (pCoord) PetscValidPointer(pCoord, 4); 643 PetscCall(DMGetCoordinateDM(dm, &cdm)); 644 PetscCall(DMGetLocalSection(cdm, &cs)); 645 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 646 PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set"); 647 PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported"); 648 PetscCall(VecGetArrayRead(coords, &arr)); 649 PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL)); 650 PetscCall(VecRestoreArrayRead(coords, &arr)); 651 if (pCoord) { 652 PetscCall(PetscSectionGetStorageSize(newcs, &n)); 653 /* set array in two steps to mimic PETSC_OWN_POINTER */ 654 PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord)); 655 PetscCall(VecReplaceArray(*pCoord, newarr)); 656 } else { 657 PetscCall(PetscFree(newarr)); 658 } 659 if (pCoordSection) { 660 *pCoordSection = newcs; 661 } else PetscCall(PetscSectionDestroy(&newcs)); 662 PetscFunctionReturn(0); 663 } 664 665 /*@ 666 DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates 667 668 Not collective 669 670 Input Parameters: 671 + dm - the DM 672 - c - coordinate vector 673 674 Notes: 675 The coordinates of ghost points can be set using DMSetCoordinates() 676 followed by DMGetCoordinatesLocal(). This is intended to enable the 677 setting of ghost coordinates outside of the domain. 678 679 The vector c should be destroyed by the caller. 680 681 Level: intermediate 682 683 .seealso: `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()` 684 @*/ 685 PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c) 686 { 687 PetscFunctionBegin; 688 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 689 if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 690 PetscCall(PetscObjectReference((PetscObject)c)); 691 PetscCall(VecDestroy(&dm->coordinates[0].xl)); 692 dm->coordinates[0].xl = c; 693 PetscCall(VecDestroy(&dm->coordinates[0].x)); 694 PetscFunctionReturn(0); 695 } 696 697 /*@ 698 DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that DMGetCellCoordinatesLocalNoncollective() can be used as non-collective afterwards. 699 700 Collective on dm 701 702 Input Parameter: 703 . dm - the DM 704 705 Level: advanced 706 707 .seealso: `DMGetCellCoordinatesLocalNoncollective()` 708 @*/ 709 PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm) 710 { 711 PetscFunctionBegin; 712 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 713 if (!dm->coordinates[1].xl && dm->coordinates[1].x) { 714 DM cdm = NULL; 715 716 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 717 PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl)); 718 PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates")); 719 PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl)); 720 PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl)); 721 } 722 PetscFunctionReturn(0); 723 } 724 725 /*@ 726 DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the DM. 727 728 Collective on dm 729 730 Input Parameter: 731 . dm - the DM 732 733 Output Parameter: 734 . c - coordinate vector 735 736 Note: 737 This is a borrowed reference, so the user should NOT destroy this vector 738 739 Each process has the local and ghost coordinates 740 741 Level: intermediate 742 743 .seealso: `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()` 744 @*/ 745 PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c) 746 { 747 PetscFunctionBegin; 748 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 749 PetscValidPointer(c, 2); 750 PetscCall(DMGetCellCoordinatesLocalSetUp(dm)); 751 *c = dm->coordinates[1].xl; 752 PetscFunctionReturn(0); 753 } 754 755 /*@ 756 DMGetCellCoordinatesLocalNoncollective - Non-collective version of DMGetCellCoordinatesLocal(). Fails if global cellwise coordinates have been set and DMGetCellCoordinatesLocalSetUp() not called. 757 758 Not collective 759 760 Input Parameter: 761 . dm - the DM 762 763 Output Parameter: 764 . c - cellwise coordinate vector 765 766 Level: advanced 767 768 .seealso: `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()` 769 @*/ 770 PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c) 771 { 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 774 PetscValidPointer(c, 2); 775 PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called"); 776 *c = dm->coordinates[1].xl; 777 PetscFunctionReturn(0); 778 } 779 780 /*@ 781 DMSetCellCoordinatesLocal - Sets into the DM a local vector that holds the cellwise coordinates 782 783 Not collective 784 785 Input Parameters: 786 + dm - the DM 787 - c - cellwise coordinate vector 788 789 Notes: 790 The coordinates of ghost points can be set using DMSetCoordinates() 791 followed by DMGetCoordinatesLocal(). This is intended to enable the 792 setting of ghost coordinates outside of the domain. 793 794 The vector c should be destroyed by the caller. 795 796 Level: intermediate 797 798 .seealso: `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()` 799 @*/ 800 PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c) 801 { 802 PetscFunctionBegin; 803 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 804 if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 805 PetscCall(PetscObjectReference((PetscObject)c)); 806 PetscCall(VecDestroy(&dm->coordinates[1].xl)); 807 dm->coordinates[1].xl = c; 808 PetscCall(VecDestroy(&dm->coordinates[1].x)); 809 PetscFunctionReturn(0); 810 } 811 812 PetscErrorCode DMGetCoordinateField(DM dm, DMField *field) 813 { 814 PetscFunctionBegin; 815 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 816 PetscValidPointer(field, 2); 817 if (!dm->coordinates[0].field) { 818 if (dm->ops->createcoordinatefield) PetscCall((*dm->ops->createcoordinatefield)(dm, &dm->coordinates[0].field)); 819 } 820 *field = dm->coordinates[0].field; 821 PetscFunctionReturn(0); 822 } 823 824 PetscErrorCode DMSetCoordinateField(DM dm, DMField field) 825 { 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 828 if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2); 829 PetscCall(PetscObjectReference((PetscObject)field)); 830 PetscCall(DMFieldDestroy(&dm->coordinates[0].field)); 831 dm->coordinates[0].field = field; 832 PetscFunctionReturn(0); 833 } 834 835 /*@ 836 DMGetLocalBoundingBox - Returns the bounding box for the piece of the DM on this process. 837 838 Not collective 839 840 Input Parameter: 841 . dm - the DM 842 843 Output Parameters: 844 + lmin - local minimum coordinates (length coord dim, optional) 845 - lmax - local maximim coordinates (length coord dim, optional) 846 847 Level: beginner 848 849 Note: If the DM is a DMDA and has no coordinates, the index bounds are returned instead. 850 851 .seealso: `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()` 852 @*/ 853 PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[]) 854 { 855 Vec coords = NULL; 856 PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 857 PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 858 PetscInt cdim, i, j; 859 860 PetscFunctionBegin; 861 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 862 PetscCall(DMGetCoordinateDim(dm, &cdim)); 863 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 864 if (coords) { 865 const PetscScalar *local_coords; 866 PetscInt N, Ni; 867 868 for (j = cdim; j < 3; ++j) { 869 min[j] = 0; 870 max[j] = 0; 871 } 872 PetscCall(VecGetArrayRead(coords, &local_coords)); 873 PetscCall(VecGetLocalSize(coords, &N)); 874 Ni = N / cdim; 875 for (i = 0; i < Ni; ++i) { 876 for (j = 0; j < cdim; ++j) { 877 min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j])); 878 max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j])); 879 } 880 } 881 PetscCall(VecRestoreArrayRead(coords, &local_coords)); 882 PetscCall(DMGetCellCoordinatesLocal(dm, &coords)); 883 if (coords) { 884 PetscCall(VecGetArrayRead(coords, &local_coords)); 885 PetscCall(VecGetLocalSize(coords, &N)); 886 Ni = N / cdim; 887 for (i = 0; i < Ni; ++i) { 888 for (j = 0; j < cdim; ++j) { 889 min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j])); 890 max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j])); 891 } 892 } 893 PetscCall(VecRestoreArrayRead(coords, &local_coords)); 894 } 895 } else { 896 PetscBool isda; 897 898 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 899 if (isda) PetscCall(DMGetLocalBoundingIndices_DMDA(dm, min, max)); 900 } 901 if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim)); 902 if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim)); 903 PetscFunctionReturn(0); 904 } 905 906 /*@ 907 DMGetBoundingBox - Returns the global bounding box for the DM. 908 909 Collective 910 911 Input Parameter: 912 . dm - the DM 913 914 Output Parameters: 915 + gmin - global minimum coordinates (length coord dim, optional) 916 - gmax - global maximim coordinates (length coord dim, optional) 917 918 Level: beginner 919 920 .seealso: `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 921 @*/ 922 PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[]) 923 { 924 PetscReal lmin[3], lmax[3]; 925 PetscInt cdim; 926 PetscMPIInt count; 927 928 PetscFunctionBegin; 929 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 930 PetscCall(DMGetCoordinateDim(dm, &cdim)); 931 PetscCall(PetscMPIIntCast(cdim, &count)); 932 PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax)); 933 if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm))); 934 if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm))); 935 PetscFunctionReturn(0); 936 } 937 938 /*@ 939 DMProjectCoordinates - Project coordinates to a different space 940 941 Input Parameters: 942 + dm - The DM object 943 - disc - The new coordinate discretization or NULL to ensure a coordinate discretization exists 944 945 Level: intermediate 946 947 .seealso: `DMGetCoordinateField()` 948 @*/ 949 PetscErrorCode DMProjectCoordinates(DM dm, PetscFE disc) 950 { 951 PetscFE discOld; 952 PetscClassId classid; 953 DM cdmOld, cdmNew; 954 Vec coordsOld, coordsNew; 955 Mat matInterp; 956 PetscBool same_space = PETSC_TRUE; 957 958 PetscFunctionBegin; 959 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 960 if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2); 961 962 PetscCall(DMGetCoordinateDM(dm, &cdmOld)); 963 /* Check current discretization is compatible */ 964 PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld)); 965 PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid)); 966 if (classid != PETSCFE_CLASSID) { 967 if (classid == PETSC_CONTAINER_CLASSID) { 968 PetscFE feLinear; 969 DMPolytopeType ct; 970 PetscInt dim, dE, cStart, cEnd; 971 PetscBool simplex; 972 973 /* Assume linear vertex coordinates */ 974 PetscCall(DMGetDimension(dm, &dim)); 975 PetscCall(DMGetCoordinateDim(dm, &dE)); 976 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 977 if (cEnd > cStart) { 978 PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 979 switch (ct) { 980 case DM_POLYTOPE_TRI_PRISM: 981 case DM_POLYTOPE_TRI_PRISM_TENSOR: 982 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot autoamtically create coordinate space for prisms"); 983 default: 984 break; 985 } 986 } 987 PetscCall(DMPlexIsSimplex(dm, &simplex)); 988 PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, simplex, 1, -1, &feLinear)); 989 PetscCall(DMSetField(cdmOld, 0, NULL, (PetscObject)feLinear)); 990 PetscCall(PetscFEDestroy(&feLinear)); 991 PetscCall(DMCreateDS(cdmOld)); 992 PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld)); 993 } else { 994 const char *discname; 995 996 PetscCall(PetscObjectGetType((PetscObject)discOld, &discname)); 997 SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname); 998 } 999 } 1000 if (!disc) PetscFunctionReturn(0); 1001 { // Check if the new space is the same as the old modulo quadrature 1002 PetscDualSpace dsOld, ds; 1003 PetscCall(PetscFEGetDualSpace(discOld, &dsOld)); 1004 PetscCall(PetscFEGetDualSpace(disc, &ds)); 1005 PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space)); 1006 } 1007 /* Make a fresh clone of the coordinate DM */ 1008 PetscCall(DMClone(cdmOld, &cdmNew)); 1009 PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc)); 1010 PetscCall(DMCreateDS(cdmNew)); 1011 PetscCall(DMGetCoordinates(dm, &coordsOld)); 1012 if (same_space) { 1013 PetscCall(PetscObjectReference((PetscObject)coordsOld)); 1014 coordsNew = coordsOld; 1015 } else { // Project the coordinate vector from old to new space 1016 PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew)); 1017 PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &matInterp, NULL)); 1018 PetscCall(MatInterpolate(matInterp, coordsOld, coordsNew)); 1019 PetscCall(MatDestroy(&matInterp)); 1020 } 1021 /* Set new coordinate structures */ 1022 PetscCall(DMSetCoordinateField(dm, NULL)); 1023 PetscCall(DMSetCoordinateDM(dm, cdmNew)); 1024 PetscCall(DMSetCoordinates(dm, coordsNew)); 1025 PetscCall(VecDestroy(&coordsNew)); 1026 PetscCall(DMDestroy(&cdmNew)); 1027 PetscFunctionReturn(0); 1028 } 1029 1030 /*@ 1031 DMLocatePoints - Locate the points in v in the mesh and return a PetscSF of the containing cells 1032 1033 Collective on v (see explanation below) 1034 1035 Input Parameters: 1036 + dm - The DM 1037 - ltype - The type of point location, e.g. DM_POINTLOCATION_NONE or DM_POINTLOCATION_NEAREST 1038 1039 Input/Output Parameters: 1040 + v - The Vec of points, on output contains the nearest mesh points to the given points if DM_POINTLOCATION_NEAREST is used 1041 - cellSF - Points to either NULL, or a PetscSF with guesses for which cells contain each point; 1042 on output, the PetscSF containing the ranks and local indices of the containing points 1043 1044 Level: developer 1045 1046 Notes: 1047 To do a search of the local cells of the mesh, v should have PETSC_COMM_SELF as its communicator. 1048 To do a search of all the cells in the distributed mesh, v should have the same communicator as dm. 1049 1050 Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions. 1051 1052 If *cellSF is NULL on input, a PetscSF will be created. 1053 If *cellSF is not NULL on input, it should point to an existing PetscSF, whose graph will be used as initial guesses. 1054 1055 An array that maps each point to its containing cell can be obtained with 1056 1057 $ const PetscSFNode *cells; 1058 $ PetscInt nFound; 1059 $ const PetscInt *found; 1060 $ 1061 $ PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells); 1062 1063 Where cells[i].rank is the rank of the cell containing point found[i] (or i if found == NULL), and cells[i].index is 1064 the index of the cell in its rank's local numbering. 1065 1066 .seealso: `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType` 1067 @*/ 1068 PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF) 1069 { 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1072 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 1073 PetscValidPointer(cellSF, 4); 1074 if (*cellSF) { 1075 PetscMPIInt result; 1076 1077 PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4); 1078 PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result)); 1079 PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's"); 1080 } else { 1081 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF)); 1082 } 1083 PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0)); 1084 PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF); 1085 PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0)); 1086 PetscFunctionReturn(0); 1087 } 1088