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