Actual source code: shellpc.c


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

 25:   char *name;
 26: } PC_Shell;

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

 31:     Not Collective

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

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

 39:     Level: advanced

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

 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.

 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 = NULL;
 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:
 76:     To use this from Fortran you must write a Fortran interface definition for this
 77:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.



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

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

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

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

107: static PetscErrorCode PCApply_Shell(PC pc,Vec x,Vec y)
108: {
109:   PC_Shell         *shell = (PC_Shell*)pc->data;
110:   PetscErrorCode   ierr;
111:   PetscObjectState instate,outstate;

114:   if (!shell->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
115:   PetscObjectStateGet((PetscObject)y, &instate);
116:   PetscStackCall("PCSHELL user function apply()",(*shell->apply)(pc,x,y);CHKERRQ(ierr));
117:   PetscObjectStateGet((PetscObject)y, &outstate);
118:   if (instate == outstate) {
119:     /* increase the state of the output vector since the user did not update its state themselve as should have been done */
120:     PetscObjectStateIncrease((PetscObject)y);
121:   }
122:   return(0);
123: }

125: static PetscErrorCode PCMatApply_Shell(PC pc,Mat X,Mat Y)
126: {
127:   PC_Shell         *shell = (PC_Shell*)pc->data;
128:   PetscErrorCode   ierr;
129:   PetscObjectState instate,outstate;

132:   if (!shell->matapply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
133:   PetscObjectStateGet((PetscObject)Y, &instate);
134:   PetscStackCall("PCSHELL user function apply()",(*shell->matapply)(pc,X,Y);CHKERRQ(ierr));
135:   PetscObjectStateGet((PetscObject)Y, &outstate);
136:   if (instate == outstate) {
137:     /* increase the state of the output vector since the user did not update its state themselve as should have been done */
138:     PetscObjectStateIncrease((PetscObject)Y);
139:   }
140:   return(0);
141: }

143: static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc,Vec x,Vec y)
144: {
145:   PC_Shell       *shell = (PC_Shell*)pc->data;

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

154: static PetscErrorCode PCApplySymmetricRight_Shell(PC pc,Vec x,Vec y)
155: {
156:   PC_Shell       *shell = (PC_Shell*)pc->data;

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

165: static PetscErrorCode PCApplyBA_Shell(PC pc,PCSide side,Vec x,Vec y,Vec w)
166: {
167:   PC_Shell         *shell = (PC_Shell*)pc->data;
168:   PetscErrorCode   ierr;
169:   PetscObjectState instate,outstate;

172:   if (!shell->applyBA) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyBA() routine provided to Shell PC");
173:   PetscObjectStateGet((PetscObject)w, &instate);
174:   PetscStackCall("PCSHELL user function applyBA()",(*shell->applyBA)(pc,side,x,y,w);CHKERRQ(ierr));
175:   PetscObjectStateGet((PetscObject)w, &outstate);
176:   if (instate == outstate) {
177:     /* increase the state of the output vector since the user did not update its state themselve as should have been done */
178:     PetscObjectStateIncrease((PetscObject)w);
179:   }
180:   return(0);
181: }

183: static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc,PetscBool* change)
184: {
186:   *change = PETSC_TRUE;
187:   return(0);
188: }

190: static PetscErrorCode PCPreSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)
191: {
192:   PC_Shell       *shell = (PC_Shell*)pc->data;

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

201: static PetscErrorCode PCPostSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)
202: {
203:   PC_Shell       *shell = (PC_Shell*)pc->data;

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

212: static PetscErrorCode PCApplyTranspose_Shell(PC pc,Vec x,Vec y)
213: {
214:   PC_Shell         *shell = (PC_Shell*)pc->data;
215:   PetscErrorCode   ierr;
216:   PetscObjectState instate,outstate;

219:   if (!shell->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applytranspose() routine provided to Shell PC");
220:   PetscObjectStateGet((PetscObject)y, &instate);
221:   PetscStackCall("PCSHELL user function applytranspose()",(*shell->applytranspose)(pc,x,y);CHKERRQ(ierr));
222:   PetscObjectStateGet((PetscObject)y, &outstate);
223:   if (instate == outstate) {
224:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
225:     PetscObjectStateIncrease((PetscObject)y);
226:   }
227:   return(0);
228: }

230: 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)
231: {
232:   PetscErrorCode   ierr;
233:   PC_Shell         *shell = (PC_Shell*)pc->data;
234:   PetscObjectState instate,outstate;

237:   if (!shell->applyrich) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyrichardson() routine provided to Shell PC");
238:   PetscObjectStateGet((PetscObject)y, &instate);
239:   PetscStackCall("PCSHELL user function applyrichardson()",(*shell->applyrich)(pc,x,y,w,rtol,abstol,dtol,it,guesszero,outits,reason);CHKERRQ(ierr));
240:   PetscObjectStateGet((PetscObject)y, &outstate);
241:   if (instate == outstate) {
242:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
243:     PetscObjectStateIncrease((PetscObject)y);
244:   }
245:   return(0);
246: }

248: static PetscErrorCode PCDestroy_Shell(PC pc)
249: {
250:   PC_Shell       *shell = (PC_Shell*)pc->data;

254:   PetscFree(shell->name);
255:   if (shell->destroy) PetscStackCall("PCSHELL user function destroy()",(*shell->destroy)(pc);CHKERRQ(ierr));
256:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",NULL);
257:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",NULL);
258:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",NULL);
259:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",NULL);
260:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",NULL);
261:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",NULL);
262:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",NULL);
263:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",NULL);
264:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",NULL);
265:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",NULL);
266:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",NULL);
267:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",NULL);
268:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",NULL);
269:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",NULL);
270:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
271:   PetscFree(pc->data);
272:   return(0);
273: }

275: static PetscErrorCode PCView_Shell(PC pc,PetscViewer viewer)
276: {
277:   PC_Shell       *shell = (PC_Shell*)pc->data;
279:   PetscBool      iascii;

282:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
283:   if (iascii) {
284:     if (shell->name) {
285:       PetscViewerASCIIPrintf(viewer,"  %s\n",shell->name);
286:     } else {
287:       PetscViewerASCIIPrintf(viewer,"  no name\n");
288:     }
289:   }
290:   if (shell->view) {
291:     PetscViewerASCIIPushTab(viewer);
292:     (*shell->view)(pc,viewer);
293:     PetscViewerASCIIPopTab(viewer);
294:   }
295:   return(0);
296: }

298: /* ------------------------------------------------------------------------------*/
299: static PetscErrorCode  PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC))
300: {
301:   PC_Shell *shell= (PC_Shell*)pc->data;

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

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

313:   shell->setup = setup;
314:   if (setup) pc->ops->setup = PCSetUp_Shell;
315:   else       pc->ops->setup = NULL;
316:   return(0);
317: }

