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) *(void**)ctx = NULL;
 60:   else      *(void**)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.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

181: static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc,PetscBool* change)
182: {
184:   *change = PETSC_TRUE;
185:   return(0);
186: }

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

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

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

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

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

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

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

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

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

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

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

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

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

302:   shell->destroy = destroy;
303:   return(0);
304: }

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

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

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

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

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

331:   shell->matapply = matapply;
332:   return(0);
333: }

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

340:   shell->applysymmetricleft = apply;
341:   return(0);
342: }

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

349:   shell->applysymmetricright = apply;
350:   return(0);
351: }

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

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

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

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

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

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

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

397:   shell->view = view;
398:   return(0);
399: }

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

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

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

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

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

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

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

439:   *name = shell->name;
440:   return(0);
441: }

443: /* -------------------------------------------------------------------------------*/

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

449:    Logically Collective on PC

451:    Input Parameters:
452: +  pc - the preconditioner context
453: -  destroy - the application-provided destroy routine

455:    Calling sequence of destroy:
456: .vb
457:    PetscErrorCode destroy (PC)
458: .ve

460: .  ptr - the application context

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

465:    Level: developer

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

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

479: /*@C
480:    PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the
481:    matrix operator is changed.

483:    Logically Collective on PC

485:    Input Parameters:
486: +  pc - the preconditioner context
487: -  setup - the application-provided setup routine

489:    Calling sequence of setup:
490: .vb
491:    PetscErrorCode setup (PC pc)
492: .ve

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

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

499:    Level: developer

501: .seealso: PCShellSetApplyRichardson(), PCShellSetApply(), PCShellSetContext()
502: @*/
503: PetscErrorCode  PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC))
504: {

509:   PetscTryMethod(pc,"PCShellSetSetUp_C",(PC,PetscErrorCode (*)(PC)),(pc,setup));
510:   return(0);
511: }

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

516:    Logically Collective on PC

518:    Input Parameters:
519: +  pc - the preconditioner context
520: -  view - the application-provided view routine

522:    Calling sequence of view:
523: .vb
524:    PetscErrorCode view(PC pc,PetscViewer v)
525: .ve

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

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

533:    Level: developer

535: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose()
536: @*/
537: PetscErrorCode  PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
538: {

543:   PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view));
544:   return(0);
545: }

547: /*@C
548:    PCShellSetApply - Sets routine to use as preconditioner.

550:    Logically Collective on PC

552:    Input Parameters:
553: +  pc - the preconditioner context
554: -  apply - the application-provided preconditioning routine

556:    Calling sequence of apply:
557: .vb
558:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
559: .ve

561: +  pc - the preconditioner, get the application context with PCShellGetContext()
562: .  xin - input vector
563: -  xout - output vector

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

568:    Level: developer

570: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApplyBA(), PCShellSetApplySymmetricRight(),PCShellSetApplySymmetricLeft()
571: @*/
572: PetscErrorCode  PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
573: {

578:   PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
579:   return(0);
580: }

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

585:    Logically Collective on PC

587:    Input Parameters:
588: +  pc - the preconditioner context
589: -  apply - the application-provided preconditioning routine

591:    Calling sequence of apply:
592: .vb
593:    PetscErrorCode apply (PC pc,Mat Xin,Mat Xout)
594: .ve

596: +  pc - the preconditioner, get the application context with PCShellGetContext()
597: .  Xin - input block of vectors
598: -  Xout - output block of vectors

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

603:    Level: developer

605: .seealso: PCShellSetApply()
606: @*/
607: PetscErrorCode  PCShellSetMatApply(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
608: {

613:   PetscTryMethod(pc,"PCShellSetMatApply_C",(PC,PetscErrorCode (*)(PC,Mat,Mat)),(pc,matapply));
614:   return(0);
615: }

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

620:    Logically Collective on PC

622:    Input Parameters:
623: +  pc - the preconditioner context
624: -  apply - the application-provided left preconditioning routine

626:    Calling sequence of apply:
627: .vb
628:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
629: .ve

631: +  pc - the preconditioner, get the application context with PCShellGetContext()
632: .  xin - input vector
633: -  xout - output vector

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

638:    Level: developer

640: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
641: @*/
642: PetscErrorCode  PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
643: {

648:   PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
649:   return(0);
650: }

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

655:    Logically Collective on PC

657:    Input Parameters:
658: +  pc - the preconditioner context
659: -  apply - the application-provided right preconditioning routine

661:    Calling sequence of apply:
662: .vb
663:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
664: .ve

666: +  pc - the preconditioner, get the application context with PCShellGetContext()
667: .  xin - input vector
668: -  xout - output vector

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

673:    Level: developer

675: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
676: @*/
677: PetscErrorCode  PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
678: {

683:   PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
684:   return(0);
685: }

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

690:    Logically Collective on PC

692:    Input Parameters:
693: +  pc - the preconditioner context
694: -  applyBA - the application-provided BA routine

696:    Calling sequence of applyBA:
697: .vb
698:    PetscErrorCode applyBA (PC pc,Vec xin,Vec xout)
699: .ve

701: +  pc - the preconditioner, get the application context with PCShellGetContext()
702: .  xin - input vector
703: -  xout - output vector

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

708:    Level: developer

710: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApply()
711: @*/
712: PetscErrorCode  PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
713: {

718:   PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA));
719:   return(0);
720: }

