Actual source code: shellpc.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: /*
  3:    This provides a simple shell for Fortran (and C programmers) to
  4:   create their own preconditioner without writing much interface code.
  5: */

  7:  #include <petsc/private/pcimpl.h>

  9: typedef struct {
 10:   void *ctx;                     /* user provided contexts for preconditioner */

 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*);

 24:   char *name;
 25: } PC_Shell;

 27: /*@C
 28:     PCShellGetContext - Returns the user-provided context associated with a shell PC

 30:     Not Collective

 32:     Input Parameter:
 33: .   pc - should have been created with PCSetType(pc,shell)

 35:     Output Parameter:
 36: .   ctx - the user provided context

 38:     Level: advanced

 40:     Notes:
 41:     This routine is intended for use within various shell routines

 43:    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
 44:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

 46: .keywords: PC, shell, get, context

 48: .seealso: PCShellSetContext()
 49: @*/
 50: PetscErrorCode  PCShellGetContext(PC pc,void **ctx)
 51: {
 53:   PetscBool      flg;

 58:   PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&flg);
 59:   if (!flg) *ctx = 0;
 60:   else      *ctx = ((PC_Shell*)(pc->data))->ctx;
 61:   return(0);
 62: }

 64: /*@
 65:     PCShellSetContext - sets the context for a shell PC

 67:    Logically Collective on PC

 69:     Input Parameters:
 70: +   pc - the shell PC
 71: -   ctx - the context

 73:    Level: advanced

 75:    Fortran Notes: 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.



 80: .seealso: PCShellGetContext(), PCSHELL
 81: @*/
 82: PetscErrorCode  PCShellSetContext(PC pc,void *ctx)
 83: {
 84:   PC_Shell       *shell = (PC_Shell*)pc->data;
 86:   PetscBool      flg;

 90:   PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&flg);
 91:   if (flg) shell->ctx = ctx;
 92:   return(0);
 93: }

 95: static PetscErrorCode PCSetUp_Shell(PC pc)
 96: {
 97:   PC_Shell       *shell = (PC_Shell*)pc->data;

101:   if (!shell->setup) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No setup() routine provided to Shell PC");
102:   PetscStackCall("PCSHELL user function setup()",(*shell->setup)(pc);CHKERRQ(ierr));
103:   return(0);
104: }

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;

113:   if (!shell->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
114:   PetscObjectStateGet((PetscObject)y, &instate);
115:   PetscStackCall("PCSHELL user function apply()",(*shell->apply)(pc,x,y);CHKERRQ(ierr));
116:   PetscObjectStateGet((PetscObject)y, &outstate);
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:     PetscObjectStateIncrease((PetscObject)y);
120:   }
121:   return(0);
122: }

124: static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc,Vec x,Vec y)
125: {
126:   PC_Shell       *shell = (PC_Shell*)pc->data;

130:   if (!shell->applysymmetricleft) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
131:   PetscStackCall("PCSHELL user function apply()",(*shell->applysymmetricleft)(pc,x,y);CHKERRQ(ierr));
132:   return(0);
133: }

135: static PetscErrorCode PCApplySymmetricRight_Shell(PC pc,Vec x,Vec y)
136: {
137:   PC_Shell       *shell = (PC_Shell*)pc->data;

141:   if (!shell->applysymmetricright) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
142:   PetscStackCall("PCSHELL user function apply()",(*shell->applysymmetricright)(pc,x,y);CHKERRQ(ierr));
143:   return(0);
144: }

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;

153:   if (!shell->applyBA) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyBA() routine provided to Shell PC");
154:   PetscObjectStateGet((PetscObject)w, &instate);
155:   PetscStackCall("PCSHELL user function applyBA()",(*shell->applyBA)(pc,side,x,y,w);CHKERRQ(ierr));
156:   PetscObjectStateGet((PetscObject)w, &outstate);
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:     PetscObjectStateIncrease((PetscObject)w);
160:   }
161:   return(0);
162: }

