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