319: static PetscErrorCode  PCShellSetApply_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
320: {
321:   PC_Shell *shell = (PC_Shell*)pc->data;

324:   shell->apply = apply;
325:   return(0);
326: }

328: static PetscErrorCode  PCShellSetMatApply_Shell(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
329: {
330:   PC_Shell *shell = (PC_Shell*)pc->data;

333:   shell->matapply = matapply;
334:   return(0);
335: }

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

342:   shell->applysymmetricleft = apply;
343:   return(0);
344: }

346: static PetscErrorCode  PCShellSetApplySymmetricRight_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
347: {
348:   PC_Shell *shell = (PC_Shell*)pc->data;

351:   shell->applysymmetricright = apply;
352:   return(0);
353: }

355: static PetscErrorCode  PCShellSetApplyBA_Shell(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
356: {
357:   PC_Shell *shell = (PC_Shell*)pc->data;

360:   shell->applyBA = applyBA;
361:   if (applyBA) pc->ops->applyBA  = PCApplyBA_Shell;
362:   else         pc->ops->applyBA  = NULL;
363:   return(0);
364: }

366: static PetscErrorCode  PCShellSetPreSolve_Shell(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
367: {
368:   PC_Shell       *shell = (PC_Shell*)pc->data;

372:   shell->presolve = presolve;
373:   if (presolve) {
374:     pc->ops->presolve = PCPreSolve_Shell;
375:     PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Shell);
376:   } else {
377:     pc->ops->presolve = NULL;
378:     PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
379:   }
380:   return(0);
381: }

383: static PetscErrorCode  PCShellSetPostSolve_Shell(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
384: {
385:   PC_Shell *shell = (PC_Shell*)pc->data;

388:   shell->postsolve = postsolve;
389:   if (postsolve) pc->ops->postsolve = PCPostSolve_Shell;
390:   else           pc->ops->postsolve = NULL;
391:   return(0);
392: }

394: static PetscErrorCode  PCShellSetView_Shell(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
395: {
396:   PC_Shell *shell = (PC_Shell*)pc->data;

399:   shell->view = view;
400:   return(0);
401: }

403: static PetscErrorCode  PCShellSetApplyTranspose_Shell(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
404: {
405:   PC_Shell *shell = (PC_Shell*)pc->data;

408:   shell->applytranspose = applytranspose;
409:   if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell;
410:   else                pc->ops->applytranspose = NULL;
411:   return(0);
412: }

414: static PetscErrorCode  PCShellSetApplyRichardson_Shell(PC pc,PetscErrorCode (*applyrich)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*))
415: {
416:   PC_Shell *shell = (PC_Shell*)pc->data;

419:   shell->applyrich = applyrich;
420:   if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell;
421:   else           pc->ops->applyrichardson = NULL;
422:   return(0);
423: }

425: static PetscErrorCode  PCShellSetName_Shell(PC pc,const char name[])
426: {
427:   PC_Shell       *shell = (PC_Shell*)pc->data;

431:   PetscFree(shell->name);
432:   PetscStrallocpy(name,&shell->name);
433:   return(0);
434: }

436: static PetscErrorCode  PCShellGetName_Shell(PC pc,const char *name[])
437: {
438:   PC_Shell *shell = (PC_Shell*)pc->data;

441:   *name = shell->name;
442:   return(0);
443: }

445: /* -------------------------------------------------------------------------------*/

447: /*@C
448:    PCShellSetDestroy - Sets routine to use to destroy the user-provided
449:    application context.

451:    Logically Collective on PC

453:    Input Parameters:
454: +  pc - the preconditioner context
455: -  destroy - the application-provided destroy routine

457:    Calling sequence of destroy:
458: .vb
459:    PetscErrorCode destroy (PC)
460: .ve

462: .  ptr - the application context

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

467:    Level: developer

469: .seealso: PCShellSetApply(), PCShellSetContext()
470: @*/
471: PetscErrorCode  PCShellSetDestroy(PC pc,PetscErrorCode (*destroy)(PC))
472: {

477:   PetscTryMethod(pc,"PCShellSetDestroy_C",(PC,PetscErrorCode (*)(PC)),(pc,destroy));
478:   return(0);
479: }


482: /*@C
483:    PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the
484:    matrix operator is changed.

486:    Logically Collective on PC

488:    Input Parameters:
489: +  pc - the preconditioner context
490: -  setup - the application-provided setup routine

492:    Calling sequence of setup:
493: .vb
494:    PetscErrorCode setup (PC pc)
495: .ve

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

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

502:    Level: developer

504: .seealso: PCShellSetApplyRichardson(), PCShellSetApply(), PCShellSetContext()
505: @*/
506: PetscErrorCode  PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC))
507: {

512:   PetscTryMethod(pc,"PCShellSetSetUp_C",(PC,PetscErrorCode (*)(PC)),(pc,setup));
513:   return(0);
514: }


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

520:    Logically Collective on PC

522:    Input Parameters:
523: +  pc - the preconditioner context
524: -  view - the application-provided view routine

526:    Calling sequence of view:
527: .vb
528:    PetscErrorCode view(PC pc,PetscViewer v)
529: .ve

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

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

537:    Level: developer

539: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose()
540: @*/
541: PetscErrorCode  PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
542: {

547:   PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view));
548:   return(0);
549: }