164: static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc,PetscBool* change)
165: {
167:   *change = PETSC_TRUE;
168:   return(0);
169: }

171: static PetscErrorCode PCPreSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)
172: {
173:   PC_Shell       *shell = (PC_Shell*)pc->data;

177:   if (!shell->presolve) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No presolve() routine provided to Shell PC");
178:   PetscStackCall("PCSHELL user function presolve()",(*shell->presolve)(pc,ksp,b,x);CHKERRQ(ierr));
179:   return(0);
180: }

182: static PetscErrorCode PCPostSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)
183: {
184:   PC_Shell       *shell = (PC_Shell*)pc->data;

188:   if (!shell->postsolve) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No postsolve() routine provided to Shell PC");
189:   PetscStackCall("PCSHELL user function postsolve()",(*shell->postsolve)(pc,ksp,b,x);CHKERRQ(ierr));
190:   return(0);
191: }

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;

200:   if (!shell->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applytranspose() routine provided to Shell PC");
201:   PetscObjectStateGet((PetscObject)y, &instate);
202:   PetscStackCall("PCSHELL user function applytranspose()",(*shell->applytranspose)(pc,x,y);CHKERRQ(ierr));
203:   PetscObjectStateGet((PetscObject)y, &outstate);
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:     PetscObjectStateIncrease((PetscObject)y);
207:   }
208:   return(0);
209: }

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;

218:   if (!shell->applyrich) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyrichardson() routine provided to Shell PC");
219:   PetscObjectStateGet((PetscObject)y, &instate);
220:   PetscStackCall("PCSHELL user function applyrichardson()",(*shell->applyrich)(pc,x,y,w,rtol,abstol,dtol,it,guesszero,outits,reason);CHKERRQ(ierr));
221:   PetscObjectStateGet((PetscObject)y, &outstate);
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:     PetscObjectStateIncrease((PetscObject)y);
225:   }
226:   return(0);
227: }

229: static PetscErrorCode PCDestroy_Shell(PC pc)
230: {
231:   PC_Shell       *shell = (PC_Shell*)pc->data;

235:   PetscFree(shell->name);
236:   if (shell->destroy) PetscStackCall("PCSHELL user function destroy()",(*shell->destroy)(pc);CHKERRQ(ierr));
237:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",NULL);
238:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",NULL);
239:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",NULL);
240:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",NULL);
241:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",NULL);
242:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",NULL);
243:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",NULL);
244:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",NULL);
245:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",NULL);
246:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",NULL);
247:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",NULL);
248:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",NULL);
249:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",NULL);
250:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
251:   PetscFree(pc->data);
252:   return(0);
253: }

255: static PetscErrorCode PCView_Shell(PC pc,PetscViewer viewer)
256: {
257:   PC_Shell       *shell = (PC_Shell*)pc->data;
259:   PetscBool      iascii;

262:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
263:   if (iascii) {
264:     if (shell->name) {
265:       PetscViewerASCIIPrintf(viewer,"  %s\n",shell->name);
266:     } else {
267:       PetscViewerASCIIPrintf(viewer,"  no name\n");
268:     }
269:   }
270:   if (shell->view) {
271:     PetscViewerASCIIPushTab(viewer);
272:     (*shell->view)(pc,viewer);
273:     PetscViewerASCIIPopTab(viewer);
274:   }
275:   return(0);
276: }

278: /* ------------------------------------------------------------------------------*/
279: static PetscErrorCode  PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC))
280: {
281:   PC_Shell *shell= (PC_Shell*)pc->data;

284:   shell->destroy = destroy;
285:   return(0);
286: }

288: static PetscErrorCode  PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC))
289: {
290:   PC_Shell *shell = (PC_Shell*)pc->data;;

293:   shell->setup = setup;
294:   if (setup) pc->ops->setup = PCSetUp_Shell;
295:   else       pc->ops->setup = 0;
296:   return(0);
297: }

