Actual source code: shellpc.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

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

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

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

 12:   PetscErrorCode (*destroy)(PC);
 13:   PetscErrorCode (*setup)(PC);
 14:   PetscErrorCode (*apply)(PC,Vec,Vec);
 15:   PetscErrorCode (*applysymmetricleft)(PC,Vec,Vec);
 16:   PetscErrorCode (*applysymmetricright)(PC,Vec,Vec);
 17:   PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec);
 18:   PetscErrorCode (*presolve)(PC,KSP,Vec,Vec);
 19:   PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec);
 20:   PetscErrorCode (*view)(PC,PetscViewer);
 21:   PetscErrorCode (*applytranspose)(PC,Vec,Vec);
 22:   PetscErrorCode (*applyrich)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);

 24:   char *name;
 25: } PC_Shell;

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

 30:     Not Collective

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

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

 38:     Level: advanced

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

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

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

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

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

 66:    Logically Collective on PC

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

 72:    Level: advanced

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



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

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

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

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

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

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

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

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

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

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

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

153:   if (!shell->applyBA) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyBA() routine provided to Shell PC");
154:   PetscObjectStateGet((PetscObject)w, &instate);
155:   PetscStackCall("PCSHELL user function applyBA()",(*shell->applyBA)(pc,side,x,y,w);CHKERRQ(ierr));
156:   PetscObjectStateGet((PetscObject)w, &outstate);
157:   if (instate == outstate) {
158:     /* increase the state of the output vector since the user did not update its state themselve as should have been done */
159:     PetscObjectStateIncrease((PetscObject)w);
160:   }
161:   return(0);
162: }

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

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

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

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

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

193: static PetscErrorCode PCApplyTranspose_Shell(PC pc,Vec x,Vec y)
194: {
195:   PC_Shell         *shell = (PC_Shell*)pc->data;
196:   PetscErrorCode   ierr;
197:   PetscObjectState instate,outstate;

200:   if (!shell->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applytranspose() routine provided to Shell PC");
201:   PetscObjectStateGet((PetscObject)y, &instate);
202:   PetscStackCall("PCSHELL user function applytranspose()",(*shell->applytranspose)(pc,x,y);CHKERRQ(ierr));
203:   PetscObjectStateGet((PetscObject)y, &outstate);
204:   if (instate == outstate) {
205:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
206:     PetscObjectStateIncrease((PetscObject)y);
207:   }
208:   return(0);
209: }

211: static PetscErrorCode PCApplyRichardson_Shell(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt it,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
212: {
213:   PetscErrorCode   ierr;
214:   PC_Shell         *shell = (PC_Shell*)pc->data;
215:   PetscObjectState instate,outstate;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

418: /*@C
419:    PCShellSetDestroy - Sets routine to use to destroy the user-provided
420:    Section 1.5 Writing Application Codes with PETSc context.

422:    Logically Collective on PC

424:    Input Parameters:
425: +  pc - the preconditioner context
426: -  destroy - the Section 1.5 Writing Application Codes with PETSc-provided destroy routine

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

433: .  ptr - the Section 1.5 Writing Application Codes with PETSc context

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

438:    Level: developer

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

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


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

457:    Logically Collective on PC

459:    Input Parameters:
460: +  pc - the preconditioner context
461: -  setup - the Section 1.5 Writing Application Codes with PETSc-provided setup routine

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

468: .  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()

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

473:    Level: developer

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

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


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

491:    Logically Collective on PC

493:    Input Parameters:
494: +  pc - the preconditioner context
495: -  view - the Section 1.5 Writing Application Codes with PETSc-provided view routine

497:    Calling sequence of view:
498: .vb
499:    PetscErrorCode view(PC pc,PetscViewer v)
500: .ve

502: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
503: -  v   - viewer

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

508:    Level: developer

510: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose()
511: @*/
512: PetscErrorCode  PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
513: {

518:   PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view));
519:   return(0);
520: }

522: /*@C
523:    PCShellSetApply - Sets routine to use as preconditioner.

525:    Logically Collective on PC

527:    Input Parameters:
528: +  pc - the preconditioner context
529: -  apply - the Section 1.5 Writing Application Codes with PETSc-provided preconditioning routine

531:    Calling sequence of apply:
532: .vb
533:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
534: .ve

536: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
537: .  xin - input vector
538: -  xout - output vector

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

543:    Level: developer

545: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApplyBA(), PCShellSetApplySymmetricRight(),PCShellSetApplySymmetricLeft()
546: @*/
547: PetscErrorCode  PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
548: {

553:   PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
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 Section 1.5 Writing Application Codes with PETSc-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 Section 1.5 Writing Application Codes with PETSc 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: {

588:   PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
589:   return(0);
590: }

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

595:    Logically Collective on PC

597:    Input Parameters:
598: +  pc - the preconditioner context
599: -  apply - the Section 1.5 Writing Application Codes with PETSc-provided right preconditioning routine

601:    Calling sequence of apply:
602: .vb
603:    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
604: .ve

606: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
607: .  xin - input vector
608: -  xout - output vector

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

613:    Level: developer

615: .seealso: PCShellSetApply(), PCShellSetApplySymmetricLeft(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext()
616: @*/
617: PetscErrorCode  PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
618: {

623:   PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
624:   return(0);
625: }

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

630:    Logically Collective on PC

632:    Input Parameters:
633: +  pc - the preconditioner context
634: -  applyBA - the Section 1.5 Writing Application Codes with PETSc-provided BA routine

636:    Calling sequence of applyBA:
637: .vb
638:    PetscErrorCode applyBA (PC pc,Vec xin,Vec xout)
639: .ve

641: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
642: .  xin - input vector
643: -  xout - output vector

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

648:    Level: developer

650: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetContext(), PCShellSetApply()
651: @*/
652: PetscErrorCode  PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
653: {

658:   PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA));
659:   return(0);
660: }

662: /*@C
663:    PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.

665:    Logically Collective on PC

667:    Input Parameters:
668: +  pc - the preconditioner context
669: -  apply - the Section 1.5 Writing Application Codes with PETSc-provided preconditioning transpose routine

671:    Calling sequence of apply:
672: .vb
673:    PetscErrorCode applytranspose (PC pc,Vec xin,Vec xout)
674: .ve

676: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
677: .  xin - input vector
678: -  xout - output vector

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

683:    Level: developer

685:    Notes:
686:    Uses the same context variable as PCShellSetApply().

688: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApply(), PCSetContext(), PCShellSetApplyBA()
689: @*/
690: PetscErrorCode  PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
691: {

696:   PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose));
697:   return(0);
698: }