722: /*@C
723:    PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.

725:    Logically Collective on PC

727:    Input Parameters:
728: +  pc - the preconditioner context
729: -  apply - the application-provided preconditioning transpose routine

731:    Calling sequence of apply:
732: .vb
733:    PetscErrorCode applytranspose (PC pc,Vec xin,Vec xout)
734: .ve

736: +  pc - the preconditioner, get the application context with PCShellGetContext()
737: .  xin - input vector
738: -  xout - output vector

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

743:    Level: developer

745:    Notes:
746:    Uses the same context variable as PCShellSetApply().

748: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApply(), PCSetContext(), PCShellSetApplyBA()
749: @*/
750: PetscErrorCode  PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
751: {

756:   PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose));
757:   return(0);
758: }

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

765:    Logically Collective on PC

767:    Input Parameters:
768: +  pc - the preconditioner context
769: -  presolve - the application-provided presolve routine

771:    Calling sequence of presolve:
772: .vb
773:    PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x)
774: .ve

776: +  pc - the preconditioner, get the application context with PCShellGetContext()
777: .  xin - input vector
778: -  xout - output vector

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

783:    Level: developer

785: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPostSolve(), PCShellSetContext()
786: @*/
787: PetscErrorCode  PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
788: {

793:   PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve));
794:   return(0);
795: }

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

802:    Logically Collective on PC

804:    Input Parameters:
805: +  pc - the preconditioner context
806: -  postsolve - the application-provided presolve routine

808:    Calling sequence of postsolve:
809: .vb
810:    PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x)
811: .ve

813: +  pc - the preconditioner, get the application context with PCShellGetContext()
814: .  xin - input vector
815: -  xout - output vector

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

820:    Level: developer

822: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPreSolve(), PCShellSetContext()
823: @*/
824: PetscErrorCode  PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
825: {

830:   PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve));
831:   return(0);
832: }

834: /*@C
835:    PCShellSetName - Sets an optional name to associate with a shell
836:    preconditioner.

838:    Not Collective

840:    Input Parameters:
841: +  pc - the preconditioner context
842: -  name - character string describing shell preconditioner

844:    Level: developer

846: .seealso: PCShellGetName()
847: @*/
848: PetscErrorCode  PCShellSetName(PC pc,const char name[])
849: {

854:   PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name));
855:   return(0);
856: }

858: /*@C
859:    PCShellGetName - Gets an optional name that the user has set for a shell
860:    preconditioner.

862:    Not Collective

864:    Input Parameter:
865: .  pc - the preconditioner context

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

870:    Level: developer

872: .seealso: PCShellSetName()
873: @*/
874: PetscErrorCode  PCShellGetName(PC pc,const char *name[])
875: {

881:   PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name));
882:   return(0);
883: }

885: /*@C
886:    PCShellSetApplyRichardson - Sets routine to use as preconditioner
887:    in Richardson iteration.

889:    Logically Collective on PC

891:    Input Parameters:
892: +  pc - the preconditioner context
893: -  apply - the application-provided preconditioning routine

895:    Calling sequence of apply:
896: .vb
897:    PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
898: .ve

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

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

912:    Level: developer

914: .seealso: PCShellSetApply(), PCShellSetContext()
915: @*/
916: PetscErrorCode  PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*))
917: {

922:   PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply));
923:   return(0);
924: }

926: /*MC
927:    PCSHELL - Creates a new preconditioner class for use with your
928:               own private data storage format.

930:    Level: advanced

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

948: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
949:            MATSHELL, PCShellSetSetUp(), PCShellSetApply(), PCShellSetView(),
950:            PCShellSetApplyTranspose(), PCShellSetName(), PCShellSetApplyRichardson(),
951:            PCShellGetName(), PCShellSetContext(), PCShellGetContext(), PCShellSetApplyBA()
952: M*/

954: PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
955: {
957:   PC_Shell       *shell;

960:   PetscNewLog(pc,&shell);
961:   pc->data = (void*)shell;

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

975:   shell->apply          = NULL;
976:   shell->applytranspose = NULL;
977:   shell->name           = NULL;
978:   shell->applyrich      = NULL;
979:   shell->presolve       = NULL;
980:   shell->postsolve      = NULL;
981:   shell->ctx            = NULL;
982:   shell->setup          = NULL;
983:   shell->view           = NULL;
984:   shell->destroy        = NULL;
985:   shell->applysymmetricleft  = NULL;
986:   shell->applysymmetricright = NULL;

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