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: {
 52:   PetscBool      flg;

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

 62: /*@
 63:     PCShellSetContext - sets the context for a shell PC

 65:    Logically Collective on PC

 67:     Input Parameters:
 68: +   pc - the shell PC
 69: -   ctx - the context

 71:    Level: advanced

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

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

 85:   PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&flg);
 86:   if (flg) shell->ctx = ctx;
 87:   return 0;
 88: }

 90: static PetscErrorCode PCSetUp_Shell(PC pc)
 91: {
 92:   PC_Shell       *shell = (PC_Shell*)pc->data;

 95:   PetscStackCall("PCSHELL user function setup()",(*shell->setup)(pc));
 96:   return 0;
 97: }

 99: static PetscErrorCode PCApply_Shell(PC pc,Vec x,Vec y)
100: {
101:   PC_Shell         *shell = (PC_Shell*)pc->data;
102:   PetscObjectState instate,outstate;

105:   PetscObjectStateGet((PetscObject)y, &instate);
106:   PetscStackCall("PCSHELL user function apply()",(*shell->apply)(pc,x,y));
107:   PetscObjectStateGet((PetscObject)y, &outstate);
108:   if (instate == outstate) {
109:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
110:     PetscObjectStateIncrease((PetscObject)y);
111:   }
112:   return 0;
113: }

115: static PetscErrorCode PCMatApply_Shell(PC pc,Mat X,Mat Y)
116: {
117:   PC_Shell         *shell = (PC_Shell*)pc->data;
118:   PetscObjectState instate,outstate;

121:   PetscObjectStateGet((PetscObject)Y, &instate);
122:   PetscStackCall("PCSHELL user function apply()",(*shell->matapply)(pc,X,Y));
123:   PetscObjectStateGet((PetscObject)Y, &outstate);
124:   if (instate == outstate) {
125:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
126:     PetscObjectStateIncrease((PetscObject)Y);
127:   }
128:   return 0;
129: }

131: static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc,Vec x,Vec y)
132: {
133:   PC_Shell       *shell = (PC_Shell*)pc->data;

136:   PetscStackCall("PCSHELL user function apply()",(*shell->applysymmetricleft)(pc,x,y));
137:   return 0;
138: }

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

145:   PetscStackCall("PCSHELL user function apply()",(*shell->applysymmetricright)(pc,x,y));
146:   return 0;
147: }

149: static PetscErrorCode PCApplyBA_Shell(PC pc,PCSide side,Vec x,Vec y,Vec w)
150: {
151:   PC_Shell         *shell = (PC_Shell*)pc->data;
152:   PetscObjectState instate,outstate;

155:   PetscObjectStateGet((PetscObject)w, &instate);
156:   PetscStackCall("PCSHELL user function applyBA()",(*shell->applyBA)(pc,side,x,y,w));
157:   PetscObjectStateGet((PetscObject)w, &outstate);
158:   if (instate == outstate) {
159:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
160:     PetscObjectStateIncrease((PetscObject)w);
161:   }
162:   return 0;
163: }

165: static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc,PetscBool* change)
166: {
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;

176:   PetscStackCall("PCSHELL user function presolve()",(*shell->presolve)(pc,ksp,b,x));
177:   return 0;
178: }

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

185:   PetscStackCall("PCSHELL user function postsolve()",(*shell->postsolve)(pc,ksp,b,x));
186:   return 0;
187: }

189: static PetscErrorCode PCApplyTranspose_Shell(PC pc,Vec x,Vec y)
190: {
191:   PC_Shell         *shell = (PC_Shell*)pc->data;
192:   PetscObjectState instate,outstate;

195:   PetscObjectStateGet((PetscObject)y, &instate);
196:   PetscStackCall("PCSHELL user function applytranspose()",(*shell->applytranspose)(pc,x,y));
197:   PetscObjectStateGet((PetscObject)y, &outstate);
198:   if (instate == outstate) {
199:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
200:     PetscObjectStateIncrease((PetscObject)y);
201:   }
202:   return 0;
203: }

205: 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)
206: {
207:   PC_Shell         *shell = (PC_Shell*)pc->data;
208:   PetscObjectState instate,outstate;

211:   PetscObjectStateGet((PetscObject)y, &instate);
212:   PetscStackCall("PCSHELL user function applyrichardson()",(*shell->applyrich)(pc,x,y,w,rtol,abstol,dtol,it,guesszero,outits,reason));
213:   PetscObjectStateGet((PetscObject)y, &outstate);
214:   if (instate == outstate) {
215:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
216:     PetscObjectStateIncrease((PetscObject)y);
217:   }
218:   return 0;
219: }

221: static PetscErrorCode PCDestroy_Shell(PC pc)
222: {
223:   PC_Shell       *shell = (PC_Shell*)pc->data;

225:   PetscFree(shell->name);
226:   if (shell->destroy) PetscStackCall("PCSHELL user function destroy()",(*shell->destroy)(pc));
227:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",NULL);
228:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",NULL);
229:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",NULL);
230:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",NULL);
231:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",NULL);
232:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",NULL);
233:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",NULL);
234:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",NULL);
235:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",NULL);
236:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",NULL);
237:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",NULL);
238:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",NULL);
239:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",NULL);
240:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",NULL);
241:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
242:   PetscFree(pc->data);
243:   return 0;
244: }

246: static PetscErrorCode PCView_Shell(PC pc,PetscViewer viewer)
247: {
248:   PC_Shell       *shell = (PC_Shell*)pc->data;
249:   PetscBool      iascii;

251:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
252:   if (iascii) {
253:     if (shell->name) {
254:       PetscViewerASCIIPrintf(viewer,"  %s\n",shell->name);
255:     } else {
256:       PetscViewerASCIIPrintf(viewer,"  no name\n");
257:     }
258:   }
259:   if (shell->view) {
260:     PetscViewerASCIIPushTab(viewer);
261:     (*shell->view)(pc,viewer);
262:     PetscViewerASCIIPopTab(viewer);
263:   }
264:   return 0;
265: }

267: /* ------------------------------------------------------------------------------*/
268: static PetscErrorCode  PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC))
269: {
270:   PC_Shell *shell= (PC_Shell*)pc->data;

272:   shell->destroy = destroy;
273:   return 0;
274: }