551: /*@C
552:    PCShellSetApply - Sets routine to use as preconditioner.

554:    Logically Collective on PC

556:    Input Parameters:
557: +  pc - the preconditioner context
558: -  apply - the application-provided preconditioning routine

560:    Calling sequence of apply:
561: .vb
562:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
563: .ve

565: +  pc - the preconditioner, get the application context with PCShellGetContext()
566: .  xin - input vector
567: -  xout - output vector

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

572:    Level: developer

574: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApplyBA(), PCShellSetApplySymmetricRight(),PCShellSetApplySymmetricLeft()
575: @*/
576: PetscErrorCode  PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
577: {

582:   PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
583:   return(0);
584: }

586: /*@C
587:    PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors.

589:    Logically Collective on PC

591:    Input Parameters:
592: +  pc - the preconditioner context
593: -  apply - the application-provided preconditioning routine

595:    Calling sequence of apply:
596: .vb
597:    PetscErrorCode apply (PC pc,Mat Xin,Mat Xout)
598: .ve

600: +  pc - the preconditioner, get the application context with PCShellGetContext()
601: .  Xin - input block of vectors
602: -  Xout - output block of vectors

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

607:    Level: developer

609: .seealso: PCShellSetApply()
610: @*/
611: PetscErrorCode  PCShellSetMatApply(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
612: {

617:   PetscTryMethod(pc,"PCShellSetMatApply_C",(PC,PetscErrorCode (*)(PC,Mat,Mat)),(pc,matapply));
618:   return(0);
619: }

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

624:    Logically Collective on PC

626:    Input Parameters:
627: +  pc - the preconditioner context
628: -  apply - the application-provided left preconditioning routine

630:    Calling sequence of apply:
631: .vb
632:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
633: .ve

635: +  pc - the preconditioner, get the application context with PCShellGetContext()
636: .  xin - input vector
637: -  xout - output vector

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

642:    Level: developer

644: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
645: @*/
646: PetscErrorCode  PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
647: {

652:   PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
653:   return(0);
654: }

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

659:    Logically Collective on PC

661:    Input Parameters:
662: +  pc - the preconditioner context
663: -  apply - the application-provided right preconditioning routine

665:    Calling sequence of apply:
666: .vb
667:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
668: .ve

670: +  pc - the preconditioner, get the application context with PCShellGetContext()
671: .  xin - input vector
672: -  xout - output vector

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

677:    Level: developer

679: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
680: @*/
681: PetscErrorCode  PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
682: {

687:   PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
688:   return(0);
689: }

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

694:    Logically Collective on PC

696:    Input Parameters:
697: +  pc - the preconditioner context
698: -  applyBA - the application-provided BA routine

700:    Calling sequence of applyBA:
701: .vb
702:    PetscErrorCode applyBA (PC pc,Vec xin,Vec xout)
703: .ve

705: +  pc - the preconditioner, get the application context with PCShellGetContext()
706: .  xin - input vector
707: -  xout - output vector

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

712:    Level: developer

714: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApply()
715: @*/
716: PetscErrorCode  PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
717: {

722:   PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA));
723:   return(0);
724: }

