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