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