700: /*@C
701:    PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is
702:       applied. This usually does something like scale the linear system in some Section 1.5 Writing Application Codes with PETSc
703:       specific way.

705:    Logically Collective on PC

707:    Input Parameters:
708: +  pc - the preconditioner context
709: -  presolve - the Section 1.5 Writing Application Codes with PETSc-provided presolve routine

711:    Calling sequence of presolve:
712: .vb
713:    PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x)
714: .ve

716: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
717: .  xin - input vector
718: -  xout - output vector

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

723:    Level: developer

725: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPostSolve(), PCShellSetContext()
726: @*/
727: PetscErrorCode  PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
728: {

733:   PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve));
734:   return(0);
735: }

737: /*@C
738:    PCShellSetPostSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is
739:       applied. This usually does something like scale the linear system in some Section 1.5 Writing Application Codes with PETSc
740:       specific way.

742:    Logically Collective on PC

744:    Input Parameters:
745: +  pc - the preconditioner context
746: -  postsolve - the Section 1.5 Writing Application Codes with PETSc-provided presolve routine

748:    Calling sequence of postsolve:
749: .vb
750:    PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x)
751: .ve

753: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
754: .  xin - input vector
755: -  xout - output vector

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

760:    Level: developer

762: .seealso: PCShellSetApplyRichardson(), PCShellSetSetUp(), PCShellSetApplyTranspose(), PCShellSetPreSolve(), PCShellSetContext()
763: @*/
764: PetscErrorCode  PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
765: {

770:   PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve));
771:   return(0);
772: }

774: /*@C
775:    PCShellSetName - Sets an optional name to associate with a shell
776:    preconditioner.

778:    Not Collective

780:    Input Parameters:
781: +  pc - the preconditioner context
782: -  name - character string describing shell preconditioner

784:    Level: developer

786: .seealso: PCShellGetName()
787: @*/
788: PetscErrorCode  PCShellSetName(PC pc,const char name[])
789: {

794:   PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name));
795:   return(0);
796: }

