1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 5ac09b921SBarry Smith 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8ac09b921SBarry Smith /* SUBMANSEC = KSP */ 9ac09b921SBarry Smith 10607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 114bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void); 121dbb0a54SBarry Smith 1328ce4d24SBarry Smith /*S 140b4b7b1cSBarry Smith KSP - Abstract PETSc object that manages the linear solves in PETSc (even those such as direct factorization-based solvers that 150b4b7b1cSBarry Smith do not use Krylov accelerators). 1628ce4d24SBarry Smith 1728ce4d24SBarry Smith Level: beginner 1828ce4d24SBarry Smith 190b4b7b1cSBarry Smith Notes: 20a4d1885cSBarry Smith When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a 2116a05f60SBarry Smith `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 230b4b7b1cSBarry Smith Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object 240b4b7b1cSBarry Smith 250b4b7b1cSBarry Smith The `PC` object is used to control preconditioners in PETSc. 260b4b7b1cSBarry Smith 27789736e1SBarry Smith `KSP` can also be used to solve some least squares problems (over or under-determined linear systems), using, for example, `KSPLSQR`, see `PETSCREGRESSORLINEAR` 28789736e1SBarry Smith for additional methods that can be used to solve least squares problems and other linear regressions). 29789736e1SBarry Smith 301cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES` 3128ce4d24SBarry Smith S*/ 32e2a1c21fSSatish Balay typedef struct _p_KSP *KSP; 332eac72dbSBarry Smith 3476bdecfbSBarry Smith /*J 350b4b7b1cSBarry Smith KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides. 3628ce4d24SBarry Smith 3728ce4d24SBarry Smith Level: beginner 3828ce4d24SBarry Smith 391cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()` 4076bdecfbSBarry Smith J*/ 4119fd82e9SBarry Smith typedef const char *KSPType; 4282bf6240SBarry Smith #define KSPRICHARDSON "richardson" 436c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 4482bf6240SBarry Smith #define KSPCG "cg" 452c8fc646SJed Brown #define KSPGROPPCG "groppcg" 462c8fc646SJed Brown #define KSPPIPECG "pipecg" 47901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 48265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 49b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 50325e8391SManasi #define KSPPIPECG2 "pipecg2" 51df2a969aSvictorle #define KSPCGNE "cgne" 5205de396fSBarry Smith #define KSPNASH "nash" 5305de396fSBarry Smith #define KSPSTCG "stcg" 5405de396fSBarry Smith #define KSPGLTR "gltr" 55edd03b47SJacob Faibussowitsch #define KSPCGNASH PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash" 56edd03b47SJacob Faibussowitsch #define KSPCGSTCG PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg" 57edd03b47SJacob Faibussowitsch #define KSPCGGLTR PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr" 58640f4f53SPatrick Sanan #define KSPFCG "fcg" 59390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 6082bf6240SBarry Smith #define KSPGMRES "gmres" 61483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 629dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 639dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 644f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 6561b468f9SJed Brown #define KSPPGMRES "pgmres" 6682bf6240SBarry Smith #define KSPTCQMR "tcqmr" 6782bf6240SBarry Smith #define KSPBCGS "bcgs" 68e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 69345ecf0bSXiangmin Jiao #define KSPQMRCGS "qmrcgs" 7018ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 71c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 72f0037002Svictorle #define KSPBCGSL "bcgsl" 73f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 7482bf6240SBarry Smith #define KSPCGS "cgs" 7582bf6240SBarry Smith #define KSPTFQMR "tfqmr" 7682bf6240SBarry Smith #define KSPCR "cr" 772c8fc646SJed Brown #define KSPPIPECR "pipecr" 7882bf6240SBarry Smith #define KSPLSQR "lsqr" 7982bf6240SBarry Smith #define KSPPREONLY "preonly" 803c2be86cSBarry Smith #define KSPNONE "none" 8182bf6240SBarry Smith #define KSPQCG "qcg" 82c9cf9b11SBarry Smith #define KSPBICG "bicg" 83b4ac9ba4SBarry Smith #define KSPMINRES "minres" 8401c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 8521356919SSatish Balay #define KSPLCD "lcd" 861d6018f0SLisandro Dalcin #define KSPPYTHON "python" 8758ad63f4SBarry Smith #define KSPGCR "gcr" 88fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 89e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 90e4d80e07Szianekhodja #define KSPCGLS "cgls" 91329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 9238cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 932eac72dbSBarry Smith 948ba1e511SMatthew Knepley /* Logging support */ 95014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 96ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 975358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 988ba1e511SMatthew Knepley 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *); 10019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType); 101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec); 10675f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool); 1075ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat); 108bbbebc2cSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat); 1093e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 110edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) 111d71ae5a4SJacob Faibussowitsch { 1129371c9d4SSatish Balay return KSPSetMatSolveBatchSize(ksp, n); 1139371c9d4SSatish Balay } 1143e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *); 115edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) 116d71ae5a4SJacob Faibussowitsch { 1179371c9d4SSatish Balay return KSPGetMatSolveBatchSize(ksp, n); 1189371c9d4SSatish Balay } 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 120ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *); 12223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool); 1237d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *); 12425c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool); 125c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec); 1262eac72dbSBarry Smith 127140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 128ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 129798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 130798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 131798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 132bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP)); 133*4d4d2bdcSBarry Smith 134*4d4d2bdcSBarry Smith /*S 135*4d4d2bdcSBarry Smith KSPMonitorRegisterFn - A function prototype for functions provided to `KSPMonitorRegister()` 136*4d4d2bdcSBarry Smith 137*4d4d2bdcSBarry Smith Calling Sequence: 138*4d4d2bdcSBarry Smith + ksp - iterative solver obtained from `KSPCreate()` 139*4d4d2bdcSBarry Smith . it - iteration number 140*4d4d2bdcSBarry Smith . rnorm - (estimated) 2-norm of (preconditioned) residual 141*4d4d2bdcSBarry Smith - ctx - `PetscViewerAndFormat` object 142*4d4d2bdcSBarry Smith 143*4d4d2bdcSBarry Smith Level: beginner 144*4d4d2bdcSBarry Smith 145*4d4d2bdcSBarry Smith Note: 146*4d4d2bdcSBarry Smith This is a `KSPMonitorFn` specialized for a context of `PetscViewerAndFormat` 147*4d4d2bdcSBarry Smith 148*4d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn`, `KSPMonitorRegisterDestroyFn` 149*4d4d2bdcSBarry Smith S*/ 150*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterFn(KSP ksp, PetscInt it, PetscReal rnorm, PetscViewerAndFormat *ctx); 151*4d4d2bdcSBarry Smith 152*4d4d2bdcSBarry Smith /*S 153*4d4d2bdcSBarry Smith KSPMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to `KSPMonitorRegister()` 154*4d4d2bdcSBarry Smith 155*4d4d2bdcSBarry Smith Calling Sequence: 156*4d4d2bdcSBarry Smith + viewer - the viewer to be used with the `KSPMonitorRegisterFn` 157*4d4d2bdcSBarry Smith . format - the format of the viewer 158*4d4d2bdcSBarry Smith . ctx - a context for the monitor 159*4d4d2bdcSBarry Smith - result - a `PetscViewerAndFormat` object 160*4d4d2bdcSBarry Smith 161*4d4d2bdcSBarry Smith Level: beginner 162*4d4d2bdcSBarry Smith 163*4d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterDestroyFn` 164*4d4d2bdcSBarry Smith S*/ 165*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterCreateFn(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **result); 166*4d4d2bdcSBarry Smith 167*4d4d2bdcSBarry Smith /*S 168*4d4d2bdcSBarry Smith KSPMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to `KSPMonitorRegister()` 169*4d4d2bdcSBarry Smith 170*4d4d2bdcSBarry Smith Calling Sequence: 171*4d4d2bdcSBarry Smith . vf - a `PetscViewerAndFormat` object to be destroyed, including any context 172*4d4d2bdcSBarry Smith 173*4d4d2bdcSBarry Smith Level: beginner 174*4d4d2bdcSBarry Smith 175*4d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn` 176*4d4d2bdcSBarry Smith S*/ 177*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterDestroyFn(PetscViewerAndFormat **result); 178*4d4d2bdcSBarry Smith 179*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, KSPMonitorRegisterFn *, KSPMonitorRegisterCreateFn *, KSPMonitorRegisterDestroyFn *); 18030de9b25SBarry Smith 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt); 184c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *); 18525c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt); 18625c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 1927d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 193c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 195c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 200734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 2012a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 202edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 203d71ae5a4SJacob Faibussowitsch { 2049371c9d4SSatish Balay return KSPCreateVecs(ksp, n, x, m, y); 2059371c9d4SSatish Balay } 2062eac72dbSBarry Smith 207*4d4d2bdcSBarry Smith /*S 208*4d4d2bdcSBarry Smith KSPPSolveFn - A function prototype for functions provided to `KSPSetPreSolve()` and `KSPSetPostSolve()` 209*4d4d2bdcSBarry Smith 210*4d4d2bdcSBarry Smith Calling Sequence: 211*4d4d2bdcSBarry Smith + ksp - the `KSP` context 212*4d4d2bdcSBarry Smith . rhs - the right-hand side vector 213*4d4d2bdcSBarry Smith . x - the solution vector 214*4d4d2bdcSBarry Smith - ctx - optional context that was provided with `KSPSetPreSolve()` or `KSPSetPostSolve()` 215*4d4d2bdcSBarry Smith 216*4d4d2bdcSBarry Smith Level: intermediate 217*4d4d2bdcSBarry Smith 218*4d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`, `PCShellPSolveFn` 219*4d4d2bdcSBarry Smith S*/ 220*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPPSolveFn(KSP ksp, Vec rhs, Vec x, void *ctx); 221*4d4d2bdcSBarry Smith 222*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, KSPPSolveFn *, void *); 223*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, KSPPSolveFn *, void *); 224d4211eb9SBarry Smith 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 2273821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt); 2283821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *); 229aabeff55SBarry Smith 230*4d4d2bdcSBarry Smith /*S 231*4d4d2bdcSBarry Smith KSPMonitorFn - A function prototype for functions provided to `KSPMonitorSet()` 232*4d4d2bdcSBarry Smith 233*4d4d2bdcSBarry Smith Calling Sequence: 234*4d4d2bdcSBarry Smith + ksp - iterative solver obtained from `KSPCreate()` 235*4d4d2bdcSBarry Smith . it - iteration number 236*4d4d2bdcSBarry Smith . rnorm - (estimated) 2-norm of (preconditioned) residual 237*4d4d2bdcSBarry Smith - ctx - optional monitoring context, as provided with `KSPMonitorSet()` 238*4d4d2bdcSBarry Smith 239*4d4d2bdcSBarry Smith Level: beginner 240*4d4d2bdcSBarry Smith 241*4d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()` 242*4d4d2bdcSBarry Smith S*/ 243*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorFn(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx); 244*4d4d2bdcSBarry Smith 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 246*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, KSPMonitorFn *, void *, PetscCtxDestroyFn *); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 2483ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *); 24994a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 2506497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool); 25194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 2526497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool); 2534b0e389bSBarry Smith 254fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 255fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 256fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 257fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 258fa0ddf94SBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 260cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 264fb80e629SPablo Brubeck PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 266285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 267b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 268b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 269b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 270b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 2724a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 273*4d4d2bdcSBarry Smith 274*4d4d2bdcSBarry Smith /*S 275*4d4d2bdcSBarry Smith PCMGCoarseSpaceConstructorFn - A function prototype for functions registered with `PCMGRegisterCoarseSpaceConstructor()` 276*4d4d2bdcSBarry Smith 277*4d4d2bdcSBarry Smith Calling Sequence: 278*4d4d2bdcSBarry Smith + pc - The `PC` object 279*4d4d2bdcSBarry Smith . l - The multigrid level, 0 is the coarse level 280*4d4d2bdcSBarry Smith . dm - The `DM` for this level 281*4d4d2bdcSBarry Smith . smooth - The level smoother 282*4d4d2bdcSBarry Smith . Nc - The size of the coarse space 283*4d4d2bdcSBarry Smith . initGuess - Basis for an initial guess for the space 284*4d4d2bdcSBarry Smith - coarseSp - A basis for the computed coarse space 285*4d4d2bdcSBarry Smith 286*4d4d2bdcSBarry Smith Level: beginner 287*4d4d2bdcSBarry Smith 288*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PCMGRegisterCoarseSpaceConstructor()`, `PCMGGetCoarseSpaceConstructor()` 289*4d4d2bdcSBarry Smith S*/ 290*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCMGCoarseSpaceConstructorFn(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp); 291*4d4d2bdcSBarry Smith 292f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 293*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn *); 294*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn **); 295f3b08a26SMatthew G. Knepley 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 2982eac72dbSBarry Smith 299f2edd1f0SMalachi Phillips /*E 30095bd0b28SBarry Smith KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV` 301f2edd1f0SMalachi Phillips 302f2edd1f0SMalachi Phillips Values: 303f2edd1f0SMalachi Phillips + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial 304f2edd1f0SMalachi Phillips . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial 305f2edd1f0SMalachi Phillips - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial 306f2edd1f0SMalachi Phillips 307f2edd1f0SMalachi Phillips Level: intermediate 308f2edd1f0SMalachi Phillips 3091cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()` 310f2edd1f0SMalachi Phillips E*/ 311f2edd1f0SMalachi Phillips typedef enum { 312f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FIRST, 313f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FOURTH, 314f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_OPT_FOURTH 315f2edd1f0SMalachi Phillips } KSPChebyshevKind; 316f2edd1f0SMalachi Phillips 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 32058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 321b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 322f2edd1f0SMalachi Phillips PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind); 32316a05f60SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *); 32458450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 326d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 327d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 3287d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 3294b0e389bSBarry Smith 330640f4f53SPatrick Sanan /*E 331640f4f53SPatrick Sanan 33206137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 333640f4f53SPatrick Sanan 33467b8a455SSatish Balay Values: 335a1cb98faSBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions 336a1cb98faSBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 337640f4f53SPatrick Sanan 3382b26979fSBarry Smith Level: intermediate 339640f4f53SPatrick Sanan 3401cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 341640f4f53SPatrick Sanan E*/ 3429371c9d4SSatish Balay typedef enum { 3439371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_STANDARD, 3449371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_NOTAY 3459371c9d4SSatish Balay } KSPFCDTruncationType; 34606137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 347640f4f53SPatrick Sanan 348640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 349640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 350640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 351640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 35206137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 35306137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 354640f4f53SPatrick Sanan 355390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 356390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 357390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 358390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 359390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 360390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 361390d8e47SPatrick Sanan 362fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 363fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 364fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 365fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 366fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 367fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 368fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 369fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 370*4d4d2bdcSBarry Smith 371*4d4d2bdcSBarry Smith /*S 372*4d4d2bdcSBarry Smith KSPFlexibleModifyPCFn - A prototype of a function used to modify the preconditioner during the use of flexible `KSP` methods, such as `KSPFGMRES` 373*4d4d2bdcSBarry Smith 374*4d4d2bdcSBarry Smith Calling Sequence: 375*4d4d2bdcSBarry Smith + ksp - the `KSP` context being used. 376*4d4d2bdcSBarry Smith . total_its - the total number of iterations that have occurred. 377*4d4d2bdcSBarry Smith . local_its - the number of iterations since last restart if applicable 378*4d4d2bdcSBarry Smith . res_norm - the current residual norm 379*4d4d2bdcSBarry Smith - ctx - optional context variable set with `KSPFlexibleSetModifyPC()`, `KSPPIPEGCRSetModifyPC()`, `KSPGCRSetModifyPC()`, `KSPFGMRESSetModifyPC()` 380*4d4d2bdcSBarry Smith 381*4d4d2bdcSBarry Smith Level: beginner 382*4d4d2bdcSBarry Smith 383*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFlexibleSetModifyPC()`, `KSPPIPEGCRSetModifyPC()`, `KSPGCRSetModifyPC()`, `KSPFGMRESSetModifyPC()` 384*4d4d2bdcSBarry Smith S*/ 385*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPFlexibleModifyPCFn(KSP ksp, PetscInt total_its, PetscInt local_its, PetscReal res_norm, void *ctx); 386*4d4d2bdcSBarry Smith 387*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPFlexibleSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *); 388*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *); 38983f0b094SBarry Smith 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 391014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 393e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 3949f236934SBarry Smith 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 4001d73ed98SBarry Smith 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 4031d73ed98SBarry Smith 404483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 405483d6965SPatrick Sanan 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 408*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *); 40958ad63f4SBarry Smith 410f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 411f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 412f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 413f87a0b54SStefano Zampini 41404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 41504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 41604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 417568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 418e82af88dSprj- 4196cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 4206cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 4216cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 422edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 423d71ae5a4SJacob Faibussowitsch { 4249371c9d4SSatish Balay return KSPHPDDMSetDeflationMat(ksp, U); 4259371c9d4SSatish Balay } 426edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 427d71ae5a4SJacob Faibussowitsch { 4289371c9d4SSatish Balay return KSPHPDDMGetDeflationMat(ksp, U); 4299371c9d4SSatish Balay } 4306cac28cbSPierre Jolivet #endif 431edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 432d71ae5a4SJacob Faibussowitsch { 4339371c9d4SSatish Balay return KSPMatSolve(ksp, B, X); 4349371c9d4SSatish Balay } 435d55ab951SPierre Jolivet /*E 43687497f52SBarry Smith KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 437d55ab951SPierre Jolivet 438d55ab951SPierre Jolivet Values: 439a4d1885cSBarry Smith + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 440a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 441a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 442a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BCG` - block CG 443a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 444a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 445a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 446a4d1885cSBarry Smith - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 447d55ab951SPierre Jolivet 44816a05f60SBarry Smith Level: intermediate 44916a05f60SBarry Smith 4501cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 451d55ab951SPierre Jolivet E*/ 452d55ab951SPierre Jolivet typedef enum { 453d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 454d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 455d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 456d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 457d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 458d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 459d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 460d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 461d55ab951SPierre Jolivet } KSPHPDDMType; 462d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 463a4d1885cSBarry Smith 4642dd49c90SPierre Jolivet /*E 46587497f52SBarry Smith KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 4662dd49c90SPierre Jolivet 4672dd49c90SPierre Jolivet Values: 468a4d1885cSBarry Smith + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16` 469a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single` 470a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double` 471a4d1885cSBarry Smith - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128` 4722dd49c90SPierre Jolivet 47316a05f60SBarry Smith Level: intermediate 47416a05f60SBarry Smith 4751cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPHPDDM` 4762dd49c90SPierre Jolivet E*/ 4772dd49c90SPierre Jolivet typedef enum { 4782dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_HALF = 0, 4792dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_SINGLE = 1, 4802dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_DOUBLE = 2, 4812dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_QUADRUPLE = 3 4822dd49c90SPierre Jolivet } KSPHPDDMPrecision; 483d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 484d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 485e82af88dSprj- 486b49fd9e1SBarry Smith /*E 48716a05f60SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers 488b49fd9e1SBarry Smith 489a4d1885cSBarry Smith Values: 490a4d1885cSBarry Smith + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 491a4d1885cSBarry Smith . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 492a4d1885cSBarry Smith - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 493b49fd9e1SBarry Smith 49416a05f60SBarry Smith Level: advanced 49516a05f60SBarry Smith 49695bd0b28SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 497a4d1885cSBarry Smith `KSPGMRESGetOrthogonalization()`, 498a4d1885cSBarry Smith `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 499b49fd9e1SBarry Smith E*/ 5009371c9d4SSatish Balay typedef enum { 5019371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_NEVER, 5029371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_IFNEEDED, 5039371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_ALWAYS 5049371c9d4SSatish Balay } KSPGMRESCGSRefinementType; 5056a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 50616a05f60SBarry Smith 5071f7e983dSSatish Balay /*MC 5081957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 5098c5b8ba0SBarry Smith 5108c5b8ba0SBarry Smith Level: advanced 5118c5b8ba0SBarry Smith 51287497f52SBarry Smith Note: 51387497f52SBarry Smith Possibly unstable, but the fastest to compute 5148c5b8ba0SBarry Smith 51595bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 516a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 517db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 518db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 5198c5b8ba0SBarry Smith M*/ 5208c5b8ba0SBarry Smith 5211f7e983dSSatish Balay /*MC 5228c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 5238c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 5248c5b8ba0SBarry Smith poor orthogonality. 5258c5b8ba0SBarry Smith 5268c5b8ba0SBarry Smith Level: advanced 5278c5b8ba0SBarry Smith 528a4d1885cSBarry Smith Note: 529a4d1885cSBarry Smith This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 5308c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 5318c5b8ba0SBarry Smith 53295bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 533a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 534db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 535db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 5368c5b8ba0SBarry Smith M*/ 5378c5b8ba0SBarry Smith 5381f7e983dSSatish Balay /*MC 53916a05f60SBarry Smith KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process. 5408c5b8ba0SBarry Smith 5418c5b8ba0SBarry Smith Level: advanced 5428c5b8ba0SBarry Smith 54387497f52SBarry Smith Notes: 54487497f52SBarry Smith This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 54587497f52SBarry Smith but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 5468c5b8ba0SBarry Smith 5478c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 5488c5b8ba0SBarry Smith 54995bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 550a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 551db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 552db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 5538c5b8ba0SBarry Smith M*/ 5548c5b8ba0SBarry Smith 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 55708480c60SBarry Smith 558*4d4d2bdcSBarry Smith PETSC_EXTERN KSPFlexibleModifyPCFn KSPFGMRESModifyPCNoChange; 559*4d4d2bdcSBarry Smith PETSC_EXTERN KSPFlexibleModifyPCFn KSPFGMRESModifyPCKSP; 560*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *); 561c38d4ed2SBarry Smith 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 565121fd945SKris Buschelman 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 569e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 570d9492815SBarry Smith 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 57287d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 5732eac72dbSBarry Smith 574798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 575*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidual; 576*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDraw; 577*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDrawLG; 578798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 579*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualShort; 580*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualRange; 581*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidual; 582*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDraw; 583*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDrawLG; 584798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 585*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualMax; 586*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorError; 587*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDraw; 588*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDrawLG; 589798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 590*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolution; 591*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDraw; 592*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDrawLG; 593798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 594*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSingularValue; 595798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 596edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 597d71ae5a4SJacob Faibussowitsch { 5989371c9d4SSatish Balay return KSPMonitorResidual(ksp, n, rnorm, vf); 5999371c9d4SSatish Balay } 600edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 601d71ae5a4SJacob Faibussowitsch { 6029371c9d4SSatish Balay return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 6039371c9d4SSatish Balay } 604edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 605d71ae5a4SJacob Faibussowitsch { 6069371c9d4SSatish Balay return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 6079371c9d4SSatish Balay } 608798534f6SMatthew G. Knepley 609798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 610341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 611341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 612341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 613341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 614e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 615e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 616e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 61784cb2905SBarry Smith 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 620c01c455dSBarry Smith 62123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 62223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 626014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 6271eb62cbbSBarry Smith 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 630014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 6321f7f0c4fSBarry Smith 633*4d4d2bdcSBarry Smith /*S 634*4d4d2bdcSBarry Smith KSPConvergedReasonViewFn - A prototype of a function used with `KSPConvergedReasonViewSet()` 635*4d4d2bdcSBarry Smith 636*4d4d2bdcSBarry Smith Calling Sequence: 637*4d4d2bdcSBarry Smith + ksp - the `KSP` object whose `KSPConvergedReason` is to be viewed 638*4d4d2bdcSBarry Smith - ctx - context used by the function, set with `KSPConvergedReasonViewSet()` 639*4d4d2bdcSBarry Smith 640*4d4d2bdcSBarry Smith Level: beginner 641*4d4d2bdcSBarry Smith 642*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`, `KSPConvergedReasonViewFromOptions()`, `KSPView()` 643*4d4d2bdcSBarry Smith S*/ 644*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergedReasonViewFn(KSP ksp, void *ctx); 645*4d4d2bdcSBarry Smith 646014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 64755849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 648fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 64919a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 650*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, KSPConvergedReasonViewFn *, void *, PetscCtxDestroyFn *); 6511b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 652c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 65394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 6541b2b9847SBarry Smith 655edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 656d71ae5a4SJacob Faibussowitsch { 6579371c9d4SSatish Balay return KSPConvergedReasonView(ksp, v); 6589371c9d4SSatish Balay } 659edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 660d71ae5a4SJacob Faibussowitsch { 6619371c9d4SSatish Balay return KSPConvergedReasonViewFromOptions(ksp); 6629371c9d4SSatish Balay } 66355849f57SBarry Smith 66455849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 6651eb62cbbSBarry Smith 666dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 6670e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 668014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 669884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 670*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidual; 671*4d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidualDrawLG; 672798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 673db9b2ab1SHong Zhang 674014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 675014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 67668ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 6779f0612e4SBarry Smith PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *); 67883ab6a24SBarry Smith 67928ce4d24SBarry Smith /*E 680a4d1885cSBarry Smith KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 6818a4b9c5cSBarry Smith test routines. 6828a4b9c5cSBarry Smith 683a4d1885cSBarry Smith Values: 684a4d1885cSBarry Smith + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 685a4d1885cSBarry Smith . `KSP_NORM_NONE` - use no norm calculation 686a4d1885cSBarry Smith . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 687a4d1885cSBarry Smith . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 688a4d1885cSBarry Smith - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 689a4d1885cSBarry Smith 6908a4b9c5cSBarry Smith Level: advanced 6918a4b9c5cSBarry Smith 692a4d1885cSBarry Smith Note: 693a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 69416a05f60SBarry Smith depending on whether left or right preconditioning is used, see `KSPSetPCSide()` 695a3f661c8SBarry Smith 6961cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 69795bd0b28SBarry Smith `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 6988a4b9c5cSBarry Smith E*/ 6999371c9d4SSatish Balay typedef enum { 7009371c9d4SSatish Balay KSP_NORM_DEFAULT = -1, 7019371c9d4SSatish Balay KSP_NORM_NONE = 0, 7029371c9d4SSatish Balay KSP_NORM_PRECONDITIONED = 1, 7039371c9d4SSatish Balay KSP_NORM_UNPRECONDITIONED = 2, 7049371c9d4SSatish Balay KSP_NORM_NATURAL = 3 7059371c9d4SSatish Balay } KSPNormType; 7069e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 707014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes; 708a21b2a99SBarry Smith 7091f7e983dSSatish Balay /*MC 7109793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 7118c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 7128c5b8ba0SBarry Smith be based on a norm of a residual etc. 7138c5b8ba0SBarry Smith 7148c5b8ba0SBarry Smith Level: advanced 7158c5b8ba0SBarry Smith 716a4d1885cSBarry Smith Note: 717a4d1885cSBarry Smith Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 7188c5b8ba0SBarry Smith 7191cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 7208c5b8ba0SBarry Smith M*/ 7218c5b8ba0SBarry Smith 7221f7e983dSSatish Balay /*MC 7231957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 7248c5b8ba0SBarry Smith convergence test routine. 7258c5b8ba0SBarry Smith 7268c5b8ba0SBarry Smith Level: advanced 7278c5b8ba0SBarry Smith 7281cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 7298c5b8ba0SBarry Smith M*/ 7308c5b8ba0SBarry Smith 7311f7e983dSSatish Balay /*MC 732ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 7338c5b8ba0SBarry Smith convergence test routine. 7348c5b8ba0SBarry Smith 7358c5b8ba0SBarry Smith Level: advanced 7368c5b8ba0SBarry Smith 7371cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 7388c5b8ba0SBarry Smith M*/ 7398c5b8ba0SBarry Smith 7401f7e983dSSatish Balay /*MC 741ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 74287497f52SBarry Smith convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 7438c5b8ba0SBarry Smith 7448c5b8ba0SBarry Smith Level: advanced 7458c5b8ba0SBarry Smith 7461cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 7478c5b8ba0SBarry Smith M*/ 7488c5b8ba0SBarry Smith 749014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 750014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 751ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt); 752014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 753014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 7548a4b9c5cSBarry Smith 755edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", ) 756edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", ) 757b5c1ed96SAlexander #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", ) 758b5c1ed96SAlexander #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", ) 7598a4b9c5cSBarry Smith /*E 760a4d1885cSBarry Smith KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 76128ce4d24SBarry Smith 762a4d1885cSBarry Smith Values: 76378daaec8SBarry Smith + `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR` 76478daaec8SBarry Smith . `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR` 765a4d1885cSBarry Smith . `KSP_CONVERGED_RTOL` - requested decrease in the residual 766a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 767a4d1885cSBarry Smith . `KSP_CONVERGED_ITS` - requested number of iterations 7684a221d59SStefano Zampini . `KSP_CONVERGED_NEG_CURVE` - see note below 769a4d1885cSBarry Smith . `KSP_CONVERGED_STEP_LENGTH` - see note below 7706e7835fbSStefano Zampini . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred). 77178daaec8SBarry Smith . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within `KSPGMRES` 772a4d1885cSBarry Smith . `KSP_DIVERGED_ITS` - requested number of iterations 77378daaec8SBarry Smith . `KSP_DIVERGED_DTOL` - large increase in the residual norm indicating the solution is diverging 774a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 7758f47816eSPierre Jolivet . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBCGS` Krylov method 776a4d1885cSBarry Smith . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 77778daaec8SBarry Smith . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG` 778a4d1885cSBarry Smith . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 77978daaec8SBarry Smith . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG` 780a4d1885cSBarry Smith - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 781a4d1885cSBarry Smith 78216a05f60SBarry Smith Level: beginner 78316a05f60SBarry Smith 784a4d1885cSBarry Smith Note: 7856e7835fbSStefano Zampini The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by 7866e7835fbSStefano Zampini the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 78728ce4d24SBarry Smith 78895bd0b28SBarry Smith Developer Note: 78987497f52SBarry Smith The string versions of these are `KSPConvergedReasons`; if you change 7904d0a8057SBarry Smith any of the values here also change them that array of names. 79186c02ca4SBarry Smith 7921cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 79328ce4d24SBarry Smith E*/ 794d15094e1SBarry Smith typedef enum { /* converged */ 79578daaec8SBarry Smith KSP_CONVERGED_RTOL_NORMAL_DEPRECATED = 1, 79678daaec8SBarry Smith KSP_CONVERGED_RTOL_NORMAL_EQUATIONS = 1, 79778daaec8SBarry Smith KSP_CONVERGED_ATOL_NORMAL_DEPRECATED = 9, 79878daaec8SBarry Smith KSP_CONVERGED_ATOL_NORMAL_EQUATIONS = 9, 799d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 800d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 801b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 8024a221d59SStefano Zampini KSP_CONVERGED_NEG_CURVE = 5, 8034a221d59SStefano Zampini KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 8044a221d59SStefano Zampini KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 8054a221d59SStefano Zampini KSP_CONVERGED_STEP_LENGTH = 6, 8064a221d59SStefano Zampini KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 807d15094e1SBarry Smith /* diverged */ 808b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 809d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 810d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 811d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 812b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 813b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 814b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 8154d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 8166aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 817c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 818aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 819d15094e1SBarry Smith 8209371c9d4SSatish Balay KSP_CONVERGED_ITERATING = 0 8219371c9d4SSatish Balay } KSPConvergedReason; 822014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons; 823d15094e1SBarry Smith 824c838673bSBarry Smith /*MC 825af27ebaaSBarry Smith KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called 826c838673bSBarry Smith 827c838673bSBarry Smith Level: beginner 828c838673bSBarry Smith 82995bd0b28SBarry Smith Notes: 83087497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 831c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 832c838673bSBarry Smith 2-norm of the residual for right preconditioning 833c838673bSBarry Smith 83487497f52SBarry Smith See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 835f9fed41fSBarry Smith 8361cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 837c838673bSBarry Smith M*/ 838c838673bSBarry Smith 839c838673bSBarry Smith /*MC 840af27ebaaSBarry Smith KSP_CONVERGED_ATOL - $||r|| \le atol$ 841c838673bSBarry Smith 842c838673bSBarry Smith Level: beginner 843c838673bSBarry Smith 84495bd0b28SBarry Smith Notes: 84587497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 846c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 847c838673bSBarry Smith 2-norm of the residual for right preconditioning 848c838673bSBarry Smith 84987497f52SBarry Smith See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 850c838673bSBarry Smith 8511cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 852c838673bSBarry Smith M*/ 853c838673bSBarry Smith 854c838673bSBarry Smith /*MC 855af27ebaaSBarry Smith KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$ 856c838673bSBarry Smith 857c838673bSBarry Smith Level: beginner 858c838673bSBarry Smith 85995bd0b28SBarry Smith Note: 86087497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 861c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 862c838673bSBarry Smith 2-norm of the residual for right preconditioning 863c838673bSBarry Smith 8640241b274SPierre Jolivet .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 865c838673bSBarry Smith M*/ 866c838673bSBarry Smith 867c838673bSBarry Smith /*MC 868c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 869c838673bSBarry Smith reached 870c838673bSBarry Smith 871c838673bSBarry Smith Level: beginner 872c838673bSBarry Smith 8731cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 874c838673bSBarry Smith M*/ 875c838673bSBarry Smith 876c838673bSBarry Smith /*MC 87787497f52SBarry Smith KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 87887497f52SBarry Smith the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 879af27ebaaSBarry Smith test routine is set in `KSP`. 880c838673bSBarry Smith 881c838673bSBarry Smith Level: beginner 882c838673bSBarry Smith 8831cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 884c838673bSBarry Smith M*/ 885c838673bSBarry Smith 886c838673bSBarry Smith /*MC 887c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 8881de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 88987497f52SBarry Smith preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 89046091a0eSPierre Jolivet are collinear. 891c838673bSBarry Smith 892c838673bSBarry Smith Level: beginner 893c838673bSBarry Smith 8941cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 895c838673bSBarry Smith M*/ 896c838673bSBarry Smith 897c838673bSBarry Smith /*MC 89887497f52SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 899c838673bSBarry Smith method could not continue to enlarge the Krylov space. 900c838673bSBarry Smith 901c838673bSBarry Smith Level: beginner 902c838673bSBarry Smith 9031cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 904c838673bSBarry Smith M*/ 905c838673bSBarry Smith 906c838673bSBarry Smith /*MC 907c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 90887497f52SBarry Smith symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 909c838673bSBarry Smith 910c838673bSBarry Smith Level: beginner 911c838673bSBarry Smith 9121cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 913c838673bSBarry Smith M*/ 914c838673bSBarry Smith 915c838673bSBarry Smith /*MC 916c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 91787497f52SBarry Smith positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 9180b4b7b1cSBarry Smith be symmetric positive definite (SPD). 919c838673bSBarry Smith 920c838673bSBarry Smith Level: beginner 921c838673bSBarry Smith 92287497f52SBarry Smith Note: 923a4d1885cSBarry Smith This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force 92487497f52SBarry Smith the `PCICC` preconditioner to generate a positive definite preconditioner 925c838673bSBarry Smith 9261cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 927c838673bSBarry Smith M*/ 928c838673bSBarry Smith 929c838673bSBarry Smith /*MC 930c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 9319fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 93287497f52SBarry Smith such as `PCFIELDSPLIT`. 9339fc87aa7SBarry Smith 9349fc87aa7SBarry Smith Level: beginner 9359fc87aa7SBarry Smith 936a4d1885cSBarry Smith Note: 937a4d1885cSBarry Smith Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details. 9389fc87aa7SBarry Smith 9391cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 9409fc87aa7SBarry Smith M*/ 9419fc87aa7SBarry Smith 9429fc87aa7SBarry Smith /*MC 943a4d1885cSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called 944a4d1885cSBarry Smith while `KSPSolve()` is still running. 945c838673bSBarry Smith 946c838673bSBarry Smith Level: beginner 947c838673bSBarry Smith 9481cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 949c838673bSBarry Smith M*/ 950c838673bSBarry Smith 951*4d4d2bdcSBarry Smith /*S 952*4d4d2bdcSBarry Smith KSPConvergenceTestFn - A prototype of a function used with `KSPSetConvergenceTest()` 953*4d4d2bdcSBarry Smith 954*4d4d2bdcSBarry Smith Calling Sequence: 955*4d4d2bdcSBarry Smith + ksp - iterative solver obtained from `KSPCreate()` 956*4d4d2bdcSBarry Smith . it - iteration number 957*4d4d2bdcSBarry Smith . rnorm - (estimated) 2-norm of (preconditioned) residual 958*4d4d2bdcSBarry Smith . reason - the reason why it has converged or diverged 959*4d4d2bdcSBarry Smith - ctx - optional convergence context, as set by `KSPSetConvergenceTest()` 960*4d4d2bdcSBarry Smith 961*4d4d2bdcSBarry Smith Level: beginner 962*4d4d2bdcSBarry Smith 963*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()` 964*4d4d2bdcSBarry Smith S*/ 965*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergenceTestFn(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx); 966*4d4d2bdcSBarry Smith 967*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, KSPConvergenceTestFn *, void *, PetscCtxDestroyFn *); 968*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, KSPConvergenceTestFn **, void **, PetscCtxDestroyFn **); 969*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, KSPConvergenceTestFn **, void **, PetscCtxDestroyFn **); 9703ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 971*4d4d2bdcSBarry Smith PETSC_EXTERN KSPConvergenceTestFn KSPConvergedDefault; 972*4d4d2bdcSBarry Smith PETSC_EXTERN KSPConvergenceTestFn KSPLSQRConvergedDefault; 973*4d4d2bdcSBarry Smith PETSC_EXTERN PetscCtxDestroyFn KSPConvergedDefaultDestroy; 9748de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 9758de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 9768de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 97754b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 9780059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 979014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 980*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char *[]); 98194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 9822a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool); 9832a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *); 984abef13c0SSatish Balay 985edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void) 986d71ae5a4SJacob Faibussowitsch { /* never called */ 9879371c9d4SSatish Balay } 9888ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 989edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void) 990d71ae5a4SJacob Faibussowitsch { /* never called */ 9919371c9d4SSatish Balay } 9928ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 993edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void) 994d71ae5a4SJacob Faibussowitsch { /* never called */ 9959371c9d4SSatish Balay } 9968ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 997edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void) 998d71ae5a4SJacob Faibussowitsch { /* never called */ 9999371c9d4SSatish Balay } 10008ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 1001edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void) 1002d71ae5a4SJacob Faibussowitsch { /* never called */ 10039371c9d4SSatish Balay } 10048ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 1005edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void) 1006d71ae5a4SJacob Faibussowitsch { /* never called */ 10079371c9d4SSatish Balay } 10088ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 10098ea1b3e6SJed Brown 10100bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 1011edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 1012d71ae5a4SJacob Faibussowitsch { 1013f22e26b7SPierre Jolivet return KSPComputeOperator(A, PETSC_NULLPTR, B); 10149371c9d4SSatish Balay } 1015d4fbbf0eSBarry Smith 101628ce4d24SBarry Smith /*E 1017a4d1885cSBarry Smith KSPCGType - Determines what type of `KSPCG` to use 101828ce4d24SBarry Smith 1019a4d1885cSBarry Smith Values: 1020a4d1885cSBarry Smith + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric 1021a4d1885cSBarry Smith - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian 1022a4d1885cSBarry Smith 102316a05f60SBarry Smith Level: beginner 102416a05f60SBarry Smith 10251cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()` 102628ce4d24SBarry Smith E*/ 10279371c9d4SSatish Balay typedef enum { 10289371c9d4SSatish Balay KSP_CG_SYMMETRIC = 0, 10299371c9d4SSatish Balay KSP_CG_HERMITIAN = 1 10309371c9d4SSatish Balay } KSPCGType; 10316a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 103228ce4d24SBarry Smith 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 10358031f4adStmunson 1036dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 1037fb01098fSStefano Zampini PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal); 1038dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 1039dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 1040fcae7a14Stmunson 104105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 104205de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 1043edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 1044d71ae5a4SJacob Faibussowitsch { 10459371c9d4SSatish Balay return KSPGLTRGetMinEig(ksp, x); 10469371c9d4SSatish Balay } 1047edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 1048d71ae5a4SJacob Faibussowitsch { 10499371c9d4SSatish Balay return KSPGLTRGetLambda(ksp, x); 10509371c9d4SSatish Balay } 10518031f4adStmunson 1052014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 1053ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 10541d6018f0SLisandro Dalcin 1055f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 1056*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCSetPostSolve(PC, PetscErrorCode (*)(PC, KSP)); 1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 10593369ce9aSBarry Smith 1060014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 10612f2e5d10SKris Buschelman 1062*4d4d2bdcSBarry Smith /*S 1063*4d4d2bdcSBarry Smith PCShellPSolveFn - A function prototype for functions provided to `PCShellSetPreSolve()` and `PCShellSetPostSolve()` 1064*4d4d2bdcSBarry Smith 1065*4d4d2bdcSBarry Smith Calling Sequence: 1066*4d4d2bdcSBarry Smith + pc - the preconditioner `PC` context 1067*4d4d2bdcSBarry Smith . ksp - the `KSP` context 1068*4d4d2bdcSBarry Smith . xin - input vector 1069*4d4d2bdcSBarry Smith - xout - output vector 1070*4d4d2bdcSBarry Smith 1071*4d4d2bdcSBarry Smith Level: intermediate 1072*4d4d2bdcSBarry Smith 1073*4d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSPPSolveFn`, `KSP`, `PCShellSetPreSolve()`, `PCShellSetPostSolve()` 1074*4d4d2bdcSBarry Smith S*/ 1075*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCShellPSolveFn(PC pc, KSP ksp, Vec xim, Vec xout); 1076*4d4d2bdcSBarry Smith 1077*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PCShellPSolveFn *); 1078*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PCShellPSolveFn *); 107903919abeSBarry Smith 1080ba36735cSStefano Zampini /*S 1081a4d1885cSBarry Smith KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods. 1082f8a50e2bSBarry Smith 1083a4d1885cSBarry Smith Level: intermediate 1084f8a50e2bSBarry Smith 10850b4b7b1cSBarry Smith Note: 10860b4b7b1cSBarry Smith These methods generate initial guesses based on a series of previous, related, linear solves. For example, 10870b4b7b1cSBarry Smith in implicit time-stepping with `TS`. 10880b4b7b1cSBarry Smith 10899168452cSPierre Jolivet .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType` 1090ba36735cSStefano Zampini S*/ 1091ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess; 109216a05f60SBarry Smith 1093ba36735cSStefano Zampini /*J 109487497f52SBarry Smith KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 1095ba36735cSStefano Zampini 1096a4d1885cSBarry Smith Values: 1097a4d1885cSBarry Smith + `KSPGUESSFISCHER` - methodology developed by Paul Fischer 10980b4b7b1cSBarry Smith - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition (POD) 1099a4d1885cSBarry Smith 110016a05f60SBarry Smith Level: intermediate 110116a05f60SBarry Smith 11021cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGuess` 1103ba36735cSStefano Zampini J*/ 1104ba36735cSStefano Zampini typedef const char *KSPGuessType; 1105ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 1106ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 1107a4d1885cSBarry Smith 11081d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 1109ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 1110ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 1111ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 1112ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 1113ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 1114ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 1115ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 11168410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 1117ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 1118ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 1119ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 1120ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 1121ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 1122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 1123ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 1124ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 1125f8a50e2bSBarry Smith 1126470b340bSDmitry Karpeev /*E 11277addb90fSBarry Smith MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines 1128470b340bSDmitry Karpeev 1129470b340bSDmitry Karpeev Level: intermediate 1130470b340bSDmitry Karpeev 1131a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, 11320b4b7b1cSBarry Smith `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()` 1133470b340bSDmitry Karpeev E*/ 11349371c9d4SSatish Balay typedef enum { 11359371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_DIAG, 11369371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_LUMP, 11379371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 11389371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_FULL 11399371c9d4SSatish Balay } MatSchurComplementAinvType; 1140470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 1141470b340bSDmitry Karpeev 1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 1144d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 1145bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 1146aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 1147bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 1148470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 1149470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 11505bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 11515a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 1152470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 1153470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 11543f22127dSBarry Smith 115578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 115678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1157bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1158bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 1159bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *); 116078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 1161864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1162864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1163864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1164864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1165864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1166cd929ea3SAlp Dener 1167cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1168b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1169cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1170cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1171e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 11720ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1173cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1174cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1175cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1176cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1177cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 11782d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1179cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 11801ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *); 1181cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1182cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1183cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 118492f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1185bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *); 1186cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1187cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1188864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1189e96a0efeSStefano Zampini 1190bbb72809SHansol Suh /*E 11911ca5963aSToby Isaac MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix 11921ca5963aSToby Isaac 11931ca5963aSToby Isaac Values: 11941ca5963aSToby Isaac + `MAT_LMVM_MULT_RECURSIVE` - Use recursive formulas for products and solves 11951ca5963aSToby Isaac . `MAT_LMVM_MULT_DENSE` - Use dense formulas for products and solves when possible 11961ca5963aSToby Isaac - `MAT_LMVM_MULT_COMPACT_DENSE` - The same as `MATLMVM_MULT_DENSE`, but go further and ensure products and solves are computed in compact low-rank update form 11971ca5963aSToby Isaac 11981ca5963aSToby Isaac Level: advanced 11991ca5963aSToby Isaac 12001ca5963aSToby Isaac Options Database Keys: 12011ca5963aSToby Isaac . -mat_lmvm_mult_algorithm - the algorithm to use for multiplication (recursive, dense, compact_dense) 12021ca5963aSToby Isaac 12031ca5963aSToby Isaac .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()` 12041ca5963aSToby Isaac E*/ 12051ca5963aSToby Isaac typedef enum { 12061ca5963aSToby Isaac MAT_LMVM_MULT_RECURSIVE, 12071ca5963aSToby Isaac MAT_LMVM_MULT_DENSE, 12081ca5963aSToby Isaac MAT_LMVM_MULT_COMPACT_DENSE, 12091ca5963aSToby Isaac } MatLMVMMultAlgorithm; 12101ca5963aSToby Isaac 12111ca5963aSToby Isaac PETSC_EXTERN const char *const MatLMVMMultAlgorithms[]; 12121ca5963aSToby Isaac 12131ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm); 12141ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *); 12151ca5963aSToby Isaac 12161ca5963aSToby Isaac /*E 12171ca5963aSToby Isaac MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix. 1218bbb72809SHansol Suh 1219bbb72809SHansol Suh Values: 1220b6982945SToby Isaac + `MAT_LMVM_SYMBROYDEN_SCALE_NONE` - no rescaling 1221b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR` - scalar rescaling 1222b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling 1223b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_USER` - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE` 1224b6982945SToby Isaac - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE` - let PETSc decide rescaling 1225bbb72809SHansol Suh 1226bbb72809SHansol Suh Level: intermediate 1227bbb72809SHansol Suh 1228bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()` 1229bbb72809SHansol Suh E*/ 1230e96a0efeSStefano Zampini typedef enum { 1231e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 1232e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 1233e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 1234b6982945SToby Isaac MAT_LMVM_SYMBROYDEN_SCALE_USER = 3, 1235b6982945SToby Isaac MAT_LMVM_SYMBROYDEN_SCALE_DECIDE = 4 1236e96a0efeSStefano Zampini } MatLMVMSymBroydenScaleType; 1237e96a0efeSStefano Zampini PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 1238e96a0efeSStefano Zampini 1239864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 12401ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *); 12411ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal); 12421ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *); 12431ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal); 1244cd929ea3SAlp Dener 1245bbb72809SHansol Suh /*E 12460b4b7b1cSBarry Smith MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`. 1247bbb72809SHansol Suh 1248bbb72809SHansol Suh Values: 1249bbb72809SHansol Suh + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch 1250bbb72809SHansol Suh - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement 1251bbb72809SHansol Suh 1252bbb72809SHansol Suh Level: intermediate 1253bbb72809SHansol Suh 1254bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()` 1255bbb72809SHansol Suh E*/ 1256bbb72809SHansol Suh typedef enum { 1257bbb72809SHansol Suh MAT_LMVM_DENSE_REORDER, 1258bbb72809SHansol Suh MAT_LMVM_DENSE_INPLACE 1259bbb72809SHansol Suh } MatLMVMDenseType; 1260bbb72809SHansol Suh PETSC_EXTERN const char *const MatLMVMDenseTypes[]; 1261bbb72809SHansol Suh 1262bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType); 1263bbb72809SHansol Suh 1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 1266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 1267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 1268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 12697b578ef6SBarry Smith 12707b578ef6SBarry Smith /*S 12718434afd1SBarry Smith KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()` 12727b578ef6SBarry Smith 12737b578ef6SBarry Smith Calling Sequence: 12747b578ef6SBarry Smith + ksp - `ksp` context 12757b578ef6SBarry Smith . b - output vector 12767b578ef6SBarry Smith - ctx - [optional] user-defined function context 12777b578ef6SBarry Smith 12787b578ef6SBarry Smith Level: beginner 12797b578ef6SBarry Smith 1280*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn` 12817b578ef6SBarry Smith S*/ 1282*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeRHSFn(KSP ksp, Vec b, void *ctx); 12837b578ef6SBarry Smith 12848434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *); 12857b578ef6SBarry Smith 12867b578ef6SBarry Smith /*S 12878434afd1SBarry Smith KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()` 12887b578ef6SBarry Smith 12897b578ef6SBarry Smith Calling Sequence: 12907b578ef6SBarry Smith + ksp - `KSP` context 12917b578ef6SBarry Smith . A - the operator that defines the linear system 12927b578ef6SBarry Smith . P - an operator from which to build the preconditioner (often the same as `A`) 12937b578ef6SBarry Smith - ctx - [optional] user-defined function context 12947b578ef6SBarry Smith 12957b578ef6SBarry Smith Level: beginner 12967b578ef6SBarry Smith 1297*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn` 12987b578ef6SBarry Smith S*/ 1299*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeOperatorsFn(KSP ksp, Mat A, Mat P, void *ctx); 13007b578ef6SBarry Smith 13018434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *); 13027b578ef6SBarry Smith 13037b578ef6SBarry Smith /*S 13048434afd1SBarry Smith KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()` 13057b578ef6SBarry Smith 13067b578ef6SBarry Smith Calling Sequence: 13077b578ef6SBarry Smith + ksp - `ksp` context 13087b578ef6SBarry Smith . x - output vector 13097b578ef6SBarry Smith - ctx - [optional] user-defined function context 13107b578ef6SBarry Smith 13117b578ef6SBarry Smith Level: beginner 13127b578ef6SBarry Smith 1313*4d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn` 13147b578ef6SBarry Smith S*/ 1315*4d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeInitialGuessFn(KSP ksp, Vec x, void *ctx); 13167b578ef6SBarry Smith 13178434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *); 13188434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *); 13198434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *); 13208434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *); 13218434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *); 13228434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *); 13238434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *); 13246c699258SBarry Smith 132502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 1326953494dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode); 1327557cf195SMatthew G. Knepley 13282b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 13292b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 13304bf303faSJacob Faibussowitsch 13314bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP); 13324bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *); 13335d83a8b1SBarry Smith 13345d83a8b1SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM); 1335*4d4d2bdcSBarry Smith 1336*4d4d2bdcSBarry Smith #include <petscdstypes.h> 1337*4d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, PetscPointFn **, InsertMode, Vec); 1338