276: static PetscErrorCode  PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC))
277: {
278:   PC_Shell *shell = (PC_Shell*)pc->data;

280:   shell->setup = setup;
281:   if (setup) pc->ops->setup = PCSetUp_Shell;
282:   else       pc->ops->setup = NULL;
283:   return 0;
284: }

286: static PetscErrorCode  PCShellSetApply_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
287: {
288:   PC_Shell *shell = (PC_Shell*)pc->data;

290:   shell->apply = apply;
291:   return 0;
292: }

294: static PetscErrorCode  PCShellSetMatApply_Shell(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
295: {
296:   PC_Shell *shell = (PC_Shell*)pc->data;

298:   shell->matapply = matapply;
299:   return 0;
300: }

302: static PetscErrorCode  PCShellSetApplySymmetricLeft_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
303: {
304:   PC_Shell *shell = (PC_Shell*)pc->data;

306:   shell->applysymmetricleft = apply;
307:   return 0;
308: }

310: static PetscErrorCode  PCShellSetApplySymmetricRight_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
311: {
312:   PC_Shell *shell = (PC_Shell*)pc->data;

314:   shell->applysymmetricright = apply;
315:   return 0;
316: }

318: static PetscErrorCode  PCShellSetApplyBA_Shell(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
319: {
320:   PC_Shell *shell = (PC_Shell*)pc->data;

322:   shell->applyBA = applyBA;
323:   if (applyBA) pc->ops->applyBA  = PCApplyBA_Shell;
324:   else         pc->ops->applyBA  = NULL;
325:   return 0;
326: }

328: static PetscErrorCode  PCShellSetPreSolve_Shell(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
329: {
330:   PC_Shell       *shell = (PC_Shell*)pc->data;

332:   shell->presolve = presolve;
333:   if (presolve) {
334:     pc->ops->presolve = PCPreSolve_Shell;
335:     PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Shell);
336:   } else {
337:     pc->ops->presolve = NULL;
338:     PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
339:   }
340:   return 0;
341: }

343: static PetscErrorCode  PCShellSetPostSolve_Shell(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
344: {
345:   PC_Shell *shell = (PC_Shell*)pc->data;

347:   shell->postsolve = postsolve;
348:   if (postsolve) pc->ops->postsolve = PCPostSolve_Shell;
349:   else           pc->ops->postsolve = NULL;
350:   return 0;
351: }

353: static PetscErrorCode  PCShellSetView_Shell(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
354: {
355:   PC_Shell *shell = (PC_Shell*)pc->data;

357:   shell->view = view;
358:   return 0;
359: }

361: static PetscErrorCode  PCShellSetApplyTranspose_Shell(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
362: {
363:   PC_Shell *shell = (PC_Shell*)pc->data;

365:   shell->applytranspose = applytranspose;
366:   if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell;
367:   else                pc->ops->applytranspose = NULL;
368:   return 0;
369: }

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

375:   shell->applyrich = applyrich;
376:   if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell;
377:   else           pc->ops->applyrichardson = NULL;
378:   return 0;
379: }

381: static PetscErrorCode  PCShellSetName_Shell(PC pc,const char name[])
382: {
383:   PC_Shell       *shell = (PC_Shell*)pc->data;

385:   PetscFree(shell->name);
386:   PetscStrallocpy(name,&shell->name);
387:   return 0;
388: }

390: static PetscErrorCode  PCShellGetName_Shell(PC pc,const char *name[])
391: {
392:   PC_Shell *shell = (PC_Shell*)pc->data;

394:   *name = shell->name;
395:   return 0;
396: }

398: /* -------------------------------------------------------------------------------*/

400: /*@C
401:    PCShellSetDestroy - Sets routine to use to destroy the user-provided
402:    application context.

404:    Logically Collective on PC

406:    Input Parameters:
407: +  pc - the preconditioner context
408: -  destroy - the application-provided destroy routine

410:    Calling sequence of destroy:
411: .vb
412:    PetscErrorCode destroy (PC)
413: .ve

415: .  ptr - the application context

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

420:    Level: developer

422: .seealso: PCShellSetApply(), PCShellSetContext()
423: @*/
424: PetscErrorCode  PCShellSetDestroy(PC pc,PetscErrorCode (*destroy)(PC))
425: {
427:   PetscTryMethod(pc,"PCShellSetDestroy_C",(PC,PetscErrorCode (*)(PC)),(pc,destroy));
428:   return 0;
429: }

431: /*@C
432:    PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the
433:    matrix operator is changed.

435:    Logically Collective on PC

437:    Input Parameters:
438: +  pc - the preconditioner context
439: -  setup - the application-provided setup routine

441:    Calling sequence of setup:
442: .vb
443:    PetscErrorCode setup (PC pc)
444: .ve

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

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

451:    Level: developer

453: .seealso: PCShellSetApplyRichardson(), PCShellSetApply(), PCShellSetContext()
454: @*/
455: PetscErrorCode  PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC))
456: {
458:   PetscTryMethod(pc,"PCShellSetSetUp_C",(PC,PetscErrorCode (*)(PC)),(pc,setup));
459:   return 0;
460: }

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

465:    Logically Collective on PC

467:    Input Parameters:
468: +  pc - the preconditioner context
469: -  view - the application-provided view routine

471:    Calling sequence of view:
472: .vb
473:    PetscErrorCode view(PC pc,PetscViewer v)
474: .ve

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

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

482:    Level: developer

484: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose()
485: @*/
486: PetscErrorCode  PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
487: {
489:   PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view));
490:   return 0;
491: }