299: static PetscErrorCode  PCShellSetApply_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
300: {
301:   PC_Shell *shell = (PC_Shell*)pc->data;

304:   shell->apply = apply;
305:   return(0);
306: }

308: static PetscErrorCode  PCShellSetApplySymmetricLeft_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
309: {
310:   PC_Shell *shell = (PC_Shell*)pc->data;

313:   shell->applysymmetricleft = apply;
314:   return(0);
315: }

317: static PetscErrorCode  PCShellSetApplySymmetricRight_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
318: {
319:   PC_Shell *shell = (PC_Shell*)pc->data;

322:   shell->applysymmetricright = apply;
323:   return(0);
324: }

326: static PetscErrorCode  PCShellSetApplyBA_Shell(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
327: {
328:   PC_Shell *shell = (PC_Shell*)pc->data;

331:   shell->applyBA = applyBA;
332:   if (applyBA) pc->ops->applyBA  = PCApplyBA_Shell;
333:   else         pc->ops->applyBA  = 0;
334:   return(0);
335: }

337: static PetscErrorCode  PCShellSetPreSolve_Shell(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
338: {
339:   PC_Shell       *shell = (PC_Shell*)pc->data;

343:   shell->presolve = presolve;
344:   if (presolve) {
345:     pc->ops->presolve = PCPreSolve_Shell;
346:     PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Shell);
347:   } else {
348:     pc->ops->presolve = 0;
349:     PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
350:   }
351:   return(0);
352: }

354: static PetscErrorCode  PCShellSetPostSolve_Shell(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
355: {
356:   PC_Shell *shell = (PC_Shell*)pc->data;

359:   shell->postsolve = postsolve;
360:   if (postsolve) pc->ops->postsolve = PCPostSolve_Shell;
361:   else           pc->ops->postsolve = 0;
362:   return(0);
363: }

365: static PetscErrorCode  PCShellSetView_Shell(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
366: {
367:   PC_Shell *shell = (PC_Shell*)pc->data;

370:   shell->view = view;
371:   return(0);
372: }

374: static PetscErrorCode  PCShellSetApplyTranspose_Shell(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
375: {
376:   PC_Shell *shell = (PC_Shell*)pc->data;

379:   shell->applytranspose = applytranspose;
380:   if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell;
381:   else                pc->ops->applytranspose = 0;
382:   return(0);
383: }

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;

390:   shell->applyrich = applyrich;
391:   if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell;
392:   else           pc->ops->applyrichardson = 0;
393:   return(0);
394: }

396: static PetscErrorCode  PCShellSetName_Shell(PC pc,const char name[])
397: {
398:   PC_Shell       *shell = (PC_Shell*)pc->data;

402:   PetscFree(shell->name);
403:   PetscStrallocpy(name,&shell->name);
404:   return(0);
405: }

407: static PetscErrorCode  PCShellGetName_Shell(PC pc,const char *name[])
408: {
409:   PC_Shell *shell = (PC_Shell*)pc->data;

412:   *name = shell->name;
413:   return(0);
414: }

416: /* -------------------------------------------------------------------------------*/

418: /*@C
419:    PCShellSetDestroy - Sets routine to use to destroy the user-provided
420:    application context.

422:    Logically Collective on PC

424:    Input Parameters:
425: +  pc - the preconditioner context
426: .  destroy - the application-provided destroy routine

428:    Calling sequence of destroy:
429: .vb
430:    PetscErrorCode destroy (PC)
431: .ve

433: .  ptr - the application context

435:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

437:    Level: developer

439: .keywords: PC, shell, set, destroy, user-provided

441: .seealso: PCShellSetApply(), PCShellSetContext()
442: @*/
443: PetscErrorCode  PCShellSetDestroy(PC pc,PetscErrorCode (*destroy)(PC))
444: {

449:   PetscTryMethod(pc,"PCShellSetDestroy_C",(PC,PetscErrorCode (*)(PC)),(pc,destroy));
450:   return(0);
451: }


454: /*@C
455:    PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the
456:    matrix operator is changed.

458:    Logically Collective on PC

460:    Input Parameters:
461: +  pc - the preconditioner context
462: .  setup - the application-provided setup routine

464:    Calling sequence of setup:
465: .vb
466:    PetscErrorCode setup (PC pc)
467: .ve

469: .  pc - the preconditioner, get the application context with PCShellGetContext()

471:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

473:    Level: developer

475: .keywords: PC, shell, set, setup, user-provided

477: .seealso: PCShellSetApplyRichardson(), PCShellSetApply(), PCShellSetContext()
478: @*/
479: PetscErrorCode  PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC))
480: {

485:   PetscTryMethod(pc,"PCShellSetSetUp_C",(PC,PetscErrorCode (*)(PC)),(pc,setup));
486:   return(0);
487: }


490: /*@C
491:    PCShellSetView - Sets routine to use as viewer of shell preconditioner

493:    Logically Collective on PC

495:    Input Parameters:
496: +  pc - the preconditioner context
497: -  view - the application-provided view routine

499:    Calling sequence of apply:
500: .vb
501:    PetscErrorCode view(PC pc,PetscViewer v)
502: .ve

504: +  pc - the preconditioner, get the application context with PCShellGetContext()
505: -  v   - viewer

507:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

509:    Level: developer

511: .keywords: PC, shell, set, apply, user-provided

513: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose()
514: @*/
515: PetscErrorCode  PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
516: {

521:   PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view));
522:   return(0);
523: }

