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 themselve 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 themselve 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 themselve 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 PetscFunctionReturn(0); 319 } 320 321 static PetscErrorCode PCShellSetApplySymmetricLeft_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec)) 322 { 323 PC_Shell *shell = (PC_Shell*)pc->data; 324 325 PetscFunctionBegin; 326 shell->applysymmetricleft = apply; 327 PetscFunctionReturn(0); 328 } 329 330 static PetscErrorCode PCShellSetApplySymmetricRight_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec)) 331 { 332 PC_Shell *shell = (PC_Shell*)pc->data; 333 334 PetscFunctionBegin; 335 shell->applysymmetricright = apply; 336 PetscFunctionReturn(0); 337 } 338 339 static PetscErrorCode PCShellSetApplyBA_Shell(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec)) 340 { 341 PC_Shell *shell = (PC_Shell*)pc->data; 342 343 PetscFunctionBegin; 344 shell->applyBA = applyBA; 345 if (applyBA) pc->ops->applyBA = PCApplyBA_Shell; 346 else pc->ops->applyBA = NULL; 347 PetscFunctionReturn(0); 348 } 349 350 static PetscErrorCode PCShellSetPreSolve_Shell(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec)) 351 { 352 PC_Shell *shell = (PC_Shell*)pc->data; 353 354 PetscFunctionBegin; 355 shell->presolve = presolve; 356 if (presolve) { 357 pc->ops->presolve = PCPreSolve_Shell; 358 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Shell)); 359 } else { 360 pc->ops->presolve = NULL; 361 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL)); 362 } 363 PetscFunctionReturn(0); 364 } 365 366 static PetscErrorCode PCShellSetPostSolve_Shell(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec)) 367 { 368 PC_Shell *shell = (PC_Shell*)pc->data; 369 370 PetscFunctionBegin; 371 shell->postsolve = postsolve; 372 if (postsolve) pc->ops->postsolve = PCPostSolve_Shell; 373 else pc->ops->postsolve = NULL; 374 PetscFunctionReturn(0); 375 } 376 377 static PetscErrorCode PCShellSetView_Shell(PC pc,PetscErrorCode (*view)(PC,PetscViewer)) 378 { 379 PC_Shell *shell = (PC_Shell*)pc->data; 380 381 PetscFunctionBegin; 382 shell->view = view; 383 PetscFunctionReturn(0); 384 } 385 386 static PetscErrorCode PCShellSetApplyTranspose_Shell(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec)) 387 { 388 PC_Shell *shell = (PC_Shell*)pc->data; 389 390 PetscFunctionBegin; 391 shell->applytranspose = applytranspose; 392 if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell; 393 else pc->ops->applytranspose = NULL; 394 PetscFunctionReturn(0); 395 } 396 397 static PetscErrorCode PCShellSetApplyRichardson_Shell(PC pc,PetscErrorCode (*applyrich)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)) 398 { 399 PC_Shell *shell = (PC_Shell*)pc->data; 400 401 PetscFunctionBegin; 402 shell->applyrich = applyrich; 403 if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell; 404 else pc->ops->applyrichardson = NULL; 405 PetscFunctionReturn(0); 406 } 407 408 static PetscErrorCode PCShellSetName_Shell(PC pc,const char name[]) 409 { 410 PC_Shell *shell = (PC_Shell*)pc->data; 411 412 PetscFunctionBegin; 413 PetscCall(PetscFree(shell->name)); 414 PetscCall(PetscStrallocpy(name,&shell->name)); 415 PetscFunctionReturn(0); 416 } 417 418 static PetscErrorCode PCShellGetName_Shell(PC pc,const char *name[]) 419 { 420 PC_Shell *shell = (PC_Shell*)pc->data; 421 422 PetscFunctionBegin; 423 *name = shell->name; 424 PetscFunctionReturn(0); 425 } 426 427 /* -------------------------------------------------------------------------------*/ 428 429 /*@C 430 PCShellSetDestroy - Sets routine to use to destroy the user-provided 431 application context. 432 433 Logically Collective on PC 434 435 Input Parameters: 436 + pc - the preconditioner context 437 - destroy - the application-provided destroy routine 438 439 Calling sequence of destroy: 440 .vb 441 PetscErrorCode destroy (PC) 442 .ve 443 444 . ptr - the application context 445 446 Notes: 447 the function MUST return an error code of 0 on success and nonzero on failure. 448 449 Level: developer 450 451 .seealso: `PCShellSetApply()`, `PCShellSetContext()` 452 @*/ 453 PetscErrorCode PCShellSetDestroy(PC pc,PetscErrorCode (*destroy)(PC)) 454 { 455 PetscFunctionBegin; 456 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 457 PetscTryMethod(pc,"PCShellSetDestroy_C",(PC,PetscErrorCode (*)(PC)),(pc,destroy)); 458 PetscFunctionReturn(0); 459 } 460 461 /*@C 462 PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the 463 matrix operator is changed. 464 465 Logically Collective on PC 466 467 Input Parameters: 468 + pc - the preconditioner context 469 - setup - the application-provided setup routine 470 471 Calling sequence of setup: 472 .vb 473 PetscErrorCode setup (PC pc) 474 .ve 475 476 . pc - the preconditioner, get the application context with PCShellGetContext() 477 478 Notes: 479 the function MUST return an error code of 0 on success and nonzero on failure. 480 481 Level: developer 482 483 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetApply()`, `PCShellSetContext()` 484 @*/ 485 PetscErrorCode PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC)) 486 { 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 489 PetscTryMethod(pc,"PCShellSetSetUp_C",(PC,PetscErrorCode (*)(PC)),(pc,setup)); 490 PetscFunctionReturn(0); 491 } 492 493 /*@C 494 PCShellSetView - Sets routine to use as viewer of shell preconditioner 495 496 Logically Collective on PC 497 498 Input Parameters: 499 + pc - the preconditioner context 500 - view - the application-provided view routine 501 502 Calling sequence of view: 503 .vb 504 PetscErrorCode view(PC pc,PetscViewer v) 505 .ve 506 507 + pc - the preconditioner, get the application context with PCShellGetContext() 508 - v - viewer 509 510 Notes: 511 the function MUST return an error code of 0 on success and nonzero on failure. 512 513 Level: developer 514 515 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()` 516 @*/ 517 PetscErrorCode PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer)) 518 { 519 PetscFunctionBegin; 520 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 521 PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view)); 522 PetscFunctionReturn(0); 523 } 524 525 /*@C 526 PCShellSetApply - Sets routine to use as preconditioner. 527 528 Logically Collective on PC 529 530 Input Parameters: 531 + pc - the preconditioner context 532 - apply - the application-provided preconditioning routine 533 534 Calling sequence of apply: 535 .vb 536 PetscErrorCode apply (PC pc,Vec xin,Vec xout) 537 .ve 538 539 + pc - the preconditioner, get the application context with PCShellGetContext() 540 . xin - input vector 541 - xout - output vector 542 543 Notes: 544 the function MUST return an error code of 0 on success and nonzero on failure. 545 546 Level: developer 547 548 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellSetApplySymmetricRight(),PCShellSetApplySymmetricLeft()` 549 @*/ 550 PetscErrorCode PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec)) 551 { 552 PetscFunctionBegin; 553 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 554 PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply)); 555 PetscFunctionReturn(0); 556 } 557 558 /*@C 559 PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors. 560 561 Logically Collective on PC 562 563 Input Parameters: 564 + pc - the preconditioner context 565 - apply - the application-provided preconditioning routine 566 567 Calling sequence of apply: 568 .vb 569 PetscErrorCode apply (PC pc,Mat Xin,Mat Xout) 570 .ve 571 572 + pc - the preconditioner, get the application context with PCShellGetContext() 573 . Xin - input block of vectors 574 - Xout - output block of vectors 575 576 Notes: 577 the function MUST return an error code of 0 on success and nonzero on failure. 578 579 Level: developer 580 581 .seealso: `PCShellSetApply()` 582 @*/ 583 PetscErrorCode PCShellSetMatApply(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat)) 584 { 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 587 PetscTryMethod(pc,"PCShellSetMatApply_C",(PC,PetscErrorCode (*)(PC,Mat,Mat)),(pc,matapply)); 588 PetscFunctionReturn(0); 589 } 590 591 /*@C 592 PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the PC_SYMMETRIC is used). 593 594 Logically Collective on PC 595 596 Input Parameters: 597 + pc - the preconditioner context 598 - apply - the application-provided left preconditioning routine 599 600 Calling sequence of apply: 601 .vb 602 PetscErrorCode apply (PC pc,Vec xin,Vec xout) 603 .ve 604 605 + pc - the preconditioner, get the application context with PCShellGetContext() 606 . xin - input vector 607 - xout - output vector 608 609 Notes: 610 the function MUST return an error code of 0 on success and nonzero on failure. 611 612 Level: developer 613 614 .seealso: `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()` 615 @*/ 616 PetscErrorCode PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec)) 617 { 618 PetscFunctionBegin; 619 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 620 PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply)); 621 PetscFunctionReturn(0); 622 } 623 624 /*@C 625 PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the PC_SYMMETRIC is used). 626 627 Logically Collective on PC 628 629 Input Parameters: 630 + pc - the preconditioner context 631 - apply - the application-provided right preconditioning routine 632 633 Calling sequence of apply: 634 .vb 635 PetscErrorCode apply (PC pc,Vec xin,Vec xout) 636 .ve 637 638 + pc - the preconditioner, get the application context with PCShellGetContext() 639 . xin - input vector 640 - xout - output vector 641 642 Notes: 643 the function MUST return an error code of 0 on success and nonzero on failure. 644 645 Level: developer 646 647 .seealso: `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()` 648 @*/ 649 PetscErrorCode PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec)) 650 { 651 PetscFunctionBegin; 652 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 653 PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply)); 654 PetscFunctionReturn(0); 655 } 656 657 /*@C 658 PCShellSetApplyBA - Sets routine to use as preconditioner times operator. 659 660 Logically Collective on PC 661 662 Input Parameters: 663 + pc - the preconditioner context 664 - applyBA - the application-provided BA routine 665 666 Calling sequence of applyBA: 667 .vb 668 PetscErrorCode applyBA (PC pc,Vec xin,Vec xout) 669 .ve 670 671 + pc - the preconditioner, get the application context with PCShellGetContext() 672 . xin - input vector 673 - xout - output vector 674 675 Notes: 676 the function MUST return an error code of 0 on success and nonzero on failure. 677 678 Level: developer 679 680 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApply()` 681 @*/ 682 PetscErrorCode PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec)) 683 { 684 PetscFunctionBegin; 685 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 686 PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA)); 687 PetscFunctionReturn(0); 688 } 689 690 /*@C 691 PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose. 692 693 Logically Collective on PC 694 695 Input Parameters: 696 + pc - the preconditioner context 697 - apply - the application-provided preconditioning transpose routine 698 699 Calling sequence of apply: 700 .vb 701 PetscErrorCode applytranspose (PC pc,Vec xin,Vec xout) 702 .ve 703 704 + pc - the preconditioner, get the application context with PCShellGetContext() 705 . xin - input vector 706 - xout - output vector 707 708 Notes: 709 the function MUST return an error code of 0 on success and nonzero on failure. 710 711 Level: developer 712 713 Notes: 714 Uses the same context variable as PCShellSetApply(). 715 716 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCSetContext()`, `PCShellSetApplyBA()` 717 @*/ 718 PetscErrorCode PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec)) 719 { 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 722 PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose)); 723 PetscFunctionReturn(0); 724 } 725 726 /*@C 727 PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is 728 applied. This usually does something like scale the linear system in some application 729 specific way. 730 731 Logically Collective on PC 732 733 Input Parameters: 734 + pc - the preconditioner context 735 - presolve - the application-provided presolve routine 736 737 Calling sequence of presolve: 738 .vb 739 PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x) 740 .ve 741 742 + pc - the preconditioner, get the application context with PCShellGetContext() 743 . xin - input vector 744 - xout - output vector 745 746 Notes: 747 the function MUST return an error code of 0 on success and nonzero on failure. 748 749 Level: developer 750 751 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPostSolve()`, `PCShellSetContext()` 752 @*/ 753 PetscErrorCode PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec)) 754 { 755 PetscFunctionBegin; 756 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 757 PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve)); 758 PetscFunctionReturn(0); 759 } 760 761 /*@C 762 PCShellSetPostSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is 763 applied. This usually does something like scale the linear system in some application 764 specific way. 765 766 Logically Collective on PC 767 768 Input Parameters: 769 + pc - the preconditioner context 770 - postsolve - the application-provided presolve routine 771 772 Calling sequence of postsolve: 773 .vb 774 PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x) 775 .ve 776 777 + pc - the preconditioner, get the application context with PCShellGetContext() 778 . xin - input vector 779 - xout - output vector 780 781 Notes: 782 the function MUST return an error code of 0 on success and nonzero on failure. 783 784 Level: developer 785 786 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPreSolve()`, `PCShellSetContext()` 787 @*/ 788 PetscErrorCode PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec)) 789 { 790 PetscFunctionBegin; 791 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 792 PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve)); 793 PetscFunctionReturn(0); 794 } 795 796 /*@C 797 PCShellSetName - Sets an optional name to associate with a shell 798 preconditioner. 799 800 Not Collective 801 802 Input Parameters: 803 + pc - the preconditioner context 804 - name - character string describing shell preconditioner 805 806 Level: developer 807 808 .seealso: `PCShellGetName()` 809 @*/ 810 PetscErrorCode PCShellSetName(PC pc,const char name[]) 811 { 812 PetscFunctionBegin; 813 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 814 PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name)); 815 PetscFunctionReturn(0); 816 } 817 818 /*@C 819 PCShellGetName - Gets an optional name that the user has set for a shell 820 preconditioner. 821 822 Not Collective 823 824 Input Parameter: 825 . pc - the preconditioner context 826 827 Output Parameter: 828 . name - character string describing shell preconditioner (you should not free this) 829 830 Level: developer 831 832 .seealso: `PCShellSetName()` 833 @*/ 834 PetscErrorCode PCShellGetName(PC pc,const char *name[]) 835 { 836 PetscFunctionBegin; 837 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 838 PetscValidPointer(name,2); 839 PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name)); 840 PetscFunctionReturn(0); 841 } 842 843 /*@C 844 PCShellSetApplyRichardson - Sets routine to use as preconditioner 845 in Richardson iteration. 846 847 Logically Collective on PC 848 849 Input Parameters: 850 + pc - the preconditioner context 851 - apply - the application-provided preconditioning routine 852 853 Calling sequence of apply: 854 .vb 855 PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits) 856 .ve 857 858 + pc - the preconditioner, get the application context with PCShellGetContext() 859 . b - right-hand-side 860 . x - current iterate 861 . r - work space 862 . rtol - relative tolerance of residual norm to stop at 863 . abstol - absolute tolerance of residual norm to stop at 864 . dtol - if residual norm increases by this factor than return 865 - maxits - number of iterations to run 866 867 Notes: 868 the function MUST return an error code of 0 on success and nonzero on failure. 869 870 Level: developer 871 872 .seealso: `PCShellSetApply()`, `PCShellSetContext()` 873 @*/ 874 PetscErrorCode PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)) 875 { 876 PetscFunctionBegin; 877 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 878 PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply)); 879 PetscFunctionReturn(0); 880 } 881 882 /*MC 883 PCSHELL - Creates a new preconditioner class for use with your 884 own private data storage format. 885 886 Level: advanced 887 888 Usage: 889 $ extern PetscErrorCode apply(PC,Vec,Vec); 890 $ extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec); 891 $ extern PetscErrorCode applytranspose(PC,Vec,Vec); 892 $ extern PetscErrorCode setup(PC); 893 $ extern PetscErrorCode destroy(PC); 894 $ 895 $ PCCreate(comm,&pc); 896 $ PCSetType(pc,PCSHELL); 897 $ PCShellSetContext(pc,ctx) 898 $ PCShellSetApply(pc,apply); 899 $ PCShellSetApplyBA(pc,applyba); (optional) 900 $ PCShellSetApplyTranspose(pc,applytranspose); (optional) 901 $ PCShellSetSetUp(pc,setup); (optional) 902 $ PCShellSetDestroy(pc,destroy); (optional) 903 904 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 905 `MATSHELL`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetView()`, 906 `PCShellSetApplyTranspose()`, `PCShellSetName()`, `PCShellSetApplyRichardson()`, 907 `PCShellGetName()`, `PCShellSetContext()`, `PCShellGetContext()`, `PCShellSetApplyBA()` 908 M*/ 909 910 PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc) 911 { 912 PC_Shell *shell; 913 914 PetscFunctionBegin; 915 PetscCall(PetscNewLog(pc,&shell)); 916 pc->data = (void*)shell; 917 918 pc->ops->destroy = PCDestroy_Shell; 919 pc->ops->view = PCView_Shell; 920 pc->ops->apply = PCApply_Shell; 921 pc->ops->matapply = PCMatApply_Shell; 922 pc->ops->applysymmetricleft = PCApplySymmetricLeft_Shell; 923 pc->ops->applysymmetricright = PCApplySymmetricRight_Shell; 924 pc->ops->applytranspose = NULL; 925 pc->ops->applyrichardson = NULL; 926 pc->ops->setup = NULL; 927 pc->ops->presolve = NULL; 928 pc->ops->postsolve = NULL; 929 930 shell->apply = NULL; 931 shell->applytranspose = NULL; 932 shell->name = NULL; 933 shell->applyrich = NULL; 934 shell->presolve = NULL; 935 shell->postsolve = NULL; 936 shell->ctx = NULL; 937 shell->setup = NULL; 938 shell->view = NULL; 939 shell->destroy = NULL; 940 shell->applysymmetricleft = NULL; 941 shell->applysymmetricright = NULL; 942 943 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",PCShellSetDestroy_Shell)); 944 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",PCShellSetSetUp_Shell)); 945 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",PCShellSetApply_Shell)); 946 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",PCShellSetMatApply_Shell)); 947 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",PCShellSetApplySymmetricLeft_Shell)); 948 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",PCShellSetApplySymmetricRight_Shell)); 949 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",PCShellSetApplyBA_Shell)); 950 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",PCShellSetPreSolve_Shell)); 951 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",PCShellSetPostSolve_Shell)); 952 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",PCShellSetView_Shell)); 953 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",PCShellSetApplyTranspose_Shell)); 954 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",PCShellSetName_Shell)); 955 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",PCShellGetName_Shell)); 956 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",PCShellSetApplyRichardson_Shell)); 957 PetscFunctionReturn(0); 958 } 959