493: /*@C
494:    PCShellSetApply - Sets routine to use as preconditioner.

496:    Logically Collective on PC

498:    Input Parameters:
499: +  pc - the preconditioner context
500: -  apply - the application-provided preconditioning routine

502:    Calling sequence of apply:
503: .vb
504:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
505: .ve

507: +  pc - the preconditioner, get the application context with PCShellGetContext()
508: .  xin - input vector
509: -  xout - output vector

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

514:    Level: developer

516: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApplyBA(), PCShellSetApplySymmetricRight(),PCShellSetApplySymmetricLeft()
517: @*/
518: PetscErrorCode  PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
519: {
521:   PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
522:   return 0;
523: }

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

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,Mat Xin,Mat Xout)
537: .ve

539: +  pc - the preconditioner, get the application context with PCShellGetContext()
540: .  Xin - input block of vectors
541: -  Xout - output block of vectors

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

546:    Level: developer

548: .seealso: PCShellSetApply()
549: @*/
550: PetscErrorCode  PCShellSetMatApply(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
551: {
553:   PetscTryMethod(pc,"PCShellSetMatApply_C",(PC,PetscErrorCode (*)(PC,Mat,Mat)),(pc,matapply));
554:   return 0;
555: }

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

560:    Logically Collective on PC

562:    Input Parameters:
563: +  pc - the preconditioner context
564: -  apply - the application-provided left preconditioning routine

566:    Calling sequence of apply:
567: .vb
568:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
569: .ve

571: +  pc - the preconditioner, get the application context with PCShellGetContext()
572: .  xin - input vector
573: -  xout - output vector

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

578:    Level: developer

580: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
581: @*/
582: PetscErrorCode  PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
583: {
585:   PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
586:   return 0;
587: }

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

592:    Logically Collective on PC

594:    Input Parameters:
595: +  pc - the preconditioner context
596: -  apply - the application-provided right preconditioning routine

598:    Calling sequence of apply:
599: .vb
600:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
601: .ve

603: +  pc - the preconditioner, get the application context with PCShellGetContext()
604: .  xin - input vector
605: -  xout - output vector

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

610:    Level: developer

612: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
613: @*/
614: PetscErrorCode  PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
615: {
617:   PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
618:   return 0;
619: }

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

624:    Logically Collective on PC

626:    Input Parameters:
627: +  pc - the preconditioner context
628: -  applyBA - the application-provided BA routine

630:    Calling sequence of applyBA:
631: .vb
632:    PetscErrorCode applyBA (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: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApply()
645: @*/
646: PetscErrorCode  PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
647: {
649:   PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA));
650:   return 0;
651: }

