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