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