653: /*@C
654:    PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.

656:    Logically Collective on PC

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

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

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

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

674:    Level: developer

676:    Notes:
677:    Uses the same context variable as PCShellSetApply().

679: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApply(), PCSetContext(), PCShellSetApplyBA()
680: @*/
681: PetscErrorCode  PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
682: {
684:   PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose));
685:   return 0;
686: }

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

693:    Logically Collective on PC

695:    Input Parameters:
696: +  pc - the preconditioner context
697: -  presolve - the application-provided presolve routine

699:    Calling sequence of presolve:
700: .vb
701:    PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x)
702: .ve

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

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

711:    Level: developer

713: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPostSolve(), PCShellSetContext()
714: @*/
715: PetscErrorCode  PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
716: {
718:   PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve));
719:   return 0;
720: }

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

727:    Logically Collective on PC

729:    Input Parameters:
730: +  pc - the preconditioner context
731: -  postsolve - the application-provided presolve routine

733:    Calling sequence of postsolve:
734: .vb
735:    PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x)
736: .ve

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

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

745:    Level: developer

747: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPreSolve(), PCShellSetContext()
748: @*/
749: PetscErrorCode  PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
750: {
752:   PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve));
753:   return 0;
754: }

756: /*@C
757:    PCShellSetName - Sets an optional name to associate with a shell
758:    preconditioner.

760:    Not Collective

762:    Input Parameters:
763: +  pc - the preconditioner context
764: -  name - character string describing shell preconditioner

766:    Level: developer

768: .seealso: PCShellGetName()
769: @*/
770: PetscErrorCode  PCShellSetName(PC pc,const char name[])
771: {
773:   PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name));
774:   return 0;
775: }

777: /*@C
778:    PCShellGetName - Gets an optional name that the user has set for a shell
779:    preconditioner.

781:    Not Collective

783:    Input Parameter:
784: .  pc - the preconditioner context

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

789:    Level: developer

791: .seealso: PCShellSetName()
792: @*/
793: PetscErrorCode  PCShellGetName(PC pc,const char *name[])
794: {
797:   PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name));
798:   return 0;
799: }