525: /*@C
526:    PCShellSetApply - Sets routine to use as preconditioner.

528:    Logically Collective on PC

530:    Input Parameters:
531: +  pc - the preconditioner context
532: -  apply - the application-provided preconditioning routine

534:    Calling sequence of apply:
535: .vb
536:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
537: .ve

539: +  pc - the preconditioner, get the application context with PCShellGetContext()
540: .  xin - input vector
541: -  xout - output vector

543:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

545:    Level: developer

547: .keywords: PC, shell, set, apply, user-provided

549: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApplyBA(), PCShellSetApplySymmetricRight(),PCShellSetApplySymmetricLeft()
550: @*/
551: PetscErrorCode  PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
552: {

557:   PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
558:   return(0);
559: }

561: /*@C
562:    PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the PC_SYMMETRIC is used).

564:    Logically Collective on PC

566:    Input Parameters:
567: +  pc - the preconditioner context
568: -  apply - the application-provided left preconditioning routine

570:    Calling sequence of apply:
571: .vb
572:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
573: .ve

575: +  pc - the preconditioner, get the application context with PCShellGetContext()
576: .  xin - input vector
577: -  xout - output vector

579:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

581:    Level: developer

583: .keywords: PC, shell, set, apply, user-provided

585: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
586: @*/
587: PetscErrorCode  PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
588: {

593:   PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
594:   return(0);
595: }

597: /*@C
598:    PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the PC_SYMMETRIC is used).

600:    Logically Collective on PC

602:    Input Parameters:
603: +  pc - the preconditioner context
604: -  apply - the application-provided right preconditioning routine

606:    Calling sequence of apply:
607: .vb
608:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
609: .ve

611: +  pc - the preconditioner, get the application context with PCShellGetContext()
612: .  xin - input vector
613: -  xout - output vector

615:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

617:    Level: developer

619: .keywords: PC, shell, set, apply, user-provided

621: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
622: @*/
623: PetscErrorCode  PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
624: {

629:   PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
630:   return(0);
631: }

