1 #include <petscdmshell.h> /*I "petscdmshell.h" I*/ 2 #include <petscmat.h> 3 #include <petsc/private/dmimpl.h> 4 5 typedef struct { 6 Vec Xglobal; 7 Vec Xlocal; 8 Mat A; 9 VecScatter gtol; 10 VecScatter ltog; 11 VecScatter ltol; 12 void *ctx; 13 PetscErrorCode (*destroyctx)(void *); 14 } DM_Shell; 15 16 /*@ 17 DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to begin a global to local scatter 18 Collective 19 20 Input Parameters: 21 + dm - shell DM 22 . g - global vector 23 . mode - InsertMode 24 - l - local vector 25 26 Level: advanced 27 28 Note: This is not normally called directly by user code, generally user code calls DMGlobalToLocalBegin() and DMGlobalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function. 29 30 .seealso: `DMGlobalToLocalEndDefaultShell()` 31 @*/ 32 PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm, Vec g, InsertMode mode, Vec l) { 33 DM_Shell *shell = (DM_Shell *)dm->data; 34 35 PetscFunctionBegin; 36 PetscCheck(shell->gtol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 37 PetscCall(VecScatterBegin(shell->gtol, g, l, mode, SCATTER_FORWARD)); 38 PetscFunctionReturn(0); 39 } 40 41 /*@ 42 DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter 43 Collective 44 45 Input Parameters: 46 + dm - shell DM 47 . g - global vector 48 . mode - InsertMode 49 - l - local vector 50 51 Level: advanced 52 53 .seealso: `DMGlobalToLocalBeginDefaultShell()` 54 @*/ 55 PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm, Vec g, InsertMode mode, Vec l) { 56 DM_Shell *shell = (DM_Shell *)dm->data; 57 58 PetscFunctionBegin; 59 PetscCheck(shell->gtol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 60 PetscCall(VecScatterEnd(shell->gtol, g, l, mode, SCATTER_FORWARD)); 61 PetscFunctionReturn(0); 62 } 63 64 /*@ 65 DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter 66 Collective 67 68 Input Parameters: 69 + dm - shell DM 70 . l - local vector 71 . mode - InsertMode 72 - g - global vector 73 74 Level: advanced 75 76 Note: This is not normally called directly by user code, generally user code calls DMLocalToGlobalBegin() and DMLocalToGlobalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function. 77 78 .seealso: `DMLocalToGlobalEndDefaultShell()` 79 @*/ 80 PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm, Vec l, InsertMode mode, Vec g) { 81 DM_Shell *shell = (DM_Shell *)dm->data; 82 83 PetscFunctionBegin; 84 PetscCheck(shell->ltog, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()"); 85 PetscCall(VecScatterBegin(shell->ltog, l, g, mode, SCATTER_FORWARD)); 86 PetscFunctionReturn(0); 87 } 88 89 /*@ 90 DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter 91 Collective 92 93 Input Parameters: 94 + dm - shell DM 95 . l - local vector 96 . mode - InsertMode 97 - g - global vector 98 99 Level: advanced 100 101 .seealso: `DMLocalToGlobalBeginDefaultShell()` 102 @*/ 103 PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm, Vec l, InsertMode mode, Vec g) { 104 DM_Shell *shell = (DM_Shell *)dm->data; 105 106 PetscFunctionBegin; 107 PetscCheck(shell->ltog, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()"); 108 PetscCall(VecScatterEnd(shell->ltog, l, g, mode, SCATTER_FORWARD)); 109 PetscFunctionReturn(0); 110 } 111 112 /*@ 113 DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter 114 Collective 115 116 Input Parameters: 117 + dm - shell DM 118 . g - the original local vector 119 - mode - InsertMode 120 121 Output Parameter: 122 . l - the local vector with correct ghost values 123 124 Level: advanced 125 126 Note: This is not normally called directly by user code, generally user code calls DMLocalToLocalBegin() and DMLocalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToLocal() then those routines might have reason to call this function. 127 128 .seealso: `DMLocalToLocalEndDefaultShell()` 129 @*/ 130 PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm, Vec g, InsertMode mode, Vec l) { 131 DM_Shell *shell = (DM_Shell *)dm->data; 132 133 PetscFunctionBegin; 134 PetscCheck(shell->ltol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()"); 135 PetscCall(VecScatterBegin(shell->ltol, g, l, mode, SCATTER_FORWARD)); 136 PetscFunctionReturn(0); 137 } 138 139 /*@ 140 DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter 141 Collective 142 143 Input Parameters: 144 + dm - shell DM 145 . g - the original local vector 146 - mode - InsertMode 147 148 Output Parameter: 149 . l - the local vector with correct ghost values 150 151 Level: advanced 152 153 .seealso: `DMLocalToLocalBeginDefaultShell()` 154 @*/ 155 PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm, Vec g, InsertMode mode, Vec l) { 156 DM_Shell *shell = (DM_Shell *)dm->data; 157 158 PetscFunctionBegin; 159 PetscCheck(shell->ltol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 160 PetscCall(VecScatterEnd(shell->ltol, g, l, mode, SCATTER_FORWARD)); 161 PetscFunctionReturn(0); 162 } 163 164 static PetscErrorCode DMCreateMatrix_Shell(DM dm, Mat *J) { 165 DM_Shell *shell = (DM_Shell *)dm->data; 166 Mat A; 167 168 PetscFunctionBegin; 169 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 170 PetscValidPointer(J, 2); 171 if (!shell->A) { 172 if (shell->Xglobal) { 173 PetscInt m, M; 174 PetscCall(PetscInfo(dm, "Naively creating matrix using global vector distribution without preallocation\n")); 175 PetscCall(VecGetSize(shell->Xglobal, &M)); 176 PetscCall(VecGetLocalSize(shell->Xglobal, &m)); 177 PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), &shell->A)); 178 PetscCall(MatSetSizes(shell->A, m, m, M, M)); 179 PetscCall(MatSetType(shell->A, dm->mattype)); 180 PetscCall(MatSetUp(shell->A)); 181 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector"); 182 } 183 A = shell->A; 184 PetscCall(MatDuplicate(A, MAT_SHARE_NONZERO_PATTERN, J)); 185 PetscCall(MatSetDM(*J, dm)); 186 PetscFunctionReturn(0); 187 } 188 189 PetscErrorCode DMCreateGlobalVector_Shell(DM dm, Vec *gvec) { 190 DM_Shell *shell = (DM_Shell *)dm->data; 191 Vec X; 192 193 PetscFunctionBegin; 194 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 195 PetscValidPointer(gvec, 2); 196 *gvec = NULL; 197 X = shell->Xglobal; 198 PetscCheck(X, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()"); 199 /* Need to create a copy in order to attach the DM to the vector */ 200 PetscCall(VecDuplicate(X, gvec)); 201 PetscCall(VecZeroEntries(*gvec)); 202 PetscCall(VecSetDM(*gvec, dm)); 203 PetscFunctionReturn(0); 204 } 205 206 PetscErrorCode DMCreateLocalVector_Shell(DM dm, Vec *gvec) { 207 DM_Shell *shell = (DM_Shell *)dm->data; 208 Vec X; 209 210 PetscFunctionBegin; 211 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 212 PetscValidPointer(gvec, 2); 213 *gvec = NULL; 214 X = shell->Xlocal; 215 PetscCheck(X, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()"); 216 /* Need to create a copy in order to attach the DM to the vector */ 217 PetscCall(VecDuplicate(X, gvec)); 218 PetscCall(VecZeroEntries(*gvec)); 219 PetscCall(VecSetDM(*gvec, dm)); 220 PetscFunctionReturn(0); 221 } 222 223 /*@ 224 DMShellSetDestroyContext - set a function that destroys the context provided with `DMShellSetContext()` 225 226 Collective 227 228 Input Parameters: 229 . ctx - the context 230 231 Level: advanced 232 233 .seealso: `DMShellSetContext()`, `DMShellGetContext()` 234 @*/ 235 PetscErrorCode DMShellSetDestroyContext(DM dm, PetscErrorCode (*destroyctx)(void *)) { 236 DM_Shell *shell = (DM_Shell *)dm->data; 237 PetscBool isshell; 238 239 PetscFunctionBegin; 240 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 241 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 242 if (!isshell) PetscFunctionReturn(0); 243 shell->destroyctx = destroyctx; 244 PetscFunctionReturn(0); 245 } 246 247 /*@ 248 DMShellSetContext - set some data to be usable by this DM 249 250 Collective 251 252 Input Parameters: 253 + dm - shell DM 254 - ctx - the context 255 256 Level: advanced 257 258 .seealso: `DMCreateMatrix()`, `DMShellGetContext()` 259 @*/ 260 PetscErrorCode DMShellSetContext(DM dm, void *ctx) { 261 DM_Shell *shell = (DM_Shell *)dm->data; 262 PetscBool isshell; 263 264 PetscFunctionBegin; 265 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 266 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 267 if (!isshell) PetscFunctionReturn(0); 268 shell->ctx = ctx; 269 PetscFunctionReturn(0); 270 } 271 272 /*@ 273 DMShellGetContext - Returns the user-provided context associated to the DM 274 275 Collective 276 277 Input Parameter: 278 . dm - shell DM 279 280 Output Parameter: 281 . ctx - the context 282 283 Level: advanced 284 285 .seealso: `DMCreateMatrix()`, `DMShellSetContext()` 286 @*/ 287 PetscErrorCode DMShellGetContext(DM dm, void *ctx) { 288 DM_Shell *shell = (DM_Shell *)dm->data; 289 PetscBool isshell; 290 291 PetscFunctionBegin; 292 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 293 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 294 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 295 *(void **)ctx = shell->ctx; 296 PetscFunctionReturn(0); 297 } 298 299 /*@ 300 DMShellSetMatrix - sets a template matrix associated with the DMShell 301 302 Collective 303 304 Input Parameters: 305 + dm - shell DM 306 - J - template matrix 307 308 Level: advanced 309 310 Developer Notes: 311 To avoid circular references, if J is already associated to the same DM, then MatDuplicate(SHARE_NONZERO_PATTERN) is called, followed by removing the DM reference from the private template. 312 313 .seealso: `DMCreateMatrix()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()` 314 @*/ 315 PetscErrorCode DMShellSetMatrix(DM dm, Mat J) { 316 DM_Shell *shell = (DM_Shell *)dm->data; 317 PetscBool isshell; 318 DM mdm; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 322 PetscValidHeaderSpecific(J, MAT_CLASSID, 2); 323 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 324 if (!isshell) PetscFunctionReturn(0); 325 if (J == shell->A) PetscFunctionReturn(0); 326 PetscCall(MatGetDM(J, &mdm)); 327 PetscCall(PetscObjectReference((PetscObject)J)); 328 PetscCall(MatDestroy(&shell->A)); 329 if (mdm == dm) { 330 PetscCall(MatDuplicate(J, MAT_SHARE_NONZERO_PATTERN, &shell->A)); 331 PetscCall(MatSetDM(shell->A, NULL)); 332 } else shell->A = J; 333 PetscFunctionReturn(0); 334 } 335 336 /*@C 337 DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM 338 339 Logically Collective on dm 340 341 Input Parameters: 342 + dm - the shell DM 343 - func - the function to create a matrix 344 345 Level: advanced 346 347 .seealso: `DMCreateMatrix()`, `DMShellSetMatrix()`, `DMShellSetContext()`, `DMShellGetContext()` 348 @*/ 349 PetscErrorCode DMShellSetCreateMatrix(DM dm, PetscErrorCode (*func)(DM, Mat *)) { 350 PetscFunctionBegin; 351 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 352 dm->ops->creatematrix = func; 353 PetscFunctionReturn(0); 354 } 355 356 /*@ 357 DMShellSetGlobalVector - sets a template global vector associated with the DMShell 358 359 Logically Collective on dm 360 361 Input Parameters: 362 + dm - shell DM 363 - X - template vector 364 365 Level: advanced 366 367 .seealso: `DMCreateGlobalVector()`, `DMShellSetMatrix()`, `DMShellSetCreateGlobalVector()` 368 @*/ 369 PetscErrorCode DMShellSetGlobalVector(DM dm, Vec X) { 370 DM_Shell *shell = (DM_Shell *)dm->data; 371 PetscBool isshell; 372 DM vdm; 373 374 PetscFunctionBegin; 375 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 376 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 377 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 378 if (!isshell) PetscFunctionReturn(0); 379 PetscCall(VecGetDM(X, &vdm)); 380 /* 381 if the vector proposed as the new base global vector for the DM is a DM vector associated 382 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 383 we get a circular dependency that prevents the DM from being destroy when it should be. 384 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 385 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 386 to set its input vector (which is associated with the DM) as the base global vector. 387 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 388 for pointing out the problem. 389 */ 390 if (vdm == dm) PetscFunctionReturn(0); 391 PetscCall(PetscObjectReference((PetscObject)X)); 392 PetscCall(VecDestroy(&shell->Xglobal)); 393 shell->Xglobal = X; 394 PetscFunctionReturn(0); 395 } 396 397 /*@ 398 DMShellGetGlobalVector - Returns the template global vector associated with the DMShell, or NULL if it was not set 399 400 Not collective 401 402 Input Parameters: 403 + dm - shell DM 404 - X - template vector 405 406 Level: advanced 407 408 .seealso: `DMShellSetGlobalVector()`, `DMShellSetCreateGlobalVector()`, `DMCreateGlobalVector()` 409 @*/ 410 PetscErrorCode DMShellGetGlobalVector(DM dm, Vec *X) { 411 DM_Shell *shell = (DM_Shell *)dm->data; 412 PetscBool isshell; 413 414 PetscFunctionBegin; 415 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 416 PetscValidPointer(X, 2); 417 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 418 if (!isshell) PetscFunctionReturn(0); 419 *X = shell->Xglobal; 420 PetscFunctionReturn(0); 421 } 422 423 /*@C 424 DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM 425 426 Logically Collective 427 428 Input Parameters: 429 + dm - the shell DM 430 - func - the creation routine 431 432 Level: advanced 433 434 .seealso: `DMShellSetGlobalVector()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()` 435 @*/ 436 PetscErrorCode DMShellSetCreateGlobalVector(DM dm, PetscErrorCode (*func)(DM, Vec *)) { 437 PetscFunctionBegin; 438 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 439 dm->ops->createglobalvector = func; 440 PetscFunctionReturn(0); 441 } 442 443 /*@ 444 DMShellSetLocalVector - sets a template local vector associated with the DMShell 445 446 Logically Collective on dm 447 448 Input Parameters: 449 + dm - shell DM 450 - X - template vector 451 452 Level: advanced 453 454 .seealso: `DMCreateLocalVector()`, `DMShellSetMatrix()`, `DMShellSetCreateLocalVector()` 455 @*/ 456 PetscErrorCode DMShellSetLocalVector(DM dm, Vec X) { 457 DM_Shell *shell = (DM_Shell *)dm->data; 458 PetscBool isshell; 459 DM vdm; 460 461 PetscFunctionBegin; 462 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 463 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 464 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 465 if (!isshell) PetscFunctionReturn(0); 466 PetscCall(VecGetDM(X, &vdm)); 467 /* 468 if the vector proposed as the new base global vector for the DM is a DM vector associated 469 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 470 we get a circular dependency that prevents the DM from being destroy when it should be. 471 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 472 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 473 to set its input vector (which is associated with the DM) as the base global vector. 474 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 475 for pointing out the problem. 476 */ 477 if (vdm == dm) PetscFunctionReturn(0); 478 PetscCall(PetscObjectReference((PetscObject)X)); 479 PetscCall(VecDestroy(&shell->Xlocal)); 480 shell->Xlocal = X; 481 PetscFunctionReturn(0); 482 } 483 484 /*@C 485 DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM 486 487 Logically Collective 488 489 Input Parameters: 490 + dm - the shell DM 491 - func - the creation routine 492 493 Level: advanced 494 495 .seealso: `DMShellSetLocalVector()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()` 496 @*/ 497 PetscErrorCode DMShellSetCreateLocalVector(DM dm, PetscErrorCode (*func)(DM, Vec *)) { 498 PetscFunctionBegin; 499 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 500 dm->ops->createlocalvector = func; 501 PetscFunctionReturn(0); 502 } 503 504 /*@C 505 DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter 506 507 Logically Collective on dm 508 509 Input Parameters 510 + dm - the shell DM 511 . begin - the routine that begins the global to local scatter 512 - end - the routine that ends the global to local scatter 513 514 Notes: 515 If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then 516 DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers 517 518 Level: advanced 519 520 .seealso: `DMShellSetLocalToGlobal()`, `DMGlobalToLocalBeginDefaultShell()`, `DMGlobalToLocalEndDefaultShell()` 521 @*/ 522 PetscErrorCode DMShellSetGlobalToLocal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec)) { 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 525 dm->ops->globaltolocalbegin = begin; 526 dm->ops->globaltolocalend = end; 527 PetscFunctionReturn(0); 528 } 529 530 /*@C 531 DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter 532 533 Logically Collective on dm 534 535 Input Parameters 536 + dm - the shell DM 537 . begin - the routine that begins the local to global scatter 538 - end - the routine that ends the local to global scatter 539 540 Notes: 541 If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then 542 DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers 543 544 Level: advanced 545 546 .seealso: `DMShellSetGlobalToLocal()` 547 @*/ 548 PetscErrorCode DMShellSetLocalToGlobal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec)) { 549 PetscFunctionBegin; 550 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 551 dm->ops->localtoglobalbegin = begin; 552 dm->ops->localtoglobalend = end; 553 PetscFunctionReturn(0); 554 } 555 556 /*@C 557 DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter 558 559 Logically Collective on dm 560 561 Input Parameters 562 + dm - the shell DM 563 . begin - the routine that begins the local to local scatter 564 - end - the routine that ends the local to local scatter 565 566 Notes: 567 If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then 568 DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers 569 570 Level: advanced 571 572 .seealso: `DMShellSetGlobalToLocal()`, `DMLocalToLocalBeginDefaultShell()`, `DMLocalToLocalEndDefaultShell()` 573 @*/ 574 PetscErrorCode DMShellSetLocalToLocal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec)) { 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 577 dm->ops->localtolocalbegin = begin; 578 dm->ops->localtolocalend = end; 579 PetscFunctionReturn(0); 580 } 581 582 /*@ 583 DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication 584 585 Logically Collective on dm 586 587 Input Parameters 588 + dm - the shell DM 589 - gtol - the global to local VecScatter context 590 591 Level: advanced 592 593 .seealso: `DMShellSetGlobalToLocal()`, `DMGlobalToLocalBeginDefaultShell()`, `DMGlobalToLocalEndDefaultShell()` 594 @*/ 595 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol) { 596 DM_Shell *shell = (DM_Shell *)dm->data; 597 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 600 PetscValidHeaderSpecific(gtol, PETSCSF_CLASSID, 2); 601 PetscCall(PetscObjectReference((PetscObject)gtol)); 602 PetscCall(VecScatterDestroy(&shell->gtol)); 603 shell->gtol = gtol; 604 PetscFunctionReturn(0); 605 } 606 607 /*@ 608 DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication 609 610 Logically Collective on dm 611 612 Input Parameters 613 + dm - the shell DM 614 - ltog - the local to global VecScatter context 615 616 Level: advanced 617 618 .seealso: `DMShellSetLocalToGlobal()`, `DMLocalToGlobalBeginDefaultShell()`, `DMLocalToGlobalEndDefaultShell()` 619 @*/ 620 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog) { 621 DM_Shell *shell = (DM_Shell *)dm->data; 622 623 PetscFunctionBegin; 624 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 625 PetscValidHeaderSpecific(ltog, PETSCSF_CLASSID, 2); 626 PetscCall(PetscObjectReference((PetscObject)ltog)); 627 PetscCall(VecScatterDestroy(&shell->ltog)); 628 shell->ltog = ltog; 629 PetscFunctionReturn(0); 630 } 631 632 /*@ 633 DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication 634 635 Logically Collective on dm 636 637 Input Parameters 638 + dm - the shell DM 639 - ltol - the local to local VecScatter context 640 641 Level: advanced 642 643 .seealso: `DMShellSetLocalToLocal()`, `DMLocalToLocalBeginDefaultShell()`, `DMLocalToLocalEndDefaultShell()` 644 @*/ 645 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol) { 646 DM_Shell *shell = (DM_Shell *)dm->data; 647 648 PetscFunctionBegin; 649 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 650 PetscValidHeaderSpecific(ltol, PETSCSF_CLASSID, 2); 651 PetscCall(PetscObjectReference((PetscObject)ltol)); 652 PetscCall(VecScatterDestroy(&shell->ltol)); 653 shell->ltol = ltol; 654 PetscFunctionReturn(0); 655 } 656 657 /*@C 658 DMShellSetCoarsen - Set the routine used to coarsen the shell DM 659 660 Logically Collective on dm 661 662 Input Parameters 663 + dm - the shell DM 664 - coarsen - the routine that coarsens the DM 665 666 Level: advanced 667 668 .seealso: `DMShellSetRefine()`, `DMCoarsen()`, `DMShellGetCoarsen()`, `DMShellSetContext()`, `DMShellGetContext()` 669 @*/ 670 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *)) { 671 PetscBool isshell; 672 673 PetscFunctionBegin; 674 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 675 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 676 if (!isshell) PetscFunctionReturn(0); 677 dm->ops->coarsen = coarsen; 678 PetscFunctionReturn(0); 679 } 680 681 /*@C 682 DMShellGetCoarsen - Get the routine used to coarsen the shell DM 683 684 Logically Collective on dm 685 686 Input Parameter: 687 . dm - the shell DM 688 689 Output Parameter: 690 . coarsen - the routine that coarsens the DM 691 692 Level: advanced 693 694 .seealso: `DMShellSetCoarsen()`, `DMCoarsen()`, `DMShellSetRefine()`, `DMRefine()` 695 @*/ 696 PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM, MPI_Comm, DM *)) { 697 PetscBool isshell; 698 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 701 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 702 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 703 *coarsen = dm->ops->coarsen; 704 PetscFunctionReturn(0); 705 } 706 707 /*@C 708 DMShellSetRefine - Set the routine used to refine the shell DM 709 710 Logically Collective on dm 711 712 Input Parameters 713 + dm - the shell DM 714 - refine - the routine that refines the DM 715 716 Level: advanced 717 718 .seealso: `DMShellSetCoarsen()`, `DMRefine()`, `DMShellGetRefine()`, `DMShellSetContext()`, `DMShellGetContext()` 719 @*/ 720 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM, MPI_Comm, DM *)) { 721 PetscBool isshell; 722 723 PetscFunctionBegin; 724 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 725 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 726 if (!isshell) PetscFunctionReturn(0); 727 dm->ops->refine = refine; 728 PetscFunctionReturn(0); 729 } 730 731 /*@C 732 DMShellGetRefine - Get the routine used to refine the shell DM 733 734 Logically Collective on dm 735 736 Input Parameter: 737 . dm - the shell DM 738 739 Output Parameter: 740 . refine - the routine that refines the DM 741 742 Level: advanced 743 744 .seealso: `DMShellSetCoarsen()`, `DMCoarsen()`, `DMShellSetRefine()`, `DMRefine()` 745 @*/ 746 PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM, MPI_Comm, DM *)) { 747 PetscBool isshell; 748 749 PetscFunctionBegin; 750 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 751 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 752 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 753 *refine = dm->ops->refine; 754 PetscFunctionReturn(0); 755 } 756 757 /*@C 758 DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator 759 760 Logically Collective on dm 761 762 Input Parameters 763 + dm - the shell DM 764 - interp - the routine to create the interpolation 765 766 Level: advanced 767 768 .seealso: `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateInterpolation()`, `DMShellSetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()` 769 @*/ 770 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM, DM, Mat *, Vec *)) { 771 PetscBool isshell; 772 773 PetscFunctionBegin; 774 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 775 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 776 if (!isshell) PetscFunctionReturn(0); 777 dm->ops->createinterpolation = interp; 778 PetscFunctionReturn(0); 779 } 780 781 /*@C 782 DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator 783 784 Logically Collective on dm 785 786 Input Parameter: 787 . dm - the shell DM 788 789 Output Parameter: 790 . interp - the routine to create the interpolation 791 792 Level: advanced 793 794 .seealso: `DMShellGetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()` 795 @*/ 796 PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM, DM, Mat *, Vec *)) { 797 PetscBool isshell; 798 799 PetscFunctionBegin; 800 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 801 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 802 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 803 *interp = dm->ops->createinterpolation; 804 PetscFunctionReturn(0); 805 } 806 807 /*@C 808 DMShellSetCreateRestriction - Set the routine used to create the restriction operator 809 810 Logically Collective on dm 811 812 Input Parameters 813 + dm - the shell DM 814 - striction- the routine to create the restriction 815 816 Level: advanced 817 818 .seealso: `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()` 819 @*/ 820 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM, DM, Mat *)) { 821 PetscBool isshell; 822 823 PetscFunctionBegin; 824 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 825 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 826 if (!isshell) PetscFunctionReturn(0); 827 dm->ops->createrestriction = restriction; 828 PetscFunctionReturn(0); 829 } 830 831 /*@C 832 DMShellGetCreateRestriction - Get the routine used to create the restriction operator 833 834 Logically Collective on dm 835 836 Input Parameter: 837 . dm - the shell DM 838 839 Output Parameter: 840 . restriction - the routine to create the restriction 841 842 Level: advanced 843 844 .seealso: `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellSetContext()`, `DMShellGetContext()` 845 @*/ 846 PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM, DM, Mat *)) { 847 PetscBool isshell; 848 849 PetscFunctionBegin; 850 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 851 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 852 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 853 *restriction = dm->ops->createrestriction; 854 PetscFunctionReturn(0); 855 } 856 857 /*@C 858 DMShellSetCreateInjection - Set the routine used to create the injection operator 859 860 Logically Collective on dm 861 862 Input Parameters: 863 + dm - the shell DM 864 - inject - the routine to create the injection 865 866 Level: advanced 867 868 .seealso: `DMShellSetCreateInterpolation()`, `DMCreateInjection()`, `DMShellGetCreateInjection()`, `DMShellSetContext()`, `DMShellGetContext()` 869 @*/ 870 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM, DM, Mat *)) { 871 PetscBool isshell; 872 873 PetscFunctionBegin; 874 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 875 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 876 if (!isshell) PetscFunctionReturn(0); 877 dm->ops->createinjection = inject; 878 PetscFunctionReturn(0); 879 } 880 881 /*@C 882 DMShellGetCreateInjection - Get the routine used to create the injection operator 883 884 Logically Collective on dm 885 886 Input Parameter: 887 . dm - the shell DM 888 889 Output Parameter: 890 . inject - the routine to create the injection 891 892 Level: advanced 893 894 .seealso: `DMShellGetCreateInterpolation()`, `DMCreateInjection()`, `DMShellSetContext()`, `DMShellGetContext()` 895 @*/ 896 PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM, DM, Mat *)) { 897 PetscBool isshell; 898 899 PetscFunctionBegin; 900 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 901 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 902 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 903 *inject = dm->ops->createinjection; 904 PetscFunctionReturn(0); 905 } 906 907 /*@C 908 DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM 909 910 Logically Collective on dm 911 912 Input Parameters: 913 + dm - the shell DM 914 - decomp - the routine to create the decomposition 915 916 Level: advanced 917 918 .seealso: `DMCreateFieldDecomposition()`, `DMShellSetContext()`, `DMShellGetContext()` 919 @*/ 920 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM, PetscInt *, char ***, IS **, DM **)) { 921 PetscBool isshell; 922 923 PetscFunctionBegin; 924 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 925 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 926 if (!isshell) PetscFunctionReturn(0); 927 dm->ops->createfielddecomposition = decomp; 928 PetscFunctionReturn(0); 929 } 930 931 /*@C 932 DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM 933 934 Logically Collective on dm 935 936 Input Parameters: 937 + dm - the shell DM 938 - decomp - the routine to create the decomposition 939 940 Level: advanced 941 942 .seealso: `DMCreateDomainDecomposition()`, `DMShellSetContext()`, `DMShellGetContext()` 943 @*/ 944 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM, PetscInt *, char ***, IS **, IS **, DM **)) { 945 PetscBool isshell; 946 947 PetscFunctionBegin; 948 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 949 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 950 if (!isshell) PetscFunctionReturn(0); 951 dm->ops->createdomaindecomposition = decomp; 952 PetscFunctionReturn(0); 953 } 954 955 /*@C 956 DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM 957 958 Logically Collective on dm 959 960 Input Parameters: 961 + dm - the shell DM 962 - scatter - the routine to create the scatters 963 964 Level: advanced 965 966 .seealso: `DMCreateDomainDecompositionScatters()`, `DMShellSetContext()`, `DMShellGetContext()` 967 @*/ 968 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **)) { 969 PetscBool isshell; 970 971 PetscFunctionBegin; 972 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 973 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 974 if (!isshell) PetscFunctionReturn(0); 975 dm->ops->createddscatters = scatter; 976 PetscFunctionReturn(0); 977 } 978 979 /*@C 980 DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM 981 982 Logically Collective on dm 983 984 Input Parameters: 985 + dm - the shell DM 986 - subdm - the routine to create the decomposition 987 988 Level: advanced 989 990 .seealso: `DMCreateSubDM()`, `DMShellGetCreateSubDM()`, `DMShellSetContext()`, `DMShellGetContext()` 991 @*/ 992 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM, PetscInt, const PetscInt[], IS *, DM *)) { 993 PetscBool isshell; 994 995 PetscFunctionBegin; 996 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 997 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 998 if (!isshell) PetscFunctionReturn(0); 999 dm->ops->createsubdm = subdm; 1000 PetscFunctionReturn(0); 1001 } 1002 1003 /*@C 1004 DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM 1005 1006 Logically Collective on dm 1007 1008 Input Parameter: 1009 . dm - the shell DM 1010 1011 Output Parameter: 1012 . subdm - the routine to create the decomposition 1013 1014 Level: advanced 1015 1016 .seealso: `DMCreateSubDM()`, `DMShellSetCreateSubDM()`, `DMShellSetContext()`, `DMShellGetContext()` 1017 @*/ 1018 PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM, PetscInt, const PetscInt[], IS *, DM *)) { 1019 PetscBool isshell; 1020 1021 PetscFunctionBegin; 1022 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1023 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 1024 PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs"); 1025 *subdm = dm->ops->createsubdm; 1026 PetscFunctionReturn(0); 1027 } 1028 1029 static PetscErrorCode DMDestroy_Shell(DM dm) { 1030 DM_Shell *shell = (DM_Shell *)dm->data; 1031 1032 PetscFunctionBegin; 1033 if (shell->destroyctx) PetscCallBack("Destroy Context", (*shell->destroyctx)(shell->ctx)); 1034 PetscCall(MatDestroy(&shell->A)); 1035 PetscCall(VecDestroy(&shell->Xglobal)); 1036 PetscCall(VecDestroy(&shell->Xlocal)); 1037 PetscCall(VecScatterDestroy(&shell->gtol)); 1038 PetscCall(VecScatterDestroy(&shell->ltog)); 1039 PetscCall(VecScatterDestroy(&shell->ltol)); 1040 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 1041 PetscCall(PetscFree(shell)); 1042 PetscFunctionReturn(0); 1043 } 1044 1045 static PetscErrorCode DMView_Shell(DM dm, PetscViewer v) { 1046 DM_Shell *shell = (DM_Shell *)dm->data; 1047 1048 PetscFunctionBegin; 1049 if (shell->Xglobal) PetscCall(VecView(shell->Xglobal, v)); 1050 PetscFunctionReturn(0); 1051 } 1052 1053 static PetscErrorCode DMLoad_Shell(DM dm, PetscViewer v) { 1054 DM_Shell *shell = (DM_Shell *)dm->data; 1055 1056 PetscFunctionBegin; 1057 PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &shell->Xglobal)); 1058 PetscCall(VecLoad(shell->Xglobal, v)); 1059 PetscFunctionReturn(0); 1060 } 1061 1062 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) { 1063 PetscFunctionBegin; 1064 if (subdm) PetscCall(DMShellCreate(PetscObjectComm((PetscObject)dm), subdm)); 1065 PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm)); 1066 PetscFunctionReturn(0); 1067 } 1068 1069 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm) { 1070 DM_Shell *shell; 1071 1072 PetscFunctionBegin; 1073 PetscCall(PetscNew(&shell)); 1074 dm->data = shell; 1075 1076 dm->ops->destroy = DMDestroy_Shell; 1077 dm->ops->createglobalvector = DMCreateGlobalVector_Shell; 1078 dm->ops->createlocalvector = DMCreateLocalVector_Shell; 1079 dm->ops->creatematrix = DMCreateMatrix_Shell; 1080 dm->ops->view = DMView_Shell; 1081 dm->ops->load = DMLoad_Shell; 1082 dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell; 1083 dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell; 1084 dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell; 1085 dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell; 1086 dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell; 1087 dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell; 1088 dm->ops->createsubdm = DMCreateSubDM_Shell; 1089 PetscCall(DMSetMatType(dm, MATDENSE)); 1090 PetscFunctionReturn(0); 1091 } 1092 1093 /*@ 1094 DMShellCreate - Creates a shell DM object, used to manage user-defined problem data 1095 1096 Collective 1097 1098 Input Parameter: 1099 . comm - the processors that will share the global vector 1100 1101 Output Parameters: 1102 . shell - the shell DM 1103 1104 Level: advanced 1105 1106 .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMShellSetContext()`, `DMShellGetContext()` 1107 @*/ 1108 PetscErrorCode DMShellCreate(MPI_Comm comm, DM *dm) { 1109 PetscFunctionBegin; 1110 PetscValidPointer(dm, 2); 1111 PetscCall(DMCreate(comm, dm)); 1112 PetscCall(DMSetType(*dm, DMSHELL)); 1113 PetscCall(DMSetUp(*dm)); 1114 PetscFunctionReturn(0); 1115 } 1116