801: /*@C
802:    PCShellSetApplyRichardson - Sets routine to use as preconditioner
803:    in Richardson iteration.

805:    Logically Collective on PC

807:    Input Parameters:
808: +  pc - the preconditioner context
809: -  apply - the application-provided preconditioning routine

811:    Calling sequence of apply:
812: .vb
813:    PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
814: .ve

816: +  pc - the preconditioner, get the application context with PCShellGetContext()
817: .  b - right-hand-side
818: .  x - current iterate
819: .  r - work space
820: .  rtol - relative tolerance of residual norm to stop at
821: .  abstol - absolute tolerance of residual norm to stop at
822: .  dtol - if residual norm increases by this factor than return
823: -  maxits - number of iterations to run

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

828:    Level: developer

830: .seealso: PCShellSetApply(), PCShellSetContext()
831: @*/
832: PetscErrorCode  PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*))
833: {
835:   PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply));
836:   return 0;
837: }

839: /*MC
840:    PCSHELL - Creates a new preconditioner class for use with your
841:               own private data storage format.

843:    Level: advanced

845:   Usage:
846: $             extern PetscErrorCode apply(PC,Vec,Vec);
847: $             extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec);
848: $             extern PetscErrorCode applytranspose(PC,Vec,Vec);
849: $             extern PetscErrorCode setup(PC);
850: $             extern PetscErrorCode destroy(PC);
851: $
852: $             PCCreate(comm,&pc);
853: $             PCSetType(pc,PCSHELL);
854: $             PCShellSetContext(pc,ctx)
855: $             PCShellSetApply(pc,apply);
856: $             PCShellSetApplyBA(pc,applyba);               (optional)
857: $             PCShellSetApplyTranspose(pc,applytranspose); (optional)
858: $             PCShellSetSetUp(pc,setup);                   (optional)
859: $             PCShellSetDestroy(pc,destroy);               (optional)

861: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
862:            MATSHELL, PCShellSetSetUp(), PCShellSetApply(), PCShellSetView(),
863:            PCShellSetApplyTranspose(), PCShellSetName(), PCShellSetApplyRichardson(),
864:            PCShellGetName(), PCShellSetContext(), PCShellGetContext(), PCShellSetApplyBA()
865: M*/

867: PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
868: {
869:   PC_Shell       *shell;

871:   PetscNewLog(pc,&shell);
872:   pc->data = (void*)shell;

874:   pc->ops->destroy         = PCDestroy_Shell;
875:   pc->ops->view            = PCView_Shell;
876:   pc->ops->apply           = PCApply_Shell;
877:   pc->ops->matapply        = PCMatApply_Shell;
878:   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Shell;
879:   pc->ops->applysymmetricright = PCApplySymmetricRight_Shell;
880:   pc->ops->applytranspose  = NULL;
881:   pc->ops->applyrichardson = NULL;
882:   pc->ops->setup           = NULL;
883:   pc->ops->presolve        = NULL;
884:   pc->ops->postsolve       = NULL;

886:   shell->apply          = NULL;
887:   shell->applytranspose = NULL;
888:   shell->name           = NULL;
889:   shell->applyrich      = NULL;
890:   shell->presolve       = NULL;
891:   shell->postsolve      = NULL;
892:   shell->ctx            = NULL;
893:   shell->setup          = NULL;
894:   shell->view           = NULL;
895:   shell->destroy        = NULL;
896:   shell->applysymmetricleft  = NULL;
897:   shell->applysymmetricright = NULL;

899:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",PCShellSetDestroy_Shell);
900:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",PCShellSetSetUp_Shell);
901:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",PCShellSetApply_Shell);
902:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",PCShellSetMatApply_Shell);
903:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",PCShellSetApplySymmetricLeft_Shell);
904:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",PCShellSetApplySymmetricRight_Shell);
905:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",PCShellSetApplyBA_Shell);
906:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",PCShellSetPreSolve_Shell);
907:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",PCShellSetPostSolve_Shell);
908:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",PCShellSetView_Shell);
909:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",PCShellSetApplyTranspose_Shell);
910:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",PCShellSetName_Shell);
911:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",PCShellGetName_Shell);
912:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",PCShellSetApplyRichardson_Shell);
913:   return 0;
914: }