726: /*@C
727:    PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.

729:    Logically Collective on PC

731:    Input Parameters:
732: +  pc - the preconditioner context
733: -  apply - the application-provided preconditioning transpose routine

735:    Calling sequence of apply:
736: .vb
737:    PetscErrorCode applytranspose (PC pc,Vec xin,Vec xout)
738: .ve

740: +  pc - the preconditioner, get the application context with PCShellGetContext()
741: .  xin - input vector
742: -  xout - output vector

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

747:    Level: developer

749:    Notes:
750:    Uses the same context variable as PCShellSetApply().

752: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApply(), PCSetContext(), PCShellSetApplyBA()
753: @*/
754: PetscErrorCode  PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
755: {

760:   PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose));
761:   return(0);
762: }

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

769:    Logically Collective on PC

771:    Input Parameters:
772: +  pc - the preconditioner context
773: -  presolve - the application-provided presolve routine

775:    Calling sequence of presolve:
776: .vb
777:    PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x)
778: .ve

780: +  pc - the preconditioner, get the application context with PCShellGetContext()
781: .  xin - input vector
782: -  xout - output vector

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

787:    Level: developer

789: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPostSolve(), PCShellSetContext()
790: @*/
791: PetscErrorCode  PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
792: {

797:   PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve));
798:   return(0);
799: }

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

806:    Logically Collective on PC

808:    Input Parameters:
809: +  pc - the preconditioner context
810: -  postsolve - the application-provided presolve routine

812:    Calling sequence of postsolve:
813: .vb
814:    PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x)
815: .ve

817: +  pc - the preconditioner, get the application context with PCShellGetContext()
818: .  xin - input vector
819: -  xout - output vector

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

824:    Level: developer

826: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPreSolve(), PCShellSetContext()
827: @*/
828: PetscErrorCode  PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
829: {

834:   PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve));
835:   return(0);
836: }

838: /*@C
839:    PCShellSetName - Sets an optional name to associate with a shell
840:    preconditioner.

842:    Not Collective

844:    Input Parameters:
845: +  pc - the preconditioner context
846: -  name - character string describing shell preconditioner

848:    Level: developer

850: .seealso: PCShellGetName()
851: @*/
852: PetscErrorCode  PCShellSetName(PC pc,const char name[])
853: {

858:   PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name));
859:   return(0);
860: }

862: /*@C
863:    PCShellGetName - Gets an optional name that the user has set for a shell
864:    preconditioner.

866:    Not Collective

868:    Input Parameter:
869: .  pc - the preconditioner context

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

874:    Level: developer

876: .seealso: PCShellSetName()
877: @*/
878: PetscErrorCode  PCShellGetName(PC pc,const char *name[])
879: {

885:   PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name));
886:   return(0);
887: }

