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