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