889: /*@C
890:    PCShellSetApplyRichardson - Sets routine to use as preconditioner
891:    in Richardson iteration.

893:    Logically Collective on PC

895:    Input Parameters:
896: +  pc - the preconditioner context
897: -  apply - the application-provided preconditioning routine

899:    Calling sequence of apply:
900: .vb
901:    PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
902: .ve

904: +  pc - the preconditioner, get the application context with PCShellGetContext()
905: .  b - right-hand-side
906: .  x - current iterate
907: .  r - work space
908: .  rtol - relative tolerance of residual norm to stop at
909: .  abstol - absolute tolerance of residual norm to stop at
910: .  dtol - if residual norm increases by this factor than return
911: -  maxits - number of iterations to run

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

916:    Level: developer

918: .seealso: PCShellSetApply(), PCShellSetContext()
919: @*/
920: PetscErrorCode  PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*))
921: {

926:   PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply));
927:   return(0);
928: }

930: /*MC
931:    PCSHELL - Creates a new preconditioner class for use with your
932:               own private data storage format.

934:    Level: advanced

936:   Usage:
937: $             extern PetscErrorCode apply(PC,Vec,Vec);
938: $             extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec);
939: $             extern PetscErrorCode applytranspose(PC,Vec,Vec);
940: $             extern PetscErrorCode setup(PC);
941: $             extern PetscErrorCode destroy(PC);
942: $
943: $             PCCreate(comm,&pc);
944: $             PCSetType(pc,PCSHELL);
945: $             PCShellSetContext(pc,ctx)
946: $             PCShellSetApply(pc,apply);
947: $             PCShellSetApplyBA(pc,applyba);               (optional)
948: $             PCShellSetApplyTranspose(pc,applytranspose); (optional)
949: $             PCShellSetSetUp(pc,setup);                   (optional)
950: $             PCShellSetDestroy(pc,destroy);               (optional)

952: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
953:            MATSHELL, PCShellSetSetUp(), PCShellSetApply(), PCShellSetView(),
954:            PCShellSetApplyTranspose(), PCShellSetName(), PCShellSetApplyRichardson(),
955:            PCShellGetName(), PCShellSetContext(), PCShellGetContext(), PCShellSetApplyBA()
956: M*/

958: PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
959: {
961:   PC_Shell       *shell;

964:   PetscNewLog(pc,&shell);
965:   pc->data = (void*)shell;

967:   pc->ops->destroy         = PCDestroy_Shell;
968:   pc->ops->view            = PCView_Shell;
969:   pc->ops->apply           = PCApply_Shell;
970:   pc->ops->matapply        = PCMatApply_Shell;
971:   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Shell;
972:   pc->ops->applysymmetricright = PCApplySymmetricRight_Shell;
973:   pc->ops->applytranspose  = NULL;
974:   pc->ops->applyrichardson = NULL;
975:   pc->ops->setup           = NULL;
976:   pc->ops->presolve        = NULL;
977:   pc->ops->postsolve       = NULL;

979:   shell->apply          = NULL;
980:   shell->applytranspose = NULL;
981:   shell->name           = NULL;
982:   shell->applyrich      = NULL;
983:   shell->presolve       = NULL;
984:   shell->postsolve      = NULL;
985:   shell->ctx            = NULL;
986:   shell->setup          = NULL;
987:   shell->view           = NULL;
988:   shell->destroy        = NULL;
989:   shell->applysymmetricleft  = NULL;
990:   shell->applysymmetricright = NULL;

992:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",PCShellSetDestroy_Shell);
993:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",PCShellSetSetUp_Shell);
994:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",PCShellSetApply_Shell);
995:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",PCShellSetMatApply_Shell);
996:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",PCShellSetApplySymmetricLeft_Shell);
997:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",PCShellSetApplySymmetricRight_Shell);
998:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",PCShellSetApplyBA_Shell);
999:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",PCShellSetPreSolve_Shell);
1000:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",PCShellSetPostSolve_Shell);
1001:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",PCShellSetView_Shell);
1002:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",PCShellSetApplyTranspose_Shell);
1003:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",PCShellSetName_Shell);
1004:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",PCShellGetName_Shell);
1005:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",PCShellSetApplyRichardson_Shell);
1006:   return(0);
1007: }