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); 22*9fa64a6bSPierre 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 4404c3f3b8SBarry Smith Fortran Note: 4595452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 46feaf08eaSBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument. 47daf670e6SBarry Smith 48562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellSetDestroy()` 49be29d3c6SBarry Smith @*/ 50d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellGetContext(PC pc, void *ctx) 51d71ae5a4SJacob Faibussowitsch { 52ace3abfcSBarry Smith PetscBool flg; 53be29d3c6SBarry Smith 54be29d3c6SBarry Smith PetscFunctionBegin; 550700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 564f572ea9SToby Isaac PetscAssertPointer(ctx, 2); 579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg)); 583ec1f749SStefano Zampini if (!flg) *(void **)ctx = NULL; 59f4f49eeaSPierre Jolivet else *(void **)ctx = ((PC_Shell *)pc->data)->ctx; 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61be29d3c6SBarry Smith } 62be29d3c6SBarry Smith 639dd1005fSJed Brown /*@ 6404c3f3b8SBarry Smith PCShellSetContext - sets the context for a shell `PC` that can be accessed with `PCShellGetContext()` 65be29d3c6SBarry Smith 66c3339decSBarry Smith Logically Collective 67be29d3c6SBarry Smith 68be29d3c6SBarry Smith Input Parameters: 69f1580f4eSBarry Smith + pc - the `PC` of type `PCSHELL` 70be29d3c6SBarry Smith - ctx - the context 71be29d3c6SBarry Smith 72be29d3c6SBarry Smith Level: advanced 73be29d3c6SBarry Smith 7404c3f3b8SBarry Smith Notes: 75c33202f7SPierre Jolivet This routine is intended for use within the various user-provided routines set with, for example, `PCShellSetApply()` 7604c3f3b8SBarry Smith 7704c3f3b8SBarry Smith One should also provide a routine to destroy the context when `pc` is destroyed with `PCShellSetDestroy()` 7804c3f3b8SBarry Smith 79feefa0e1SJacob Faibussowitsch Fortran Notes: 8095452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 81feaf08eaSBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument. 82daf670e6SBarry Smith 83562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCShellGetContext()`, `PCSHELL`, `PCShellSetApply()`, `PCShellSetDestroy()` 84be29d3c6SBarry Smith @*/ 85d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellSetContext(PC pc, void *ctx) 86d71ae5a4SJacob Faibussowitsch { 87c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 88ace3abfcSBarry Smith PetscBool flg; 89be29d3c6SBarry Smith 90be29d3c6SBarry Smith PetscFunctionBegin; 910700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &flg)); 932fa5cd67SKarl Rupp if (flg) shell->ctx = ctx; 943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95be29d3c6SBarry Smith } 96be29d3c6SBarry Smith 97d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Shell(PC pc) 98d71ae5a4SJacob Faibussowitsch { 99c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith PetscFunctionBegin; 10228b400f6SJacob Faibussowitsch PetscCheck(shell->setup, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No setup() routine provided to Shell PC"); 10325e27a38SBarry Smith PetscCallBack("PCSHELL callback setup", (*shell->setup)(pc)); 1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1054b9ad928SBarry Smith } 1064b9ad928SBarry Smith 107d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Shell(PC pc, Vec x, Vec y) 108d71ae5a4SJacob Faibussowitsch { 109c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 110e3f487b0SBarry Smith PetscObjectState instate, outstate; 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith PetscFunctionBegin; 11328b400f6SJacob Faibussowitsch PetscCheck(shell->apply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 11525e27a38SBarry Smith PetscCallBack("PCSHELL callback apply", (*shell->apply)(pc, x, y)); 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 1171e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 1181e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 122d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Shell(PC pc, Mat X, Mat Y) 123d71ae5a4SJacob Faibussowitsch { 1247b6e2003SPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 1257b6e2003SPierre Jolivet PetscObjectState instate, outstate; 1267b6e2003SPierre Jolivet 1277b6e2003SPierre Jolivet PetscFunctionBegin; 12828b400f6SJacob Faibussowitsch PetscCheck(shell->matapply, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)Y, &instate)); 13025e27a38SBarry Smith PetscCallBack("PCSHELL callback apply", (*shell->matapply)(pc, X, Y)); 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)Y, &outstate)); 1321e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 1331e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1357b6e2003SPierre Jolivet } 1367b6e2003SPierre Jolivet 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Shell(PC pc, Vec x, Vec y) 138d71ae5a4SJacob Faibussowitsch { 1391b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1401b581b66SBarry Smith 1411b581b66SBarry Smith PetscFunctionBegin; 14228b400f6SJacob Faibussowitsch PetscCheck(shell->applysymmetricleft, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 14325e27a38SBarry Smith PetscCallBack("PCSHELL callback apply symmetric left", (*shell->applysymmetricleft)(pc, x, y)); 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1451b581b66SBarry Smith } 1461b581b66SBarry Smith 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Shell(PC pc, Vec x, Vec y) 148d71ae5a4SJacob Faibussowitsch { 1491b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1501b581b66SBarry Smith 1511b581b66SBarry Smith PetscFunctionBegin; 15228b400f6SJacob Faibussowitsch PetscCheck(shell->applysymmetricright, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No apply() routine provided to Shell PC"); 15325e27a38SBarry Smith PetscCallBack("PCSHELL callback apply symmetric right", (*shell->applysymmetricright)(pc, x, y)); 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1551b581b66SBarry Smith } 1561b581b66SBarry Smith 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyBA_Shell(PC pc, PCSide side, Vec x, Vec y, Vec w) 158d71ae5a4SJacob Faibussowitsch { 159c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 160e3f487b0SBarry Smith PetscObjectState instate, outstate; 1612bb17772SBarry Smith 1622bb17772SBarry Smith PetscFunctionBegin; 16328b400f6SJacob Faibussowitsch PetscCheck(shell->applyBA, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyBA() routine provided to Shell PC"); 1649566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)w, &instate)); 16525e27a38SBarry Smith PetscCallBack("PCSHELL callback applyBA", (*shell->applyBA)(pc, side, x, y, w)); 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)w, &outstate)); 1671e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 1681e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)w)); 1693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1702bb17772SBarry Smith } 1712bb17772SBarry Smith 172d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolveChangeRHS_Shell(PC pc, PetscBool *change) 173d71ae5a4SJacob Faibussowitsch { 174a06fd7f2SStefano Zampini PetscFunctionBegin; 175a06fd7f2SStefano Zampini *change = PETSC_TRUE; 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177a06fd7f2SStefano Zampini } 178a06fd7f2SStefano Zampini 179d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolve_Shell(PC pc, KSP ksp, Vec b, Vec x) 180d71ae5a4SJacob Faibussowitsch { 181c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1827cdd61b2SBarry Smith 1837cdd61b2SBarry Smith PetscFunctionBegin; 18428b400f6SJacob Faibussowitsch PetscCheck(shell->presolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No presolve() routine provided to Shell PC"); 18525e27a38SBarry Smith PetscCallBack("PCSHELL callback presolve", (*shell->presolve)(pc, ksp, b, x)); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1877cdd61b2SBarry Smith } 1887cdd61b2SBarry Smith 189d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPostSolve_Shell(PC pc, KSP ksp, Vec b, Vec x) 190d71ae5a4SJacob Faibussowitsch { 191c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 1927cdd61b2SBarry Smith 1937cdd61b2SBarry Smith PetscFunctionBegin; 19428b400f6SJacob Faibussowitsch PetscCheck(shell->postsolve, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No postsolve() routine provided to Shell PC"); 19525e27a38SBarry Smith PetscCallBack("PCSHELL callback postsolve()", (*shell->postsolve)(pc, ksp, b, x)); 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1977cdd61b2SBarry Smith } 1987cdd61b2SBarry Smith 199d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Shell(PC pc, Vec x, Vec y) 200d71ae5a4SJacob Faibussowitsch { 201c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 202e3f487b0SBarry Smith PetscObjectState instate, outstate; 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith PetscFunctionBegin; 20528b400f6SJacob Faibussowitsch PetscCheck(shell->applytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applytranspose() routine provided to Shell PC"); 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 20725e27a38SBarry Smith PetscCallBack("PCSHELL callback applytranspose", (*shell->applytranspose)(pc, x, y)); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 2091e66621cSBarry Smith /* increase the state of the output vector if the user did not update its state themself as should have been done */ 2101e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2124b9ad928SBarry Smith } 2134b9ad928SBarry Smith 214*9fa64a6bSPierre Jolivet static PetscErrorCode PCMatApplyTranspose_Shell(PC pc, Mat x, Mat y) 215*9fa64a6bSPierre Jolivet { 216*9fa64a6bSPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 217*9fa64a6bSPierre Jolivet PetscObjectState instate, outstate; 218*9fa64a6bSPierre Jolivet 219*9fa64a6bSPierre Jolivet PetscFunctionBegin; 220*9fa64a6bSPierre Jolivet PetscCheck(shell->matapplytranspose, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No matapplytranspose() routine provided to Shell PC"); 221*9fa64a6bSPierre Jolivet PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 222*9fa64a6bSPierre Jolivet PetscCallBack("PCSHELL callback matapplytranspose", (*shell->matapplytranspose)(pc, x, y)); 223*9fa64a6bSPierre Jolivet PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 224*9fa64a6bSPierre Jolivet /* increase the state of the output matrix if the user did not update its state themself as should have been done */ 225*9fa64a6bSPierre Jolivet if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 226*9fa64a6bSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 227*9fa64a6bSPierre Jolivet } 228*9fa64a6bSPierre Jolivet 229d71ae5a4SJacob 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) 230d71ae5a4SJacob Faibussowitsch { 231c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 232e3f487b0SBarry Smith PetscObjectState instate, outstate; 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith PetscFunctionBegin; 23528b400f6SJacob Faibussowitsch PetscCheck(shell->applyrich, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No applyrichardson() routine provided to Shell PC"); 2369566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 23725e27a38SBarry Smith PetscCallBack("PCSHELL callback applyrichardson", (*shell->applyrich)(pc, x, y, w, rtol, abstol, dtol, it, guesszero, outits, reason)); 2389566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 239e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 2401e66621cSBarry Smith if (instate == outstate) PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2424b9ad928SBarry Smith } 2434b9ad928SBarry Smith 244d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Shell(PC pc) 245d71ae5a4SJacob Faibussowitsch { 2464b9ad928SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 2474b9ad928SBarry Smith 2484b9ad928SBarry Smith PetscFunctionBegin; 2499566063dSJacob Faibussowitsch PetscCall(PetscFree(shell->name)); 25025e27a38SBarry Smith if (shell->destroy) PetscCallBack("PCSHELL callback destroy", (*shell->destroy)(pc)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", NULL)); 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", NULL)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", NULL)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", NULL)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", NULL)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", NULL)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", NULL)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", NULL)); 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", NULL)); 262*9fa64a6bSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApplyTranspose_C", NULL)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", NULL)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 2679566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2694b9ad928SBarry Smith } 2704b9ad928SBarry Smith 271d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Shell(PC pc, PetscViewer viewer) 272d71ae5a4SJacob Faibussowitsch { 2734b9ad928SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 274ace3abfcSBarry Smith PetscBool iascii; 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith PetscFunctionBegin; 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 27832077d6dSBarry Smith if (iascii) { 2791e66621cSBarry Smith if (shell->name) PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", shell->name)); 2801e66621cSBarry Smith else PetscCall(PetscViewerASCIIPrintf(viewer, " no name\n")); 2814b9ad928SBarry Smith } 2824b9ad928SBarry Smith if (shell->view) { 2839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2849566063dSJacob Faibussowitsch PetscCall((*shell->view)(pc, viewer)); 2859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2864b9ad928SBarry Smith } 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2884b9ad928SBarry Smith } 2894b9ad928SBarry Smith 290d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetDestroy_Shell(PC pc, PetscErrorCode (*destroy)(PC)) 291d71ae5a4SJacob Faibussowitsch { 292c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 29318be62a5SSatish Balay 29418be62a5SSatish Balay PetscFunctionBegin; 29518be62a5SSatish Balay shell->destroy = destroy; 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29718be62a5SSatish Balay } 29818be62a5SSatish Balay 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetSetUp_Shell(PC pc, PetscErrorCode (*setup)(PC)) 300d71ae5a4SJacob Faibussowitsch { 301feb237baSPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith PetscFunctionBegin; 3044b9ad928SBarry Smith shell->setup = setup; 305d01c8aa3SLisandro Dalcin if (setup) pc->ops->setup = PCSetUp_Shell; 3060a545947SLisandro Dalcin else pc->ops->setup = NULL; 3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3084b9ad928SBarry Smith } 3094b9ad928SBarry Smith 310d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApply_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 311d71ae5a4SJacob Faibussowitsch { 312c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3134b9ad928SBarry Smith 3144b9ad928SBarry Smith PetscFunctionBegin; 3154b9ad928SBarry Smith shell->apply = apply; 3163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3174b9ad928SBarry Smith } 3184b9ad928SBarry Smith 319d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetMatApply_Shell(PC pc, PetscErrorCode (*matapply)(PC, Mat, Mat)) 320d71ae5a4SJacob Faibussowitsch { 3217b6e2003SPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 3227b6e2003SPierre Jolivet 3237b6e2003SPierre Jolivet PetscFunctionBegin; 3247b6e2003SPierre Jolivet shell->matapply = matapply; 3250e0fe96bSStefano Zampini if (matapply) pc->ops->matapply = PCMatApply_Shell; 3260e0fe96bSStefano Zampini else pc->ops->matapply = NULL; 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3287b6e2003SPierre Jolivet } 3297b6e2003SPierre Jolivet 330d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplySymmetricLeft_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 331d71ae5a4SJacob Faibussowitsch { 3321b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3331b581b66SBarry Smith 3341b581b66SBarry Smith PetscFunctionBegin; 3351b581b66SBarry Smith shell->applysymmetricleft = apply; 3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3371b581b66SBarry Smith } 3381b581b66SBarry Smith 339d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplySymmetricRight_Shell(PC pc, PetscErrorCode (*apply)(PC, Vec, Vec)) 340d71ae5a4SJacob Faibussowitsch { 3411b581b66SBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3421b581b66SBarry Smith 3431b581b66SBarry Smith PetscFunctionBegin; 3441b581b66SBarry Smith shell->applysymmetricright = apply; 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3461b581b66SBarry Smith } 3471b581b66SBarry Smith 348d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyBA_Shell(PC pc, PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec)) 349d71ae5a4SJacob Faibussowitsch { 350c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3512bb17772SBarry Smith 3522bb17772SBarry Smith PetscFunctionBegin; 353d01c8aa3SLisandro Dalcin shell->applyBA = applyBA; 354d01c8aa3SLisandro Dalcin if (applyBA) pc->ops->applyBA = PCApplyBA_Shell; 3550a545947SLisandro Dalcin else pc->ops->applyBA = NULL; 3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3572bb17772SBarry Smith } 3582bb17772SBarry Smith 3594d4d2bdcSBarry Smith static PetscErrorCode PCShellSetPreSolve_Shell(PC pc, PCShellPSolveFn *presolve) 360d71ae5a4SJacob Faibussowitsch { 361c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3627cdd61b2SBarry Smith 3637cdd61b2SBarry Smith PetscFunctionBegin; 3647cdd61b2SBarry Smith shell->presolve = presolve; 365a06fd7f2SStefano Zampini if (presolve) { 366a06fd7f2SStefano Zampini pc->ops->presolve = PCPreSolve_Shell; 3679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Shell)); 368a06fd7f2SStefano Zampini } else { 3690a545947SLisandro Dalcin pc->ops->presolve = NULL; 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 371a06fd7f2SStefano Zampini } 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3737cdd61b2SBarry Smith } 3747cdd61b2SBarry Smith 3754d4d2bdcSBarry Smith static PetscErrorCode PCShellSetPostSolve_Shell(PC pc, PCShellPSolveFn *postsolve) 376d71ae5a4SJacob Faibussowitsch { 377c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3787cdd61b2SBarry Smith 3797cdd61b2SBarry Smith PetscFunctionBegin; 3807cdd61b2SBarry Smith shell->postsolve = postsolve; 381d01c8aa3SLisandro Dalcin if (postsolve) pc->ops->postsolve = PCPostSolve_Shell; 3820a545947SLisandro Dalcin else pc->ops->postsolve = NULL; 3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3847cdd61b2SBarry Smith } 3857cdd61b2SBarry Smith 386d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetView_Shell(PC pc, PetscErrorCode (*view)(PC, PetscViewer)) 387d71ae5a4SJacob Faibussowitsch { 388c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith PetscFunctionBegin; 3914b9ad928SBarry Smith shell->view = view; 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3934b9ad928SBarry Smith } 3944b9ad928SBarry Smith 395d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyTranspose_Shell(PC pc, PetscErrorCode (*applytranspose)(PC, Vec, Vec)) 396d71ae5a4SJacob Faibussowitsch { 397c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith PetscFunctionBegin; 4004b9ad928SBarry Smith shell->applytranspose = applytranspose; 401d01c8aa3SLisandro Dalcin if (applytranspose) pc->ops->applytranspose = PCApplyTranspose_Shell; 4020a545947SLisandro Dalcin else pc->ops->applytranspose = NULL; 4033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 404d01c8aa3SLisandro Dalcin } 405d01c8aa3SLisandro Dalcin 406*9fa64a6bSPierre Jolivet static PetscErrorCode PCShellSetMatApplyTranspose_Shell(PC pc, PetscErrorCode (*matapplytranspose)(PC, Mat, Mat)) 407*9fa64a6bSPierre Jolivet { 408*9fa64a6bSPierre Jolivet PC_Shell *shell = (PC_Shell *)pc->data; 409*9fa64a6bSPierre Jolivet 410*9fa64a6bSPierre Jolivet PetscFunctionBegin; 411*9fa64a6bSPierre Jolivet shell->matapplytranspose = matapplytranspose; 412*9fa64a6bSPierre Jolivet if (matapplytranspose) pc->ops->matapplytranspose = PCMatApplyTranspose_Shell; 413*9fa64a6bSPierre Jolivet else pc->ops->matapplytranspose = NULL; 414*9fa64a6bSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 415*9fa64a6bSPierre Jolivet } 416*9fa64a6bSPierre Jolivet 417d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetApplyRichardson_Shell(PC pc, PetscErrorCode (*applyrich)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)) 418d71ae5a4SJacob Faibussowitsch { 419c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 420d01c8aa3SLisandro Dalcin 421d01c8aa3SLisandro Dalcin PetscFunctionBegin; 422d01c8aa3SLisandro Dalcin shell->applyrich = applyrich; 423d01c8aa3SLisandro Dalcin if (applyrich) pc->ops->applyrichardson = PCApplyRichardson_Shell; 4240a545947SLisandro Dalcin else pc->ops->applyrichardson = NULL; 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4264b9ad928SBarry Smith } 4274b9ad928SBarry Smith 428d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellSetName_Shell(PC pc, const char name[]) 429d71ae5a4SJacob Faibussowitsch { 430c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith PetscFunctionBegin; 4339566063dSJacob Faibussowitsch PetscCall(PetscFree(shell->name)); 4349566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &shell->name)); 4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4364b9ad928SBarry Smith } 4374b9ad928SBarry Smith 438d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCShellGetName_Shell(PC pc, const char *name[]) 439d71ae5a4SJacob Faibussowitsch { 440c5ae4b9aSBarry Smith PC_Shell *shell = (PC_Shell *)pc->data; 4414b9ad928SBarry Smith 4424b9ad928SBarry Smith PetscFunctionBegin; 4434b9ad928SBarry Smith *name = shell->name; 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4454b9ad928SBarry Smith } 4464b9ad928SBarry Smith 44718be62a5SSatish Balay /*@C 44804c3f3b8SBarry Smith PCShellSetDestroy - Sets routine to use to destroy the user-provided application context that was provided with `PCShellSetContext()` 44918be62a5SSatish Balay 450c3339decSBarry Smith Logically Collective 45118be62a5SSatish Balay 45218be62a5SSatish Balay Input Parameters: 45318be62a5SSatish Balay + pc - the preconditioner context 454a2b725a8SWilliam Gropp - destroy - the application-provided destroy routine 45518be62a5SSatish Balay 45620f4b53cSBarry Smith Calling sequence of `destroy`: 45704c3f3b8SBarry Smith . pc - the preconditioner 4584aa34b0aSBarry Smith 459f1580f4eSBarry Smith Level: intermediate 46018be62a5SSatish Balay 461562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()` 46218be62a5SSatish Balay @*/ 46304c3f3b8SBarry Smith PetscErrorCode PCShellSetDestroy(PC pc, PetscErrorCode (*destroy)(PC pc)) 464d71ae5a4SJacob Faibussowitsch { 46518be62a5SSatish Balay PetscFunctionBegin; 4660700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 467cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetDestroy_C", (PC, PetscErrorCode (*)(PC)), (pc, destroy)); 4683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46918be62a5SSatish Balay } 47018be62a5SSatish Balay 4714b9ad928SBarry Smith /*@C 4724b9ad928SBarry Smith PCShellSetSetUp - Sets routine to use to "setup" the preconditioner whenever the 4734b9ad928SBarry Smith matrix operator is changed. 4744b9ad928SBarry Smith 475c3339decSBarry Smith Logically Collective 4764b9ad928SBarry Smith 4774b9ad928SBarry Smith Input Parameters: 4784b9ad928SBarry Smith + pc - the preconditioner context 479a2b725a8SWilliam Gropp - setup - the application-provided setup routine 4804b9ad928SBarry Smith 48120f4b53cSBarry Smith Calling sequence of `setup`: 48204c3f3b8SBarry Smith . pc - the preconditioner 4834aa34b0aSBarry Smith 484f1580f4eSBarry Smith Level: intermediate 4854b9ad928SBarry Smith 48604c3f3b8SBarry Smith Note: 48704c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `setup`. 48804c3f3b8SBarry Smith 489562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetApply()`, `PCShellSetContext()`, , `PCShellGetContext()` 4904b9ad928SBarry Smith @*/ 49104c3f3b8SBarry Smith PetscErrorCode PCShellSetSetUp(PC pc, PetscErrorCode (*setup)(PC pc)) 492d71ae5a4SJacob Faibussowitsch { 4934b9ad928SBarry Smith PetscFunctionBegin; 4940700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 495cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetSetUp_C", (PC, PetscErrorCode (*)(PC)), (pc, setup)); 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4974b9ad928SBarry Smith } 4984b9ad928SBarry Smith 4994b9ad928SBarry Smith /*@C 50004c3f3b8SBarry Smith PCShellSetView - Sets routine to use as viewer of a `PCSHELL` shell preconditioner 5014b9ad928SBarry Smith 502c3339decSBarry Smith Logically Collective 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith Input Parameters: 5054b9ad928SBarry Smith + pc - the preconditioner context 5064b9ad928SBarry Smith - view - the application-provided view routine 5074b9ad928SBarry Smith 50820f4b53cSBarry Smith Calling sequence of `view`: 50904c3f3b8SBarry Smith + pc - the preconditioner 5104b9ad928SBarry Smith - v - viewer 5114b9ad928SBarry Smith 512f1580f4eSBarry Smith Level: advanced 5134b9ad928SBarry Smith 51404c3f3b8SBarry Smith Note: 51504c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `view`. 51604c3f3b8SBarry Smith 517562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()` 5184b9ad928SBarry Smith @*/ 51904c3f3b8SBarry Smith PetscErrorCode PCShellSetView(PC pc, PetscErrorCode (*view)(PC pc, PetscViewer v)) 520d71ae5a4SJacob Faibussowitsch { 5214b9ad928SBarry Smith PetscFunctionBegin; 5220700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 523cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetView_C", (PC, PetscErrorCode (*)(PC, PetscViewer)), (pc, view)); 5243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5254b9ad928SBarry Smith } 5264b9ad928SBarry Smith 5274b9ad928SBarry Smith /*@C 5284b9ad928SBarry Smith PCShellSetApply - Sets routine to use as preconditioner. 5294b9ad928SBarry Smith 530c3339decSBarry Smith Logically Collective 5314b9ad928SBarry Smith 5324b9ad928SBarry Smith Input Parameters: 5334b9ad928SBarry Smith + pc - the preconditioner context 534be29d3c6SBarry Smith - apply - the application-provided preconditioning routine 5354b9ad928SBarry Smith 53620f4b53cSBarry Smith Calling sequence of `apply`: 53720f4b53cSBarry Smith + pc - the preconditioner, get the application context with `PCShellGetContext()` 5384b9ad928SBarry Smith . xin - input vector 5394b9ad928SBarry Smith - xout - output vector 5404b9ad928SBarry Smith 541f1580f4eSBarry Smith Level: intermediate 5424b9ad928SBarry Smith 54304c3f3b8SBarry Smith Note: 54404c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 54504c3f3b8SBarry Smith 546562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellSetApplySymmetricRight()`, `PCShellSetApplySymmetricLeft()`, `PCShellGetContext()` 5474b9ad928SBarry Smith @*/ 54804c3f3b8SBarry Smith PetscErrorCode PCShellSetApply(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 549d71ae5a4SJacob Faibussowitsch { 5504b9ad928SBarry Smith PetscFunctionBegin; 5510700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 552cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApply_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply)); 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5544b9ad928SBarry Smith } 5554b9ad928SBarry Smith 5561b581b66SBarry Smith /*@C 5576437efc7SEric Chamberland PCShellSetMatApply - Sets routine to use as preconditioner on a block of vectors. 5587b6e2003SPierre Jolivet 559c3339decSBarry Smith Logically Collective 5607b6e2003SPierre Jolivet 5617b6e2003SPierre Jolivet Input Parameters: 5627b6e2003SPierre Jolivet + pc - the preconditioner context 56304c3f3b8SBarry Smith - matapply - the application-provided preconditioning routine 5647b6e2003SPierre Jolivet 56504c3f3b8SBarry Smith Calling sequence of `matapply`: 56604c3f3b8SBarry Smith + pc - the preconditioner 56704c3f3b8SBarry Smith . Xin - input block of vectors represented as a dense `Mat` 56804c3f3b8SBarry Smith - Xout - output block of vectors represented as a dense `Mat` 5697b6e2003SPierre Jolivet 570f1580f4eSBarry Smith Level: advanced 5717b6e2003SPierre Jolivet 57204c3f3b8SBarry Smith Note: 57304c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapply`. 57404c3f3b8SBarry Smith 575562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellGetContext()` 5767b6e2003SPierre Jolivet @*/ 57704c3f3b8SBarry Smith PetscErrorCode PCShellSetMatApply(PC pc, PetscErrorCode (*matapply)(PC pc, Mat Xin, Mat Xout)) 578d71ae5a4SJacob Faibussowitsch { 5797b6e2003SPierre Jolivet PetscFunctionBegin; 5807b6e2003SPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 581cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetMatApply_C", (PC, PetscErrorCode (*)(PC, Mat, Mat)), (pc, matapply)); 5823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5837b6e2003SPierre Jolivet } 5847b6e2003SPierre Jolivet 5857b6e2003SPierre Jolivet /*@C 586f1580f4eSBarry Smith PCShellSetApplySymmetricLeft - Sets routine to use as left preconditioner (when the `PC_SYMMETRIC` is used). 5871b581b66SBarry Smith 588c3339decSBarry Smith Logically Collective 5891b581b66SBarry Smith 5901b581b66SBarry Smith Input Parameters: 5911b581b66SBarry Smith + pc - the preconditioner context 5921b581b66SBarry Smith - apply - the application-provided left preconditioning routine 5931b581b66SBarry Smith 59420f4b53cSBarry Smith Calling sequence of `apply`: 59504c3f3b8SBarry Smith + pc - the preconditioner 5961b581b66SBarry Smith . xin - input vector 5971b581b66SBarry Smith - xout - output vector 5981b581b66SBarry Smith 599f1580f4eSBarry Smith Level: advanced 6001b581b66SBarry Smith 60104c3f3b8SBarry Smith Note: 60204c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 60304c3f3b8SBarry Smith 604562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()` 6051b581b66SBarry Smith @*/ 60604c3f3b8SBarry Smith PetscErrorCode PCShellSetApplySymmetricLeft(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 607d71ae5a4SJacob Faibussowitsch { 6081b581b66SBarry Smith PetscFunctionBegin; 6091b581b66SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 610cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplySymmetricLeft_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply)); 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6121b581b66SBarry Smith } 6131b581b66SBarry Smith 6141b581b66SBarry Smith /*@C 615f1580f4eSBarry Smith PCShellSetApplySymmetricRight - Sets routine to use as right preconditioner (when the `PC_SYMMETRIC` is used). 6161b581b66SBarry Smith 617c3339decSBarry Smith Logically Collective 6181b581b66SBarry Smith 6191b581b66SBarry Smith Input Parameters: 6201b581b66SBarry Smith + pc - the preconditioner context 6211b581b66SBarry Smith - apply - the application-provided right preconditioning routine 6221b581b66SBarry Smith 62320f4b53cSBarry Smith Calling sequence of `apply`: 62404c3f3b8SBarry Smith + pc - the preconditioner 6251b581b66SBarry Smith . xin - input vector 6261b581b66SBarry Smith - xout - output vector 6271b581b66SBarry Smith 628f1580f4eSBarry Smith Level: advanced 6291b581b66SBarry Smith 63004c3f3b8SBarry Smith Note: 63104c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 63204c3f3b8SBarry Smith 633562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetApplySymmetricLeft()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellGetContext()` 6341b581b66SBarry Smith @*/ 63504c3f3b8SBarry Smith PetscErrorCode PCShellSetApplySymmetricRight(PC pc, PetscErrorCode (*apply)(PC pc, Vec xin, Vec xout)) 636d71ae5a4SJacob Faibussowitsch { 6371b581b66SBarry Smith PetscFunctionBegin; 6381b581b66SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 639cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplySymmetricRight_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, apply)); 6403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6411b581b66SBarry Smith } 6421b581b66SBarry Smith 6432bb17772SBarry Smith /*@C 64404c3f3b8SBarry Smith PCShellSetApplyBA - Sets routine to use as the preconditioner times the operator. 6452bb17772SBarry Smith 646c3339decSBarry Smith Logically Collective 6472bb17772SBarry Smith 6482bb17772SBarry Smith Input Parameters: 6492bb17772SBarry Smith + pc - the preconditioner context 6502bb17772SBarry Smith - applyBA - the application-provided BA routine 6512bb17772SBarry Smith 65220f4b53cSBarry Smith Calling sequence of `applyBA`: 65304c3f3b8SBarry Smith + pc - the preconditioner 65404c3f3b8SBarry Smith . side - `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 6552bb17772SBarry Smith . xin - input vector 65604c3f3b8SBarry Smith . xout - output vector 65704c3f3b8SBarry Smith - w - work vector 6582bb17772SBarry Smith 659f1580f4eSBarry Smith Level: intermediate 6602bb17772SBarry Smith 66104c3f3b8SBarry Smith Note: 66204c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applyBA`. 66304c3f3b8SBarry Smith 664562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCShellGetContext()`, `PCSide` 6652bb17772SBarry Smith @*/ 66604c3f3b8SBarry Smith PetscErrorCode PCShellSetApplyBA(PC pc, PetscErrorCode (*applyBA)(PC pc, PCSide side, Vec xin, Vec xout, Vec w)) 667d71ae5a4SJacob Faibussowitsch { 6682bb17772SBarry Smith PetscFunctionBegin; 6690700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 670cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplyBA_C", (PC, PetscErrorCode (*)(PC, PCSide, Vec, Vec, Vec)), (pc, applyBA)); 6713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6722bb17772SBarry Smith } 6732bb17772SBarry Smith 6744b9ad928SBarry Smith /*@C 6754b9ad928SBarry Smith PCShellSetApplyTranspose - Sets routine to use as preconditioner transpose. 6764b9ad928SBarry Smith 677c3339decSBarry Smith Logically Collective 6784b9ad928SBarry Smith 6794b9ad928SBarry Smith Input Parameters: 6804b9ad928SBarry Smith + pc - the preconditioner context 68120f4b53cSBarry Smith - applytranspose - the application-provided preconditioning transpose routine 6824b9ad928SBarry Smith 68320f4b53cSBarry Smith Calling sequence of `applytranspose`: 68404c3f3b8SBarry Smith + pc - the preconditioner 6854b9ad928SBarry Smith . xin - input vector 6864b9ad928SBarry Smith - xout - output vector 6874b9ad928SBarry Smith 688f1580f4eSBarry Smith Level: intermediate 6894b9ad928SBarry Smith 690f1580f4eSBarry Smith Note: 69104c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `applytranspose`. 6924b9ad928SBarry Smith 69354c05997SPierre Jolivet .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellGetContext()` 6944b9ad928SBarry Smith @*/ 69504c3f3b8SBarry Smith PetscErrorCode PCShellSetApplyTranspose(PC pc, PetscErrorCode (*applytranspose)(PC pc, Vec xin, Vec xout)) 696d71ae5a4SJacob Faibussowitsch { 6974b9ad928SBarry Smith PetscFunctionBegin; 6980700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 699cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplyTranspose_C", (PC, PetscErrorCode (*)(PC, Vec, Vec)), (pc, applytranspose)); 7003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7014b9ad928SBarry Smith } 7024b9ad928SBarry Smith 7037cdd61b2SBarry Smith /*@C 704*9fa64a6bSPierre Jolivet PCShellSetMatApplyTranspose - Sets routine to use as preconditioner transpose. 705*9fa64a6bSPierre Jolivet 706*9fa64a6bSPierre Jolivet Logically Collective 707*9fa64a6bSPierre Jolivet 708*9fa64a6bSPierre Jolivet Input Parameters: 709*9fa64a6bSPierre Jolivet + pc - the preconditioner context 710*9fa64a6bSPierre Jolivet - matapplytranspose - the application-provided preconditioning transpose routine 711*9fa64a6bSPierre Jolivet 712*9fa64a6bSPierre Jolivet Calling sequence of `matapplytranspose`: 713*9fa64a6bSPierre Jolivet + pc - the preconditioner 714*9fa64a6bSPierre Jolivet . xin - input matrix 715*9fa64a6bSPierre Jolivet - xout - output matrix 716*9fa64a6bSPierre Jolivet 717*9fa64a6bSPierre Jolivet Level: intermediate 718*9fa64a6bSPierre Jolivet 719*9fa64a6bSPierre Jolivet Note: 720*9fa64a6bSPierre Jolivet You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `matapplytranspose`. 721*9fa64a6bSPierre Jolivet 722*9fa64a6bSPierre Jolivet .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetContext()`, `PCShellSetApplyBA()`, `PCShellGetContext()` 723*9fa64a6bSPierre Jolivet @*/ 724*9fa64a6bSPierre Jolivet PetscErrorCode PCShellSetMatApplyTranspose(PC pc, PetscErrorCode (*matapplytranspose)(PC pc, Mat xin, Mat xout)) 725*9fa64a6bSPierre Jolivet { 726*9fa64a6bSPierre Jolivet PetscFunctionBegin; 727*9fa64a6bSPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 728*9fa64a6bSPierre Jolivet PetscTryMethod(pc, "PCShellSetMatApplyTranspose_C", (PC, PetscErrorCode (*)(PC, Mat, Mat)), (pc, matapplytranspose)); 729*9fa64a6bSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 730*9fa64a6bSPierre Jolivet } 731*9fa64a6bSPierre Jolivet 732*9fa64a6bSPierre Jolivet /*@C 733f1580f4eSBarry Smith PCShellSetPreSolve - Sets routine to apply to the operators/vectors before a `KSPSolve()` is 7347cdd61b2SBarry Smith applied. This usually does something like scale the linear system in some application 7357cdd61b2SBarry Smith specific way. 7367cdd61b2SBarry Smith 737c3339decSBarry Smith Logically Collective 7387cdd61b2SBarry Smith 7397cdd61b2SBarry Smith Input Parameters: 7407cdd61b2SBarry Smith + pc - the preconditioner context 7414d4d2bdcSBarry Smith - presolve - the application-provided presolve routine, see `PCShellPSolveFn` 7427cdd61b2SBarry Smith 743f1580f4eSBarry Smith Level: advanced 7447cdd61b2SBarry Smith 74504c3f3b8SBarry Smith Note: 74604c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `presolve`. 74704c3f3b8SBarry Smith 7484d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellPSolveFn`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPostSolve()`, `PCShellSetContext()`, `PCShellGetContext()` 7497cdd61b2SBarry Smith @*/ 7504d4d2bdcSBarry Smith PetscErrorCode PCShellSetPreSolve(PC pc, PCShellPSolveFn *presolve) 751d71ae5a4SJacob Faibussowitsch { 7527cdd61b2SBarry Smith PetscFunctionBegin; 7530700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 7544d4d2bdcSBarry Smith PetscTryMethod(pc, "PCShellSetPreSolve_C", (PC, PCShellPSolveFn *), (pc, presolve)); 7553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7567cdd61b2SBarry Smith } 7577cdd61b2SBarry Smith 7587cdd61b2SBarry Smith /*@C 75904c3f3b8SBarry Smith PCShellSetPostSolve - Sets routine to apply to the operators/vectors after a `KSPSolve()` is 7607cdd61b2SBarry Smith applied. This usually does something like scale the linear system in some application 7617cdd61b2SBarry Smith specific way. 7627cdd61b2SBarry Smith 763c3339decSBarry Smith Logically Collective 7647cdd61b2SBarry Smith 7657cdd61b2SBarry Smith Input Parameters: 7667cdd61b2SBarry Smith + pc - the preconditioner context 7674d4d2bdcSBarry Smith - postsolve - the application-provided postsolve routine, see `PCShellPSolveFn` 7687cdd61b2SBarry Smith 769f1580f4eSBarry Smith Level: advanced 7707cdd61b2SBarry Smith 77104c3f3b8SBarry Smith Note: 77204c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `postsolve`. 77304c3f3b8SBarry Smith 7744d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellPSolveFn`, `PCShellSetApplyRichardson()`, `PCShellSetSetUp()`, `PCShellSetApplyTranspose()`, `PCShellSetPreSolve()`, `PCShellSetContext()`, `PCShellGetContext()` 7757cdd61b2SBarry Smith @*/ 7764d4d2bdcSBarry Smith PetscErrorCode PCShellSetPostSolve(PC pc, PCShellPSolveFn *postsolve) 777d71ae5a4SJacob Faibussowitsch { 7787cdd61b2SBarry Smith PetscFunctionBegin; 7790700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 7804d4d2bdcSBarry Smith PetscTryMethod(pc, "PCShellSetPostSolve_C", (PC, PCShellPSolveFn *), (pc, postsolve)); 7813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7827cdd61b2SBarry Smith } 7837cdd61b2SBarry Smith 784cc4c1da9SBarry Smith /*@ 785f1580f4eSBarry Smith PCShellSetName - Sets an optional name to associate with a `PCSHELL` 7864b9ad928SBarry Smith preconditioner. 7874b9ad928SBarry Smith 7884b9ad928SBarry Smith Not Collective 7894b9ad928SBarry Smith 7904b9ad928SBarry Smith Input Parameters: 7914b9ad928SBarry Smith + pc - the preconditioner context 7924b9ad928SBarry Smith - name - character string describing shell preconditioner 7934b9ad928SBarry Smith 794f1580f4eSBarry Smith Level: intermediate 7954b9ad928SBarry Smith 79604c3f3b8SBarry Smith Note: 797baca6076SPierre Jolivet This is separate from the name you can provide with `PetscObjectSetName()` 79804c3f3b8SBarry Smith 799562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellGetName()`, `PetscObjectSetName()`, `PetscObjectGetName()` 8004b9ad928SBarry Smith @*/ 801d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellSetName(PC pc, const char name[]) 802d71ae5a4SJacob Faibussowitsch { 8034b9ad928SBarry Smith PetscFunctionBegin; 8040700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 805cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetName_C", (PC, const char[]), (pc, name)); 8063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8074b9ad928SBarry Smith } 8084b9ad928SBarry Smith 809cc4c1da9SBarry Smith /*@ 81004c3f3b8SBarry Smith PCShellGetName - Gets an optional name that the user has set for a `PCSHELL` with `PCShellSetName()` 8114b9ad928SBarry Smith preconditioner. 8124b9ad928SBarry Smith 8134b9ad928SBarry Smith Not Collective 8144b9ad928SBarry Smith 8154b9ad928SBarry Smith Input Parameter: 8164b9ad928SBarry Smith . pc - the preconditioner context 8174b9ad928SBarry Smith 8184b9ad928SBarry Smith Output Parameter: 8194b9ad928SBarry Smith . name - character string describing shell preconditioner (you should not free this) 8204b9ad928SBarry Smith 821f1580f4eSBarry Smith Level: intermediate 8224b9ad928SBarry Smith 823562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetName()`, `PetscObjectSetName()`, `PetscObjectGetName()` 8244b9ad928SBarry Smith @*/ 825d71ae5a4SJacob Faibussowitsch PetscErrorCode PCShellGetName(PC pc, const char *name[]) 826d71ae5a4SJacob Faibussowitsch { 8274b9ad928SBarry Smith PetscFunctionBegin; 8280700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8294f572ea9SToby Isaac PetscAssertPointer(name, 2); 830cac4c232SBarry Smith PetscUseMethod(pc, "PCShellGetName_C", (PC, const char *[]), (pc, name)); 8313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8324b9ad928SBarry Smith } 8334b9ad928SBarry Smith 8344b9ad928SBarry Smith /*@C 8354b9ad928SBarry Smith PCShellSetApplyRichardson - Sets routine to use as preconditioner 8364b9ad928SBarry Smith in Richardson iteration. 8374b9ad928SBarry Smith 838c3339decSBarry Smith Logically Collective 8394b9ad928SBarry Smith 8404b9ad928SBarry Smith Input Parameters: 8414b9ad928SBarry Smith + pc - the preconditioner context 842be29d3c6SBarry Smith - apply - the application-provided preconditioning routine 8434b9ad928SBarry Smith 84420f4b53cSBarry Smith Calling sequence of `apply`: 84504c3f3b8SBarry Smith + pc - the preconditioner 846dd8e379bSPierre Jolivet . b - right-hand side 8474b9ad928SBarry Smith . x - current iterate 8484b9ad928SBarry Smith . r - work space 8494b9ad928SBarry Smith . rtol - relative tolerance of residual norm to stop at 85070441072SBarry Smith . abstol - absolute tolerance of residual norm to stop at 8514b9ad928SBarry Smith . dtol - if residual norm increases by this factor than return 85204c3f3b8SBarry Smith . maxits - number of iterations to run 85304c3f3b8SBarry Smith . zeroinitialguess - `PETSC_TRUE` if `x` is known to be initially zero 85404c3f3b8SBarry Smith . its - returns the number of iterations used 85504c3f3b8SBarry Smith - reason - returns the reason the iteration has converged 8564b9ad928SBarry Smith 857f1580f4eSBarry Smith Level: advanced 8584b9ad928SBarry Smith 8590b4b7b1cSBarry Smith Notes: 86004c3f3b8SBarry Smith You can get the `PCSHELL` context set with `PCShellSetContext()` using `PCShellGetContext()` if needed by `apply`. 86104c3f3b8SBarry Smith 8620b4b7b1cSBarry 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, 8630b4b7b1cSBarry Smith recomputing the residual via $ r = b - A x $, and then computing the next step. SOR is an algorithm for which this is true. 8640b4b7b1cSBarry Smith 8650b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PCSHELL`, `PCShellSetApply()`, `PCShellSetContext()`, `PCRichardsonConvergedReason()`, `PCShellGetContext()`, `KSPRICHARDSON` 8664b9ad928SBarry Smith @*/ 86704c3f3b8SBarry 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)) 868d71ae5a4SJacob Faibussowitsch { 8694b9ad928SBarry Smith PetscFunctionBegin; 8700700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 871cac4c232SBarry Smith PetscTryMethod(pc, "PCShellSetApplyRichardson_C", (PC, PetscErrorCode (*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)), (pc, apply)); 8723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8734b9ad928SBarry Smith } 8744b9ad928SBarry Smith 8754b9ad928SBarry Smith /*MC 876f1580f4eSBarry Smith PCSHELL - Creates a new preconditioner class for use with a users 877f1580f4eSBarry Smith own private data storage format and preconditioner application code 8784b9ad928SBarry Smith 8794b9ad928SBarry Smith Level: advanced 880e0bb08deSStefano Zampini 8814b9ad928SBarry Smith Usage: 882f1580f4eSBarry Smith .vb 883f1580f4eSBarry Smith extern PetscErrorCode apply(PC,Vec,Vec); 884f1580f4eSBarry Smith extern PetscErrorCode applyba(PC,PCSide,Vec,Vec,Vec); 885f1580f4eSBarry Smith extern PetscErrorCode applytranspose(PC,Vec,Vec); 886f1580f4eSBarry Smith extern PetscErrorCode setup(PC); 887f1580f4eSBarry Smith extern PetscErrorCode destroy(PC); 888f1580f4eSBarry Smith 889f1580f4eSBarry Smith PCCreate(comm,&pc); 890f1580f4eSBarry Smith PCSetType(pc,PCSHELL); 891f1580f4eSBarry Smith PCShellSetContext(pc,ctx) 892f1580f4eSBarry Smith PCShellSetApply(pc,apply); 893f1580f4eSBarry Smith PCShellSetApplyBA(pc,applyba); (optional) 894f1580f4eSBarry Smith PCShellSetApplyTranspose(pc,applytranspose); (optional) 895f1580f4eSBarry Smith PCShellSetSetUp(pc,setup); (optional) 896f1580f4eSBarry Smith PCShellSetDestroy(pc,destroy); (optional) 897f1580f4eSBarry Smith .ve 8984b9ad928SBarry Smith 8990b4b7b1cSBarry Smith Notes: 90004c3f3b8SBarry Smith Information required for the preconditioner and its internal datastructures can be set with `PCShellSetContext()` and then accessed 9010b4b7b1cSBarry Smith with `PCShellGetContext()` inside the routines provided above. 9020b4b7b1cSBarry Smith 9030b4b7b1cSBarry Smith When using `MATSHELL`, where the explicit entries of matrix are not available to build the preconditioner, `PCSHELL` can be used 9040b4b7b1cSBarry Smith to construct a custom preconditioner for the `MATSHELL`, assuming the user knows enough about their problem to provide a 9050b4b7b1cSBarry Smith custom preconditioner. 90604c3f3b8SBarry Smith 907562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 908f1580f4eSBarry Smith `MATSHELL`, `PCShellSetSetUp()`, `PCShellSetApply()`, `PCShellSetView()`, `PCShellSetDestroy()`, `PCShellSetPostSolve()`, 909f1580f4eSBarry Smith `PCShellSetApplyTranspose()`, `PCShellSetName()`, `PCShellSetApplyRichardson()`, `PCShellSetPreSolve()`, `PCShellSetView()`, 910f1580f4eSBarry Smith `PCShellGetName()`, `PCShellSetContext()`, `PCShellGetContext()`, `PCShellSetApplyBA()`, `MATSHELL`, `PCShellSetMatApply()`, 9114b9ad928SBarry Smith M*/ 9124b9ad928SBarry Smith 913d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Shell(PC pc) 914d71ae5a4SJacob Faibussowitsch { 9154b9ad928SBarry Smith PC_Shell *shell; 9164b9ad928SBarry Smith 9174b9ad928SBarry Smith PetscFunctionBegin; 9184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&shell)); 9194b9ad928SBarry Smith pc->data = (void *)shell; 9204b9ad928SBarry Smith 921d01c8aa3SLisandro Dalcin pc->ops->destroy = PCDestroy_Shell; 9224b9ad928SBarry Smith pc->ops->view = PCView_Shell; 923d01c8aa3SLisandro Dalcin pc->ops->apply = PCApply_Shell; 9241b581b66SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_Shell; 9251b581b66SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_Shell; 9260e0fe96bSStefano Zampini pc->ops->matapply = NULL; 9270a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 9280a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 9290a545947SLisandro Dalcin pc->ops->setup = NULL; 9300a545947SLisandro Dalcin pc->ops->presolve = NULL; 9310a545947SLisandro Dalcin pc->ops->postsolve = NULL; 9324b9ad928SBarry Smith 9330a545947SLisandro Dalcin shell->apply = NULL; 9340a545947SLisandro Dalcin shell->applytranspose = NULL; 9350a545947SLisandro Dalcin shell->name = NULL; 9360a545947SLisandro Dalcin shell->applyrich = NULL; 9370a545947SLisandro Dalcin shell->presolve = NULL; 9380a545947SLisandro Dalcin shell->postsolve = NULL; 9390a545947SLisandro Dalcin shell->ctx = NULL; 9400a545947SLisandro Dalcin shell->setup = NULL; 9410a545947SLisandro Dalcin shell->view = NULL; 9420a545947SLisandro Dalcin shell->destroy = NULL; 9430a545947SLisandro Dalcin shell->applysymmetricleft = NULL; 9440a545947SLisandro Dalcin shell->applysymmetricright = NULL; 9454b9ad928SBarry Smith 9469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetDestroy_C", PCShellSetDestroy_Shell)); 9479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetSetUp_C", PCShellSetSetUp_Shell)); 9489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApply_C", PCShellSetApply_Shell)); 9499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApply_C", PCShellSetMatApply_Shell)); 9509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricLeft_C", PCShellSetApplySymmetricLeft_Shell)); 9519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplySymmetricRight_C", PCShellSetApplySymmetricRight_Shell)); 9529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyBA_C", PCShellSetApplyBA_Shell)); 9539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPreSolve_C", PCShellSetPreSolve_Shell)); 9549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetPostSolve_C", PCShellSetPostSolve_Shell)); 9559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetView_C", PCShellSetView_Shell)); 9569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyTranspose_C", PCShellSetApplyTranspose_Shell)); 957*9fa64a6bSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetMatApplyTranspose_C", PCShellSetMatApplyTranspose_Shell)); 9589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetName_C", PCShellSetName_Shell)); 9599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellGetName_C", PCShellGetName_Shell)); 9609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCShellSetApplyRichardson_C", PCShellSetApplyRichardson_Shell)); 9613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9624b9ad928SBarry Smith } 963