xref: /petsc/src/ksp/pc/impls/shell/shellpc.c (revision f312b9440342594fb778b0f6acedd27284658bfc)
1 
2 /*
3    This provides a simple shell for Fortran (and C programmers) to
4   create their own preconditioner without writing much interface code.
5 */
6 
7 #include <petsc/private/pcimpl.h>        /*I "petscpc.h" I*/
8 
9 typedef struct {
10   void *ctx;                     /* user provided contexts for preconditioner */
11 
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*);
24 
25   char *name;
26 } PC_Shell;
27 
28 /*@C
29     PCShellGetContext - Returns the user-provided context associated with a shell PC
30 
31     Not Collective
32 
33     Input Parameter:
34 .   pc - should have been created with PCSetType(pc,shell)
35 
36     Output Parameter:
37 .   ctx - the user provided context
38 
39     Level: advanced
40 
41     Notes:
42     This routine is intended for use within various shell routines
43 
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.
47 
48 .seealso: `PCShellSetContext()`
49 @*/
50 PetscErrorCode  PCShellGetContext(PC pc,void *ctx)
51 {
52   PetscBool      flg;
53 
54   PetscFunctionBegin;
55   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
56   PetscValidPointer(ctx,2);
57   PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&flg));
58   if (!flg) *(void**)ctx = NULL;
59   else      *(void**)ctx = ((PC_Shell*)(pc->data))->ctx;
60   PetscFunctionReturn(0);
61 }
62 
63 /*@
64     PCShellSetContext - sets the context for a shell PC
65 
66    Logically Collective on PC
67 
68     Input Parameters:
69 +   pc - the shell PC
70 -   ctx - the context
71 
72    Level: advanced
73 
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.
77 
78 .seealso: `PCShellGetContext()`, `PCSHELL`
79 @*/
80 PetscErrorCode  PCShellSetContext(PC pc,void *ctx)
81 {
82   PC_Shell       *shell = (PC_Shell*)pc->data;
83   PetscBool      flg;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
87   PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&flg));
88   if (flg) shell->ctx = ctx;
89   PetscFunctionReturn(0);
90 }
91 
92 static PetscErrorCode PCSetUp_Shell(PC pc)
93 {
94   PC_Shell       *shell = (PC_Shell*)pc->data;
95 
96   PetscFunctionBegin;
97   PetscCheck(shell->setup,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No setup() routine provided to Shell PC");
98   PetscCallBack("PCSHELL callback setup",(*shell->setup)(pc));
99   PetscFunctionReturn(0);
100 }
101 
102 static PetscErrorCode PCApply_Shell(PC pc,Vec x,Vec y)
103 {
104   PC_Shell         *shell = (PC_Shell*)pc->data;
105   PetscObjectState instate,outstate;
106 
107   PetscFunctionBegin;
108   PetscCheck(shell->apply,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
109   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
110   PetscCallBack("PCSHELL callback apply",(*shell->apply)(pc,x,y));
111   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
112   /* increase the state of the output vector if the user did not update its state themself as should have been done */
113   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
114   PetscFunctionReturn(0);
115 }
116 
117 static PetscErrorCode PCMatApply_Shell(PC pc,Mat X,Mat Y)
118 {
119   PC_Shell         *shell = (PC_Shell*)pc->data;
120   PetscObjectState instate,outstate;
121 
122   PetscFunctionBegin;
123   PetscCheck(shell->matapply,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
124   PetscCall(PetscObjectStateGet((PetscObject)Y, &instate));
125   PetscCallBack("PCSHELL callback apply",(*shell->matapply)(pc,X,Y));
126   PetscCall(PetscObjectStateGet((PetscObject)Y, &outstate));
127   /* increase the state of the output vector if the user did not update its state themself as should have been done */
128   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)Y));
129   PetscFunctionReturn(0);
130 }
131 
132 static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc,Vec x,Vec y)
133 {
134   PC_Shell       *shell = (PC_Shell*)pc->data;
135 
136   PetscFunctionBegin;
137   PetscCheck(shell->applysymmetricleft,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
138   PetscCallBack("PCSHELL callback apply symmetric left",(*shell->applysymmetricleft)(pc,x,y));
139   PetscFunctionReturn(0);
140 }
141 
142 static PetscErrorCode PCApplySymmetricRight_Shell(PC pc,Vec x,Vec y)
143 {
144   PC_Shell       *shell = (PC_Shell*)pc->data;
145 
146   PetscFunctionBegin;
147   PetscCheck(shell->applysymmetricright,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No apply() routine provided to Shell PC");
148   PetscCallBack("PCSHELL callback apply symmetric right",(*shell->applysymmetricright)(pc,x,y));
149   PetscFunctionReturn(0);
150 }
151 
152 static PetscErrorCode PCApplyBA_Shell(PC pc,PCSide side,Vec x,Vec y,Vec w)
153 {
154   PC_Shell         *shell = (PC_Shell*)pc->data;
155   PetscObjectState instate,outstate;
156 
157   PetscFunctionBegin;
158   PetscCheck(shell->applyBA,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyBA() routine provided to Shell PC");
159   PetscCall(PetscObjectStateGet((PetscObject)w, &instate));
160   PetscCallBack("PCSHELL callback applyBA",(*shell->applyBA)(pc,side,x,y,w));
161   PetscCall(PetscObjectStateGet((PetscObject)w, &outstate));
162   /* increase the state of the output vector if the user did not update its state themself as should have been done */
163   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)w));
164   PetscFunctionReturn(0);
165 }
166 
167 static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc,PetscBool* change)
168 {
169   PetscFunctionBegin;
170   *change = PETSC_TRUE;
171   PetscFunctionReturn(0);
172 }
173 
174 static PetscErrorCode PCPreSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)
175 {
176   PC_Shell       *shell = (PC_Shell*)pc->data;
177 
178   PetscFunctionBegin;
179   PetscCheck(shell->presolve,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No presolve() routine provided to Shell PC");
180   PetscCallBack("PCSHELL callback presolve",(*shell->presolve)(pc,ksp,b,x));
181   PetscFunctionReturn(0);
182 }
183 
184 static PetscErrorCode PCPostSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)
185 {
186   PC_Shell       *shell = (PC_Shell*)pc->data;
187 
188   PetscFunctionBegin;
189   PetscCheck(shell->postsolve,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No postsolve() routine provided to Shell PC");
190   PetscCallBack("PCSHELL callback postsolve()",(*shell->postsolve)(pc,ksp,b,x));
191   PetscFunctionReturn(0);
192 }
193 
194 static PetscErrorCode PCApplyTranspose_Shell(PC pc,Vec x,Vec y)
195 {
196   PC_Shell         *shell = (PC_Shell*)pc->data;
197   PetscObjectState instate,outstate;
198 
199   PetscFunctionBegin;
200   PetscCheck(shell->applytranspose,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applytranspose() routine provided to Shell PC");
201   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
202   PetscCallBack("PCSHELL callback applytranspose",(*shell->applytranspose)(pc,x,y));
203   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
204   /* increase the state of the output vector if the user did not update its state themself as should have been done */
205   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
206   PetscFunctionReturn(0);
207 }
208 
209 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)
210 {
211   PC_Shell         *shell = (PC_Shell*)pc->data;
212   PetscObjectState instate,outstate;
213 
214   PetscFunctionBegin;
215   PetscCheck(shell->applyrich,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No applyrichardson() routine provided to Shell PC");
216   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
217   PetscCallBack("PCSHELL callback applyrichardson",(*shell->applyrich)(pc,x,y,w,rtol,abstol,dtol,it,guesszero,outits,reason));
218   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
219   /* increase the state of the output vector since the user did not update its state themself as should have been done */
220   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
221   PetscFunctionReturn(0);
222 }
223 
224 static PetscErrorCode PCDestroy_Shell(PC pc)
225 {
226   PC_Shell       *shell = (PC_Shell*)pc->data;
227 
228   PetscFunctionBegin;
229   PetscCall(PetscFree(shell->name));
230   if (shell->destroy) PetscCallBack("PCSHELL callback destroy",(*shell->destroy)(pc));
231   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",NULL));
232   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",NULL));
233   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",NULL));
234   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",NULL));
235   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",NULL));
236   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",NULL));
237   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",NULL));
238   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",NULL));
239   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",NULL));
240   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",NULL));
241   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",NULL));
242   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",NULL));
243   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",NULL));
244   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",NULL));
245   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL));
246   PetscCall(PetscFree(pc->data));
247   PetscFunctionReturn(0);
248 }
249 
250 static PetscErrorCode PCView_Shell(PC pc,PetscViewer viewer)
251 {
252   PC_Shell       *shell = (PC_Shell*)pc->data;
253   PetscBool      iascii;
254 
255   PetscFunctionBegin;
256   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
257   if (iascii) {
258     if (shell->name) PetscCall(PetscViewerASCIIPrintf(viewer,"  %s\n",shell->name));
259     else             PetscCall(PetscViewerASCIIPrintf(viewer,"  no name\n"));
260   }
261   if (shell->view) {
262     PetscCall(PetscViewerASCIIPushTab(viewer));
263     PetscCall((*shell->view)(pc,viewer));
264     PetscCall(PetscViewerASCIIPopTab(viewer));
265   }
266   PetscFunctionReturn(0);
267 }
268 
269 /* ------------------------------------------------------------------------------*/
270 static PetscErrorCode  PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC))
271 {
272   PC_Shell *shell= (PC_Shell*)pc->data;
273 
274   PetscFunctionBegin;
275   shell->destroy = destroy;
276   PetscFunctionReturn(0);
277 }
278 
279 static PetscErrorCode  PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC))
280 {
281   PC_Shell *shell = (PC_Shell*)pc->data;
282 
283   PetscFunctionBegin;
284   shell->setup = setup;
285   if (setup) pc->ops->setup = PCSetUp_Shell;
286   else       pc->ops->setup = NULL;
287   PetscFunctionReturn(0);
288 }
289 
290 static PetscErrorCode  PCShellSetApply_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
291 {
292   PC_Shell *shell = (PC_Shell*)pc->data;
293 
294   PetscFunctionBegin;
295   shell->apply = apply;
296   PetscFunctionReturn(0);
297 }
298 
299 static PetscErrorCode  PCShellSetMatApply_Shell(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
300 {
301   PC_Shell *shell = (PC_Shell*)pc->data;
302 
303   PetscFunctionBegin;
304   shell->matapply = matapply;
305   if (matapply) pc->ops->matapply = PCMatApply_Shell;
306   else          pc->ops->matapply = NULL;
307   PetscFunctionReturn(0);
308 }
309 
310 static PetscErrorCode  PCShellSetApplySymmetricLeft_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
311 {
312   PC_Shell *shell = (PC_Shell*)pc->data;
313 
314   PetscFunctionBegin;
315   shell->applysymmetricleft = apply;
316   PetscFunctionReturn(0);
317 }
318 
319 static PetscErrorCode  PCShellSetApplySymmetricRight_Shell(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
320 {
321   PC_Shell *shell = (PC_Shell*)pc->data;
322 
323   PetscFunctionBegin;
324   shell->applysymmetricright = apply;
325   PetscFunctionReturn(0);
326 }
327 
328 static PetscErrorCode  PCShellSetApplyBA_Shell(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
329 {
330   PC_Shell *shell = (PC_Shell*)pc->data;
331 
332   PetscFunctionBegin;
333   shell->applyBA = applyBA;
334   if (applyBA) pc->ops->applyBA  = PCApplyBA_Shell;
335   else         pc->ops->applyBA  = NULL;
336   PetscFunctionReturn(0);
337 }
338 
339 static PetscErrorCode  PCShellSetPreSolve_Shell(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
340 {
341   PC_Shell       *shell = (PC_Shell*)pc->data;
342 
343   PetscFunctionBegin;
344   shell->presolve = presolve;
345   if (presolve) {
346     pc->ops->presolve = PCPreSolve_Shell;
347     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Shell));
348   } else {
349     pc->ops->presolve = NULL;
350     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL));
351   }
352   PetscFunctionReturn(0);
353 }
354 
355 static PetscErrorCode  PCShellSetPostSolve_Shell(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
356 {
357   PC_Shell *shell = (PC_Shell*)pc->data;
358 
359   PetscFunctionBegin;
360   shell->postsolve = postsolve;
361   if (postsolve) pc->ops->postsolve = PCPostSolve_Shell;
362   else           pc->ops->postsolve = NULL;
363   PetscFunctionReturn(0);
364 }
365 
366 static PetscErrorCode  PCShellSetView_Shell(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
367 {
368   PC_Shell *shell = (PC_Shell*)pc->data;
369 
370   PetscFunctionBegin;
371   shell->view = view;
372   PetscFunctionReturn(0);
373 }
374 
375 static PetscErrorCode  PCShellSetApplyTranspose_Shell(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
376 {
377   PC_Shell *shell = (PC_Shell*)pc->data;
378 
379   PetscFunctionBegin;
380   shell->applytranspose = applytranspose;
381   if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell;
382   else                pc->ops->applytranspose = NULL;
383   PetscFunctionReturn(0);
384 }
385 
386 static PetscErrorCode  PCShellSetApplyRichardson_Shell(PC pc,PetscErrorCode (*applyrich)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*))
387 {
388   PC_Shell *shell = (PC_Shell*)pc->data;
389 
390   PetscFunctionBegin;
391   shell->applyrich = applyrich;
392   if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell;
393   else           pc->ops->applyrichardson = NULL;
394   PetscFunctionReturn(0);
395 }
396 
397 static PetscErrorCode  PCShellSetName_Shell(PC pc,const char name[])
398 {
399   PC_Shell       *shell = (PC_Shell*)pc->data;
400 
401   PetscFunctionBegin;
402   PetscCall(PetscFree(shell->name));
403   PetscCall(PetscStrallocpy(name,&shell->name));
404   PetscFunctionReturn(0);
405 }
406 
407 static PetscErrorCode  PCShellGetName_Shell(PC pc,const char *name[])
408 {
409   PC_Shell *shell = (PC_Shell*)pc->data;
410 
411   PetscFunctionBegin;
412   *name = shell->name;
413   PetscFunctionReturn(0);
414 }
415 
416 /* -------------------------------------------------------------------------------*/
417 
418 /*@C
419    PCShellSetDestroy - Sets routine to use to destroy the user-provided
420    application context.
421 
422    Logically Collective on PC
423 
424    Input Parameters:
425 +  pc - the preconditioner context
426 -  destroy - the application-provided destroy routine
427 
428    Calling sequence of destroy:
429 .vb
430    PetscErrorCode destroy (PC)
431 .ve
432 
433 .  ptr - the application context
434 
435    Notes:
436     the function MUST return an error code of 0 on success and nonzero on failure.
437 
438    Level: developer
439 
440 .seealso: `PCShellSetApply()`, `PCShellSetContext()`
441 @*/
442 PetscErrorCode  PCShellSetDestroy(PC pc,PetscErrorCode (*destroy)(PC))
443 {
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
446   PetscTryMethod(pc,"PCShellSetDestroy_C",(PC,PetscErrorCode (*)(PC)),(pc,destroy));
447   PetscFunctionReturn(0);
448 }
449 
450 /*@C
451    PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the
452    matrix operator is changed.
453 
454    Logically Collective on PC
455 
456    Input Parameters:
457 +  pc - the preconditioner context
458 -  setup - the application-provided setup routine
459 
460    Calling sequence of setup:
461 .vb
462    PetscErrorCode setup (PC pc)
463 .ve
464 
465 .  pc - the preconditioner, get the application context with PCShellGetContext()
466 
467    Notes:
468     the function MUST return an error code of 0 on success and nonzero on failure.
469 
470    Level: developer
471 
472 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetApply()`, `PCShellSetContext()`
473 @*/
474 PetscErrorCode  PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC))
475 {
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
478   PetscTryMethod(pc,"PCShellSetSetUp_C",(PC,PetscErrorCode (*)(PC)),(pc,setup));
479   PetscFunctionReturn(0);
480 }
481 
482 /*@C
483    PCShellSetView - Sets routine to use as viewer of shell preconditioner
484 
485    Logically Collective on PC
486 
487    Input Parameters:
488 +  pc - the preconditioner context
489 -  view - the application-provided view routine
490 
491    Calling sequence of view:
492 .vb
493    PetscErrorCode view(PC pc,PetscViewer v)
494 .ve
495 
496 +  pc - the preconditioner, get the application context with PCShellGetContext()
497 -  v   - viewer
498 
499    Notes:
500     the function MUST return an error code of 0 on success and nonzero on failure.
501 
502    Level: developer
503 
504 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`
505 @*/
506 PetscErrorCode  PCShellSetView(PC pc,PetscErrorCode (*view)(PC,PetscViewer))
507 {
508   PetscFunctionBegin;
509   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
510   PetscTryMethod(pc,"PCShellSetView_C",(PC,PetscErrorCode (*)(PC,PetscViewer)),(pc,view));
511   PetscFunctionReturn(0);
512 }
513 
514 /*@C
515    PCShellSetApply - Sets routine to use as preconditioner.
516 
517    Logically Collective on PC
518 
519    Input Parameters:
520 +  pc - the preconditioner context
521 -  apply - the application-provided preconditioning routine
522 
523    Calling sequence of apply:
524 .vb
525    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
526 .ve
527 
528 +  pc - the preconditioner, get the application context with PCShellGetContext()
529 .  xin - input vector
530 -  xout - output vector
531 
532    Notes:
533     the function MUST return an error code of 0 on success and nonzero on failure.
534 
535    Level: developer
536 
537 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellSetApplySymmetricRight()`, `PCShellSetApplySymmetricLeft()`
538 @*/
539 PetscErrorCode  PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
540 {
541   PetscFunctionBegin;
542   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
543   PetscTryMethod(pc,"PCShellSetApply_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
544   PetscFunctionReturn(0);
545 }
546 
547 /*@C
548    PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors.
549 
550    Logically Collective on PC
551 
552    Input Parameters:
553 +  pc - the preconditioner context
554 -  apply - the application-provided preconditioning routine
555 
556    Calling sequence of apply:
557 .vb
558    PetscErrorCode apply (PC pc,Mat Xin,Mat Xout)
559 .ve
560 
561 +  pc - the preconditioner, get the application context with PCShellGetContext()
562 .  Xin - input block of vectors
563 -  Xout - output block of vectors
564 
565    Notes:
566     the function MUST return an error code of 0 on success and nonzero on failure.
567 
568    Level: developer
569 
570 .seealso: `PCShellSetApply()`
571 @*/
572 PetscErrorCode  PCShellSetMatApply(PC pc,PetscErrorCode (*matapply)(PC,Mat,Mat))
573 {
574   PetscFunctionBegin;
575   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
576   PetscTryMethod(pc,"PCShellSetMatApply_C",(PC,PetscErrorCode (*)(PC,Mat,Mat)),(pc,matapply));
577   PetscFunctionReturn(0);
578 }
579 
580 /*@C
581    PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the PC_SYMMETRIC is used).
582 
583    Logically Collective on PC
584 
585    Input Parameters:
586 +  pc - the preconditioner context
587 -  apply - the application-provided left preconditioning routine
588 
589    Calling sequence of apply:
590 .vb
591    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
592 .ve
593 
594 +  pc - the preconditioner, get the application context with PCShellGetContext()
595 .  xin - input vector
596 -  xout - output vector
597 
598    Notes:
599     the function MUST return an error code of 0 on success and nonzero on failure.
600 
601    Level: developer
602 
603 .seealso: `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`
604 @*/
605 PetscErrorCode  PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
606 {
607   PetscFunctionBegin;
608   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
609   PetscTryMethod(pc,"PCShellSetApplySymmetricLeft_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
610   PetscFunctionReturn(0);
611 }
612 
613 /*@C
614    PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the PC_SYMMETRIC is used).
615 
616    Logically Collective on PC
617 
618    Input Parameters:
619 +  pc - the preconditioner context
620 -  apply - the application-provided right preconditioning routine
621 
622    Calling sequence of apply:
623 .vb
624    PetscErrorCode apply (PC pc,Vec xin,Vec xout)
625 .ve
626 
627 +  pc - the preconditioner, get the application context with PCShellGetContext()
628 .  xin - input vector
629 -  xout - output vector
630 
631    Notes:
632     the function MUST return an error code of 0 on success and nonzero on failure.
633 
634    Level: developer
635 
636 .seealso: `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`
637 @*/
638 PetscErrorCode  PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec))
639 {
640   PetscFunctionBegin;
641   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
642   PetscTryMethod(pc,"PCShellSetApplySymmetricRight_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,apply));
643   PetscFunctionReturn(0);
644 }
645 
646 /*@C
647    PCShellSetApplyBA - Sets routine to use as preconditioner times operator.
648 
649    Logically Collective on PC
650 
651    Input Parameters:
652 +  pc - the preconditioner context
653 -  applyBA - the application-provided BA routine
654 
655    Calling sequence of applyBA:
656 .vb
657    PetscErrorCode applyBA (PC pc,Vec xin,Vec xout)
658 .ve
659 
660 +  pc - the preconditioner, get the application context with PCShellGetContext()
661 .  xin - input vector
662 -  xout - output vector
663 
664    Notes:
665     the function MUST return an error code of 0 on success and nonzero on failure.
666 
667    Level: developer
668 
669 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApply()`
670 @*/
671 PetscErrorCode  PCShellSetApplyBA(PC pc,PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec))
672 {
673   PetscFunctionBegin;
674   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
675   PetscTryMethod(pc,"PCShellSetApplyBA_C",(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)),(pc,applyBA));
676   PetscFunctionReturn(0);
677 }
678 
679 /*@C
680    PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.
681 
682    Logically Collective on PC
683 
684    Input Parameters:
685 +  pc - the preconditioner context
686 -  apply - the application-provided preconditioning transpose routine
687 
688    Calling sequence of apply:
689 .vb
690    PetscErrorCode applytranspose (PC pc,Vec xin,Vec xout)
691 .ve
692 
693 +  pc - the preconditioner, get the application context with PCShellGetContext()
694 .  xin - input vector
695 -  xout - output vector
696 
697    Notes:
698     the function MUST return an error code of 0 on success and nonzero on failure.
699 
700    Level: developer
701 
702    Notes:
703    Uses the same context variable as PCShellSetApply().
704 
705 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCSetContext()`, `PCShellSetApplyBA()`
706 @*/
707 PetscErrorCode  PCShellSetApplyTranspose(PC pc,PetscErrorCode (*applytranspose)(PC,Vec,Vec))
708 {
709   PetscFunctionBegin;
710   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
711   PetscTryMethod(pc,"PCShellSetApplyTranspose_C",(PC,PetscErrorCode (*)(PC,Vec,Vec)),(pc,applytranspose));
712   PetscFunctionReturn(0);
713 }
714 
715 /*@C
716    PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is
717       applied. This usually does something like scale the linear system in some application
718       specific way.
719 
720    Logically Collective on PC
721 
722    Input Parameters:
723 +  pc - the preconditioner context
724 -  presolve - the application-provided presolve routine
725 
726    Calling sequence of presolve:
727 .vb
728    PetscErrorCode presolve (PC,KSP ksp,Vec b,Vec x)
729 .ve
730 
731 +  pc - the preconditioner, get the application context with PCShellGetContext()
732 .  xin - input vector
733 -  xout - output vector
734 
735    Notes:
736     the function MUST return an error code of 0 on success and nonzero on failure.
737 
738    Level: developer
739 
740 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPostSolve()`, `PCShellSetContext()`
741 @*/
742 PetscErrorCode  PCShellSetPreSolve(PC pc,PetscErrorCode (*presolve)(PC,KSP,Vec,Vec))
743 {
744   PetscFunctionBegin;
745   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
746   PetscTryMethod(pc,"PCShellSetPreSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,presolve));
747   PetscFunctionReturn(0);
748 }
749 
750 /*@C
751    PCShellSetPostSolve - Sets routine to apply to the operators/vectors before a KSPSolve() is
752       applied. This usually does something like scale the linear system in some application
753       specific way.
754 
755    Logically Collective on PC
756 
757    Input Parameters:
758 +  pc - the preconditioner context
759 -  postsolve - the application-provided presolve routine
760 
761    Calling sequence of postsolve:
762 .vb
763    PetscErrorCode postsolve(PC,KSP ksp,Vec b,Vec x)
764 .ve
765 
766 +  pc - the preconditioner, get the application context with PCShellGetContext()
767 .  xin - input vector
768 -  xout - output vector
769 
770    Notes:
771     the function MUST return an error code of 0 on success and nonzero on failure.
772 
773    Level: developer
774 
775 .seealso: `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPreSolve()`, `PCShellSetContext()`
776 @*/
777 PetscErrorCode  PCShellSetPostSolve(PC pc,PetscErrorCode (*postsolve)(PC,KSP,Vec,Vec))
778 {
779   PetscFunctionBegin;
780   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
781   PetscTryMethod(pc,"PCShellSetPostSolve_C",(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)),(pc,postsolve));
782   PetscFunctionReturn(0);
783 }
784 
785 /*@C
786    PCShellSetName - Sets an optional name to associate with a shell
787    preconditioner.
788 
789    Not Collective
790 
791    Input Parameters:
792 +  pc - the preconditioner context
793 -  name - character string describing shell preconditioner
794 
795    Level: developer
796 
797 .seealso: `PCShellGetName()`
798 @*/
799 PetscErrorCode  PCShellSetName(PC pc,const char name[])
800 {
801   PetscFunctionBegin;
802   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
803   PetscTryMethod(pc,"PCShellSetName_C",(PC,const char []),(pc,name));
804   PetscFunctionReturn(0);
805 }
806 
807 /*@C
808    PCShellGetName - Gets an optional name that the user has set for a shell
809    preconditioner.
810 
811    Not Collective
812 
813    Input Parameter:
814 .  pc - the preconditioner context
815 
816    Output Parameter:
817 .  name - character string describing shell preconditioner (you should not free this)
818 
819    Level: developer
820 
821 .seealso: `PCShellSetName()`
822 @*/
823 PetscErrorCode  PCShellGetName(PC pc,const char *name[])
824 {
825   PetscFunctionBegin;
826   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
827   PetscValidPointer(name,2);
828   PetscUseMethod(pc,"PCShellGetName_C",(PC,const char*[]),(pc,name));
829   PetscFunctionReturn(0);
830 }
831 
832 /*@C
833    PCShellSetApplyRichardson - Sets routine to use as preconditioner
834    in Richardson iteration.
835 
836    Logically Collective on PC
837 
838    Input Parameters:
839 +  pc - the preconditioner context
840 -  apply - the application-provided preconditioning routine
841 
842    Calling sequence of apply:
843 .vb
844    PetscErrorCode apply (PC pc,Vec b,Vec x,Vec r,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
845 .ve
846 
847 +  pc - the preconditioner, get the application context with PCShellGetContext()
848 .  b - right-hand-side
849 .  x - current iterate
850 .  r - work space
851 .  rtol - relative tolerance of residual norm to stop at
852 .  abstol - absolute tolerance of residual norm to stop at
853 .  dtol - if residual norm increases by this factor than return
854 -  maxits - number of iterations to run
855 
856    Notes:
857     the function MUST return an error code of 0 on success and nonzero on failure.
858 
859    Level: developer
860 
861 .seealso: `PCShellSetApply()`, `PCShellSetContext()`
862 @*/
863 PetscErrorCode  PCShellSetApplyRichardson(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*))
864 {
865   PetscFunctionBegin;
866   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
867   PetscTryMethod(pc,"PCShellSetApplyRichardson_C",(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*)),(pc,apply));
868   PetscFunctionReturn(0);
869 }
870 
871 /*MC
872    PCSHELL - Creates a new preconditioner class for use with your
873               own private data storage format.
874 
875    Level: advanced
876 
877   Usage:
878 $             extern PetscErrorCode apply(PC,Vec,Vec);
879 $             extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec);
880 $             extern PetscErrorCode applytranspose(PC,Vec,Vec);
881 $             extern PetscErrorCode setup(PC);
882 $             extern PetscErrorCode destroy(PC);
883 $
884 $             PCCreate(comm,&pc);
885 $             PCSetType(pc,PCSHELL);
886 $             PCShellSetContext(pc,ctx)
887 $             PCShellSetApply(pc,apply);
888 $             PCShellSetApplyBA(pc,applyba);               (optional)
889 $             PCShellSetApplyTranspose(pc,applytranspose); (optional)
890 $             PCShellSetSetUp(pc,setup);                   (optional)
891 $             PCShellSetDestroy(pc,destroy);               (optional)
892 
893 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
894           `MATSHELL`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetView()`,
895           `PCShellSetApplyTranspose()`, `PCShellSetName()`, `PCShellSetApplyRichardson()`,
896           `PCShellGetName()`, `PCShellSetContext()`, `PCShellGetContext()`, `PCShellSetApplyBA()`
897 M*/
898 
899 PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
900 {
901   PC_Shell       *shell;
902 
903   PetscFunctionBegin;
904   PetscCall(PetscNewLog(pc,&shell));
905   pc->data = (void*)shell;
906 
907   pc->ops->destroy         = PCDestroy_Shell;
908   pc->ops->view            = PCView_Shell;
909   pc->ops->apply           = PCApply_Shell;
910   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Shell;
911   pc->ops->applysymmetricright = PCApplySymmetricRight_Shell;
912   pc->ops->matapply        = NULL;
913   pc->ops->applytranspose  = NULL;
914   pc->ops->applyrichardson = NULL;
915   pc->ops->setup           = NULL;
916   pc->ops->presolve        = NULL;
917   pc->ops->postsolve       = NULL;
918 
919   shell->apply          = NULL;
920   shell->applytranspose = NULL;
921   shell->name           = NULL;
922   shell->applyrich      = NULL;
923   shell->presolve       = NULL;
924   shell->postsolve      = NULL;
925   shell->ctx            = NULL;
926   shell->setup          = NULL;
927   shell->view           = NULL;
928   shell->destroy        = NULL;
929   shell->applysymmetricleft  = NULL;
930   shell->applysymmetricright = NULL;
931 
932   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetDestroy_C",PCShellSetDestroy_Shell));
933   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetSetUp_C",PCShellSetSetUp_Shell));
934   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApply_C",PCShellSetApply_Shell));
935   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetMatApply_C",PCShellSetMatApply_Shell));
936   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricLeft_C",PCShellSetApplySymmetricLeft_Shell));
937   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplySymmetricRight_C",PCShellSetApplySymmetricRight_Shell));
938   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyBA_C",PCShellSetApplyBA_Shell));
939   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPreSolve_C",PCShellSetPreSolve_Shell));
940   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetPostSolve_C",PCShellSetPostSolve_Shell));
941   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetView_C",PCShellSetView_Shell));
942   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyTranspose_C",PCShellSetApplyTranspose_Shell));
943   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetName_C",PCShellSetName_Shell));
944   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellGetName_C",PCShellGetName_Shell));
945   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCShellSetApplyRichardson_C",PCShellSetApplyRichardson_Shell));
946   PetscFunctionReturn(0);
947 }
948