798: /*@C
799:    PCShellGetName - Gets an optional name that the user has set for a shell
800:    preconditioner.

802:    Not Collective

804:    Input Parameter:
805: .  pc - the preconditioner context

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

810:    Level: developer

812: .seealso: PCShellSetName()
813: @*/
814: PetscErrorCode  PCShellGetName(PC pc,const char *name[])
815: {

821:   PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name));
822:   return(0);
823: }

825: /*@C
826:    PCShellSetApplyRichardson - Sets routine to use as preconditioner
827:    in Richardson iteration.

829:    Logically Collective on PC

831:    Input Parameters:
832: +  pc - the preconditioner context
833: -  apply - the Section 1.5 Writing Application Codes with PETSc-provided preconditioning routine

835:    Calling sequence of apply:
836: .vb
837:    PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
838: .ve

840: +  pc - the preconditioner, get the Section 1.5 Writing Application Codes with PETSc context with PCShellGetContext()
841: .  b - right-hand-side
842: .  x - current iterate
843: .  r - work space
844: .  rtol - relative tolerance of residual norm to stop at
845: .  abstol - absolute tolerance of residual norm to stop at
846: .  dtol - if residual norm increases by this factor than return
847: -  maxits - number of iterations to run

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

852:    Level: developer

854: .seealso: PCShellSetApply(), PCShellSetContext()
855: @*/
856: PetscErrorCode  PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*))
857: {

862:   PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply));
863:   return(0);
864: }

866: /*MC
867:    PCSHELL - Creates a new preconditioner class for use with your
868:               own private data storage format.

870:    Level: advanced

872:   Usage:
873: $             extern PetscErrorCode apply(PC,Vec,Vec);
874: $             extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec);
875: $             extern PetscErrorCode applytranspose(PC,Vec,Vec);
876: $             extern PetscErrorCode setup(PC);
877: $             extern PetscErrorCode destroy(PC);
878: $
879: $             PCCreate(comm,&pc);
880: $             PCSetType(pc,PCSHELL);
881: $             PCShellSetContext(pc,ctx)
882: $             PCShellSetApply(pc,apply);
883: $             PCShellSetApplyBA(pc,applyba);               (optional)
884: $             PCShellSetApplyTranspose(pc,applytranspose); (optional)
885: $             PCShellSetSetUp(pc,setup);                   (optional)
886: $             PCShellSetDestroy(pc,destroy);               (optional)

888: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
889:            MATSHELL, PCShellSetSetUp(), PCShellSetApply(), PCShellSetView(),
890:            PCShellSetApplyTranspose(), PCShellSetName(), PCShellSetApplyRichardson(),
891:            PCShellGetName(), PCShellSetContext(), PCShellGetContext(), PCShellSetApplyBA()
892: M*/

894: PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
895: {
897:   PC_Shell       *shell;

900:   PetscNewLog(pc,&shell);
901:   pc->data = (void*)shell;

903:   pc->ops->destroy         = PCDestroy_Shell;
904:   pc->ops->view            = PCView_Shell;
905:   pc->ops->apply           = PCApply_Shell;
906:   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Shell;
907:   pc->ops->applysymmetricright = PCApplySymmetricRight_Shell;
908:   pc->ops->applytranspose  = 0;
909:   pc->ops->applyrichardson = 0;
910:   pc->ops->setup           = 0;
911:   pc->ops->presolve        = 0;
912:   pc->ops->postsolve       = 0;

914:   shell->apply          = 0;
915:   shell->applytranspose = 0;
916:   shell->name           = 0;
917:   shell->applyrich      = 0;
918:   shell->presolve       = 0;
919:   shell->postsolve      = 0;
920:   shell->ctx            = 0;
921:   shell->setup          = 0;
922:   shell->view           = 0;
923:   shell->destroy        = 0;
924:   shell->applysymmetricleft  = 0;
925:   shell->applysymmetricright = 0;

927:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",PCShellSetDestroy_Shell);
928:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",PCShellSetSetUp_Shell);
929:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",PCShellSetApply_Shell);
930:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",PCShellSetApplySymmetricLeft_Shell);
931:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",PCShellSetApplySymmetricRight_Shell);
932:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",PCShellSetApplyBA_Shell);
933:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",PCShellSetPreSolve_Shell);
934:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",PCShellSetPostSolve_Shell);
935:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",PCShellSetView_Shell);
936:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",PCShellSetApplyTranspose_Shell);
937:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",PCShellSetName_Shell);
938:   PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",PCShellGetName_Shell);
939:   PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",PCShellSetApplyRichardson_Shell);
940:   return(0);
941: }