633: /*@C
634:    PCShellSetApplyBA - Sets routine to use as preconditioner times operator.

636:    Logically Collective on PC

638:    Input Parameters:
639: +  pc - the preconditioner context
640: -  applyBA - the application-provided BA routine

642:    Calling sequence of apply:
643: .vb
644:    PetscErrorCode applyBA (PC pc,Vec xin,Vec xout)
645: .ve

647: +  pc - the preconditioner, get the application context with PCShellGetContext()
648: .  xin - input vector
649: -  xout - output vector

651:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

653:    Level: developer

655: .keywords: PC, shell, set, apply, user-provided

657: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApply()
658: @*/
659: PetscErrorCode  PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
660: {

665:   PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA));
666:   return(0);
667: }

669: /*@C
670:    PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.

672:    Logically Collective on PC

674:    Input Parameters:
675: +  pc - the preconditioner context
676: -  apply - the application-provided preconditioning transpose routine

678:    Calling sequence of apply:
679: .vb
680:    PetscErrorCode applytranspose (PC pc,Vec xin,Vec xout)
681: .ve

683: +  pc - the preconditioner, get the application context with PCShellGetContext()
684: .  xin - input vector
685: -  xout - output vector

687:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

689:    Level: developer

691:    Notes:
692:    Uses the same context variable as PCShellSetApply().

694: .keywords: PC, shell, set, apply, user-provided

696: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApply(), PCSetContext(), PCShellSetApplyBA()
697: @*/
698: PetscErrorCode  PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
699: {

704:   PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose));
705:   return(0);
706: }

708: /*@C
709:    PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is
710:       applied. This usually does something like scale the linear system in some application
711:       specific way.

713:    Logically Collective on PC

715:    Input Parameters:
716: +  pc - the preconditioner context
717: -  presolve - the application-provided presolve routine

719:    Calling sequence of presolve:
720: .vb
721:    PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x)
722: .ve

724: +  pc - the preconditioner, get the application context with PCShellGetContext()
725: .  xin - input vector
726: -  xout - output vector

728:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

730:    Level: developer

732: .keywords: PC, shell, set, apply, user-provided

734: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPostSolve(), PCShellSetContext()
735: @*/
736: PetscErrorCode  PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
737: {

742:   PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve));
743:   return(0);
744: }

746: /*@C
747:    PCShellSetPostSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is
748:       applied. This usually does something like scale the linear system in some application
749:       specific way.

751:    Logically Collective on PC

753:    Input Parameters:
754: +  pc - the preconditioner context
755: -  postsolve - the application-provided presolve routine

757:    Calling sequence of postsolve:
758: .vb
759:    PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x)
760: .ve

762: +  pc - the preconditioner, get the application context with PCShellGetContext()
763: .  xin - input vector
764: -  xout - output vector

766:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

768:    Level: developer

770: .keywords: PC, shell, set, apply, user-provided

772: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPreSolve(), PCShellSetContext()
773: @*/
774: PetscErrorCode  PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
775: {

780:   PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve));
781:   return(0);
782: }

784: /*@C
785:    PCShellSetName - Sets an optional name to associate with a shell
786:    preconditioner.

788:    Not Collective

790:    Input Parameters:
791: +  pc - the preconditioner context
792: -  name - character string describing shell preconditioner

794:    Level: developer

796: .keywords: PC, shell, set, name, user-provided

798: .seealso: PCShellGetName()
799: @*/
800: PetscErrorCode  PCShellSetName(PC pc,const char name[])
801: {

806:   PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name));
807:   return(0);
808: }

810: /*@C
811:    PCShellGetName - Gets an optional name that the user has set for a shell
812:    preconditioner.

814:    Not Collective

816:    Input Parameter:
817: .  pc - the preconditioner context

819:    Output Parameter:
820: .  name - character string describing shell preconditioner (you should not free this)

822:    Level: developer

824: .keywords: PC, shell, get, name, user-provided

826: .seealso: PCShellSetName()
827: @*/
828: PetscErrorCode  PCShellGetName(PC pc,const char *name[])
829: {

835:   PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name));
836:   return(0);
837: }

