xref: /petsc/src/ksp/pc/impls/shell/shellpc.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith    This provides a simple shell for Fortran (and C programmers) to
34b9ad928SBarry Smith   create their own preconditioner without writing much interface code.
44b9ad928SBarry Smith */
54b9ad928SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith typedef struct {
9*2a8381b2SBarry Smith   PetscCtx ctx; /* user provided contexts for preconditioner */
102fa5cd67SKarl Rupp 
116891c3e4SJed Brown   PetscErrorCode (*destroy)(PC);
126891c3e4SJed Brown   PetscErrorCode (*setup)(PC);
136891c3e4SJed Brown   PetscErrorCode (*apply)(PC, Vec, Vec);
147b6e2003SPierre Jolivet   PetscErrorCode (*matapply)(PC, Mat, Mat);
151b581b66SBarry Smith   PetscErrorCode (*applysymmetricleft)(PC, Vec, Vec);
161b581b66SBarry Smith   PetscErrorCode (*applysymmetricright)(PC, Vec, Vec);
176891c3e4SJed Brown   PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec);
186891c3e4SJed Brown   PetscErrorCode (*presolve)(PC, KSP, Vec, Vec);
196891c3e4SJed Brown   PetscErrorCode (*postsolve)(PC, KSP, Vec, Vec);
206891c3e4SJed Brown   PetscErrorCode (*view)(PC, PetscViewer);
216891c3e4SJed Brown   PetscErrorCode (*applytranspose)(PC, Vec, Vec);
229fa64a6bSPierre Jolivet   PetscErrorCode (*matapplytranspose)(PC, Mat, Mat);
23ace3abfcSBarry Smith   PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
242fa5cd67SKarl Rupp 
254b9ad928SBarry Smith   char *name;
264b9ad928SBarry Smith } PC_Shell;
274b9ad928SBarry Smith 
28b29801fcSSatish Balay /*@C
2904c3f3b8SBarry Smith   PCShellGetContext - Returns the user-provided context associated with a shell `PC` that was provided with `PCShellSetContext()`
30be29d3c6SBarry Smith 
31be29d3c6SBarry Smith   Not Collective
32be29d3c6SBarry Smith 
33be29d3c6SBarry Smith   Input Parameter:
3404c3f3b8SBarry Smith . pc - of type `PCSHELL`
35be29d3c6SBarry Smith 
36be29d3c6SBarry Smith   Output Parameter:
37be29d3c6SBarry Smith . ctx - the user provided context
38be29d3c6SBarry Smith 
39be29d3c6SBarry Smith   Level: advanced
40be29d3c6SBarry Smith 
41f1580f4eSBarry Smith   Note:
42c33202f7SPierre Jolivet   This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()`
43be29d3c6SBarry Smith 
44*2a8381b2SBarry Smith   Fortran Notes:
45*2a8381b2SBarry Smith   This only works when the context is a Fortran derived type or a `PetscObject`. Define `ctx` with
46*2a8381b2SBarry Smith .vb
47*2a8381b2SBarry Smith   type(tUsertype), pointer :: ctx
48*2a8381b2SBarry Smith .ve
49daf670e6SBarry Smith 
50562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellSetDestroy()`
51be29d3c6SBarry Smith @*/
PCShellGetContext(PC pc,PetscCtxRt ctx)52*2a8381b2SBarry Smith PetscErrorCode PCShellGetContext(PC pc, PetscCtxRt ctx)
53d71ae5a4SJacob Faibussowitsch {
54ace3abfcSBarry Smith   PetscBool flg;
55be29d3c6SBarry Smith 
56be29d3c6SBarry Smith   PetscFunctionBegin;
570700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
584f572ea9SToby Isaac   PetscAssertPointer(ctx, 2);
599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg));
603ec1f749SStefano Zampini   if (!flg) *(void **)ctx = NULL;
61f4f49eeaSPierre Jolivet   else *(void **)ctx = ((PC_Shell *)pc->data)->ctx;
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63be29d3c6SBarry Smith }
64be29d3c6SBarry Smith 
659dd1005fSJed Brown /*@
6604c3f3b8SBarry Smith   PCShellSetContext - sets the context for a shell `PC` that can be accessed with `PCShellGetContext()`
67be29d3c6SBarry Smith 
68c3339decSBarry Smith   Logically Collective
69be29d3c6SBarry Smith 
70be29d3c6SBarry Smith   Input Parameters:
71f1580f4eSBarry Smith + pc  - the `PC` of type `PCSHELL`
72be29d3c6SBarry Smith - ctx - the context
73be29d3c6SBarry Smith 
74be29d3c6SBarry Smith   Level: advanced
75be29d3c6SBarry Smith 
7604c3f3b8SBarry Smith   Notes:
77c33202f7SPierre Jolivet   This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()`
7804c3f3b8SBarry Smith 
7904c3f3b8SBarry Smith   One should also provide a routine to destroy the context when `pc` is destroyed with `PCShellSetDestroy()`
8004c3f3b8SBarry Smith 
81feefa0e1SJacob Faibussowitsch   Fortran Notes:
8295452b02SPatrick Sanan   To use this from Fortran you must write a Fortran interface definition for this
83feaf08eaSBarry Smith   function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument.
84daf670e6SBarry Smith 
85562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCShellGetContext()`, `PCSHELL`, `PCShellSetApply()`, `PCShellSetDestroy()`
86be29d3c6SBarry Smith @*/
PCShellSetContext(PC pc,PetscCtx ctx)87*2a8381b2SBarry Smith PetscErrorCode PCShellSetContext(PC pc, PetscCtx ctx)
88d71ae5a4SJacob Faibussowitsch {
89c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
90ace3abfcSBarry Smith   PetscBool flg;
91be29d3c6SBarry Smith 
92be29d3c6SBarry Smith   PetscFunctionBegin;
930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg));
952fa5cd67SKarl Rupp   if (flg) shell->ctx = ctx;
963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97be29d3c6SBarry Smith }
98be29d3c6SBarry Smith 
PCSetUp_Shell(PC pc)99d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Shell(PC pc)
100d71ae5a4SJacob Faibussowitsch {
101c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith   PetscFunctionBegin;
10428b400f6SJacob Faibussowitsch   PetscCheck(shell->setup, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No setup() routine provided to Shell PC");
10525e27a38SBarry Smith   PetscCallBack("PCSHELL callback setup", (*shell->setup)(pc));
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1074b9ad928SBarry Smith }
1084b9ad928SBarry Smith 
PCApply_Shell(PC pc,Vec x,Vec y)109d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Shell(PC pc, Vec x, Vec y)
110d71ae5a4SJacob Faibussowitsch {
111c5ae4b9aSBarry Smith   PC_Shell        *shell = (PC_Shell *)pc->data;
112e3f487b0SBarry Smith   PetscObjectState instate, outstate;
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith   PetscFunctionBegin;
11528b400f6SJacob Faibussowitsch   PetscCheck(shell->apply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC");
1169566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
11725e27a38SBarry Smith   PetscCallBack("PCSHELL callback apply", (*shell->apply)(pc, x, y));
1189566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
1191e66621cSBarry Smith   /* increase the state of the output vector if the user did not update its state themself as should have been done */
1201e66621cSBarry Smith   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1224b9ad928SBarry Smith }
1234b9ad928SBarry Smith 
PCMatApply_Shell(PC pc,Mat X,Mat Y)124d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Shell(PC pc, Mat X, Mat Y)
125d71ae5a4SJacob Faibussowitsch {
1267b6e2003SPierre Jolivet   PC_Shell        *shell = (PC_Shell *)pc->data;
1277b6e2003SPierre Jolivet   PetscObjectState instate, outstate;
1287b6e2003SPierre Jolivet 
1297b6e2003SPierre Jolivet   PetscFunctionBegin;
13028b400f6SJacob Faibussowitsch   PetscCheck(shell->matapply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC");
1319566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)Y, &instate));
13225e27a38SBarry Smith   PetscCallBack("PCSHELL callback apply", (*shell->matapply)(pc, X, Y));
1339566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)Y, &outstate));
1341e66621cSBarry Smith   /* increase the state of the output vector if the user did not update its state themself as should have been done */
1351e66621cSBarry Smith   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)Y));
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1377b6e2003SPierre Jolivet }
1387b6e2003SPierre Jolivet 
PCApplySymmetricLeft_Shell(PC pc,Vec x,Vec y)139d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc, Vec x, Vec y)
140d71ae5a4SJacob Faibussowitsch {
1411b581b66SBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
1421b581b66SBarry Smith 
1431b581b66SBarry Smith   PetscFunctionBegin;
14428b400f6SJacob Faibussowitsch   PetscCheck(shell->applysymmetricleft, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC");
14525e27a38SBarry Smith   PetscCallBack("PCSHELL callback apply symmetric left", (*shell->applysymmetricleft)(pc, x, y));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1471b581b66SBarry Smith }
1481b581b66SBarry Smith 
PCApplySymmetricRight_Shell(PC pc,Vec x,Vec y)149d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Shell(PC pc, Vec x, Vec y)
150d71ae5a4SJacob Faibussowitsch {
1511b581b66SBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
1521b581b66SBarry Smith 
1531b581b66SBarry Smith   PetscFunctionBegin;
15428b400f6SJacob Faibussowitsch   PetscCheck(shell->applysymmetricright, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC");
15525e27a38SBarry Smith   PetscCallBack("PCSHELL callback apply symmetric right", (*shell->applysymmetricright)(pc, x, y));
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1571b581b66SBarry Smith }
1581b581b66SBarry Smith 
PCApplyBA_Shell(PC pc,PCSide side,Vec x,Vec y,Vec w)159d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyBA_Shell(PC pc, PCSide side, Vec x, Vec y, Vec w)
160d71ae5a4SJacob Faibussowitsch {
161c5ae4b9aSBarry Smith   PC_Shell        *shell = (PC_Shell *)pc->data;
162e3f487b0SBarry Smith   PetscObjectState instate, outstate;
1632bb17772SBarry Smith 
1642bb17772SBarry Smith   PetscFunctionBegin;
16528b400f6SJacob Faibussowitsch   PetscCheck(shell->applyBA, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyBA() routine provided to Shell PC");
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)w, &instate));
16725e27a38SBarry Smith   PetscCallBack("PCSHELL callback applyBA", (*shell->applyBA)(pc, side, x, y, w));
1689566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)w, &outstate));
1691e66621cSBarry Smith   /* increase the state of the output vector if the user did not update its state themself as should have been done */
1701e66621cSBarry Smith   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)w));
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1722bb17772SBarry Smith }
1732bb17772SBarry Smith 
PCPreSolveChangeRHS_Shell(PC pc,PetscBool * change)174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc, PetscBool *change)
175d71ae5a4SJacob Faibussowitsch {
176a06fd7f2SStefano Zampini   PetscFunctionBegin;
177a06fd7f2SStefano Zampini   *change = PETSC_TRUE;
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179a06fd7f2SStefano Zampini }
180a06fd7f2SStefano Zampini 
PCPreSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)181d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolve_Shell(PC pc, KSP ksp, Vec b, Vec x)
182d71ae5a4SJacob Faibussowitsch {
183c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
1847cdd61b2SBarry Smith 
1857cdd61b2SBarry Smith   PetscFunctionBegin;
18628b400f6SJacob Faibussowitsch   PetscCheck(shell->presolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No presolve() routine provided to Shell PC");
18725e27a38SBarry Smith   PetscCallBack("PCSHELL callback presolve", (*shell->presolve)(pc, ksp, b, x));
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1897cdd61b2SBarry Smith }
1907cdd61b2SBarry Smith 
PCPostSolve_Shell(PC pc,KSP ksp,Vec b,Vec x)191d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPostSolve_Shell(PC pc, KSP ksp, Vec b, Vec x)
192d71ae5a4SJacob Faibussowitsch {
193c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
1947cdd61b2SBarry Smith 
1957cdd61b2SBarry Smith   PetscFunctionBegin;
19628b400f6SJacob Faibussowitsch   PetscCheck(shell->postsolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No postsolve() routine provided to Shell PC");
19725e27a38SBarry Smith   PetscCallBack("PCSHELL callback postsolve()", (*shell->postsolve)(pc, ksp, b, x));
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1997cdd61b2SBarry Smith }
2007cdd61b2SBarry Smith 
PCApplyTranspose_Shell(PC pc,Vec x,Vec y)201d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Shell(PC pc, Vec x, Vec y)
202d71ae5a4SJacob Faibussowitsch {
203c5ae4b9aSBarry Smith   PC_Shell        *shell = (PC_Shell *)pc->data;
204e3f487b0SBarry Smith   PetscObjectState instate, outstate;
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith   PetscFunctionBegin;
20728b400f6SJacob Faibussowitsch   PetscCheck(shell->applytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applytranspose() routine provided to Shell PC");
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
20925e27a38SBarry Smith   PetscCallBack("PCSHELL callback applytranspose", (*shell->applytranspose)(pc, x, y));
2109566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
2111e66621cSBarry Smith   /* increase the state of the output vector if the user did not update its state themself as should have been done */
2121e66621cSBarry Smith   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2144b9ad928SBarry Smith }
2154b9ad928SBarry Smith 
PCMatApplyTranspose_Shell(PC pc,Mat x,Mat y)2169fa64a6bSPierre Jolivet static PetscErrorCode PCMatApplyTranspose_Shell(PC pc, Mat x, Mat y)
2179fa64a6bSPierre Jolivet {
2189fa64a6bSPierre Jolivet   PC_Shell        *shell = (PC_Shell *)pc->data;
2199fa64a6bSPierre Jolivet   PetscObjectState instate, outstate;
2209fa64a6bSPierre Jolivet 
2219fa64a6bSPierre Jolivet   PetscFunctionBegin;
2229fa64a6bSPierre Jolivet   PetscCheck(shell->matapplytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No matapplytranspose() routine provided to Shell PC");
2239fa64a6bSPierre Jolivet   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
2249fa64a6bSPierre Jolivet   PetscCallBack("PCSHELL callback matapplytranspose", (*shell->matapplytranspose)(pc, x, y));
2259fa64a6bSPierre Jolivet   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
2269fa64a6bSPierre Jolivet   /* increase the state of the output matrix if the user did not update its state themself as should have been done */
2279fa64a6bSPierre Jolivet   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
2289fa64a6bSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
2299fa64a6bSPierre Jolivet }
2309fa64a6bSPierre Jolivet 
PCApplyRichardson_Shell(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt it,PetscBool guesszero,PetscInt * outits,PCRichardsonConvergedReason * reason)231d71ae5a4SJacob Faibussowitsch 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)
232d71ae5a4SJacob Faibussowitsch {
233c5ae4b9aSBarry Smith   PC_Shell        *shell = (PC_Shell *)pc->data;
234e3f487b0SBarry Smith   PetscObjectState instate, outstate;
2354b9ad928SBarry Smith 
2364b9ad928SBarry Smith   PetscFunctionBegin;
23728b400f6SJacob Faibussowitsch   PetscCheck(shell->applyrich, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyrichardson() routine provided to Shell PC");
2389566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)y, &instate));
23925e27a38SBarry Smith   PetscCallBack("PCSHELL callback applyrichardson", (*shell->applyrich)(pc, x, y, w, rtol, abstol, dtol, it, guesszero, outits, reason));
2409566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)y, &outstate));
241e3f487b0SBarry Smith   /* increase the state of the output vector since the user did not update its state themself as should have been done */
2421e66621cSBarry Smith   if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y));
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2444b9ad928SBarry Smith }
2454b9ad928SBarry Smith 
PCDestroy_Shell(PC pc)246d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Shell(PC pc)
247d71ae5a4SJacob Faibussowitsch {
2484b9ad928SBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith   PetscFunctionBegin;
2519566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell->name));
25225e27a38SBarry Smith   if (shell->destroy) PetscCallBack("PCSHELL callback destroy", (*shell->destroy)(pc));
2539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", NULL));
2549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", NULL));
2559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", NULL));
2569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", NULL));
2579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", NULL));
2589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", NULL));
2599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", NULL));
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", NULL));
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", NULL));
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", NULL));
2639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", NULL));
2649fa64a6bSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApplyTranspose_C", NULL));
2659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", NULL));
2669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", NULL));
2679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", NULL));
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2714b9ad928SBarry Smith }
2724b9ad928SBarry Smith 
PCView_Shell(PC pc,PetscViewer viewer)273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Shell(PC pc, PetscViewer viewer)
274d71ae5a4SJacob Faibussowitsch {
2754b9ad928SBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
2769f196a02SMartin Diehl   PetscBool isascii;
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith   PetscFunctionBegin;
2799f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2809f196a02SMartin Diehl   if (isascii) {
2811e66621cSBarry Smith     if (shell->name) PetscCall(PetscViewerASCIIPrintf(viewer, "  %s\n", shell->name));
2821e66621cSBarry Smith     else PetscCall(PetscViewerASCIIPrintf(viewer, "  no name\n"));
2834b9ad928SBarry Smith   }
2844b9ad928SBarry Smith   if (shell->view) {
2859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2869566063dSJacob Faibussowitsch     PetscCall((*shell->view)(pc, viewer));
2879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
2884b9ad928SBarry Smith   }
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2904b9ad928SBarry Smith }
2914b9ad928SBarry Smith 
PCShellSetDestroy_Shell(PC pc,PetscErrorCode (* destroy)(PC))292d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC))
293d71ae5a4SJacob Faibussowitsch {
294c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
29518be62a5SSatish Balay 
29618be62a5SSatish Balay   PetscFunctionBegin;
29718be62a5SSatish Balay   shell->destroy = destroy;
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29918be62a5SSatish Balay }
30018be62a5SSatish Balay 
PCShellSetSetUp_Shell(PC pc,PetscErrorCode (* setup)(PC))301d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC))
302d71ae5a4SJacob Faibussowitsch {
303feb237baSPierre Jolivet   PC_Shell *shell = (PC_Shell *)pc->data;
3044b9ad928SBarry Smith 
3054b9ad928SBarry Smith   PetscFunctionBegin;
3064b9ad928SBarry Smith   shell->setup = setup;
307d01c8aa3SLisandro Dalcin   if (setup) pc->ops->setup = PCSetUp_Shell;
3080a545947SLisandro Dalcin   else pc->ops->setup = NULL;
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3104b9ad928SBarry Smith }
3114b9ad928SBarry Smith 
PCShellSetApply_Shell(PC pc,PetscErrorCode (* apply)(PC,Vec,Vec))312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApply_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec))
313d71ae5a4SJacob Faibussowitsch {
314c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith   PetscFunctionBegin;
3174b9ad928SBarry Smith   shell->apply = apply;
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3194b9ad928SBarry Smith }
3204b9ad928SBarry Smith 
PCShellSetMatApply_Shell(PC pc,PetscErrorCode (* matapply)(PC,Mat,Mat))321d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetMatApply_Shell(PC pc, PetscErrorCode (*matapply)(PC, Mat, Mat))
322d71ae5a4SJacob Faibussowitsch {
3237b6e2003SPierre Jolivet   PC_Shell *shell = (PC_Shell *)pc->data;
3247b6e2003SPierre Jolivet 
3257b6e2003SPierre Jolivet   PetscFunctionBegin;
3267b6e2003SPierre Jolivet   shell->matapply = matapply;
3270e0fe96bSStefano Zampini   if (matapply) pc->ops->matapply = PCMatApply_Shell;
3280e0fe96bSStefano Zampini   else pc->ops->matapply = NULL;
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3307b6e2003SPierre Jolivet }
3317b6e2003SPierre Jolivet 
PCShellSetApplySymmetricLeft_Shell(PC pc,PetscErrorCode (* apply)(PC,Vec,Vec))332d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplySymmetricLeft_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec))
333d71ae5a4SJacob Faibussowitsch {
3341b581b66SBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3351b581b66SBarry Smith 
3361b581b66SBarry Smith   PetscFunctionBegin;
3371b581b66SBarry Smith   shell->applysymmetricleft = apply;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3391b581b66SBarry Smith }
3401b581b66SBarry Smith 
PCShellSetApplySymmetricRight_Shell(PC pc,PetscErrorCode (* apply)(PC,Vec,Vec))341d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplySymmetricRight_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec))
342d71ae5a4SJacob Faibussowitsch {
3431b581b66SBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3441b581b66SBarry Smith 
3451b581b66SBarry Smith   PetscFunctionBegin;
3461b581b66SBarry Smith   shell->applysymmetricright = apply;
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3481b581b66SBarry Smith }
3491b581b66SBarry Smith 
PCShellSetApplyBA_Shell(PC pc,PetscErrorCode (* applyBA)(PC,PCSide,Vec,Vec,Vec))350d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyBA_Shell(PC pc, PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec))
351d71ae5a4SJacob Faibussowitsch {
352c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3532bb17772SBarry Smith 
3542bb17772SBarry Smith   PetscFunctionBegin;
355d01c8aa3SLisandro Dalcin   shell->applyBA = applyBA;
356d01c8aa3SLisandro Dalcin   if (applyBA) pc->ops->applyBA = PCApplyBA_Shell;
3570a545947SLisandro Dalcin   else pc->ops->applyBA = NULL;
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3592bb17772SBarry Smith }
3602bb17772SBarry Smith 
PCShellSetPreSolve_Shell(PC pc,PCShellPSolveFn * presolve)3614d4d2bdcSBarry Smith static PetscErrorCode PCShellSetPreSolve_Shell(PC pc, PCShellPSolveFn *presolve)
362d71ae5a4SJacob Faibussowitsch {
363c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3647cdd61b2SBarry Smith 
3657cdd61b2SBarry Smith   PetscFunctionBegin;
3667cdd61b2SBarry Smith   shell->presolve = presolve;
367a06fd7f2SStefano Zampini   if (presolve) {
368a06fd7f2SStefano Zampini     pc->ops->presolve = PCPreSolve_Shell;
3699566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Shell));
370a06fd7f2SStefano Zampini   } else {
3710a545947SLisandro Dalcin     pc->ops->presolve = NULL;
3729566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL));
373a06fd7f2SStefano Zampini   }
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3757cdd61b2SBarry Smith }
3767cdd61b2SBarry Smith 
PCShellSetPostSolve_Shell(PC pc,PCShellPSolveFn * postsolve)3774d4d2bdcSBarry Smith static PetscErrorCode PCShellSetPostSolve_Shell(PC pc, PCShellPSolveFn *postsolve)
378d71ae5a4SJacob Faibussowitsch {
379c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3807cdd61b2SBarry Smith 
3817cdd61b2SBarry Smith   PetscFunctionBegin;
3827cdd61b2SBarry Smith   shell->postsolve = postsolve;
383d01c8aa3SLisandro Dalcin   if (postsolve) pc->ops->postsolve = PCPostSolve_Shell;
3840a545947SLisandro Dalcin   else pc->ops->postsolve = NULL;
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3867cdd61b2SBarry Smith }
3877cdd61b2SBarry Smith 
PCShellSetView_Shell(PC pc,PetscErrorCode (* view)(PC,PetscViewer))388d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetView_Shell(PC pc, PetscErrorCode (*view)(PC, PetscViewer))
389d71ae5a4SJacob Faibussowitsch {
390c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith   PetscFunctionBegin;
3934b9ad928SBarry Smith   shell->view = view;
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3954b9ad928SBarry Smith }
3964b9ad928SBarry Smith 
PCShellSetApplyTranspose_Shell(PC pc,PetscErrorCode (* applytranspose)(PC,Vec,Vec))397d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyTranspose_Shell(PC pc, PetscErrorCode (*applytranspose)(PC, Vec, Vec))
398d71ae5a4SJacob Faibussowitsch {
399c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith   PetscFunctionBegin;
4024b9ad928SBarry Smith   shell->applytranspose = applytranspose;
403d01c8aa3SLisandro Dalcin   if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell;
4040a545947SLisandro Dalcin   else pc->ops->applytranspose = NULL;
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406d01c8aa3SLisandro Dalcin }
407d01c8aa3SLisandro Dalcin 
PCShellSetMatApplyTranspose_Shell(PC pc,PetscErrorCode (* matapplytranspose)(PC,Mat,Mat))4089fa64a6bSPierre Jolivet static PetscErrorCode PCShellSetMatApplyTranspose_Shell(PC pc, PetscErrorCode (*matapplytranspose)(PC, Mat, Mat))
4099fa64a6bSPierre Jolivet {
4109fa64a6bSPierre Jolivet   PC_Shell *shell = (PC_Shell *)pc->data;
4119fa64a6bSPierre Jolivet 
4129fa64a6bSPierre Jolivet   PetscFunctionBegin;
4139fa64a6bSPierre Jolivet   shell->matapplytranspose = matapplytranspose;
4149fa64a6bSPierre Jolivet   if (matapplytranspose) pc->ops->matapplytranspose = PCMatApplyTranspose_Shell;
4159fa64a6bSPierre Jolivet   else pc->ops->matapplytranspose = NULL;
4169fa64a6bSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
4179fa64a6bSPierre Jolivet }
4189fa64a6bSPierre Jolivet 
PCShellSetApplyRichardson_Shell(PC pc,PetscErrorCode (* applyrich)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt *,PCRichardsonConvergedReason *))419d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyRichardson_Shell(PC pc, PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *))
420d71ae5a4SJacob Faibussowitsch {
421c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
422d01c8aa3SLisandro Dalcin 
423d01c8aa3SLisandro Dalcin   PetscFunctionBegin;
424d01c8aa3SLisandro Dalcin   shell->applyrich = applyrich;
425d01c8aa3SLisandro Dalcin   if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell;
4260a545947SLisandro Dalcin   else pc->ops->applyrichardson = NULL;
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4284b9ad928SBarry Smith }
4294b9ad928SBarry Smith 
PCShellSetName_Shell(PC pc,const char name[])430d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetName_Shell(PC pc, const char name[])
431d71ae5a4SJacob Faibussowitsch {
432c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith   PetscFunctionBegin;
4359566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell->name));
4369566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &shell->name));
4373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4384b9ad928SBarry Smith }
4394b9ad928SBarry Smith 
PCShellGetName_Shell(PC pc,const char * name[])440d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellGetName_Shell(PC pc, const char *name[])
441d71ae5a4SJacob Faibussowitsch {
442c5ae4b9aSBarry Smith   PC_Shell *shell = (PC_Shell *)pc->data;
4434b9ad928SBarry Smith 
4444b9ad928SBarry Smith   PetscFunctionBegin;
4454b9ad928SBarry Smith   *name = shell->name;
4463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4474b9ad928SBarry Smith }
4484b9ad928SBarry Smith 
44918be62a5SSatish Balay /*@C
45004c3f3b8SBarry Smith   PCShellSetDestroy - Sets routine to use to destroy the user-provided application context that was provided with `PCShellSetContext()`
45118be62a5SSatish Balay 
452c3339decSBarry Smith   Logically Collective
45318be62a5SSatish Balay 
45418be62a5SSatish Balay   Input Parameters:
45518be62a5SSatish Balay + pc      - the preconditioner context
456a2b725a8SWilliam Gropp - destroy - the application-provided destroy routine
45718be62a5SSatish Balay 
45820f4b53cSBarry Smith   Calling sequence of `destroy`:
45904c3f3b8SBarry Smith . pc - the preconditioner
4604aa34b0aSBarry Smith 
461f1580f4eSBarry Smith   Level: intermediate
46218be62a5SSatish Balay 
463562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()`
46418be62a5SSatish Balay @*/
PCShellSetDestroy(PC pc,PetscErrorCode (* destroy)(PC pc))46504c3f3b8SBarry Smith PetscErrorCode PCShellSetDestroy(PC pc, PetscErrorCode (*destroy)(PC pc))
466d71ae5a4SJacob Faibussowitsch {
46718be62a5SSatish Balay   PetscFunctionBegin;
4680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
469cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetDestroy_C", (PC, PetscErrorCode (*)(PC)), (pc, destroy));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47118be62a5SSatish Balay }
47218be62a5SSatish Balay 
4734b9ad928SBarry Smith /*@C
4744b9ad928SBarry Smith   PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the
4754b9ad928SBarry Smith   matrix operator is changed.
4764b9ad928SBarry Smith 
477c3339decSBarry Smith   Logically Collective
4784b9ad928SBarry Smith 
4794b9ad928SBarry Smith   Input Parameters:
4804b9ad928SBarry Smith + pc    - the preconditioner context
481a2b725a8SWilliam Gropp - setup - the application-provided setup routine
4824b9ad928SBarry Smith 
48320f4b53cSBarry Smith   Calling sequence of `setup`:
48404c3f3b8SBarry Smith . pc - the preconditioner
4854aa34b0aSBarry Smith 
486f1580f4eSBarry Smith   Level: intermediate
4874b9ad928SBarry Smith 
48804c3f3b8SBarry Smith   Note:
48904c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `setup`.
49004c3f3b8SBarry Smith 
491562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetApply()`, `PCShellSetContext()`, , `PCShellGetContext()`
4924b9ad928SBarry Smith @*/
PCShellSetSetUp(PC pc,PetscErrorCode (* setup)(PC pc))49304c3f3b8SBarry Smith PetscErrorCode PCShellSetSetUp(PC pc, PetscErrorCode (*setup)(PC pc))
494d71ae5a4SJacob Faibussowitsch {
4954b9ad928SBarry Smith   PetscFunctionBegin;
4960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
497cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetSetUp_C", (PC, PetscErrorCode (*)(PC)), (pc, setup));
4983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4994b9ad928SBarry Smith }
5004b9ad928SBarry Smith 
5014b9ad928SBarry Smith /*@C
50204c3f3b8SBarry Smith   PCShellSetView - Sets routine to use as viewer of a `PCSHELL` shell preconditioner
5034b9ad928SBarry Smith 
504c3339decSBarry Smith   Logically Collective
5054b9ad928SBarry Smith 
5064b9ad928SBarry Smith   Input Parameters:
5074b9ad928SBarry Smith + pc   - the preconditioner context
5084b9ad928SBarry Smith - view - the application-provided view routine
5094b9ad928SBarry Smith 
51020f4b53cSBarry Smith   Calling sequence of `view`:
51104c3f3b8SBarry Smith + pc - the preconditioner
5124b9ad928SBarry Smith - v  - viewer
5134b9ad928SBarry Smith 
514f1580f4eSBarry Smith   Level: advanced
5154b9ad928SBarry Smith 
51604c3f3b8SBarry Smith   Note:
51704c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `view`.
51804c3f3b8SBarry Smith 
519562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()`
5204b9ad928SBarry Smith @*/
PCShellSetView(PC pc,PetscErrorCode (* view)(PC pc,PetscViewer v))52104c3f3b8SBarry Smith PetscErrorCode PCShellSetView(PC pc, PetscErrorCode (*view)(PC pc, PetscViewer v))
522d71ae5a4SJacob Faibussowitsch {
5234b9ad928SBarry Smith   PetscFunctionBegin;
5240700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
525cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetView_C", (PC, PetscErrorCode (*)(PC, PetscViewer)), (pc, view));
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5274b9ad928SBarry Smith }
5284b9ad928SBarry Smith 
5294b9ad928SBarry Smith /*@C
5304b9ad928SBarry Smith   PCShellSetApply - Sets routine to use as preconditioner.
5314b9ad928SBarry Smith 
532c3339decSBarry Smith   Logically Collective
5334b9ad928SBarry Smith 
5344b9ad928SBarry Smith   Input Parameters:
5354b9ad928SBarry Smith + pc    - the preconditioner context
536be29d3c6SBarry Smith - apply - the application-provided preconditioning routine
5374b9ad928SBarry Smith 
53820f4b53cSBarry Smith   Calling sequence of `apply`:
53920f4b53cSBarry Smith + pc   - the preconditioner, get the application context with `PCShellGetContext()`
5404b9ad928SBarry Smith . xin  - input vector
5414b9ad928SBarry Smith - xout - output vector
5424b9ad928SBarry Smith 
543f1580f4eSBarry Smith   Level: intermediate
5444b9ad928SBarry Smith 
54504c3f3b8SBarry Smith   Note:
54604c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`.
54704c3f3b8SBarry Smith 
548562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellSetApplySymmetricRight()`, `PCShellSetApplySymmetricLeft()`, `PCShellGetContext()`
5494b9ad928SBarry Smith @*/
PCShellSetApply(PC pc,PetscErrorCode (* apply)(PC pc,Vec xin,Vec xout))55004c3f3b8SBarry Smith PetscErrorCode PCShellSetApply(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout))
551d71ae5a4SJacob Faibussowitsch {
5524b9ad928SBarry Smith   PetscFunctionBegin;
5530700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
554cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetApply_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply));
5553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5564b9ad928SBarry Smith }
5574b9ad928SBarry Smith 
5581b581b66SBarry Smith /*@C
5596437efc7SEric Chamberland   PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors.
5607b6e2003SPierre Jolivet 
561c3339decSBarry Smith   Logically Collective
5627b6e2003SPierre Jolivet 
5637b6e2003SPierre Jolivet   Input Parameters:
5647b6e2003SPierre Jolivet + pc       - the preconditioner context
56504c3f3b8SBarry Smith - matapply - the application-provided preconditioning routine
5667b6e2003SPierre Jolivet 
56704c3f3b8SBarry Smith   Calling sequence of `matapply`:
56804c3f3b8SBarry Smith + pc   - the preconditioner
56904c3f3b8SBarry Smith . Xin  - input block of vectors represented as a dense `Mat`
57004c3f3b8SBarry Smith - Xout - output block of vectors represented as a dense `Mat`
5717b6e2003SPierre Jolivet 
572f1580f4eSBarry Smith   Level: advanced
5737b6e2003SPierre Jolivet 
57404c3f3b8SBarry Smith   Note:
57504c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapply`.
57604c3f3b8SBarry Smith 
577562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()`
5787b6e2003SPierre Jolivet @*/
PCShellSetMatApply(PC pc,PetscErrorCode (* matapply)(PC pc,Mat Xin,Mat Xout))57904c3f3b8SBarry Smith PetscErrorCode PCShellSetMatApply(PC pc, PetscErrorCode (*matapply)(PC pc, Mat Xin, Mat Xout))
580d71ae5a4SJacob Faibussowitsch {
5817b6e2003SPierre Jolivet   PetscFunctionBegin;
5827b6e2003SPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
583cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetMatApply_C", (PC, PetscErrorCode (*)(PC, Mat, Mat)), (pc, matapply));
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5857b6e2003SPierre Jolivet }
5867b6e2003SPierre Jolivet 
5877b6e2003SPierre Jolivet /*@C
588f1580f4eSBarry Smith   PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the `PC_SYMMETRIC` is used).
5891b581b66SBarry Smith 
590c3339decSBarry Smith   Logically Collective
5911b581b66SBarry Smith 
5921b581b66SBarry Smith   Input Parameters:
5931b581b66SBarry Smith + pc    - the preconditioner context
5941b581b66SBarry Smith - apply - the application-provided left preconditioning routine
5951b581b66SBarry Smith 
59620f4b53cSBarry Smith   Calling sequence of `apply`:
59704c3f3b8SBarry Smith + pc   - the preconditioner
5981b581b66SBarry Smith . xin  - input vector
5991b581b66SBarry Smith - xout - output vector
6001b581b66SBarry Smith 
601f1580f4eSBarry Smith   Level: advanced
6021b581b66SBarry Smith 
60304c3f3b8SBarry Smith   Note:
60404c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`.
60504c3f3b8SBarry Smith 
606562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`
6071b581b66SBarry Smith @*/
PCShellSetApplySymmetricLeft(PC pc,PetscErrorCode (* apply)(PC pc,Vec xin,Vec xout))60804c3f3b8SBarry Smith PetscErrorCode PCShellSetApplySymmetricLeft(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout))
609d71ae5a4SJacob Faibussowitsch {
6101b581b66SBarry Smith   PetscFunctionBegin;
6111b581b66SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
612cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetApplySymmetricLeft_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply));
6133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6141b581b66SBarry Smith }
6151b581b66SBarry Smith 
6161b581b66SBarry Smith /*@C
617f1580f4eSBarry Smith   PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the `PC_SYMMETRIC` is used).
6181b581b66SBarry Smith 
619c3339decSBarry Smith   Logically Collective
6201b581b66SBarry Smith 
6211b581b66SBarry Smith   Input Parameters:
6221b581b66SBarry Smith + pc    - the preconditioner context
6231b581b66SBarry Smith - apply - the application-provided right preconditioning routine
6241b581b66SBarry Smith 
62520f4b53cSBarry Smith   Calling sequence of `apply`:
62604c3f3b8SBarry Smith + pc   - the preconditioner
6271b581b66SBarry Smith . xin  - input vector
6281b581b66SBarry Smith - xout - output vector
6291b581b66SBarry Smith 
630f1580f4eSBarry Smith   Level: advanced
6311b581b66SBarry Smith 
63204c3f3b8SBarry Smith   Note:
63304c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`.
63404c3f3b8SBarry Smith 
635562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()`
6361b581b66SBarry Smith @*/
PCShellSetApplySymmetricRight(PC pc,PetscErrorCode (* apply)(PC pc,Vec xin,Vec xout))63704c3f3b8SBarry Smith PetscErrorCode PCShellSetApplySymmetricRight(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout))
638d71ae5a4SJacob Faibussowitsch {
6391b581b66SBarry Smith   PetscFunctionBegin;
6401b581b66SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
641cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetApplySymmetricRight_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply));
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6431b581b66SBarry Smith }
6441b581b66SBarry Smith 
6452bb17772SBarry Smith /*@C
64604c3f3b8SBarry Smith   PCShellSetApplyBA - Sets routine to use as the preconditioner times the operator.
6472bb17772SBarry Smith 
648c3339decSBarry Smith   Logically Collective
6492bb17772SBarry Smith 
6502bb17772SBarry Smith   Input Parameters:
6512bb17772SBarry Smith + pc      - the preconditioner context
6522bb17772SBarry Smith - applyBA - the application-provided BA routine
6532bb17772SBarry Smith 
65420f4b53cSBarry Smith   Calling sequence of `applyBA`:
65504c3f3b8SBarry Smith + pc   - the preconditioner
65604c3f3b8SBarry Smith . side - `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
6572bb17772SBarry Smith . xin  - input vector
65804c3f3b8SBarry Smith . xout - output vector
65904c3f3b8SBarry Smith - w    - work vector
6602bb17772SBarry Smith 
661f1580f4eSBarry Smith   Level: intermediate
6622bb17772SBarry Smith 
66304c3f3b8SBarry Smith   Note:
66404c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applyBA`.
66504c3f3b8SBarry Smith 
666562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellGetContext()`, `PCSide`
6672bb17772SBarry Smith @*/
PCShellSetApplyBA(PC pc,PetscErrorCode (* applyBA)(PC pc,PCSide side,Vec xin,Vec xout,Vec w))66804c3f3b8SBarry Smith PetscErrorCode PCShellSetApplyBA(PC pc, PetscErrorCode (*applyBA)(PC pc, PCSide side, Vec xin, Vec xout, Vec w))
669d71ae5a4SJacob Faibussowitsch {
6702bb17772SBarry Smith   PetscFunctionBegin;
6710700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
672cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetApplyBA_C", (PC, PetscErrorCode (*)(PC, PCSide, Vec, Vec, Vec)), (pc, applyBA));
6733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6742bb17772SBarry Smith }
6752bb17772SBarry Smith 
6764b9ad928SBarry Smith /*@C
6774b9ad928SBarry Smith   PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose.
6784b9ad928SBarry Smith 
679c3339decSBarry Smith   Logically Collective
6804b9ad928SBarry Smith 
6814b9ad928SBarry Smith   Input Parameters:
6824b9ad928SBarry Smith + pc             - the preconditioner context
68320f4b53cSBarry Smith - applytranspose - the application-provided preconditioning transpose routine
6844b9ad928SBarry Smith 
68520f4b53cSBarry Smith   Calling sequence of `applytranspose`:
68604c3f3b8SBarry Smith + pc   - the preconditioner
6874b9ad928SBarry Smith . xin  - input vector
6884b9ad928SBarry Smith - xout - output vector
6894b9ad928SBarry Smith 
690f1580f4eSBarry Smith   Level: intermediate
6914b9ad928SBarry Smith 
692f1580f4eSBarry Smith   Note:
69304c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applytranspose`.
6944b9ad928SBarry Smith 
69554c05997SPierre Jolivet .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellGetContext()`
6964b9ad928SBarry Smith @*/
PCShellSetApplyTranspose(PC pc,PetscErrorCode (* applytranspose)(PC pc,Vec xin,Vec xout))69704c3f3b8SBarry Smith PetscErrorCode PCShellSetApplyTranspose(PC pc, PetscErrorCode (*applytranspose)(PC pc, Vec xin, Vec xout))
698d71ae5a4SJacob Faibussowitsch {
6994b9ad928SBarry Smith   PetscFunctionBegin;
7000700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
701cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetApplyTranspose_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, applytranspose));
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7034b9ad928SBarry Smith }
7044b9ad928SBarry Smith 
7057cdd61b2SBarry Smith /*@C
7069fa64a6bSPierre Jolivet   PCShellSetMatApplyTranspose - Sets routine to use as preconditioner transpose.
7079fa64a6bSPierre Jolivet 
7089fa64a6bSPierre Jolivet   Logically Collective
7099fa64a6bSPierre Jolivet 
7109fa64a6bSPierre Jolivet   Input Parameters:
7119fa64a6bSPierre Jolivet + pc                - the preconditioner context
7129fa64a6bSPierre Jolivet - matapplytranspose - the application-provided preconditioning transpose routine
7139fa64a6bSPierre Jolivet 
7149fa64a6bSPierre Jolivet   Calling sequence of `matapplytranspose`:
7159fa64a6bSPierre Jolivet + pc   - the preconditioner
7169fa64a6bSPierre Jolivet . xin  - input matrix
7179fa64a6bSPierre Jolivet - xout - output matrix
7189fa64a6bSPierre Jolivet 
7199fa64a6bSPierre Jolivet   Level: intermediate
7209fa64a6bSPierre Jolivet 
7219fa64a6bSPierre Jolivet   Note:
7229fa64a6bSPierre Jolivet   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapplytranspose`.
7239fa64a6bSPierre Jolivet 
7249fa64a6bSPierre Jolivet .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellGetContext()`
7259fa64a6bSPierre Jolivet @*/
PCShellSetMatApplyTranspose(PC pc,PetscErrorCode (* matapplytranspose)(PC pc,Mat xin,Mat xout))7269fa64a6bSPierre Jolivet PetscErrorCode PCShellSetMatApplyTranspose(PC pc, PetscErrorCode (*matapplytranspose)(PC pc, Mat xin, Mat xout))
7279fa64a6bSPierre Jolivet {
7289fa64a6bSPierre Jolivet   PetscFunctionBegin;
7299fa64a6bSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7309fa64a6bSPierre Jolivet   PetscTryMethod(pc, "PCShellSetMatApplyTranspose_C", (PC, PetscErrorCode (*)(PC, Mat, Mat)), (pc, matapplytranspose));
7319fa64a6bSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
7329fa64a6bSPierre Jolivet }
7339fa64a6bSPierre Jolivet 
7349fa64a6bSPierre Jolivet /*@C
735f1580f4eSBarry Smith   PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a `KSPSolve()` is
7367cdd61b2SBarry Smith   applied. This usually does something like scale the linear system in some application
7377cdd61b2SBarry Smith   specific way.
7387cdd61b2SBarry Smith 
739c3339decSBarry Smith   Logically Collective
7407cdd61b2SBarry Smith 
7417cdd61b2SBarry Smith   Input Parameters:
7427cdd61b2SBarry Smith + pc       - the preconditioner context
7434d4d2bdcSBarry Smith - presolve - the application-provided presolve routine, see `PCShellPSolveFn`
7447cdd61b2SBarry Smith 
745f1580f4eSBarry Smith   Level: advanced
7467cdd61b2SBarry Smith 
74704c3f3b8SBarry Smith   Note:
74804c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `presolve`.
74904c3f3b8SBarry Smith 
7504d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellPSolveFn`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPostSolve()`, `PCShellSetContext()`, `PCShellGetContext()`
7517cdd61b2SBarry Smith @*/
PCShellSetPreSolve(PC pc,PCShellPSolveFn * presolve)7524d4d2bdcSBarry Smith PetscErrorCode PCShellSetPreSolve(PC pc, PCShellPSolveFn *presolve)
753d71ae5a4SJacob Faibussowitsch {
7547cdd61b2SBarry Smith   PetscFunctionBegin;
7550700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7564d4d2bdcSBarry Smith   PetscTryMethod(pc, "PCShellSetPreSolve_C", (PC, PCShellPSolveFn *), (pc, presolve));
7573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7587cdd61b2SBarry Smith }
7597cdd61b2SBarry Smith 
7607cdd61b2SBarry Smith /*@C
76104c3f3b8SBarry Smith   PCShellSetPostSolve - Sets routine to apply to the operators/vectors after a `KSPSolve()` is
7627cdd61b2SBarry Smith   applied. This usually does something like scale the linear system in some application
7637cdd61b2SBarry Smith   specific way.
7647cdd61b2SBarry Smith 
765c3339decSBarry Smith   Logically Collective
7667cdd61b2SBarry Smith 
7677cdd61b2SBarry Smith   Input Parameters:
7687cdd61b2SBarry Smith + pc        - the preconditioner context
7694d4d2bdcSBarry Smith - postsolve - the application-provided postsolve routine, see `PCShellPSolveFn`
7707cdd61b2SBarry Smith 
771f1580f4eSBarry Smith   Level: advanced
7727cdd61b2SBarry Smith 
77304c3f3b8SBarry Smith   Note:
77404c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `postsolve`.
77504c3f3b8SBarry Smith 
7764d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellPSolveFn`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPreSolve()`, `PCShellSetContext()`, `PCShellGetContext()`
7777cdd61b2SBarry Smith @*/
PCShellSetPostSolve(PC pc,PCShellPSolveFn * postsolve)7784d4d2bdcSBarry Smith PetscErrorCode PCShellSetPostSolve(PC pc, PCShellPSolveFn *postsolve)
779d71ae5a4SJacob Faibussowitsch {
7807cdd61b2SBarry Smith   PetscFunctionBegin;
7810700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7824d4d2bdcSBarry Smith   PetscTryMethod(pc, "PCShellSetPostSolve_C", (PC, PCShellPSolveFn *), (pc, postsolve));
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7847cdd61b2SBarry Smith }
7857cdd61b2SBarry Smith 
786cc4c1da9SBarry Smith /*@
787f1580f4eSBarry Smith   PCShellSetName - Sets an optional name to associate with a `PCSHELL`
7884b9ad928SBarry Smith   preconditioner.
7894b9ad928SBarry Smith 
7904b9ad928SBarry Smith   Not Collective
7914b9ad928SBarry Smith 
7924b9ad928SBarry Smith   Input Parameters:
7934b9ad928SBarry Smith + pc   - the preconditioner context
7944b9ad928SBarry Smith - name - character string describing shell preconditioner
7954b9ad928SBarry Smith 
796f1580f4eSBarry Smith   Level: intermediate
7974b9ad928SBarry Smith 
79804c3f3b8SBarry Smith   Note:
799baca6076SPierre Jolivet   This is separate from the name you can provide with `PetscObjectSetName()`
80004c3f3b8SBarry Smith 
801562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellGetName()`, `PetscObjectSetName()`, `PetscObjectGetName()`
8024b9ad928SBarry Smith @*/
PCShellSetName(PC pc,const char name[])803d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellSetName(PC pc, const char name[])
804d71ae5a4SJacob Faibussowitsch {
8054b9ad928SBarry Smith   PetscFunctionBegin;
8060700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
807cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetName_C", (PC, const char[]), (pc, name));
8083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8094b9ad928SBarry Smith }
8104b9ad928SBarry Smith 
811cc4c1da9SBarry Smith /*@
81204c3f3b8SBarry Smith   PCShellGetName - Gets an optional name that the user has set for a `PCSHELL` with `PCShellSetName()`
8134b9ad928SBarry Smith   preconditioner.
8144b9ad928SBarry Smith 
8154b9ad928SBarry Smith   Not Collective
8164b9ad928SBarry Smith 
8174b9ad928SBarry Smith   Input Parameter:
8184b9ad928SBarry Smith . pc - the preconditioner context
8194b9ad928SBarry Smith 
8204b9ad928SBarry Smith   Output Parameter:
8214b9ad928SBarry Smith . name - character string describing shell preconditioner (you should not free this)
8224b9ad928SBarry Smith 
823f1580f4eSBarry Smith   Level: intermediate
8244b9ad928SBarry Smith 
825562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetName()`, `PetscObjectSetName()`, `PetscObjectGetName()`
8264b9ad928SBarry Smith @*/
PCShellGetName(PC pc,const char * name[])827d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellGetName(PC pc, const char *name[])
828d71ae5a4SJacob Faibussowitsch {
8294b9ad928SBarry Smith   PetscFunctionBegin;
8300700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8314f572ea9SToby Isaac   PetscAssertPointer(name, 2);
832cac4c232SBarry Smith   PetscUseMethod(pc, "PCShellGetName_C", (PC, const char *[]), (pc, name));
8333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8344b9ad928SBarry Smith }
8354b9ad928SBarry Smith 
8364b9ad928SBarry Smith /*@C
8374b9ad928SBarry Smith   PCShellSetApplyRichardson - Sets routine to use as preconditioner
8384b9ad928SBarry Smith   in Richardson iteration.
8394b9ad928SBarry Smith 
840c3339decSBarry Smith   Logically Collective
8414b9ad928SBarry Smith 
8424b9ad928SBarry Smith   Input Parameters:
8434b9ad928SBarry Smith + pc    - the preconditioner context
844be29d3c6SBarry Smith - apply - the application-provided preconditioning routine
8454b9ad928SBarry Smith 
84620f4b53cSBarry Smith   Calling sequence of `apply`:
84704c3f3b8SBarry Smith + pc               - the preconditioner
848dd8e379bSPierre Jolivet . b                - right-hand side
8494b9ad928SBarry Smith . x                - current iterate
8504b9ad928SBarry Smith . r                - work space
8514b9ad928SBarry Smith . rtol             - relative tolerance of residual norm to stop at
85270441072SBarry Smith . abstol           - absolute tolerance of residual norm to stop at
8534b9ad928SBarry Smith . dtol             - if residual norm increases by this factor than return
85404c3f3b8SBarry Smith . maxits           - number of iterations to run
85504c3f3b8SBarry Smith . zeroinitialguess - `PETSC_TRUE` if `x` is known to be initially zero
85604c3f3b8SBarry Smith . its              - returns the number of iterations used
85704c3f3b8SBarry Smith - reason           - returns the reason the iteration has converged
8584b9ad928SBarry Smith 
859f1580f4eSBarry Smith   Level: advanced
8604b9ad928SBarry Smith 
8610b4b7b1cSBarry Smith   Notes:
86204c3f3b8SBarry Smith   You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`.
86304c3f3b8SBarry Smith 
8640b4b7b1cSBarry Smith   This is used when one can provide code for multiple steps of Richardson's method that is more efficient than computing a single step,
8650b4b7b1cSBarry Smith   recomputing the residual via $ r = b - A x $, and then computing the next step. SOR is an algorithm for which this is true.
8660b4b7b1cSBarry Smith 
8670b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCRichardsonConvergedReason()`, `PCShellGetContext()`, `KSPRICHARDSON`
8684b9ad928SBarry Smith @*/
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))86904c3f3b8SBarry Smith 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))
870d71ae5a4SJacob Faibussowitsch {
8714b9ad928SBarry Smith   PetscFunctionBegin;
8720700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
873cac4c232SBarry Smith   PetscTryMethod(pc, "PCShellSetApplyRichardson_C", (PC, PetscErrorCode (*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)), (pc, apply));
8743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8754b9ad928SBarry Smith }
8764b9ad928SBarry Smith 
8774b9ad928SBarry Smith /*MC
878f1580f4eSBarry Smith   PCSHELL - Creates a new preconditioner class for use with a users
879f1580f4eSBarry Smith             own private data storage format and preconditioner application code
8804b9ad928SBarry Smith 
8814b9ad928SBarry Smith   Level: advanced
882e0bb08deSStefano Zampini 
8834b9ad928SBarry Smith   Usage:
884f1580f4eSBarry Smith .vb
885f1580f4eSBarry Smith   extern PetscErrorCode apply(PC,Vec,Vec);
886f1580f4eSBarry Smith   extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec);
887f1580f4eSBarry Smith   extern PetscErrorCode applytranspose(PC,Vec,Vec);
888f1580f4eSBarry Smith   extern PetscErrorCode setup(PC);
889f1580f4eSBarry Smith   extern PetscErrorCode destroy(PC);
890f1580f4eSBarry Smith 
891f1580f4eSBarry Smith   PCCreate(comm,&pc);
892f1580f4eSBarry Smith   PCSetType(pc,PCSHELL);
893f1580f4eSBarry Smith   PCShellSetContext(pc,ctx)
894f1580f4eSBarry Smith   PCShellSetApply(pc,apply);
895f1580f4eSBarry Smith   PCShellSetApplyBA(pc,applyba);               (optional)
896f1580f4eSBarry Smith   PCShellSetApplyTranspose(pc,applytranspose); (optional)
897f1580f4eSBarry Smith   PCShellSetSetUp(pc,setup);                   (optional)
898f1580f4eSBarry Smith   PCShellSetDestroy(pc,destroy);               (optional)
899f1580f4eSBarry Smith .ve
9004b9ad928SBarry Smith 
9010b4b7b1cSBarry Smith   Notes:
90204c3f3b8SBarry Smith   Information required for the preconditioner and its internal datastructures can be set with `PCShellSetContext()` and then accessed
9030b4b7b1cSBarry Smith   with `PCShellGetContext()` inside the routines provided above.
9040b4b7b1cSBarry Smith 
9050b4b7b1cSBarry Smith   When using `MATSHELL`, where the explicit entries of matrix are not available to build the preconditioner, `PCSHELL` can be used
9060b4b7b1cSBarry Smith   to construct a custom preconditioner for the `MATSHELL`, assuming the user knows enough about their problem to provide a
9070b4b7b1cSBarry Smith   custom preconditioner.
90804c3f3b8SBarry Smith 
909562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
910f1580f4eSBarry Smith           `MATSHELL`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetView()`, `PCShellSetDestroy()`, `PCShellSetPostSolve()`,
911f1580f4eSBarry Smith           `PCShellSetApplyTranspose()`, `PCShellSetName()`, `PCShellSetApplyRichardson()`, `PCShellSetPreSolve()`, `PCShellSetView()`,
912f1580f4eSBarry Smith           `PCShellGetName()`, `PCShellSetContext()`, `PCShellGetContext()`, `PCShellSetApplyBA()`, `MATSHELL`, `PCShellSetMatApply()`,
9134b9ad928SBarry Smith M*/
9144b9ad928SBarry Smith 
PCCreate_Shell(PC pc)915d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc)
916d71ae5a4SJacob Faibussowitsch {
9174b9ad928SBarry Smith   PC_Shell *shell;
9184b9ad928SBarry Smith 
9194b9ad928SBarry Smith   PetscFunctionBegin;
9204dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&shell));
9214b9ad928SBarry Smith   pc->data = (void *)shell;
9224b9ad928SBarry Smith 
923d01c8aa3SLisandro Dalcin   pc->ops->destroy             = PCDestroy_Shell;
9244b9ad928SBarry Smith   pc->ops->view                = PCView_Shell;
925d01c8aa3SLisandro Dalcin   pc->ops->apply               = PCApply_Shell;
9261b581b66SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Shell;
9271b581b66SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricRight_Shell;
9280e0fe96bSStefano Zampini   pc->ops->matapply            = NULL;
9290a545947SLisandro Dalcin   pc->ops->applytranspose      = NULL;
9300a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
9310a545947SLisandro Dalcin   pc->ops->setup               = NULL;
9320a545947SLisandro Dalcin   pc->ops->presolve            = NULL;
9330a545947SLisandro Dalcin   pc->ops->postsolve           = NULL;
9344b9ad928SBarry Smith 
9350a545947SLisandro Dalcin   shell->apply               = NULL;
9360a545947SLisandro Dalcin   shell->applytranspose      = NULL;
9370a545947SLisandro Dalcin   shell->name                = NULL;
9380a545947SLisandro Dalcin   shell->applyrich           = NULL;
9390a545947SLisandro Dalcin   shell->presolve            = NULL;
9400a545947SLisandro Dalcin   shell->postsolve           = NULL;
9410a545947SLisandro Dalcin   shell->ctx                 = NULL;
9420a545947SLisandro Dalcin   shell->setup               = NULL;
9430a545947SLisandro Dalcin   shell->view                = NULL;
9440a545947SLisandro Dalcin   shell->destroy             = NULL;
9450a545947SLisandro Dalcin   shell->applysymmetricleft  = NULL;
9460a545947SLisandro Dalcin   shell->applysymmetricright = NULL;
9474b9ad928SBarry Smith 
9489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", PCShellSetDestroy_Shell));
9499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", PCShellSetSetUp_Shell));
9509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", PCShellSetApply_Shell));
9519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", PCShellSetMatApply_Shell));
9529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", PCShellSetApplySymmetricLeft_Shell));
9539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", PCShellSetApplySymmetricRight_Shell));
9549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", PCShellSetApplyBA_Shell));
9559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", PCShellSetPreSolve_Shell));
9569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", PCShellSetPostSolve_Shell));
9579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", PCShellSetView_Shell));
9589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", PCShellSetApplyTranspose_Shell));
9599fa64a6bSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApplyTranspose_C", PCShellSetMatApplyTranspose_Shell));
9609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", PCShellSetName_Shell));
9619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", PCShellGetName_Shell));
9629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", PCShellSetApplyRichardson_Shell));
9633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9644b9ad928SBarry Smith }
965