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 { 9be29d3c6SBarry Smith void *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); 22ace3abfcSBarry Smith PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *); 232fa5cd67SKarl Rupp 244b9ad928SBarry Smith char *name; 254b9ad928SBarry Smith } PC_Shell; 264b9ad928SBarry Smith 27b29801fcSSatish Balay /*@C 2804c3f3b8SBarry Smith PCShellGetContext - Returns the user-provided context associated with a shell `PC` that was provided with `PCShellSetContext()` 29be29d3c6SBarry Smith 30be29d3c6SBarry Smith Not Collective 31be29d3c6SBarry Smith 32be29d3c6SBarry Smith Input Parameter: 3304c3f3b8SBarry Smith . pc - of type `PCSHELL` 34be29d3c6SBarry Smith 35be29d3c6SBarry Smith Output Parameter: 36be29d3c6SBarry Smith . ctx - the user provided context 37be29d3c6SBarry Smith 38be29d3c6SBarry Smith Level: advanced 39be29d3c6SBarry Smith 40f1580f4eSBarry Smith Note: 41c33202f7SPierre Jolivet This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()` 42be29d3c6SBarry Smith 4304c3f3b8SBarry Smith Fortran Note: 4495452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 45daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 46daf670e6SBarry Smith 47*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellSetDestroy()` 48be29d3c6SBarry Smith @*/ 49d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellGetContext(PC pc, void *ctx) 50d71ae5a4SJacob Faibussowitsch { 51ace3abfcSBarry Smith PetscBool flg; 52be29d3c6SBarry Smith 53be29d3c6SBarry Smith PetscFunctionBegin; 540700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 554f572ea9SToby Isaac PetscAssertPointer(ctx, 2); 569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg)); 573ec1f749SStefano Zampini if (!flg) *(void **)ctx = NULL; 583ec1f749SStefano Zampini else *(void **)ctx = ((PC_Shell *)(pc->data))->ctx; 593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60be29d3c6SBarry Smith } 61be29d3c6SBarry Smith 629dd1005fSJed Brown /*@ 6304c3f3b8SBarry Smith PCShellSetContext - sets the context for a shell `PC` that can be accessed with `PCShellGetContext()` 64be29d3c6SBarry Smith 65c3339decSBarry Smith Logically Collective 66be29d3c6SBarry Smith 67be29d3c6SBarry Smith Input Parameters: 68f1580f4eSBarry Smith + pc - the `PC` of type `PCSHELL` 69be29d3c6SBarry Smith - ctx - the context 70be29d3c6SBarry Smith 71be29d3c6SBarry Smith Level: advanced 72be29d3c6SBarry Smith 7304c3f3b8SBarry Smith Notes: 74c33202f7SPierre Jolivet This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()` 7504c3f3b8SBarry Smith 7604c3f3b8SBarry Smith One should also provide a routine to destroy the context when `pc` is destroyed with `PCShellSetDestroy()` 7704c3f3b8SBarry Smith 78feefa0e1SJacob Faibussowitsch Fortran Notes: 7995452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 80daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 81daf670e6SBarry Smith 82*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCShellGetContext()`, `PCSHELL`, `PCShellSetApply()`, `PCShellSetDestroy()` 83be29d3c6SBarry Smith @*/ 84d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellSetContext(PC pc, void *ctx) 85d71ae5a4SJacob Faibussowitsch { 86c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 87ace3abfcSBarry Smith PetscBool flg; 88be29d3c6SBarry Smith 89be29d3c6SBarry Smith PetscFunctionBegin; 900700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg)); 922fa5cd67SKarl Rupp if (flg) shell->ctx = ctx; 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 94be29d3c6SBarry Smith } 95be29d3c6SBarry Smith 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Shell(PC pc) 97d71ae5a4SJacob Faibussowitsch { 98c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 994b9ad928SBarry Smith 1004b9ad928SBarry Smith PetscFunctionBegin; 10128b400f6SJacob Faibussowitsch PetscCheck(shell->setup, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No setup() routine provided to Shell PC"); 10225e27a38SBarry Smith PetscCallBack("PCSHELL callback setup", (*shell->setup)(pc)); 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1044b9ad928SBarry Smith } 1054b9ad928SBarry Smith 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Shell(PC pc, Vec x, Vec y) 107d71ae5a4SJacob Faibussowitsch { 108c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 109e3f487b0SBarry Smith PetscObjectState instate, outstate; 1104b9ad928SBarry Smith 1114b9ad928SBarry Smith PetscFunctionBegin; 11228b400f6SJacob Faibussowitsch PetscCheck(shell->apply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 11425e27a38SBarry Smith PetscCallBack("PCSHELL callback apply", (*shell->apply)(pc, x, y)); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 1161e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 1171e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1194b9ad928SBarry Smith } 1204b9ad928SBarry Smith 121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Shell(PC pc, Mat X, Mat Y) 122d71ae5a4SJacob Faibussowitsch { 1237b6e2003SPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 1247b6e2003SPierre Jolivet PetscObjectState instate, outstate; 1257b6e2003SPierre Jolivet 1267b6e2003SPierre Jolivet PetscFunctionBegin; 12728b400f6SJacob Faibussowitsch PetscCheck(shell->matapply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)Y, &instate)); 12925e27a38SBarry Smith PetscCallBack("PCSHELL callback apply", (*shell->matapply)(pc, X, Y)); 1309566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)Y, &outstate)); 1311e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 1321e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1347b6e2003SPierre Jolivet } 1357b6e2003SPierre Jolivet 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc, Vec x, Vec y) 137d71ae5a4SJacob Faibussowitsch { 1381b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1391b581b66SBarry Smith 1401b581b66SBarry Smith PetscFunctionBegin; 14128b400f6SJacob Faibussowitsch PetscCheck(shell->applysymmetricleft, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 14225e27a38SBarry Smith PetscCallBack("PCSHELL callback apply symmetric left", (*shell->applysymmetricleft)(pc, x, y)); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1441b581b66SBarry Smith } 1451b581b66SBarry Smith 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Shell(PC pc, Vec x, Vec y) 147d71ae5a4SJacob Faibussowitsch { 1481b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1491b581b66SBarry Smith 1501b581b66SBarry Smith PetscFunctionBegin; 15128b400f6SJacob Faibussowitsch PetscCheck(shell->applysymmetricright, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 15225e27a38SBarry Smith PetscCallBack("PCSHELL callback apply symmetric right", (*shell->applysymmetricright)(pc, x, y)); 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1541b581b66SBarry Smith } 1551b581b66SBarry Smith 156d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyBA_Shell(PC pc, PCSide side, Vec x, Vec y, Vec w) 157d71ae5a4SJacob Faibussowitsch { 158c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 159e3f487b0SBarry Smith PetscObjectState instate, outstate; 1602bb17772SBarry Smith 1612bb17772SBarry Smith PetscFunctionBegin; 16228b400f6SJacob Faibussowitsch PetscCheck(shell->applyBA, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyBA() routine provided to Shell PC"); 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)w, &instate)); 16425e27a38SBarry Smith PetscCallBack("PCSHELL callback applyBA", (*shell->applyBA)(pc, side, x, y, w)); 1659566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)w, &outstate)); 1661e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 1671e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)w)); 1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1692bb17772SBarry Smith } 1702bb17772SBarry Smith 171d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc, PetscBool *change) 172d71ae5a4SJacob Faibussowitsch { 173a06fd7f2SStefano Zampini PetscFunctionBegin; 174a06fd7f2SStefano Zampini *change = PETSC_TRUE; 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176a06fd7f2SStefano Zampini } 177a06fd7f2SStefano Zampini 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolve_Shell(PC pc, KSP ksp, Vec b, Vec x) 179d71ae5a4SJacob Faibussowitsch { 180c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1817cdd61b2SBarry Smith 1827cdd61b2SBarry Smith PetscFunctionBegin; 18328b400f6SJacob Faibussowitsch PetscCheck(shell->presolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No presolve() routine provided to Shell PC"); 18425e27a38SBarry Smith PetscCallBack("PCSHELL callback presolve", (*shell->presolve)(pc, ksp, b, x)); 1853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1867cdd61b2SBarry Smith } 1877cdd61b2SBarry Smith 188d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPostSolve_Shell(PC pc, KSP ksp, Vec b, Vec x) 189d71ae5a4SJacob Faibussowitsch { 190c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1917cdd61b2SBarry Smith 1927cdd61b2SBarry Smith PetscFunctionBegin; 19328b400f6SJacob Faibussowitsch PetscCheck(shell->postsolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No postsolve() routine provided to Shell PC"); 19425e27a38SBarry Smith PetscCallBack("PCSHELL callback postsolve()", (*shell->postsolve)(pc, ksp, b, x)); 1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1967cdd61b2SBarry Smith } 1977cdd61b2SBarry Smith 198d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Shell(PC pc, Vec x, Vec y) 199d71ae5a4SJacob Faibussowitsch { 200c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 201e3f487b0SBarry Smith PetscObjectState instate, outstate; 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith PetscFunctionBegin; 20428b400f6SJacob Faibussowitsch PetscCheck(shell->applytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applytranspose() routine provided to Shell PC"); 2059566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 20625e27a38SBarry Smith PetscCallBack("PCSHELL callback applytranspose", (*shell->applytranspose)(pc, x, y)); 2079566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 2081e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 2091e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2114b9ad928SBarry Smith } 2124b9ad928SBarry Smith 213d71ae5a4SJacob 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) 214d71ae5a4SJacob Faibussowitsch { 215c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 216e3f487b0SBarry Smith PetscObjectState instate, outstate; 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith PetscFunctionBegin; 21928b400f6SJacob Faibussowitsch PetscCheck(shell->applyrich, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyrichardson() routine provided to Shell PC"); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 22125e27a38SBarry Smith PetscCallBack("PCSHELL callback applyrichardson", (*shell->applyrich)(pc, x, y, w, rtol, abstol, dtol, it, guesszero, outits, reason)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 223e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 2241e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2264b9ad928SBarry Smith } 2274b9ad928SBarry Smith 228d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Shell(PC pc) 229d71ae5a4SJacob Faibussowitsch { 2304b9ad928SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 2314b9ad928SBarry Smith 2324b9ad928SBarry Smith PetscFunctionBegin; 2339566063dSJacob Faibussowitsch PetscCall(PetscFree(shell->name)); 23425e27a38SBarry Smith if (shell->destroy) PetscCallBack("PCSHELL callback destroy", (*shell->destroy)(pc)); 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", NULL)); 2369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", NULL)); 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", NULL)); 2389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", NULL)); 2399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", NULL)); 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", NULL)); 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", NULL)); 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", NULL)); 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", NULL)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", NULL)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", NULL)); 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", NULL)); 2479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", NULL)); 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", NULL)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 2509566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2524b9ad928SBarry Smith } 2534b9ad928SBarry Smith 254d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Shell(PC pc, PetscViewer viewer) 255d71ae5a4SJacob Faibussowitsch { 2564b9ad928SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 257ace3abfcSBarry Smith PetscBool iascii; 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith PetscFunctionBegin; 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 26132077d6dSBarry Smith if (iascii) { 2621e66621cSBarry Smith if (shell->name) PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", shell->name)); 2631e66621cSBarry Smith else PetscCall(PetscViewerASCIIPrintf(viewer, " no name\n")); 2644b9ad928SBarry Smith } 2654b9ad928SBarry Smith if (shell->view) { 2669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2679566063dSJacob Faibussowitsch PetscCall((*shell->view)(pc, viewer)); 2689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2694b9ad928SBarry Smith } 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2714b9ad928SBarry Smith } 2724b9ad928SBarry Smith 273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC)) 274d71ae5a4SJacob Faibussowitsch { 275c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 27618be62a5SSatish Balay 27718be62a5SSatish Balay PetscFunctionBegin; 27818be62a5SSatish Balay shell->destroy = destroy; 2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28018be62a5SSatish Balay } 28118be62a5SSatish Balay 282d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC)) 283d71ae5a4SJacob Faibussowitsch { 284feb237baSPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 2854b9ad928SBarry Smith 2864b9ad928SBarry Smith PetscFunctionBegin; 2874b9ad928SBarry Smith shell->setup = setup; 288d01c8aa3SLisandro Dalcin if (setup) pc->ops->setup = PCSetUp_Shell; 2890a545947SLisandro Dalcin else pc->ops->setup = NULL; 2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2914b9ad928SBarry Smith } 2924b9ad928SBarry Smith 293d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApply_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 294d71ae5a4SJacob Faibussowitsch { 295c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 2964b9ad928SBarry Smith 2974b9ad928SBarry Smith PetscFunctionBegin; 2984b9ad928SBarry Smith shell->apply = apply; 2993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3004b9ad928SBarry Smith } 3014b9ad928SBarry Smith 302d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetMatApply_Shell(PC pc, PetscErrorCode (*matapply)(PC, Mat, Mat)) 303d71ae5a4SJacob Faibussowitsch { 3047b6e2003SPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 3057b6e2003SPierre Jolivet 3067b6e2003SPierre Jolivet PetscFunctionBegin; 3077b6e2003SPierre Jolivet shell->matapply = matapply; 3080e0fe96bSStefano Zampini if (matapply) pc->ops->matapply = PCMatApply_Shell; 3090e0fe96bSStefano Zampini else pc->ops->matapply = NULL; 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3117b6e2003SPierre Jolivet } 3127b6e2003SPierre Jolivet 313d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplySymmetricLeft_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 314d71ae5a4SJacob Faibussowitsch { 3151b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3161b581b66SBarry Smith 3171b581b66SBarry Smith PetscFunctionBegin; 3181b581b66SBarry Smith shell->applysymmetricleft = apply; 3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3201b581b66SBarry Smith } 3211b581b66SBarry Smith 322d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplySymmetricRight_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 323d71ae5a4SJacob Faibussowitsch { 3241b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3251b581b66SBarry Smith 3261b581b66SBarry Smith PetscFunctionBegin; 3271b581b66SBarry Smith shell->applysymmetricright = apply; 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3291b581b66SBarry Smith } 3301b581b66SBarry Smith 331d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyBA_Shell(PC pc, PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec)) 332d71ae5a4SJacob Faibussowitsch { 333c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3342bb17772SBarry Smith 3352bb17772SBarry Smith PetscFunctionBegin; 336d01c8aa3SLisandro Dalcin shell->applyBA = applyBA; 337d01c8aa3SLisandro Dalcin if (applyBA) pc->ops->applyBA = PCApplyBA_Shell; 3380a545947SLisandro Dalcin else pc->ops->applyBA = NULL; 3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3402bb17772SBarry Smith } 3412bb17772SBarry Smith 342d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetPreSolve_Shell(PC pc, PetscErrorCode (*presolve)(PC, KSP, Vec, Vec)) 343d71ae5a4SJacob Faibussowitsch { 344c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3457cdd61b2SBarry Smith 3467cdd61b2SBarry Smith PetscFunctionBegin; 3477cdd61b2SBarry Smith shell->presolve = presolve; 348a06fd7f2SStefano Zampini if (presolve) { 349a06fd7f2SStefano Zampini pc->ops->presolve = PCPreSolve_Shell; 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Shell)); 351a06fd7f2SStefano Zampini } else { 3520a545947SLisandro Dalcin pc->ops->presolve = NULL; 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 354a06fd7f2SStefano Zampini } 3553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3567cdd61b2SBarry Smith } 3577cdd61b2SBarry Smith 358d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetPostSolve_Shell(PC pc, PetscErrorCode (*postsolve)(PC, KSP, Vec, Vec)) 359d71ae5a4SJacob Faibussowitsch { 360c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3617cdd61b2SBarry Smith 3627cdd61b2SBarry Smith PetscFunctionBegin; 3637cdd61b2SBarry Smith shell->postsolve = postsolve; 364d01c8aa3SLisandro Dalcin if (postsolve) pc->ops->postsolve = PCPostSolve_Shell; 3650a545947SLisandro Dalcin else pc->ops->postsolve = NULL; 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3677cdd61b2SBarry Smith } 3687cdd61b2SBarry Smith 369d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetView_Shell(PC pc, PetscErrorCode (*view)(PC, PetscViewer)) 370d71ae5a4SJacob Faibussowitsch { 371c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3724b9ad928SBarry Smith 3734b9ad928SBarry Smith PetscFunctionBegin; 3744b9ad928SBarry Smith shell->view = view; 3753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3764b9ad928SBarry Smith } 3774b9ad928SBarry Smith 378d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyTranspose_Shell(PC pc, PetscErrorCode (*applytranspose)(PC, Vec, Vec)) 379d71ae5a4SJacob Faibussowitsch { 380c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith PetscFunctionBegin; 3834b9ad928SBarry Smith shell->applytranspose = applytranspose; 384d01c8aa3SLisandro Dalcin if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell; 3850a545947SLisandro Dalcin else pc->ops->applytranspose = NULL; 3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 387d01c8aa3SLisandro Dalcin } 388d01c8aa3SLisandro Dalcin 389d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyRichardson_Shell(PC pc, PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)) 390d71ae5a4SJacob Faibussowitsch { 391c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 392d01c8aa3SLisandro Dalcin 393d01c8aa3SLisandro Dalcin PetscFunctionBegin; 394d01c8aa3SLisandro Dalcin shell->applyrich = applyrich; 395d01c8aa3SLisandro Dalcin if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell; 3960a545947SLisandro Dalcin else pc->ops->applyrichardson = NULL; 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3984b9ad928SBarry Smith } 3994b9ad928SBarry Smith 400d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetName_Shell(PC pc, const char name[]) 401d71ae5a4SJacob Faibussowitsch { 402c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 4034b9ad928SBarry Smith 4044b9ad928SBarry Smith PetscFunctionBegin; 4059566063dSJacob Faibussowitsch PetscCall(PetscFree(shell->name)); 4069566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &shell->name)); 4073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4084b9ad928SBarry Smith } 4094b9ad928SBarry Smith 410d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellGetName_Shell(PC pc, const char *name[]) 411d71ae5a4SJacob Faibussowitsch { 412c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith PetscFunctionBegin; 4154b9ad928SBarry Smith *name = shell->name; 4163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4174b9ad928SBarry Smith } 4184b9ad928SBarry Smith 41918be62a5SSatish Balay /*@C 42004c3f3b8SBarry Smith PCShellSetDestroy - Sets routine to use to destroy the user-provided application context that was provided with `PCShellSetContext()` 42118be62a5SSatish Balay 422c3339decSBarry Smith Logically Collective 42318be62a5SSatish Balay 42418be62a5SSatish Balay Input Parameters: 42518be62a5SSatish Balay + pc - the preconditioner context 426a2b725a8SWilliam Gropp - destroy - the application-provided destroy routine 42718be62a5SSatish Balay 42820f4b53cSBarry Smith Calling sequence of `destroy`: 42904c3f3b8SBarry Smith . pc - the preconditioner 4304aa34b0aSBarry Smith 431f1580f4eSBarry Smith Level: intermediate 43218be62a5SSatish Balay 433*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()` 43418be62a5SSatish Balay @*/ 43504c3f3b8SBarry Smith PetscErrorCode PCShellSetDestroy(PC pc, PetscErrorCode (*destroy)(PC pc)) 436d71ae5a4SJacob Faibussowitsch { 43718be62a5SSatish Balay PetscFunctionBegin; 4380700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 439cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetDestroy_C", (PC, PetscErrorCode(*)(PC)), (pc, destroy)); 4403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44118be62a5SSatish Balay } 44218be62a5SSatish Balay 4434b9ad928SBarry Smith /*@C 4444b9ad928SBarry Smith PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the 4454b9ad928SBarry Smith matrix operator is changed. 4464b9ad928SBarry Smith 447c3339decSBarry Smith Logically Collective 4484b9ad928SBarry Smith 4494b9ad928SBarry Smith Input Parameters: 4504b9ad928SBarry Smith + pc - the preconditioner context 451a2b725a8SWilliam Gropp - setup - the application-provided setup routine 4524b9ad928SBarry Smith 45320f4b53cSBarry Smith Calling sequence of `setup`: 45404c3f3b8SBarry Smith . pc - the preconditioner 4554aa34b0aSBarry Smith 456f1580f4eSBarry Smith Level: intermediate 4574b9ad928SBarry Smith 45804c3f3b8SBarry Smith Note: 45904c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `setup`. 46004c3f3b8SBarry Smith 461*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetApply()`, `PCShellSetContext()`, , `PCShellGetContext()` 4624b9ad928SBarry Smith @*/ 46304c3f3b8SBarry Smith PetscErrorCode PCShellSetSetUp(PC pc, PetscErrorCode (*setup)(PC pc)) 464d71ae5a4SJacob Faibussowitsch { 4654b9ad928SBarry Smith PetscFunctionBegin; 4660700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 467cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetSetUp_C", (PC, PetscErrorCode(*)(PC)), (pc, setup)); 4683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4694b9ad928SBarry Smith } 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith /*@C 47204c3f3b8SBarry Smith PCShellSetView - Sets routine to use as viewer of a `PCSHELL` shell preconditioner 4734b9ad928SBarry Smith 474c3339decSBarry Smith Logically Collective 4754b9ad928SBarry Smith 4764b9ad928SBarry Smith Input Parameters: 4774b9ad928SBarry Smith + pc - the preconditioner context 4784b9ad928SBarry Smith - view - the application-provided view routine 4794b9ad928SBarry Smith 48020f4b53cSBarry Smith Calling sequence of `view`: 48104c3f3b8SBarry Smith + pc - the preconditioner 4824b9ad928SBarry Smith - v - viewer 4834b9ad928SBarry Smith 484f1580f4eSBarry Smith Level: advanced 4854b9ad928SBarry Smith 48604c3f3b8SBarry Smith Note: 48704c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `view`. 48804c3f3b8SBarry Smith 489*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()` 4904b9ad928SBarry Smith @*/ 49104c3f3b8SBarry Smith PetscErrorCode PCShellSetView(PC pc, PetscErrorCode (*view)(PC pc, PetscViewer v)) 492d71ae5a4SJacob Faibussowitsch { 4934b9ad928SBarry Smith PetscFunctionBegin; 4940700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 495cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetView_C", (PC, PetscErrorCode(*)(PC, PetscViewer)), (pc, view)); 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4974b9ad928SBarry Smith } 4984b9ad928SBarry Smith 4994b9ad928SBarry Smith /*@C 5004b9ad928SBarry Smith PCShellSetApply - Sets routine to use as preconditioner. 5014b9ad928SBarry Smith 502c3339decSBarry Smith Logically Collective 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith Input Parameters: 5054b9ad928SBarry Smith + pc - the preconditioner context 506be29d3c6SBarry Smith - apply - the application-provided preconditioning routine 5074b9ad928SBarry Smith 50820f4b53cSBarry Smith Calling sequence of `apply`: 50920f4b53cSBarry Smith + pc - the preconditioner, get the application context with `PCShellGetContext()` 5104b9ad928SBarry Smith . xin - input vector 5114b9ad928SBarry Smith - xout - output vector 5124b9ad928SBarry Smith 513f1580f4eSBarry Smith Level: intermediate 5144b9ad928SBarry Smith 51504c3f3b8SBarry Smith Note: 51604c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 51704c3f3b8SBarry Smith 518*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellSetApplySymmetricRight()`, `PCShellSetApplySymmetricLeft()`, `PCShellGetContext()` 5194b9ad928SBarry Smith @*/ 52004c3f3b8SBarry Smith PetscErrorCode PCShellSetApply(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 521d71ae5a4SJacob Faibussowitsch { 5224b9ad928SBarry Smith PetscFunctionBegin; 5230700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 524cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApply_C", (PC, PetscErrorCode(*)(PC, Vec, Vec)), (pc, apply)); 5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5264b9ad928SBarry Smith } 5274b9ad928SBarry Smith 5281b581b66SBarry Smith /*@C 5296437efc7SEric Chamberland PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors. 5307b6e2003SPierre Jolivet 531c3339decSBarry Smith Logically Collective 5327b6e2003SPierre Jolivet 5337b6e2003SPierre Jolivet Input Parameters: 5347b6e2003SPierre Jolivet + pc - the preconditioner context 53504c3f3b8SBarry Smith - matapply - the application-provided preconditioning routine 5367b6e2003SPierre Jolivet 53704c3f3b8SBarry Smith Calling sequence of `matapply`: 53804c3f3b8SBarry Smith + pc - the preconditioner 53904c3f3b8SBarry Smith . Xin - input block of vectors represented as a dense `Mat` 54004c3f3b8SBarry Smith - Xout - output block of vectors represented as a dense `Mat` 5417b6e2003SPierre Jolivet 542f1580f4eSBarry Smith Level: advanced 5437b6e2003SPierre Jolivet 54404c3f3b8SBarry Smith Note: 54504c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapply`. 54604c3f3b8SBarry Smith 547*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()` 5487b6e2003SPierre Jolivet @*/ 54904c3f3b8SBarry Smith PetscErrorCode PCShellSetMatApply(PC pc, PetscErrorCode (*matapply)(PC pc, Mat Xin, Mat Xout)) 550d71ae5a4SJacob Faibussowitsch { 5517b6e2003SPierre Jolivet PetscFunctionBegin; 5527b6e2003SPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 553cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetMatApply_C", (PC, PetscErrorCode(*)(PC, Mat, Mat)), (pc, matapply)); 5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5557b6e2003SPierre Jolivet } 5567b6e2003SPierre Jolivet 5577b6e2003SPierre Jolivet /*@C 558f1580f4eSBarry Smith PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the `PC_SYMMETRIC` is used). 5591b581b66SBarry Smith 560c3339decSBarry Smith Logically Collective 5611b581b66SBarry Smith 5621b581b66SBarry Smith Input Parameters: 5631b581b66SBarry Smith + pc - the preconditioner context 5641b581b66SBarry Smith - apply - the application-provided left preconditioning routine 5651b581b66SBarry Smith 56620f4b53cSBarry Smith Calling sequence of `apply`: 56704c3f3b8SBarry Smith + pc - the preconditioner 5681b581b66SBarry Smith . xin - input vector 5691b581b66SBarry Smith - xout - output vector 5701b581b66SBarry Smith 571f1580f4eSBarry Smith Level: advanced 5721b581b66SBarry Smith 57304c3f3b8SBarry Smith Note: 57404c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 57504c3f3b8SBarry Smith 576*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()` 5771b581b66SBarry Smith @*/ 57804c3f3b8SBarry Smith PetscErrorCode PCShellSetApplySymmetricLeft(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 579d71ae5a4SJacob Faibussowitsch { 5801b581b66SBarry Smith PetscFunctionBegin; 5811b581b66SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 582cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplySymmetricLeft_C", (PC, PetscErrorCode(*)(PC, Vec, Vec)), (pc, apply)); 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5841b581b66SBarry Smith } 5851b581b66SBarry Smith 5861b581b66SBarry Smith /*@C 587f1580f4eSBarry Smith PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the `PC_SYMMETRIC` is used). 5881b581b66SBarry Smith 589c3339decSBarry Smith Logically Collective 5901b581b66SBarry Smith 5911b581b66SBarry Smith Input Parameters: 5921b581b66SBarry Smith + pc - the preconditioner context 5931b581b66SBarry Smith - apply - the application-provided right preconditioning routine 5941b581b66SBarry Smith 59520f4b53cSBarry Smith Calling sequence of `apply`: 59604c3f3b8SBarry Smith + pc - the preconditioner 5971b581b66SBarry Smith . xin - input vector 5981b581b66SBarry Smith - xout - output vector 5991b581b66SBarry Smith 600f1580f4eSBarry Smith Level: advanced 6011b581b66SBarry Smith 60204c3f3b8SBarry Smith Note: 60304c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 60404c3f3b8SBarry Smith 605*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()` 6061b581b66SBarry Smith @*/ 60704c3f3b8SBarry Smith PetscErrorCode PCShellSetApplySymmetricRight(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 608d71ae5a4SJacob Faibussowitsch { 6091b581b66SBarry Smith PetscFunctionBegin; 6101b581b66SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 611cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplySymmetricRight_C", (PC, PetscErrorCode(*)(PC, Vec, Vec)), (pc, apply)); 6123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6131b581b66SBarry Smith } 6141b581b66SBarry Smith 6152bb17772SBarry Smith /*@C 61604c3f3b8SBarry Smith PCShellSetApplyBA - Sets routine to use as the preconditioner times the operator. 6172bb17772SBarry Smith 618c3339decSBarry Smith Logically Collective 6192bb17772SBarry Smith 6202bb17772SBarry Smith Input Parameters: 6212bb17772SBarry Smith + pc - the preconditioner context 6222bb17772SBarry Smith - applyBA - the application-provided BA routine 6232bb17772SBarry Smith 62420f4b53cSBarry Smith Calling sequence of `applyBA`: 62504c3f3b8SBarry Smith + pc - the preconditioner 62604c3f3b8SBarry Smith . side - `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 6272bb17772SBarry Smith . xin - input vector 62804c3f3b8SBarry Smith . xout - output vector 62904c3f3b8SBarry Smith - w - work vector 6302bb17772SBarry Smith 631f1580f4eSBarry Smith Level: intermediate 6322bb17772SBarry Smith 63304c3f3b8SBarry Smith Note: 63404c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applyBA`. 63504c3f3b8SBarry Smith 636*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellGetContext()`, `PCSide` 6372bb17772SBarry Smith @*/ 63804c3f3b8SBarry Smith PetscErrorCode PCShellSetApplyBA(PC pc, PetscErrorCode (*applyBA)(PC pc, PCSide side, Vec xin, Vec xout, Vec w)) 639d71ae5a4SJacob Faibussowitsch { 6402bb17772SBarry Smith PetscFunctionBegin; 6410700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 642cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplyBA_C", (PC, PetscErrorCode(*)(PC, PCSide, Vec, Vec, Vec)), (pc, applyBA)); 6433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6442bb17772SBarry Smith } 6452bb17772SBarry Smith 6464b9ad928SBarry Smith /*@C 6474b9ad928SBarry Smith PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose. 6484b9ad928SBarry Smith 649c3339decSBarry Smith Logically Collective 6504b9ad928SBarry Smith 6514b9ad928SBarry Smith Input Parameters: 6524b9ad928SBarry Smith + pc - the preconditioner context 65320f4b53cSBarry Smith - applytranspose - the application-provided preconditioning transpose routine 6544b9ad928SBarry Smith 65520f4b53cSBarry Smith Calling sequence of `applytranspose`: 65604c3f3b8SBarry Smith + pc - the preconditioner 6574b9ad928SBarry Smith . xin - input vector 6584b9ad928SBarry Smith - xout - output vector 6594b9ad928SBarry Smith 660f1580f4eSBarry Smith Level: intermediate 6614b9ad928SBarry Smith 662f1580f4eSBarry Smith Note: 66304c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applytranspose`. 6644b9ad928SBarry Smith 665*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCSetContext()`, `PCShellSetApplyBA()`, `PCGetContext()` 6664b9ad928SBarry Smith @*/ 66704c3f3b8SBarry Smith PetscErrorCode PCShellSetApplyTranspose(PC pc, PetscErrorCode (*applytranspose)(PC pc, Vec xin, Vec xout)) 668d71ae5a4SJacob Faibussowitsch { 6694b9ad928SBarry Smith PetscFunctionBegin; 6700700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 671cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplyTranspose_C", (PC, PetscErrorCode(*)(PC, Vec, Vec)), (pc, applytranspose)); 6723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6734b9ad928SBarry Smith } 6744b9ad928SBarry Smith 6757cdd61b2SBarry Smith /*@C 676f1580f4eSBarry Smith PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a `KSPSolve()` is 6777cdd61b2SBarry Smith applied. This usually does something like scale the linear system in some application 6787cdd61b2SBarry Smith specific way. 6797cdd61b2SBarry Smith 680c3339decSBarry Smith Logically Collective 6817cdd61b2SBarry Smith 6827cdd61b2SBarry Smith Input Parameters: 6837cdd61b2SBarry Smith + pc - the preconditioner context 6847cdd61b2SBarry Smith - presolve - the application-provided presolve routine 6857cdd61b2SBarry Smith 68620f4b53cSBarry Smith Calling sequence of `presolve`: 68704c3f3b8SBarry Smith + pc - the preconditioner 68804c3f3b8SBarry Smith . ksp - the `KSP` that contains `pc` 6897cdd61b2SBarry Smith . xin - input vector 6907cdd61b2SBarry Smith - xout - output vector 6917cdd61b2SBarry Smith 692f1580f4eSBarry Smith Level: advanced 6937cdd61b2SBarry Smith 69404c3f3b8SBarry Smith Note: 69504c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `presolve`. 69604c3f3b8SBarry Smith 697*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPostSolve()`, `PCShellSetContext()`, `PCGetContext()` 6987cdd61b2SBarry Smith @*/ 69904c3f3b8SBarry Smith PetscErrorCode PCShellSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp, Vec xin, Vec xout)) 700d71ae5a4SJacob Faibussowitsch { 7017cdd61b2SBarry Smith PetscFunctionBegin; 7020700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 703cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetPreSolve_C", (PC, PetscErrorCode(*)(PC, KSP, Vec, Vec)), (pc, presolve)); 7043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7057cdd61b2SBarry Smith } 7067cdd61b2SBarry Smith 7077cdd61b2SBarry Smith /*@C 70804c3f3b8SBarry Smith PCShellSetPostSolve - Sets routine to apply to the operators/vectors after a `KSPSolve()` is 7097cdd61b2SBarry Smith applied. This usually does something like scale the linear system in some application 7107cdd61b2SBarry Smith specific way. 7117cdd61b2SBarry Smith 712c3339decSBarry Smith Logically Collective 7137cdd61b2SBarry Smith 7147cdd61b2SBarry Smith Input Parameters: 7157cdd61b2SBarry Smith + pc - the preconditioner context 7167cdd61b2SBarry Smith - postsolve - the application-provided presolve routine 7177cdd61b2SBarry Smith 71820f4b53cSBarry Smith Calling sequence of `postsolve`: 71904c3f3b8SBarry Smith + pc - the preconditioner 72004c3f3b8SBarry Smith . ksp - the `KSP` that contains `pc` 7217cdd61b2SBarry Smith . xin - input vector 7227cdd61b2SBarry Smith - xout - output vector 7237cdd61b2SBarry Smith 724f1580f4eSBarry Smith Level: advanced 7257cdd61b2SBarry Smith 72604c3f3b8SBarry Smith Note: 72704c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `postsolve`. 72804c3f3b8SBarry Smith 729*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPreSolve()`, `PCShellSetContext()`, `PCGetContext()` 7307cdd61b2SBarry Smith @*/ 73104c3f3b8SBarry Smith PetscErrorCode PCShellSetPostSolve(PC pc, PetscErrorCode (*postsolve)(PC pc, KSP ksp, Vec xin, Vec xout)) 732d71ae5a4SJacob Faibussowitsch { 7337cdd61b2SBarry Smith PetscFunctionBegin; 7340700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 735cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetPostSolve_C", (PC, PetscErrorCode(*)(PC, KSP, Vec, Vec)), (pc, postsolve)); 7363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7377cdd61b2SBarry Smith } 7387cdd61b2SBarry Smith 7394b9ad928SBarry Smith /*@C 740f1580f4eSBarry Smith PCShellSetName - Sets an optional name to associate with a `PCSHELL` 7414b9ad928SBarry Smith preconditioner. 7424b9ad928SBarry Smith 7434b9ad928SBarry Smith Not Collective 7444b9ad928SBarry Smith 7454b9ad928SBarry Smith Input Parameters: 7464b9ad928SBarry Smith + pc - the preconditioner context 7474b9ad928SBarry Smith - name - character string describing shell preconditioner 7484b9ad928SBarry Smith 749f1580f4eSBarry Smith Level: intermediate 7504b9ad928SBarry Smith 75104c3f3b8SBarry Smith Note: 752baca6076SPierre Jolivet This is separate from the name you can provide with `PetscObjectSetName()` 75304c3f3b8SBarry Smith 754*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellGetName()`, `PetscObjectSetName()`, `PetscObjectGetName()` 7554b9ad928SBarry Smith @*/ 756d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellSetName(PC pc, const char name[]) 757d71ae5a4SJacob Faibussowitsch { 7584b9ad928SBarry Smith PetscFunctionBegin; 7590700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 760cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetName_C", (PC, const char[]), (pc, name)); 7613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7624b9ad928SBarry Smith } 7634b9ad928SBarry Smith 7644b9ad928SBarry Smith /*@C 76504c3f3b8SBarry Smith PCShellGetName - Gets an optional name that the user has set for a `PCSHELL` with `PCShellSetName()` 7664b9ad928SBarry Smith preconditioner. 7674b9ad928SBarry Smith 7684b9ad928SBarry Smith Not Collective 7694b9ad928SBarry Smith 7704b9ad928SBarry Smith Input Parameter: 7714b9ad928SBarry Smith . pc - the preconditioner context 7724b9ad928SBarry Smith 7734b9ad928SBarry Smith Output Parameter: 7744b9ad928SBarry Smith . name - character string describing shell preconditioner (you should not free this) 7754b9ad928SBarry Smith 776f1580f4eSBarry Smith Level: intermediate 7774b9ad928SBarry Smith 778*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetName()`, `PetscObjectSetName()`, `PetscObjectGetName()` 7794b9ad928SBarry Smith @*/ 780d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellGetName(PC pc, const char *name[]) 781d71ae5a4SJacob Faibussowitsch { 7824b9ad928SBarry Smith PetscFunctionBegin; 7830700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 7844f572ea9SToby Isaac PetscAssertPointer(name, 2); 785cac4c232SBarry Smith PetscUseMethod(pc, "PCShellGetName_C", (PC, const char *[]), (pc, name)); 7863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7874b9ad928SBarry Smith } 7884b9ad928SBarry Smith 7894b9ad928SBarry Smith /*@C 7904b9ad928SBarry Smith PCShellSetApplyRichardson - Sets routine to use as preconditioner 7914b9ad928SBarry Smith in Richardson iteration. 7924b9ad928SBarry Smith 793c3339decSBarry Smith Logically Collective 7944b9ad928SBarry Smith 7954b9ad928SBarry Smith Input Parameters: 7964b9ad928SBarry Smith + pc - the preconditioner context 797be29d3c6SBarry Smith - apply - the application-provided preconditioning routine 7984b9ad928SBarry Smith 79920f4b53cSBarry Smith Calling sequence of `apply`: 80004c3f3b8SBarry Smith + pc - the preconditioner 8014b9ad928SBarry Smith . b - right-hand-side 8024b9ad928SBarry Smith . x - current iterate 8034b9ad928SBarry Smith . r - work space 8044b9ad928SBarry Smith . rtol - relative tolerance of residual norm to stop at 80570441072SBarry Smith . abstol - absolute tolerance of residual norm to stop at 8064b9ad928SBarry Smith . dtol - if residual norm increases by this factor than return 80704c3f3b8SBarry Smith . maxits - number of iterations to run 80804c3f3b8SBarry Smith . zeroinitialguess - `PETSC_TRUE` if `x` is known to be initially zero 80904c3f3b8SBarry Smith . its - returns the number of iterations used 81004c3f3b8SBarry Smith - reason - returns the reason the iteration has converged 8114b9ad928SBarry Smith 812f1580f4eSBarry Smith Level: advanced 8134b9ad928SBarry Smith 81404c3f3b8SBarry Smith Note: 81504c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 81604c3f3b8SBarry Smith 817*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCRichardsonConvergedReason()`, `PCShellGetContext()` 8184b9ad928SBarry Smith @*/ 81904c3f3b8SBarry 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)) 820d71ae5a4SJacob Faibussowitsch { 8214b9ad928SBarry Smith PetscFunctionBegin; 8220700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 823cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplyRichardson_C", (PC, PetscErrorCode(*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)), (pc, apply)); 8243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8254b9ad928SBarry Smith } 8264b9ad928SBarry Smith 8274b9ad928SBarry Smith /*MC 828f1580f4eSBarry Smith PCSHELL - Creates a new preconditioner class for use with a users 829f1580f4eSBarry Smith own private data storage format and preconditioner application code 8304b9ad928SBarry Smith 8314b9ad928SBarry Smith Level: advanced 832e0bb08deSStefano Zampini 8334b9ad928SBarry Smith Usage: 834f1580f4eSBarry Smith .vb 835f1580f4eSBarry Smith extern PetscErrorCode apply(PC,Vec,Vec); 836f1580f4eSBarry Smith extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec); 837f1580f4eSBarry Smith extern PetscErrorCode applytranspose(PC,Vec,Vec); 838f1580f4eSBarry Smith extern PetscErrorCode setup(PC); 839f1580f4eSBarry Smith extern PetscErrorCode destroy(PC); 840f1580f4eSBarry Smith 841f1580f4eSBarry Smith PCCreate(comm,&pc); 842f1580f4eSBarry Smith PCSetType(pc,PCSHELL); 843f1580f4eSBarry Smith PCShellSetContext(pc,ctx) 844f1580f4eSBarry Smith PCShellSetApply(pc,apply); 845f1580f4eSBarry Smith PCShellSetApplyBA(pc,applyba); (optional) 846f1580f4eSBarry Smith PCShellSetApplyTranspose(pc,applytranspose); (optional) 847f1580f4eSBarry Smith PCShellSetSetUp(pc,setup); (optional) 848f1580f4eSBarry Smith PCShellSetDestroy(pc,destroy); (optional) 849f1580f4eSBarry Smith .ve 8504b9ad928SBarry Smith 85104c3f3b8SBarry Smith Note: 85204c3f3b8SBarry Smith Information required for the preconditioner and its internal datastructures can be set with `PCShellSetContext()` and then accessed 85304c3f3b8SBarry Smith with `PCShellGetContext()` inside the routines provided above 85404c3f3b8SBarry Smith 855*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 856f1580f4eSBarry Smith `MATSHELL`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetView()`, `PCShellSetDestroy()`, `PCShellSetPostSolve()`, 857f1580f4eSBarry Smith `PCShellSetApplyTranspose()`, `PCShellSetName()`, `PCShellSetApplyRichardson()`, `PCShellSetPreSolve()`, `PCShellSetView()`, 858f1580f4eSBarry Smith `PCShellGetName()`, `PCShellSetContext()`, `PCShellGetContext()`, `PCShellSetApplyBA()`, `MATSHELL`, `PCShellSetMatApply()`, 8594b9ad928SBarry Smith M*/ 8604b9ad928SBarry Smith 861d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc) 862d71ae5a4SJacob Faibussowitsch { 8634b9ad928SBarry Smith PC_Shell *shell; 8644b9ad928SBarry Smith 8654b9ad928SBarry Smith PetscFunctionBegin; 8664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&shell)); 8674b9ad928SBarry Smith pc->data = (void *)shell; 8684b9ad928SBarry Smith 869d01c8aa3SLisandro Dalcin pc->ops->destroy = PCDestroy_Shell; 8704b9ad928SBarry Smith pc->ops->view = PCView_Shell; 871d01c8aa3SLisandro Dalcin pc->ops->apply = PCApply_Shell; 8721b581b66SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_Shell; 8731b581b66SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_Shell; 8740e0fe96bSStefano Zampini pc->ops->matapply = NULL; 8750a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 8760a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 8770a545947SLisandro Dalcin pc->ops->setup = NULL; 8780a545947SLisandro Dalcin pc->ops->presolve = NULL; 8790a545947SLisandro Dalcin pc->ops->postsolve = NULL; 8804b9ad928SBarry Smith 8810a545947SLisandro Dalcin shell->apply = NULL; 8820a545947SLisandro Dalcin shell->applytranspose = NULL; 8830a545947SLisandro Dalcin shell->name = NULL; 8840a545947SLisandro Dalcin shell->applyrich = NULL; 8850a545947SLisandro Dalcin shell->presolve = NULL; 8860a545947SLisandro Dalcin shell->postsolve = NULL; 8870a545947SLisandro Dalcin shell->ctx = NULL; 8880a545947SLisandro Dalcin shell->setup = NULL; 8890a545947SLisandro Dalcin shell->view = NULL; 8900a545947SLisandro Dalcin shell->destroy = NULL; 8910a545947SLisandro Dalcin shell->applysymmetricleft = NULL; 8920a545947SLisandro Dalcin shell->applysymmetricright = NULL; 8934b9ad928SBarry Smith 8949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", PCShellSetDestroy_Shell)); 8959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", PCShellSetSetUp_Shell)); 8969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", PCShellSetApply_Shell)); 8979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", PCShellSetMatApply_Shell)); 8989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", PCShellSetApplySymmetricLeft_Shell)); 8999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", PCShellSetApplySymmetricRight_Shell)); 9009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", PCShellSetApplyBA_Shell)); 9019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", PCShellSetPreSolve_Shell)); 9029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", PCShellSetPostSolve_Shell)); 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", PCShellSetView_Shell)); 9049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", PCShellSetApplyTranspose_Shell)); 9059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", PCShellSetName_Shell)); 9069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", PCShellGetName_Shell)); 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", PCShellSetApplyRichardson_Shell)); 9083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9094b9ad928SBarry Smith } 910