1 /* 2 This provides a simple shell for Fortran (and C programmers) to 3 create their own preconditioner without writing much interface code. 4 */ 5 6 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 7 8 typedef struct { 9 void *ctx; /* user provided contexts for preconditioner */ 10 11 PetscErrorCode (*destroy)(PC); 12 PetscErrorCode (*setup)(PC); 13 PetscErrorCode (*apply)(PC, Vec, Vec); 14 PetscErrorCode (*matapply)(PC, Mat, Mat); 15 PetscErrorCode (*applysymmetricleft)(PC, Vec, Vec); 16 PetscErrorCode (*applysymmetricright)(PC, Vec, Vec); 17 PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec); 18 PetscErrorCode (*presolve)(PC, KSP, Vec, Vec); 19 PetscErrorCode (*postsolve)(PC, KSP, Vec, Vec); 20 PetscErrorCode (*view)(PC, PetscViewer); 21 PetscErrorCode (*applytranspose)(PC, Vec, Vec); 22 PetscErrorCode (*matapplytranspose)(PC, Mat, Mat); 23 PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *); 24 25 char *name; 26 } PC_Shell; 27 28 /*@C 29 PCShellGetContext - Returns the user-provided context associated with a shell `PC` that was provided with `PCShellSetContext()` 30 31 Not Collective 32 33 Input Parameter: 34 . pc - of type `PCSHELL` 35 36 Output Parameter: 37 . ctx - the user provided context 38 39 Level: advanced 40 41 Note: 42 This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()` 43 44 Fortran Note: 45 To use this from Fortran you must write a Fortran interface definition for this 46 function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument. 47 48 .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellSetDestroy()` 49 @*/ 50 PetscErrorCode PCShellGetContext(PC pc, void *ctx) 51 { 52 PetscBool flg; 53 54 PetscFunctionBegin; 55 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 56 PetscAssertPointer(ctx, 2); 57 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg)); 58 if (!flg) *(void **)ctx = NULL; 59 else *(void **)ctx = ((PC_Shell *)pc->data)->ctx; 60 PetscFunctionReturn(PETSC_SUCCESS); 61 } 62 63 /*@ 64 PCShellSetContext - sets the context for a shell `PC` that can be accessed with `PCShellGetContext()` 65 66 Logically Collective 67 68 Input Parameters: 69 + pc - the `PC` of type `PCSHELL` 70 - ctx - the context 71 72 Level: advanced 73 74 Notes: 75 This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()` 76 77 One should also provide a routine to destroy the context when `pc` is destroyed with `PCShellSetDestroy()` 78 79 Fortran Notes: 80 To use this from Fortran you must write a Fortran interface definition for this 81 function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument. 82 83 .seealso: [](ch_ksp), `PC`, `PCShellGetContext()`, `PCSHELL`, `PCShellSetApply()`, `PCShellSetDestroy()` 84 @*/ 85 PetscErrorCode PCShellSetContext(PC pc, void *ctx) 86 { 87 PC_Shell *shell = (PC_Shell *)pc->data; 88 PetscBool flg; 89 90 PetscFunctionBegin; 91 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 92 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg)); 93 if (flg) shell->ctx = ctx; 94 PetscFunctionReturn(PETSC_SUCCESS); 95 } 96 97 static PetscErrorCode PCSetUp_Shell(PC pc) 98 { 99 PC_Shell *shell = (PC_Shell *)pc->data; 100 101 PetscFunctionBegin; 102 PetscCheck(shell->setup, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No setup() routine provided to Shell PC"); 103 PetscCallBack("PCSHELL callback setup", (*shell->setup)(pc)); 104 PetscFunctionReturn(PETSC_SUCCESS); 105 } 106 107 static PetscErrorCode PCApply_Shell(PC pc, Vec x, Vec y) 108 { 109 PC_Shell *shell = (PC_Shell *)pc->data; 110 PetscObjectState instate, outstate; 111 112 PetscFunctionBegin; 113 PetscCheck(shell->apply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 114 PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 115 PetscCallBack("PCSHELL callback apply", (*shell->apply)(pc, x, y)); 116 PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 117 /* increase the state of the output vector if the user did not update its state themself as should have been done */ 118 if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 119 PetscFunctionReturn(PETSC_SUCCESS); 120 } 121 122 static PetscErrorCode PCMatApply_Shell(PC pc, Mat X, Mat Y) 123 { 124 PC_Shell *shell = (PC_Shell *)pc->data; 125 PetscObjectState instate, outstate; 126 127 PetscFunctionBegin; 128 PetscCheck(shell->matapply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 129 PetscCall(PetscObjectStateGet((PetscObject)Y, &instate)); 130 PetscCallBack("PCSHELL callback apply", (*shell->matapply)(pc, X, Y)); 131 PetscCall(PetscObjectStateGet((PetscObject)Y, &outstate)); 132 /* increase the state of the output vector if the user did not update its state themself as should have been done */ 133 if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 134 PetscFunctionReturn(PETSC_SUCCESS); 135 } 136 137 static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc, Vec x, Vec y) 138 { 139 PC_Shell *shell = (PC_Shell *)pc->data; 140 141 PetscFunctionBegin; 142 PetscCheck(shell->applysymmetricleft, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 143 PetscCallBack("PCSHELL callback apply symmetric left", (*shell->applysymmetricleft)(pc, x, y)); 144 PetscFunctionReturn(PETSC_SUCCESS); 145 } 146 147 static PetscErrorCode PCApplySymmetricRight_Shell(PC pc, Vec x, Vec y) 148 { 149 PC_Shell *shell = (PC_Shell *)pc->data; 150 151 PetscFunctionBegin; 152 PetscCheck(shell->applysymmetricright, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 153 PetscCallBack("PCSHELL callback apply symmetric right", (*shell->applysymmetricright)(pc, x, y)); 154 PetscFunctionReturn(PETSC_SUCCESS); 155 } 156 157 static PetscErrorCode PCApplyBA_Shell(PC pc, PCSide side, Vec x, Vec y, Vec w) 158 { 159 PC_Shell *shell = (PC_Shell *)pc->data; 160 PetscObjectState instate, outstate; 161 162 PetscFunctionBegin; 163 PetscCheck(shell->applyBA, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyBA() routine provided to Shell PC"); 164 PetscCall(PetscObjectStateGet((PetscObject)w, &instate)); 165 PetscCallBack("PCSHELL callback applyBA", (*shell->applyBA)(pc, side, x, y, w)); 166 PetscCall(PetscObjectStateGet((PetscObject)w, &outstate)); 167 /* increase the state of the output vector if the user did not update its state themself as should have been done */ 168 if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)w)); 169 PetscFunctionReturn(PETSC_SUCCESS); 170 } 171 172 static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc, PetscBool *change) 173 { 174 PetscFunctionBegin; 175 *change = PETSC_TRUE; 176 PetscFunctionReturn(PETSC_SUCCESS); 177 } 178 179 static PetscErrorCode PCPreSolve_Shell(PC pc, KSP ksp, Vec b, Vec x) 180 { 181 PC_Shell *shell = (PC_Shell *)pc->data; 182 183 PetscFunctionBegin; 184 PetscCheck(shell->presolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No presolve() routine provided to Shell PC"); 185 PetscCallBack("PCSHELL callback presolve", (*shell->presolve)(pc, ksp, b, x)); 186 PetscFunctionReturn(PETSC_SUCCESS); 187 } 188 189 static PetscErrorCode PCPostSolve_Shell(PC pc, KSP ksp, Vec b, Vec x) 190 { 191 PC_Shell *shell = (PC_Shell *)pc->data; 192 193 PetscFunctionBegin; 194 PetscCheck(shell->postsolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No postsolve() routine provided to Shell PC"); 195 PetscCallBack("PCSHELL callback postsolve()", (*shell->postsolve)(pc, ksp, b, x)); 196 PetscFunctionReturn(PETSC_SUCCESS); 197 } 198 199 static PetscErrorCode PCApplyTranspose_Shell(PC pc, Vec x, Vec y) 200 { 201 PC_Shell *shell = (PC_Shell *)pc->data; 202 PetscObjectState instate, outstate; 203 204 PetscFunctionBegin; 205 PetscCheck(shell->applytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applytranspose() routine provided to Shell PC"); 206 PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 207 PetscCallBack("PCSHELL callback applytranspose", (*shell->applytranspose)(pc, x, y)); 208 PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 209 /* increase the state of the output vector if the user did not update its state themself as should have been done */ 210 if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 211 PetscFunctionReturn(PETSC_SUCCESS); 212 } 213 214 static PetscErrorCode PCMatApplyTranspose_Shell(PC pc, Mat x, Mat y) 215 { 216 PC_Shell *shell = (PC_Shell *)pc->data; 217 PetscObjectState instate, outstate; 218 219 PetscFunctionBegin; 220 PetscCheck(shell->matapplytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No matapplytranspose() routine provided to Shell PC"); 221 PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 222 PetscCallBack("PCSHELL callback matapplytranspose", (*shell->matapplytranspose)(pc, x, y)); 223 PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 224 /* increase the state of the output matrix if the user did not update its state themself as should have been done */ 225 if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 226 PetscFunctionReturn(PETSC_SUCCESS); 227 } 228 229 static PetscErrorCode PCApplyRichardson_Shell(PC pc, Vec x, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt it, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason) 230 { 231 PC_Shell *shell = (PC_Shell *)pc->data; 232 PetscObjectState instate, outstate; 233 234 PetscFunctionBegin; 235 PetscCheck(shell->applyrich, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyrichardson() routine provided to Shell PC"); 236 PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 237 PetscCallBack("PCSHELL callback applyrichardson", (*shell->applyrich)(pc, x, y, w, rtol, abstol, dtol, it, guesszero, outits, reason)); 238 PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 239 /* increase the state of the output vector since the user did not update its state themself as should have been done */ 240 if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 241 PetscFunctionReturn(PETSC_SUCCESS); 242 } 243 244 static PetscErrorCode PCDestroy_Shell(PC pc) 245 { 246 PC_Shell *shell = (PC_Shell *)pc->data; 247 248 PetscFunctionBegin; 249 PetscCall(PetscFree(shell->name)); 250 if (shell->destroy) PetscCallBack("PCSHELL callback destroy", (*shell->destroy)(pc)); 251 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", NULL)); 252 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", NULL)); 253 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", NULL)); 254 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", NULL)); 255 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", NULL)); 256 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", NULL)); 257 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", NULL)); 258 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", NULL)); 259 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", NULL)); 260 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", NULL)); 261 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", NULL)); 262 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApplyTranspose_C", NULL)); 263 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", NULL)); 264 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", NULL)); 265 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", NULL)); 266 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 267 PetscCall(PetscFree(pc->data)); 268 PetscFunctionReturn(PETSC_SUCCESS); 269 } 270 271 static PetscErrorCode PCView_Shell(PC pc, PetscViewer viewer) 272 { 273 PC_Shell *shell = (PC_Shell *)pc->data; 274 PetscBool isascii; 275 276 PetscFunctionBegin; 277 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 278 if (isascii) { 279 if (shell->name) PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", shell->name)); 280 else PetscCall(PetscViewerASCIIPrintf(viewer, " no name\n")); 281 } 282 if (shell->view) { 283 PetscCall(PetscViewerASCIIPushTab(viewer)); 284 PetscCall((*shell->view)(pc, viewer)); 285 PetscCall(PetscViewerASCIIPopTab(viewer)); 286 } 287 PetscFunctionReturn(PETSC_SUCCESS); 288 } 289 290 static PetscErrorCode PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC)) 291 { 292 PC_Shell *shell = (PC_Shell *)pc->data; 293 294 PetscFunctionBegin; 295 shell->destroy = destroy; 296 PetscFunctionReturn(PETSC_SUCCESS); 297 } 298 299 static PetscErrorCode PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC)) 300 { 301 PC_Shell *shell = (PC_Shell *)pc->data; 302 303 PetscFunctionBegin; 304 shell->setup = setup; 305 if (setup) pc->ops->setup = PCSetUp_Shell; 306 else pc->ops->setup = NULL; 307 PetscFunctionReturn(PETSC_SUCCESS); 308 } 309 310 static PetscErrorCode PCShellSetApply_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 311 { 312 PC_Shell *shell = (PC_Shell *)pc->data; 313 314 PetscFunctionBegin; 315 shell->apply = apply; 316 PetscFunctionReturn(PETSC_SUCCESS); 317 } 318 319 static PetscErrorCode PCShellSetMatApply_Shell(PC pc, PetscErrorCode (*matapply)(PC, Mat, Mat)) 320 { 321 PC_Shell *shell = (PC_Shell *)pc->data; 322 323 PetscFunctionBegin; 324 shell->matapply = matapply; 325 if (matapply) pc->ops->matapply = PCMatApply_Shell; 326 else pc->ops->matapply = NULL; 327 PetscFunctionReturn(PETSC_SUCCESS); 328 } 329 330 static PetscErrorCode PCShellSetApplySymmetricLeft_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 331 { 332 PC_Shell *shell = (PC_Shell *)pc->data; 333 334 PetscFunctionBegin; 335 shell->applysymmetricleft = apply; 336 PetscFunctionReturn(PETSC_SUCCESS); 337 } 338 339 static PetscErrorCode PCShellSetApplySymmetricRight_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 340 { 341 PC_Shell *shell = (PC_Shell *)pc->data; 342 343 PetscFunctionBegin; 344 shell->applysymmetricright = apply; 345 PetscFunctionReturn(PETSC_SUCCESS); 346 } 347 348 static PetscErrorCode PCShellSetApplyBA_Shell(PC pc, PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec)) 349 { 350 PC_Shell *shell = (PC_Shell *)pc->data; 351 352 PetscFunctionBegin; 353 shell->applyBA = applyBA; 354 if (applyBA) pc->ops->applyBA = PCApplyBA_Shell; 355 else pc->ops->applyBA = NULL; 356 PetscFunctionReturn(PETSC_SUCCESS); 357 } 358 359 static PetscErrorCode PCShellSetPreSolve_Shell(PC pc, PCShellPSolveFn *presolve) 360 { 361 PC_Shell *shell = (PC_Shell *)pc->data; 362 363 PetscFunctionBegin; 364 shell->presolve = presolve; 365 if (presolve) { 366 pc->ops->presolve = PCPreSolve_Shell; 367 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Shell)); 368 } else { 369 pc->ops->presolve = NULL; 370 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 371 } 372 PetscFunctionReturn(PETSC_SUCCESS); 373 } 374 375 static PetscErrorCode PCShellSetPostSolve_Shell(PC pc, PCShellPSolveFn *postsolve) 376 { 377 PC_Shell *shell = (PC_Shell *)pc->data; 378 379 PetscFunctionBegin; 380 shell->postsolve = postsolve; 381 if (postsolve) pc->ops->postsolve = PCPostSolve_Shell; 382 else pc->ops->postsolve = NULL; 383 PetscFunctionReturn(PETSC_SUCCESS); 384 } 385 386 static PetscErrorCode PCShellSetView_Shell(PC pc, PetscErrorCode (*view)(PC, PetscViewer)) 387 { 388 PC_Shell *shell = (PC_Shell *)pc->data; 389 390 PetscFunctionBegin; 391 shell->view = view; 392 PetscFunctionReturn(PETSC_SUCCESS); 393 } 394 395 static PetscErrorCode PCShellSetApplyTranspose_Shell(PC pc, PetscErrorCode (*applytranspose)(PC, Vec, Vec)) 396 { 397 PC_Shell *shell = (PC_Shell *)pc->data; 398 399 PetscFunctionBegin; 400 shell->applytranspose = applytranspose; 401 if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell; 402 else pc->ops->applytranspose = NULL; 403 PetscFunctionReturn(PETSC_SUCCESS); 404 } 405 406 static PetscErrorCode PCShellSetMatApplyTranspose_Shell(PC pc, PetscErrorCode (*matapplytranspose)(PC, Mat, Mat)) 407 { 408 PC_Shell *shell = (PC_Shell *)pc->data; 409 410 PetscFunctionBegin; 411 shell->matapplytranspose = matapplytranspose; 412 if (matapplytranspose) pc->ops->matapplytranspose = PCMatApplyTranspose_Shell; 413 else pc->ops->matapplytranspose = NULL; 414 PetscFunctionReturn(PETSC_SUCCESS); 415 } 416 417 static PetscErrorCode PCShellSetApplyRichardson_Shell(PC pc, PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)) 418 { 419 PC_Shell *shell = (PC_Shell *)pc->data; 420 421 PetscFunctionBegin; 422 shell->applyrich = applyrich; 423 if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell; 424 else pc->ops->applyrichardson = NULL; 425 PetscFunctionReturn(PETSC_SUCCESS); 426 } 427 428 static PetscErrorCode PCShellSetName_Shell(PC pc, const char name[]) 429 { 430 PC_Shell *shell = (PC_Shell *)pc->data; 431 432 PetscFunctionBegin; 433 PetscCall(PetscFree(shell->name)); 434 PetscCall(PetscStrallocpy(name, &shell->name)); 435 PetscFunctionReturn(PETSC_SUCCESS); 436 } 437 438 static PetscErrorCode PCShellGetName_Shell(PC pc, const char *name[]) 439 { 440 PC_Shell *shell = (PC_Shell *)pc->data; 441 442 PetscFunctionBegin; 443 *name = shell->name; 444 PetscFunctionReturn(PETSC_SUCCESS); 445 } 446 447 /*@C 448 PCShellSetDestroy - Sets routine to use to destroy the user-provided application context that was provided with `PCShellSetContext()` 449 450 Logically Collective 451 452 Input Parameters: 453 + pc - the preconditioner context 454 - destroy - the application-provided destroy routine 455 456 Calling sequence of `destroy`: 457 . pc - the preconditioner 458 459 Level: intermediate 460 461 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()` 462 @*/ 463 PetscErrorCode PCShellSetDestroy(PC pc, PetscErrorCode (*destroy)(PC pc)) 464 { 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 467 PetscTryMethod(pc, "PCShellSetDestroy_C", (PC, PetscErrorCode (*)(PC)), (pc, destroy)); 468 PetscFunctionReturn(PETSC_SUCCESS); 469 } 470 471 /*@C 472 PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the 473 matrix operator is changed. 474 475 Logically Collective 476 477 Input Parameters: 478 + pc - the preconditioner context 479 - setup - the application-provided setup routine 480 481 Calling sequence of `setup`: 482 . pc - the preconditioner 483 484 Level: intermediate 485 486 Note: 487 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `setup`. 488 489 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetApply()`, `PCShellSetContext()`, , `PCShellGetContext()` 490 @*/ 491 PetscErrorCode PCShellSetSetUp(PC pc, PetscErrorCode (*setup)(PC pc)) 492 { 493 PetscFunctionBegin; 494 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 495 PetscTryMethod(pc, "PCShellSetSetUp_C", (PC, PetscErrorCode (*)(PC)), (pc, setup)); 496 PetscFunctionReturn(PETSC_SUCCESS); 497 } 498 499 /*@C 500 PCShellSetView - Sets routine to use as viewer of a `PCSHELL` shell preconditioner 501 502 Logically Collective 503 504 Input Parameters: 505 + pc - the preconditioner context 506 - view - the application-provided view routine 507 508 Calling sequence of `view`: 509 + pc - the preconditioner 510 - v - viewer 511 512 Level: advanced 513 514 Note: 515 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `view`. 516 517 .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()` 518 @*/ 519 PetscErrorCode PCShellSetView(PC pc, PetscErrorCode (*view)(PC pc, PetscViewer v)) 520 { 521 PetscFunctionBegin; 522 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 523 PetscTryMethod(pc, "PCShellSetView_C", (PC, PetscErrorCode (*)(PC, PetscViewer)), (pc, view)); 524 PetscFunctionReturn(PETSC_SUCCESS); 525 } 526 527 /*@C 528 PCShellSetApply - Sets routine to use as preconditioner. 529 530 Logically Collective 531 532 Input Parameters: 533 + pc - the preconditioner context 534 - apply - the application-provided preconditioning routine 535 536 Calling sequence of `apply`: 537 + pc - the preconditioner, get the application context with `PCShellGetContext()` 538 . xin - input vector 539 - xout - output vector 540 541 Level: intermediate 542 543 Note: 544 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 545 546 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellSetApplySymmetricRight()`, `PCShellSetApplySymmetricLeft()`, `PCShellGetContext()` 547 @*/ 548 PetscErrorCode PCShellSetApply(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 549 { 550 PetscFunctionBegin; 551 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 552 PetscTryMethod(pc, "PCShellSetApply_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply)); 553 PetscFunctionReturn(PETSC_SUCCESS); 554 } 555 556 /*@C 557 PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors. 558 559 Logically Collective 560 561 Input Parameters: 562 + pc - the preconditioner context 563 - matapply - the application-provided preconditioning routine 564 565 Calling sequence of `matapply`: 566 + pc - the preconditioner 567 . Xin - input block of vectors represented as a dense `Mat` 568 - Xout - output block of vectors represented as a dense `Mat` 569 570 Level: advanced 571 572 Note: 573 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapply`. 574 575 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()` 576 @*/ 577 PetscErrorCode PCShellSetMatApply(PC pc, PetscErrorCode (*matapply)(PC pc, Mat Xin, Mat Xout)) 578 { 579 PetscFunctionBegin; 580 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 581 PetscTryMethod(pc, "PCShellSetMatApply_C", (PC, PetscErrorCode (*)(PC, Mat, Mat)), (pc, matapply)); 582 PetscFunctionReturn(PETSC_SUCCESS); 583 } 584 585 /*@C 586 PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the `PC_SYMMETRIC` is used). 587 588 Logically Collective 589 590 Input Parameters: 591 + pc - the preconditioner context 592 - apply - the application-provided left preconditioning routine 593 594 Calling sequence of `apply`: 595 + pc - the preconditioner 596 . xin - input vector 597 - xout - output vector 598 599 Level: advanced 600 601 Note: 602 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 603 604 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()` 605 @*/ 606 PetscErrorCode PCShellSetApplySymmetricLeft(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 607 { 608 PetscFunctionBegin; 609 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 610 PetscTryMethod(pc, "PCShellSetApplySymmetricLeft_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply)); 611 PetscFunctionReturn(PETSC_SUCCESS); 612 } 613 614 /*@C 615 PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the `PC_SYMMETRIC` is used). 616 617 Logically Collective 618 619 Input Parameters: 620 + pc - the preconditioner context 621 - apply - the application-provided right preconditioning routine 622 623 Calling sequence of `apply`: 624 + pc - the preconditioner 625 . xin - input vector 626 - xout - output vector 627 628 Level: advanced 629 630 Note: 631 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 632 633 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()` 634 @*/ 635 PetscErrorCode PCShellSetApplySymmetricRight(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 636 { 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 639 PetscTryMethod(pc, "PCShellSetApplySymmetricRight_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply)); 640 PetscFunctionReturn(PETSC_SUCCESS); 641 } 642 643 /*@C 644 PCShellSetApplyBA - Sets routine to use as the preconditioner times the operator. 645 646 Logically Collective 647 648 Input Parameters: 649 + pc - the preconditioner context 650 - applyBA - the application-provided BA routine 651 652 Calling sequence of `applyBA`: 653 + pc - the preconditioner 654 . side - `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 655 . xin - input vector 656 . xout - output vector 657 - w - work vector 658 659 Level: intermediate 660 661 Note: 662 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applyBA`. 663 664 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellGetContext()`, `PCSide` 665 @*/ 666 PetscErrorCode PCShellSetApplyBA(PC pc, PetscErrorCode (*applyBA)(PC pc, PCSide side, Vec xin, Vec xout, Vec w)) 667 { 668 PetscFunctionBegin; 669 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 670 PetscTryMethod(pc, "PCShellSetApplyBA_C", (PC, PetscErrorCode (*)(PC, PCSide, Vec, Vec, Vec)), (pc, applyBA)); 671 PetscFunctionReturn(PETSC_SUCCESS); 672 } 673 674 /*@C 675 PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose. 676 677 Logically Collective 678 679 Input Parameters: 680 + pc - the preconditioner context 681 - applytranspose - the application-provided preconditioning transpose routine 682 683 Calling sequence of `applytranspose`: 684 + pc - the preconditioner 685 . xin - input vector 686 - xout - output vector 687 688 Level: intermediate 689 690 Note: 691 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applytranspose`. 692 693 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellGetContext()` 694 @*/ 695 PetscErrorCode PCShellSetApplyTranspose(PC pc, PetscErrorCode (*applytranspose)(PC pc, Vec xin, Vec xout)) 696 { 697 PetscFunctionBegin; 698 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 699 PetscTryMethod(pc, "PCShellSetApplyTranspose_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, applytranspose)); 700 PetscFunctionReturn(PETSC_SUCCESS); 701 } 702 703 /*@C 704 PCShellSetMatApplyTranspose - Sets routine to use as preconditioner transpose. 705 706 Logically Collective 707 708 Input Parameters: 709 + pc - the preconditioner context 710 - matapplytranspose - the application-provided preconditioning transpose routine 711 712 Calling sequence of `matapplytranspose`: 713 + pc - the preconditioner 714 . xin - input matrix 715 - xout - output matrix 716 717 Level: intermediate 718 719 Note: 720 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapplytranspose`. 721 722 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellGetContext()` 723 @*/ 724 PetscErrorCode PCShellSetMatApplyTranspose(PC pc, PetscErrorCode (*matapplytranspose)(PC pc, Mat xin, Mat xout)) 725 { 726 PetscFunctionBegin; 727 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 728 PetscTryMethod(pc, "PCShellSetMatApplyTranspose_C", (PC, PetscErrorCode (*)(PC, Mat, Mat)), (pc, matapplytranspose)); 729 PetscFunctionReturn(PETSC_SUCCESS); 730 } 731 732 /*@C 733 PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a `KSPSolve()` is 734 applied. This usually does something like scale the linear system in some application 735 specific way. 736 737 Logically Collective 738 739 Input Parameters: 740 + pc - the preconditioner context 741 - presolve - the application-provided presolve routine, see `PCShellPSolveFn` 742 743 Level: advanced 744 745 Note: 746 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `presolve`. 747 748 .seealso: [](ch_ksp), `PCSHELL`, `PCShellPSolveFn`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPostSolve()`, `PCShellSetContext()`, `PCShellGetContext()` 749 @*/ 750 PetscErrorCode PCShellSetPreSolve(PC pc, PCShellPSolveFn *presolve) 751 { 752 PetscFunctionBegin; 753 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 754 PetscTryMethod(pc, "PCShellSetPreSolve_C", (PC, PCShellPSolveFn *), (pc, presolve)); 755 PetscFunctionReturn(PETSC_SUCCESS); 756 } 757 758 /*@C 759 PCShellSetPostSolve - Sets routine to apply to the operators/vectors after a `KSPSolve()` is 760 applied. This usually does something like scale the linear system in some application 761 specific way. 762 763 Logically Collective 764 765 Input Parameters: 766 + pc - the preconditioner context 767 - postsolve - the application-provided postsolve routine, see `PCShellPSolveFn` 768 769 Level: advanced 770 771 Note: 772 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `postsolve`. 773 774 .seealso: [](ch_ksp), `PCSHELL`, `PCShellPSolveFn`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPreSolve()`, `PCShellSetContext()`, `PCShellGetContext()` 775 @*/ 776 PetscErrorCode PCShellSetPostSolve(PC pc, PCShellPSolveFn *postsolve) 777 { 778 PetscFunctionBegin; 779 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 780 PetscTryMethod(pc, "PCShellSetPostSolve_C", (PC, PCShellPSolveFn *), (pc, postsolve)); 781 PetscFunctionReturn(PETSC_SUCCESS); 782 } 783 784 /*@ 785 PCShellSetName - Sets an optional name to associate with a `PCSHELL` 786 preconditioner. 787 788 Not Collective 789 790 Input Parameters: 791 + pc - the preconditioner context 792 - name - character string describing shell preconditioner 793 794 Level: intermediate 795 796 Note: 797 This is separate from the name you can provide with `PetscObjectSetName()` 798 799 .seealso: [](ch_ksp), `PCSHELL`, `PCShellGetName()`, `PetscObjectSetName()`, `PetscObjectGetName()` 800 @*/ 801 PetscErrorCode PCShellSetName(PC pc, const char name[]) 802 { 803 PetscFunctionBegin; 804 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 805 PetscTryMethod(pc, "PCShellSetName_C", (PC, const char[]), (pc, name)); 806 PetscFunctionReturn(PETSC_SUCCESS); 807 } 808 809 /*@ 810 PCShellGetName - Gets an optional name that the user has set for a `PCSHELL` with `PCShellSetName()` 811 preconditioner. 812 813 Not Collective 814 815 Input Parameter: 816 . pc - the preconditioner context 817 818 Output Parameter: 819 . name - character string describing shell preconditioner (you should not free this) 820 821 Level: intermediate 822 823 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetName()`, `PetscObjectSetName()`, `PetscObjectGetName()` 824 @*/ 825 PetscErrorCode PCShellGetName(PC pc, const char *name[]) 826 { 827 PetscFunctionBegin; 828 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 829 PetscAssertPointer(name, 2); 830 PetscUseMethod(pc, "PCShellGetName_C", (PC, const char *[]), (pc, name)); 831 PetscFunctionReturn(PETSC_SUCCESS); 832 } 833 834 /*@C 835 PCShellSetApplyRichardson - Sets routine to use as preconditioner 836 in Richardson iteration. 837 838 Logically Collective 839 840 Input Parameters: 841 + pc - the preconditioner context 842 - apply - the application-provided preconditioning routine 843 844 Calling sequence of `apply`: 845 + pc - the preconditioner 846 . b - right-hand side 847 . x - current iterate 848 . r - work space 849 . rtol - relative tolerance of residual norm to stop at 850 . abstol - absolute tolerance of residual norm to stop at 851 . dtol - if residual norm increases by this factor than return 852 . maxits - number of iterations to run 853 . zeroinitialguess - `PETSC_TRUE` if `x` is known to be initially zero 854 . its - returns the number of iterations used 855 - reason - returns the reason the iteration has converged 856 857 Level: advanced 858 859 Notes: 860 You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 861 862 This is used when one can provide code for multiple steps of Richardson's method that is more efficient than computing a single step, 863 recomputing the residual via $ r = b - A x $, and then computing the next step. SOR is an algorithm for which this is true. 864 865 .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCRichardsonConvergedReason()`, `PCShellGetContext()`, `KSPRICHARDSON` 866 @*/ 867 PetscErrorCode PCShellSetApplyRichardson(PC pc, PetscErrorCode (*apply)(PC pc, Vec b, Vec x, Vec r, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits, PetscBool zeroinitialguess, PetscInt *its, PCRichardsonConvergedReason *reason)) 868 { 869 PetscFunctionBegin; 870 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 871 PetscTryMethod(pc, "PCShellSetApplyRichardson_C", (PC, PetscErrorCode (*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)), (pc, apply)); 872 PetscFunctionReturn(PETSC_SUCCESS); 873 } 874 875 /*MC 876 PCSHELL - Creates a new preconditioner class for use with a users 877 own private data storage format and preconditioner application code 878 879 Level: advanced 880 881 Usage: 882 .vb 883 extern PetscErrorCode apply(PC,Vec,Vec); 884 extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec); 885 extern PetscErrorCode applytranspose(PC,Vec,Vec); 886 extern PetscErrorCode setup(PC); 887 extern PetscErrorCode destroy(PC); 888 889 PCCreate(comm,&pc); 890 PCSetType(pc,PCSHELL); 891 PCShellSetContext(pc,ctx) 892 PCShellSetApply(pc,apply); 893 PCShellSetApplyBA(pc,applyba); (optional) 894 PCShellSetApplyTranspose(pc,applytranspose); (optional) 895 PCShellSetSetUp(pc,setup); (optional) 896 PCShellSetDestroy(pc,destroy); (optional) 897 .ve 898 899 Notes: 900 Information required for the preconditioner and its internal datastructures can be set with `PCShellSetContext()` and then accessed 901 with `PCShellGetContext()` inside the routines provided above. 902 903 When using `MATSHELL`, where the explicit entries of matrix are not available to build the preconditioner, `PCSHELL` can be used 904 to construct a custom preconditioner for the `MATSHELL`, assuming the user knows enough about their problem to provide a 905 custom preconditioner. 906 907 .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 908 `MATSHELL`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetView()`, `PCShellSetDestroy()`, `PCShellSetPostSolve()`, 909 `PCShellSetApplyTranspose()`, `PCShellSetName()`, `PCShellSetApplyRichardson()`, `PCShellSetPreSolve()`, `PCShellSetView()`, 910 `PCShellGetName()`, `PCShellSetContext()`, `PCShellGetContext()`, `PCShellSetApplyBA()`, `MATSHELL`, `PCShellSetMatApply()`, 911 M*/ 912 913 PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc) 914 { 915 PC_Shell *shell; 916 917 PetscFunctionBegin; 918 PetscCall(PetscNew(&shell)); 919 pc->data = (void *)shell; 920 921 pc->ops->destroy = PCDestroy_Shell; 922 pc->ops->view = PCView_Shell; 923 pc->ops->apply = PCApply_Shell; 924 pc->ops->applysymmetricleft = PCApplySymmetricLeft_Shell; 925 pc->ops->applysymmetricright = PCApplySymmetricRight_Shell; 926 pc->ops->matapply = NULL; 927 pc->ops->applytranspose = NULL; 928 pc->ops->applyrichardson = NULL; 929 pc->ops->setup = NULL; 930 pc->ops->presolve = NULL; 931 pc->ops->postsolve = NULL; 932 933 shell->apply = NULL; 934 shell->applytranspose = NULL; 935 shell->name = NULL; 936 shell->applyrich = NULL; 937 shell->presolve = NULL; 938 shell->postsolve = NULL; 939 shell->ctx = NULL; 940 shell->setup = NULL; 941 shell->view = NULL; 942 shell->destroy = NULL; 943 shell->applysymmetricleft = NULL; 944 shell->applysymmetricright = NULL; 945 946 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", PCShellSetDestroy_Shell)); 947 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", PCShellSetSetUp_Shell)); 948 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", PCShellSetApply_Shell)); 949 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", PCShellSetMatApply_Shell)); 950 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", PCShellSetApplySymmetricLeft_Shell)); 951 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", PCShellSetApplySymmetricRight_Shell)); 952 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", PCShellSetApplyBA_Shell)); 953 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", PCShellSetPreSolve_Shell)); 954 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", PCShellSetPostSolve_Shell)); 955 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", PCShellSetView_Shell)); 956 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", PCShellSetApplyTranspose_Shell)); 957 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApplyTranspose_C", PCShellSetMatApplyTranspose_Shell)); 958 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", PCShellSetName_Shell)); 959 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", PCShellGetName_Shell)); 960 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", PCShellSetApplyRichardson_Shell)); 961 PetscFunctionReturn(PETSC_SUCCESS); 962 } 963