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