839: /*@C
840:    PCShellSetApplyRichardson - Sets routine to use as preconditioner
841:    in Richardson iteration.

843:    Logically Collective on PC

845:    Input Parameters:
846: +  pc - the preconditioner context
847: -  apply - the application-provided preconditioning routine

849:    Calling sequence of apply:
850: .vb
851:    PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
852: .ve

854: +  pc - the preconditioner, get the application context with PCShellGetContext()
855: .  b - right-hand-side
856: .  x - current iterate
857: .  r - work space
858: .  rtol - relative tolerance of residual norm to stop at
859: .  abstol - absolute tolerance of residual norm to stop at
860: .  dtol - if residual norm increases by this factor than return
861: -  maxits - number of iterations to run

863:    Notes: the function MUST return an error code of 0 on success and nonzero on failure.

865:    Level: developer

867: .keywords: PC, shell, set, apply, Richardson, user-provided

869: .seealso: PCShellSetApply(), PCShellSetContext()
870: @*/
871: PetscErrorCode  PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*))
872: {

877:   PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply));
878:   return(0);
879: }

881: /*MC
882:    PCSHELL - Creates a new preconditioner class for use with your
883:               own private data storage format.

885:    Level: advanced
886: >
887:    Concepts: providing your own preconditioner

889:   Usage:
890: $             extern PetscErrorCode apply(PC,Vec,Vec);
891: $             extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec);
892: $             extern PetscErrorCode applytranspose(PC,Vec,Vec);
893: $             extern PetscErrorCode setup(PC);
894: $             extern PetscErrorCode destroy(PC);
895: $
896: $             PCCreate(comm,&pc);
897: $             PCSetType(pc,PCSHELL);
898: $             PCShellSetContext(pc,ctx)
899: $             PCShellSetApply(pc,apply);
900: $             PCShellSetApplyBA(pc,applyba);               (optional)
901: $             PCShellSetApplyTranspose(pc,applytranspose); (optional)
902: $             PCShellSetSetUp(pc,setup);                   (optional)
903: $             PCShellSetDestroy(pc,destroy);               (optional)

905: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
906:            MATSHELL, PCShellSetSetUp(), PCShellSetApply(), PCShellSetView(),
907:            PCShellSetApplyTranspose(), PCShellSetName(), PCShellSetApplyRichardson(),
908:            PCShellGetName(), PCShellSetContext(), PCShellGetContext(), PCShellSetApplyBA()
909: M*/

911: PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
912: {
914:   PC_Shell       *shell;

917:   PetscNewLog(pc,&shell);
918:   pc->data = (void*)shell;

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->applytranspose  = 0;
926:   pc->ops->applyrichardson = 0;
927:   pc->ops->setup           = 0;
928:   pc->ops->presolve        = 0;
929:   pc->ops->postsolve       = 0;

931:   shell->apply          = 0;
932:   shell->applytranspose = 0;
933:   shell->name           = 0;
934:   shell->applyrich      = 0;
935:   shell->presolve       = 0;
936:   shell->postsolve      = 0;
937:   shell->ctx            = 0;
938:   shell->setup          = 0;
939:   shell->view           = 0;
940:   shell->destroy        = 0;
941:   shell->applysymmetricleft  = 0;
942:   shell->applysymmetricright = 0;

944:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",PCShellSetDestroy_Shell);
945:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",PCShellSetSetUp_Shell);
946:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",PCShellSetApply_Shell);
947:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",PCShellSetApplySymmetricLeft_Shell);
948:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",PCShellSetApplySymmetricRight_Shell);
949:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",PCShellSetApplyBA_Shell);
950:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",PCShellSetPreSolve_Shell);
951:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",PCShellSetPostSolve_Shell);
952:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",PCShellSetView_Shell);
953:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",PCShellSetApplyTranspose_Shell);
954:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",PCShellSetName_Shell);
955:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",PCShellGetName_Shell);
956:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",PCShellSetApplyRichardson_Shell);
957:   return(0);
958: }