xref: /petsc/include/petscksp.h (revision 7addb90f52a7936ba144cdab1bb2cc37152af090)
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 
271cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
2828ce4d24SBarry Smith S*/
29e2a1c21fSSatish Balay typedef struct _p_KSP *KSP;
302eac72dbSBarry Smith 
3176bdecfbSBarry Smith /*J
320b4b7b1cSBarry Smith    KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides.
3328ce4d24SBarry Smith 
3428ce4d24SBarry Smith    Level: beginner
3528ce4d24SBarry Smith 
361cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
3776bdecfbSBarry Smith J*/
3819fd82e9SBarry Smith typedef const char *KSPType;
3982bf6240SBarry Smith #define KSPRICHARDSON "richardson"
406c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
4182bf6240SBarry Smith #define KSPCG         "cg"
422c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
432c8fc646SJed Brown #define KSPPIPECG     "pipecg"
44901ccb91SSiegfried Cools #define KSPPIPECGRR   "pipecgrr"
45265663fdSSiegfried Cools #define KSPPIPELCG    "pipelcg"
46b21a8899STyler Chen #define KSPPIPEPRCG   "pipeprcg"
47325e8391SManasi #define KSPPIPECG2    "pipecg2"
48df2a969aSvictorle #define KSPCGNE       "cgne"
4905de396fSBarry Smith #define KSPNASH       "nash"
5005de396fSBarry Smith #define KSPSTCG       "stcg"
5105de396fSBarry Smith #define KSPGLTR       "gltr"
52edd03b47SJacob Faibussowitsch #define KSPCGNASH     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash"
53edd03b47SJacob Faibussowitsch #define KSPCGSTCG     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg"
54edd03b47SJacob Faibussowitsch #define KSPCGGLTR     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr"
55640f4f53SPatrick Sanan #define KSPFCG        "fcg"
56390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
5782bf6240SBarry Smith #define KSPGMRES      "gmres"
58483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
599dcbbd2bSBarry Smith #define KSPFGMRES     "fgmres"
609dcbbd2bSBarry Smith #define KSPLGMRES     "lgmres"
614f8e6cd9SSatish Balay #define KSPDGMRES     "dgmres"
6261b468f9SJed Brown #define KSPPGMRES     "pgmres"
6382bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
6482bf6240SBarry Smith #define KSPBCGS       "bcgs"
65e1c61ce8SBarry Smith #define KSPIBCGS      "ibcgs"
66345ecf0bSXiangmin Jiao #define KSPQMRCGS     "qmrcgs"
6718ac38e6SHong Zhang #define KSPFBCGS      "fbcgs"
68c2b71004SHong Zhang #define KSPFBCGSR     "fbcgsr"
69f0037002Svictorle #define KSPBCGSL      "bcgsl"
70f154af2dSSiegfried Cools #define KSPPIPEBCGS   "pipebcgs"
7182bf6240SBarry Smith #define KSPCGS        "cgs"
7282bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
7382bf6240SBarry Smith #define KSPCR         "cr"
742c8fc646SJed Brown #define KSPPIPECR     "pipecr"
7582bf6240SBarry Smith #define KSPLSQR       "lsqr"
7682bf6240SBarry Smith #define KSPPREONLY    "preonly"
773c2be86cSBarry Smith #define KSPNONE       "none"
7882bf6240SBarry Smith #define KSPQCG        "qcg"
79c9cf9b11SBarry Smith #define KSPBICG       "bicg"
80b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
8101c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
8221356919SSatish Balay #define KSPLCD        "lcd"
831d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
8458ad63f4SBarry Smith #define KSPGCR        "gcr"
85fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
86e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
87e4d80e07Szianekhodja #define KSPCGLS       "cgls"
88329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
8938cfc46eSPierre Jolivet #define KSPHPDDM      "hpddm"
902eac72dbSBarry Smith 
918ba1e511SMatthew Knepley /* Logging support */
92014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
93ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
945358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
958ba1e511SMatthew Knepley 
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
9719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
98c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
10375f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
1045ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
105bbbebc2cSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat);
1063e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
107edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
108d71ae5a4SJacob Faibussowitsch {
1099371c9d4SSatish Balay   return KSPSetMatSolveBatchSize(ksp, n);
1109371c9d4SSatish Balay }
1113e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
112edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
113d71ae5a4SJacob Faibussowitsch {
1149371c9d4SSatish Balay   return KSPGetMatSolveBatchSize(ksp, n);
1159371c9d4SSatish Balay }
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
117ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
11923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
1207d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
12125c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
122c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
1232eac72dbSBarry Smith 
124140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
125ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
126798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList;
127798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
128798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
129bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
130798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
13130de9b25SBarry Smith 
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
135c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
13625c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt);
13725c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
1437d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
144c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
146c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
151734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
1522a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
153edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
154d71ae5a4SJacob Faibussowitsch {
1559371c9d4SSatish Balay   return KSPCreateVecs(ksp, n, x, m, y);
1569371c9d4SSatish Balay }
1572eac72dbSBarry Smith 
158d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
159d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
160d4211eb9SBarry Smith 
161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
1633821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt);
1643821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *);
165aabeff55SBarry Smith 
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
16749abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscCtxDestroyFn *);
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
1693ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
17094a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
1716497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool);
17294a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
1736497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool);
1744b0e389bSBarry Smith 
175fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
176fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
177fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
178fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
179fa0ddf94SBarry Smith 
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
181cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
185fb80e629SPablo Brubeck PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
187285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
188b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
189b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
190b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
191b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
1934a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
194f3b08a26SMatthew G. Knepley /*
195f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
196f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
197f3b08a26SMatthew G. Knepley */
198f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
1992b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
2002b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
201f3b08a26SMatthew G. Knepley 
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
2042eac72dbSBarry Smith 
205f2edd1f0SMalachi Phillips /*E
20695bd0b28SBarry Smith   KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV`
207f2edd1f0SMalachi Phillips 
208f2edd1f0SMalachi Phillips   Values:
209f2edd1f0SMalachi Phillips + `KSP_CHEBYSHEV_FIRST`      - "classic" first-kind Chebyshev polynomial
210f2edd1f0SMalachi Phillips . `KSP_CHEBYSHEV_FOURTH`     - fourth-kind Chebyshev polynomial
211f2edd1f0SMalachi Phillips - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial
212f2edd1f0SMalachi Phillips 
213f2edd1f0SMalachi Phillips   Level: intermediate
214f2edd1f0SMalachi Phillips 
2151cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()`
216f2edd1f0SMalachi Phillips E*/
217f2edd1f0SMalachi Phillips typedef enum {
218f2edd1f0SMalachi Phillips   KSP_CHEBYSHEV_FIRST,
219f2edd1f0SMalachi Phillips   KSP_CHEBYSHEV_FOURTH,
220f2edd1f0SMalachi Phillips   KSP_CHEBYSHEV_OPT_FOURTH
221f2edd1f0SMalachi Phillips } KSPChebyshevKind;
222f2edd1f0SMalachi Phillips 
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
22658450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
227b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
228f2edd1f0SMalachi Phillips PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind);
22916a05f60SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *);
23058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
232d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
233d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
2347d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
2354b0e389bSBarry Smith 
236640f4f53SPatrick Sanan /*E
237640f4f53SPatrick Sanan 
23806137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
239640f4f53SPatrick Sanan 
24067b8a455SSatish Balay   Values:
241a1cb98faSBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
242a1cb98faSBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY`    - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
243640f4f53SPatrick Sanan 
2442b26979fSBarry Smith    Level: intermediate
245640f4f53SPatrick Sanan 
2461cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
247640f4f53SPatrick Sanan E*/
2489371c9d4SSatish Balay typedef enum {
2499371c9d4SSatish Balay   KSP_FCD_TRUNC_TYPE_STANDARD,
2509371c9d4SSatish Balay   KSP_FCD_TRUNC_TYPE_NOTAY
2519371c9d4SSatish Balay } KSPFCDTruncationType;
25206137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
253640f4f53SPatrick Sanan 
254640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
255640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
256640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
257640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
25806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
25906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
260640f4f53SPatrick Sanan 
261390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
262390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
263390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
264390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
265390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
266390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
267390d8e47SPatrick Sanan 
268fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
269fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
270fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
271fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
272fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
273fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
274fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
275fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
2764bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
27783f0b094SBarry Smith 
278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
281e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
2829f236934SBarry Smith 
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
2881d73ed98SBarry Smith 
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2911d73ed98SBarry Smith 
292483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
293483d6965SPatrick Sanan 
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
29758ad63f4SBarry Smith 
298f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
299f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
300f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);
301f87a0b54SStefano Zampini 
30204d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
30304d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
30404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
305568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
306e82af88dSprj- 
3076cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
3086cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
3096cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM)
310edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
311d71ae5a4SJacob Faibussowitsch {
3129371c9d4SSatish Balay   return KSPHPDDMSetDeflationMat(ksp, U);
3139371c9d4SSatish Balay }
314edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
315d71ae5a4SJacob Faibussowitsch {
3169371c9d4SSatish Balay   return KSPHPDDMGetDeflationMat(ksp, U);
3179371c9d4SSatish Balay }
3186cac28cbSPierre Jolivet #endif
319edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
320d71ae5a4SJacob Faibussowitsch {
3219371c9d4SSatish Balay   return KSPMatSolve(ksp, B, X);
3229371c9d4SSatish Balay }
323d55ab951SPierre Jolivet /*E
32487497f52SBarry Smith     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
325d55ab951SPierre Jolivet 
326d55ab951SPierre Jolivet     Values:
327a4d1885cSBarry Smith +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
328a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BGMRES`          - block GMRES
329a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_CG`              - Conjugate Gradient
330a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BCG`             - block CG
331a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_GCRODR`          - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
332a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BGCRODR`         - block GCRODR
333a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BFBCG`           - breakdown-free BCG
334a4d1885cSBarry Smith -   `KSP_HPDDM_TYPE_PREONLY`         - apply the preconditioner only
335d55ab951SPierre Jolivet 
33616a05f60SBarry Smith     Level: intermediate
33716a05f60SBarry Smith 
3381cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
339d55ab951SPierre Jolivet E*/
340d55ab951SPierre Jolivet typedef enum {
341d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES   = 0,
342d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES  = 1,
343d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG      = 2,
344d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG     = 3,
345d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR  = 4,
346d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
347d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG   = 6,
348d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
349d55ab951SPierre Jolivet } KSPHPDDMType;
350d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
351a4d1885cSBarry Smith 
3522dd49c90SPierre Jolivet /*E
35387497f52SBarry Smith     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
3542dd49c90SPierre Jolivet 
3552dd49c90SPierre Jolivet     Values:
356a4d1885cSBarry Smith +   `KSP_HPDDM_PRECISION_HALF`      - default when PETSc is configured `--with-precision=__fp16`
357a4d1885cSBarry Smith .   `KSP_HPDDM_PRECISION_SINGLE`    - default when PETSc is configured `--with-precision=single`
358a4d1885cSBarry Smith .   `KSP_HPDDM_PRECISION_DOUBLE`    - default when PETSc is configured `--with-precision=double`
359a4d1885cSBarry Smith -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
3602dd49c90SPierre Jolivet 
36116a05f60SBarry Smith     Level: intermediate
36216a05f60SBarry Smith 
3631cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPHPDDM`
3642dd49c90SPierre Jolivet E*/
3652dd49c90SPierre Jolivet typedef enum {
3662dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_HALF      = 0,
3672dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_SINGLE    = 1,
3682dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_DOUBLE    = 2,
3692dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_QUADRUPLE = 3
3702dd49c90SPierre Jolivet } KSPHPDDMPrecision;
371d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
372d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
373e82af88dSprj- 
374b49fd9e1SBarry Smith /*E
37516a05f60SBarry Smith    KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers
376b49fd9e1SBarry Smith 
377a4d1885cSBarry Smith    Values:
378a4d1885cSBarry Smith +  `KSP_GMRES_CGS_REFINE_NEVER`    - one step of classical Gram-Schmidt
379a4d1885cSBarry Smith .  `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
380a4d1885cSBarry Smith -  `KSP_GMRES_CGS_REFINE_ALWAYS`   - always perform two steps
381b49fd9e1SBarry Smith 
38216a05f60SBarry Smith    Level: advanced
38316a05f60SBarry Smith 
38495bd0b28SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
385a4d1885cSBarry Smith           `KSPGMRESGetOrthogonalization()`,
386a4d1885cSBarry Smith           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
387b49fd9e1SBarry Smith E*/
3889371c9d4SSatish Balay typedef enum {
3899371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_NEVER,
3909371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_IFNEEDED,
3919371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_ALWAYS
3929371c9d4SSatish Balay } KSPGMRESCGSRefinementType;
3936a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
39416a05f60SBarry Smith 
3951f7e983dSSatish Balay /*MC
3961957e957SBarry Smith    KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
3978c5b8ba0SBarry Smith 
3988c5b8ba0SBarry Smith    Level: advanced
3998c5b8ba0SBarry Smith 
40087497f52SBarry Smith    Note:
40187497f52SBarry Smith    Possibly unstable, but the fastest to compute
4028c5b8ba0SBarry Smith 
40395bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
404a4d1885cSBarry Smith           `KSP`, `KSPGMRESGetOrthogonalization()`,
405db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
406db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
4078c5b8ba0SBarry Smith M*/
4088c5b8ba0SBarry Smith 
4091f7e983dSSatish Balay /*MC
4108c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
4118c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
4128c5b8ba0SBarry Smith           poor orthogonality.
4138c5b8ba0SBarry Smith 
4148c5b8ba0SBarry Smith    Level: advanced
4158c5b8ba0SBarry Smith 
416a4d1885cSBarry Smith    Note:
417a4d1885cSBarry Smith    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
4188c5b8ba0SBarry Smith    estimate the orthogonality but is more stable.
4198c5b8ba0SBarry Smith 
42095bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
421a4d1885cSBarry Smith           `KSP`, `KSPGMRESGetOrthogonalization()`,
422db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
423db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
4248c5b8ba0SBarry Smith M*/
4258c5b8ba0SBarry Smith 
4261f7e983dSSatish Balay /*MC
42716a05f60SBarry Smith    KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process.
4288c5b8ba0SBarry Smith 
4298c5b8ba0SBarry Smith    Level: advanced
4308c5b8ba0SBarry Smith 
43187497f52SBarry Smith    Notes:
43287497f52SBarry Smith    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
43387497f52SBarry Smith    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
4348c5b8ba0SBarry Smith 
4358c5b8ba0SBarry Smith    You should only use this if you absolutely know that the iterative refinement is needed.
4368c5b8ba0SBarry Smith 
43795bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
438a4d1885cSBarry Smith           `KSP`, `KSPGMRESGetOrthogonalization()`,
439db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
440db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
4418c5b8ba0SBarry Smith M*/
4428c5b8ba0SBarry Smith 
443014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
444014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
44508480c60SBarry Smith 
446014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
448014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
449c38d4ed2SBarry Smith 
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
453121fd945SKris Buschelman 
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
457e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
458d9492815SBarry Smith 
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
46087d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
4612eac72dbSBarry Smith 
462798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
463798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
464798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
465798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
466798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
467798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
468798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
469798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
470798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
471798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
472798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
473798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
474798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
475798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
476798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
477798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
478798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
479798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
480798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
481798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
482fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
483798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
484edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
485d71ae5a4SJacob Faibussowitsch {
4869371c9d4SSatish Balay   return KSPMonitorResidual(ksp, n, rnorm, vf);
4879371c9d4SSatish Balay }
488edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
489d71ae5a4SJacob Faibussowitsch {
4909371c9d4SSatish Balay   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
4919371c9d4SSatish Balay }
492edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
493d71ae5a4SJacob Faibussowitsch {
4949371c9d4SSatish Balay   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
4959371c9d4SSatish Balay }
496798534f6SMatthew G. Knepley 
497798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
498341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
499341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
500341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
501341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
502e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
503e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
504e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
50584cb2905SBarry Smith 
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
508c01c455dSBarry Smith 
50923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
51023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
5151eb62cbbSBarry Smith 
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
5201f7f0c4fSBarry Smith 
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
52255849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
523fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
52419a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
52549abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *, PetscCtxDestroyFn *);
5261b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
527c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
52894a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
5291b2b9847SBarry Smith 
530edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
531d71ae5a4SJacob Faibussowitsch {
5329371c9d4SSatish Balay   return KSPConvergedReasonView(ksp, v);
5339371c9d4SSatish Balay }
534edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
535d71ae5a4SJacob Faibussowitsch {
5369371c9d4SSatish Balay   return KSPConvergedReasonViewFromOptions(ksp);
5379371c9d4SSatish Balay }
53855849f57SBarry Smith 
53955849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
5401eb62cbbSBarry Smith 
541dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
5420e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
544884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
545798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
546798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
547798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
548db9b2ab1SHong Zhang 
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
55168ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
5529f0612e4SBarry Smith PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);
55383ab6a24SBarry Smith 
55428ce4d24SBarry Smith /*E
555a4d1885cSBarry Smith    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
5568a4b9c5cSBarry Smith        test routines.
5578a4b9c5cSBarry Smith 
558a4d1885cSBarry Smith    Values:
559a4d1885cSBarry Smith +  `KSP_NORM_DEFAULT`          - use the default for the current `KSPType`
560a4d1885cSBarry Smith .  `KSP_NORM_NONE`             - use no norm calculation
561a4d1885cSBarry Smith .  `KSP_NORM_PRECONDITIONED`   - use the preconditioned residual norm
562a4d1885cSBarry Smith .  `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
563a4d1885cSBarry Smith -  `KSP_NORM_NATURAL`          - use the natural norm (the norm induced by the linear operator)
564a4d1885cSBarry Smith 
5658a4b9c5cSBarry Smith    Level: advanced
5668a4b9c5cSBarry Smith 
567a4d1885cSBarry Smith    Note:
568a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
56916a05f60SBarry Smith    depending on whether left or right preconditioning is used, see `KSPSetPCSide()`
570a3f661c8SBarry Smith 
5711cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
57295bd0b28SBarry Smith           `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
5738a4b9c5cSBarry Smith E*/
5749371c9d4SSatish Balay typedef enum {
5759371c9d4SSatish Balay   KSP_NORM_DEFAULT          = -1,
5769371c9d4SSatish Balay   KSP_NORM_NONE             = 0,
5779371c9d4SSatish Balay   KSP_NORM_PRECONDITIONED   = 1,
5789371c9d4SSatish Balay   KSP_NORM_UNPRECONDITIONED = 2,
5799371c9d4SSatish Balay   KSP_NORM_NATURAL          = 3
5809371c9d4SSatish Balay } KSPNormType;
5819e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
582014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes;
583a21b2a99SBarry Smith 
5841f7e983dSSatish Balay /*MC
5859793b452SBarry Smith    KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
5868c5b8ba0SBarry Smith    possibly save some computation but means the convergence test cannot
5878c5b8ba0SBarry Smith    be based on a norm of a residual etc.
5888c5b8ba0SBarry Smith 
5898c5b8ba0SBarry Smith    Level: advanced
5908c5b8ba0SBarry Smith 
591a4d1885cSBarry Smith    Note:
592a4d1885cSBarry Smith    Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
5938c5b8ba0SBarry Smith 
5941cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
5958c5b8ba0SBarry Smith M*/
5968c5b8ba0SBarry Smith 
5971f7e983dSSatish Balay /*MC
5981957e957SBarry Smith    KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
5998c5b8ba0SBarry Smith    convergence test routine.
6008c5b8ba0SBarry Smith 
6018c5b8ba0SBarry Smith    Level: advanced
6028c5b8ba0SBarry Smith 
6031cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
6048c5b8ba0SBarry Smith M*/
6058c5b8ba0SBarry Smith 
6061f7e983dSSatish Balay /*MC
607ce9499c7SBarry Smith    KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
6088c5b8ba0SBarry Smith    convergence test routine.
6098c5b8ba0SBarry Smith 
6108c5b8ba0SBarry Smith    Level: advanced
6118c5b8ba0SBarry Smith 
6121cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
6138c5b8ba0SBarry Smith M*/
6148c5b8ba0SBarry Smith 
6151f7e983dSSatish Balay /*MC
616ce9499c7SBarry Smith    KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
61787497f52SBarry Smith    convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
6188c5b8ba0SBarry Smith 
6198c5b8ba0SBarry Smith    Level: advanced
6208c5b8ba0SBarry Smith 
6211cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
6228c5b8ba0SBarry Smith M*/
6238c5b8ba0SBarry Smith 
624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
626ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt);
627014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
6298a4b9c5cSBarry Smith 
630edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
631edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
632edd03b47SJacob Faibussowitsch #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM(3, 11, 0, "KSP_DIVERGED_PC_FAILED", )
6338a4b9c5cSBarry Smith /*E
634a4d1885cSBarry Smith    KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
63528ce4d24SBarry Smith 
636a4d1885cSBarry Smith    Values:
637a4d1885cSBarry Smith +  `KSP_CONVERGED_RTOL_NORMAL`     - requested decrease in the residual for the normal equations
638a4d1885cSBarry Smith .  `KSP_CONVERGED_ATOL_NORMAL`     - requested absolute value in the residual for the normal equations
639a4d1885cSBarry Smith .  `KSP_CONVERGED_RTOL`            - requested decrease in the residual
640a4d1885cSBarry Smith .  `KSP_CONVERGED_ATOL`            - requested absolute value in the residual
641a4d1885cSBarry Smith .  `KSP_CONVERGED_ITS`             - requested number of iterations
6424a221d59SStefano Zampini .  `KSP_CONVERGED_NEG_CURVE`       - see note below
643a4d1885cSBarry Smith .  `KSP_CONVERGED_STEP_LENGTH`     - see note below
6446e7835fbSStefano Zampini .  `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred).
645a4d1885cSBarry Smith .  `KSP_DIVERGED_NULL`             - breakdown when solving the Hessenberg system within GMRES
646a4d1885cSBarry Smith .  `KSP_DIVERGED_ITS`              - requested number of iterations
647a4d1885cSBarry Smith .  `KSP_DIVERGED_DTOL`             - large increase in the residual norm
648a4d1885cSBarry Smith .  `KSP_DIVERGED_BREAKDOWN`        - breakdown in the Krylov method
6498f47816eSPierre Jolivet .  `KSP_DIVERGED_BREAKDOWN_BICG`   - breakdown in the `KSPBCGS` Krylov method
650a4d1885cSBarry Smith .  `KSP_DIVERGED_NONSYMMETRIC`     - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
651a4d1885cSBarry Smith .  `KSP_DIVERGED_INDEFINITE_PC`    - the preconditioner was indefinite for a `KSPType` that requires it be definite
652a4d1885cSBarry Smith .  `KSP_DIVERGED_NANORINF`         - a not a number of infinity was detected in a vector during the computation
653a4d1885cSBarry Smith .  `KSP_DIVERGED_INDEFINITE_MAT`   - the operator was indefinite for a `KSPType` that requires it be definite
654a4d1885cSBarry Smith -  `KSP_DIVERGED_PC_FAILED`        - the action of the preconditioner failed for some reason
655a4d1885cSBarry Smith 
65616a05f60SBarry Smith    Level: beginner
65716a05f60SBarry Smith 
658a4d1885cSBarry Smith    Note:
6596e7835fbSStefano Zampini    The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
6606e7835fbSStefano Zampini    the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
66128ce4d24SBarry Smith 
66295bd0b28SBarry Smith    Developer Note:
66387497f52SBarry Smith    The string versions of these are `KSPConvergedReasons`; if you change
6644d0a8057SBarry Smith    any of the values here also change them that array of names.
66586c02ca4SBarry Smith 
6661cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
66728ce4d24SBarry Smith E*/
668d15094e1SBarry Smith typedef enum { /* converged */
6699ecb1286SBarry Smith   KSP_CONVERGED_RTOL_NORMAL               = 1,
6709ecb1286SBarry Smith   KSP_CONVERGED_ATOL_NORMAL               = 9,
671d15094e1SBarry Smith   KSP_CONVERGED_RTOL                      = 2,
672d15094e1SBarry Smith   KSP_CONVERGED_ATOL                      = 3,
673b335793eSSatish Balay   KSP_CONVERGED_ITS                       = 4,
6744a221d59SStefano Zampini   KSP_CONVERGED_NEG_CURVE                 = 5,
6754a221d59SStefano Zampini   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
6764a221d59SStefano Zampini   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
6774a221d59SStefano Zampini   KSP_CONVERGED_STEP_LENGTH               = 6,
6784a221d59SStefano Zampini   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
679d15094e1SBarry Smith   /* diverged */
680b3cc6726SBarry Smith   KSP_DIVERGED_NULL                      = -2,
681d15094e1SBarry Smith   KSP_DIVERGED_ITS                       = -3,
682d15094e1SBarry Smith   KSP_DIVERGED_DTOL                      = -4,
683d15094e1SBarry Smith   KSP_DIVERGED_BREAKDOWN                 = -5,
684b4ac9ba4SBarry Smith   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
685b4ac9ba4SBarry Smith   KSP_DIVERGED_NONSYMMETRIC              = -7,
686b4ac9ba4SBarry Smith   KSP_DIVERGED_INDEFINITE_PC             = -8,
6874d51c080SBarry Smith   KSP_DIVERGED_NANORINF                  = -9,
6886aee1118SBarry Smith   KSP_DIVERGED_INDEFINITE_MAT            = -10,
689c0decd05SBarry Smith   KSP_DIVERGED_PC_FAILED                 = -11,
690aa4d2078SSatish Balay   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
691d15094e1SBarry Smith 
6929371c9d4SSatish Balay   KSP_CONVERGED_ITERATING = 0
6939371c9d4SSatish Balay } KSPConvergedReason;
694014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons;
695d15094e1SBarry Smith 
696c838673bSBarry Smith /*MC
697af27ebaaSBarry Smith    KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called
698c838673bSBarry Smith 
699c838673bSBarry Smith    Level: beginner
700c838673bSBarry Smith 
70195bd0b28SBarry Smith    Notes:
70287497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
703c838673bSBarry Smith    for left preconditioning it is the 2-norm of the preconditioned residual, and the
704c838673bSBarry Smith    2-norm of the residual for right preconditioning
705c838673bSBarry Smith 
70687497f52SBarry Smith    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
707f9fed41fSBarry Smith 
7081cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
709c838673bSBarry Smith M*/
710c838673bSBarry Smith 
711c838673bSBarry Smith /*MC
712af27ebaaSBarry Smith    KSP_CONVERGED_ATOL - $||r|| \le atol$
713c838673bSBarry Smith 
714c838673bSBarry Smith    Level: beginner
715c838673bSBarry Smith 
71695bd0b28SBarry Smith    Notes:
71787497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
718c838673bSBarry Smith    for left preconditioning it is the 2-norm of the preconditioned residual, and the
719c838673bSBarry Smith    2-norm of the residual for right preconditioning
720c838673bSBarry Smith 
72187497f52SBarry Smith    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
722c838673bSBarry Smith 
7231cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
724c838673bSBarry Smith M*/
725c838673bSBarry Smith 
726c838673bSBarry Smith /*MC
727af27ebaaSBarry Smith    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$
728c838673bSBarry Smith 
729c838673bSBarry Smith    Level: beginner
730c838673bSBarry Smith 
73195bd0b28SBarry Smith    Note:
73287497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
733c838673bSBarry Smith    for left preconditioning it is the 2-norm of the preconditioned residual, and the
734c838673bSBarry Smith    2-norm of the residual for right preconditioning
735c838673bSBarry Smith 
7360241b274SPierre Jolivet .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
737c838673bSBarry Smith M*/
738c838673bSBarry Smith 
739c838673bSBarry Smith /*MC
740c838673bSBarry Smith    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
741c838673bSBarry Smith    reached
742c838673bSBarry Smith 
743c838673bSBarry Smith    Level: beginner
744c838673bSBarry Smith 
7451cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
746c838673bSBarry Smith M*/
747c838673bSBarry Smith 
748c838673bSBarry Smith /*MC
74987497f52SBarry Smith    KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
75087497f52SBarry Smith    the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
751af27ebaaSBarry Smith    test routine is set in `KSP`.
752c838673bSBarry Smith 
753c838673bSBarry Smith    Level: beginner
754c838673bSBarry Smith 
7551cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
756c838673bSBarry Smith M*/
757c838673bSBarry Smith 
758c838673bSBarry Smith /*MC
759c838673bSBarry Smith    KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
7601de96524SPierre Jolivet    method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
76187497f52SBarry Smith    preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
76246091a0eSPierre Jolivet    are collinear.
763c838673bSBarry Smith 
764c838673bSBarry Smith    Level: beginner
765c838673bSBarry Smith 
7661cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
767c838673bSBarry Smith M*/
768c838673bSBarry Smith 
769c838673bSBarry Smith /*MC
77087497f52SBarry Smith    KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
771c838673bSBarry Smith    method could not continue to enlarge the Krylov space.
772c838673bSBarry Smith 
773c838673bSBarry Smith    Level: beginner
774c838673bSBarry Smith 
7751cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
776c838673bSBarry Smith M*/
777c838673bSBarry Smith 
778c838673bSBarry Smith /*MC
779c838673bSBarry Smith    KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
78087497f52SBarry Smith    symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
781c838673bSBarry Smith 
782c838673bSBarry Smith    Level: beginner
783c838673bSBarry Smith 
7841cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
785c838673bSBarry Smith M*/
786c838673bSBarry Smith 
787c838673bSBarry Smith /*MC
788c838673bSBarry Smith    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
78987497f52SBarry Smith    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
7900b4b7b1cSBarry Smith    be symmetric positive definite (SPD).
791c838673bSBarry Smith 
792c838673bSBarry Smith    Level: beginner
793c838673bSBarry Smith 
79487497f52SBarry Smith    Note:
795a4d1885cSBarry Smith    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
79687497f52SBarry Smith    the `PCICC` preconditioner to generate a positive definite preconditioner
797c838673bSBarry Smith 
7981cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
799c838673bSBarry Smith M*/
800c838673bSBarry Smith 
801c838673bSBarry Smith /*MC
802c0decd05SBarry Smith    KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
8039fc87aa7SBarry Smith    zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
80487497f52SBarry Smith    such as `PCFIELDSPLIT`.
8059fc87aa7SBarry Smith 
8069fc87aa7SBarry Smith    Level: beginner
8079fc87aa7SBarry Smith 
808a4d1885cSBarry Smith    Note:
809a4d1885cSBarry Smith    Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
8109fc87aa7SBarry Smith 
8111cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
8129fc87aa7SBarry Smith M*/
8139fc87aa7SBarry Smith 
8149fc87aa7SBarry Smith /*MC
815a4d1885cSBarry Smith    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
816a4d1885cSBarry Smith    while `KSPSolve()` is still running.
817c838673bSBarry Smith 
818c838673bSBarry Smith    Level: beginner
819c838673bSBarry Smith 
8201cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
821c838673bSBarry Smith M*/
822c838673bSBarry Smith 
823014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
824df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
825df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
8263ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
8278de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
8283eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
8298de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
8308de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
8318de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
8328de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
83354b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
8340059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
835014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
836c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
83794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
8382a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
8392a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
840abef13c0SSatish Balay 
841edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
842d71ae5a4SJacob Faibussowitsch { /* never called */
8439371c9d4SSatish Balay }
8448ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
845edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
846d71ae5a4SJacob Faibussowitsch { /* never called */
8479371c9d4SSatish Balay }
8488ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
849edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
850d71ae5a4SJacob Faibussowitsch { /* never called */
8519371c9d4SSatish Balay }
8528ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
853edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
854d71ae5a4SJacob Faibussowitsch { /* never called */
8559371c9d4SSatish Balay }
8568ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
857edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
858d71ae5a4SJacob Faibussowitsch { /* never called */
8599371c9d4SSatish Balay }
8608ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
861edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
862d71ae5a4SJacob Faibussowitsch { /* never called */
8639371c9d4SSatish Balay }
8648ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
8658ea1b3e6SJed Brown 
8660bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
867edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
868d71ae5a4SJacob Faibussowitsch {
869f22e26b7SPierre Jolivet   return KSPComputeOperator(A, PETSC_NULLPTR, B);
8709371c9d4SSatish Balay }
871d4fbbf0eSBarry Smith 
87228ce4d24SBarry Smith /*E
873a4d1885cSBarry Smith    KSPCGType - Determines what type of `KSPCG` to use
87428ce4d24SBarry Smith 
875a4d1885cSBarry Smith    Values:
876a4d1885cSBarry Smith  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
877a4d1885cSBarry Smith  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
878a4d1885cSBarry Smith 
87916a05f60SBarry Smith    Level: beginner
88016a05f60SBarry Smith 
8811cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
88228ce4d24SBarry Smith E*/
8839371c9d4SSatish Balay typedef enum {
8849371c9d4SSatish Balay   KSP_CG_SYMMETRIC = 0,
8859371c9d4SSatish Balay   KSP_CG_HERMITIAN = 1
8869371c9d4SSatish Balay } KSPCGType;
8876a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
88828ce4d24SBarry Smith 
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
890014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
8918031f4adStmunson 
892dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
893fb01098fSStefano Zampini PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
894dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
895dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
896fcae7a14Stmunson 
89705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
89805de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
899edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
900d71ae5a4SJacob Faibussowitsch {
9019371c9d4SSatish Balay   return KSPGLTRGetMinEig(ksp, x);
9029371c9d4SSatish Balay }
903edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
904d71ae5a4SJacob Faibussowitsch {
9059371c9d4SSatish Balay   return KSPGLTRGetLambda(ksp, x);
9069371c9d4SSatish Balay }
9078031f4adStmunson 
908014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
909ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
9101d6018f0SLisandro Dalcin 
911f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
912014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
913014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
9143369ce9aSBarry Smith 
915014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
9162f2e5d10SKris Buschelman 
917014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
918014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
91903919abeSBarry Smith 
920ba36735cSStefano Zampini /*S
921a4d1885cSBarry Smith    KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
922f8a50e2bSBarry Smith 
923a4d1885cSBarry Smith    Level: intermediate
924f8a50e2bSBarry Smith 
9250b4b7b1cSBarry Smith    Note:
9260b4b7b1cSBarry Smith    These methods generate initial guesses based on a series of previous, related, linear solves. For example,
9270b4b7b1cSBarry Smith    in implicit time-stepping with `TS`.
9280b4b7b1cSBarry Smith 
9299168452cSPierre Jolivet .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
930ba36735cSStefano Zampini S*/
931ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess;
93216a05f60SBarry Smith 
933ba36735cSStefano Zampini /*J
93487497f52SBarry Smith    KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
935ba36735cSStefano Zampini 
936a4d1885cSBarry Smith    Values:
937a4d1885cSBarry Smith  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
9380b4b7b1cSBarry Smith  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition (POD)
939a4d1885cSBarry Smith 
94016a05f60SBarry Smith    Level: intermediate
94116a05f60SBarry Smith 
9421cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGuess`
943ba36735cSStefano Zampini J*/
944ba36735cSStefano Zampini typedef const char *KSPGuessType;
945ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
946ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
947a4d1885cSBarry Smith 
9481d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
949ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
950ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
951ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
952ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
953ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
954ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
955ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
9568410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
957ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
958ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
959ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
960ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
961ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
962014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
963ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
964ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
965f8a50e2bSBarry Smith 
966470b340bSDmitry Karpeev /*E
967*7addb90fSBarry Smith     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines
968470b340bSDmitry Karpeev 
969470b340bSDmitry Karpeev     Level: intermediate
970470b340bSDmitry Karpeev 
971a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
9720b4b7b1cSBarry Smith           `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
973470b340bSDmitry Karpeev E*/
9749371c9d4SSatish Balay typedef enum {
9759371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
9769371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
9779371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
9789371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_FULL
9799371c9d4SSatish Balay } MatSchurComplementAinvType;
980470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
981470b340bSDmitry Karpeev 
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
984d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
985bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
986aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
987bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
988470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
989470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
9905bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
9915a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
992470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
993470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
9943f22127dSBarry Smith 
99578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
99678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
997bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
998bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
999bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
100078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1001864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1002864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1003864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1004864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1005864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1006cd929ea3SAlp Dener 
1007cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1008b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1009cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1010cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1011e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
10120ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1013cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1014cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1015cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1016cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1017cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
10182d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1019cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1020cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1021cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1022cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
102392f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1024bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1025cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1026cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1027864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1028e96a0efeSStefano Zampini 
1029bbb72809SHansol Suh /*E
1030bbb72809SHansol Suh   MatLMVMSymBroydenScaleType - Scaling type for symmetric Broyden.
1031bbb72809SHansol Suh 
1032bbb72809SHansol Suh   Values:
1033bbb72809SHansol Suh + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - No scaling
1034bbb72809SHansol Suh . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar scaling
1035bbb72809SHansol Suh . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal scaling
1036bbb72809SHansol Suh - `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - user-provided scale option
1037bbb72809SHansol Suh 
1038bbb72809SHansol Suh   Level: intermediate
1039bbb72809SHansol Suh 
1040bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1041bbb72809SHansol Suh E*/
1042e96a0efeSStefano Zampini typedef enum {
1043e96a0efeSStefano Zampini   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1044e96a0efeSStefano Zampini   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1045e96a0efeSStefano Zampini   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1046e96a0efeSStefano Zampini   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
1047e96a0efeSStefano Zampini } MatLMVMSymBroydenScaleType;
1048e96a0efeSStefano Zampini PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
1049e96a0efeSStefano Zampini 
1050864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1051cd929ea3SAlp Dener 
1052bbb72809SHansol Suh /*E
10530b4b7b1cSBarry Smith   MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.
1054bbb72809SHansol Suh 
1055bbb72809SHansol Suh   Values:
1056bbb72809SHansol Suh + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1057bbb72809SHansol Suh - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement
1058bbb72809SHansol Suh 
1059bbb72809SHansol Suh   Level: intermediate
1060bbb72809SHansol Suh 
1061bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1062bbb72809SHansol Suh E*/
1063bbb72809SHansol Suh typedef enum {
1064bbb72809SHansol Suh   MAT_LMVM_DENSE_REORDER,
1065bbb72809SHansol Suh   MAT_LMVM_DENSE_INPLACE
1066bbb72809SHansol Suh } MatLMVMDenseType;
1067bbb72809SHansol Suh PETSC_EXTERN const char *const MatLMVMDenseTypes[];
1068bbb72809SHansol Suh 
1069bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);
1070bbb72809SHansol Suh 
1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1073014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1074014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1075014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
10767b578ef6SBarry Smith 
10777b578ef6SBarry Smith /*S
10788434afd1SBarry Smith   KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`
10797b578ef6SBarry Smith 
10807b578ef6SBarry Smith   Calling Sequence:
10817b578ef6SBarry Smith + ksp  - `ksp` context
10827b578ef6SBarry Smith . b    - output vector
10837b578ef6SBarry Smith - ctx - [optional] user-defined function context
10847b578ef6SBarry Smith 
10857b578ef6SBarry Smith   Level: beginner
10867b578ef6SBarry Smith 
10878434afd1SBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
10887b578ef6SBarry Smith S*/
10898434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeRHSFn)(KSP ksp, Vec b, void *ctx);
10907b578ef6SBarry Smith 
10918434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);
10927b578ef6SBarry Smith 
10937b578ef6SBarry Smith /*S
10948434afd1SBarry Smith   KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`
10957b578ef6SBarry Smith 
10967b578ef6SBarry Smith   Calling Sequence:
10977b578ef6SBarry Smith + ksp - `KSP` context
10987b578ef6SBarry Smith . A   - the operator that defines the linear system
10997b578ef6SBarry Smith . P   - an operator from which to build the preconditioner (often the same as `A`)
11007b578ef6SBarry Smith - ctx - [optional] user-defined function context
11017b578ef6SBarry Smith 
11027b578ef6SBarry Smith   Level: beginner
11037b578ef6SBarry Smith 
11048434afd1SBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
11057b578ef6SBarry Smith S*/
11068434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeOperatorsFn)(KSP ksp, Mat A, Mat P, void *ctx);
11077b578ef6SBarry Smith 
11088434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);
11097b578ef6SBarry Smith 
11107b578ef6SBarry Smith /*S
11118434afd1SBarry Smith   KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`
11127b578ef6SBarry Smith 
11137b578ef6SBarry Smith   Calling Sequence:
11147b578ef6SBarry Smith + ksp  - `ksp` context
11157b578ef6SBarry Smith . x    - output vector
11167b578ef6SBarry Smith - ctx - [optional] user-defined function context
11177b578ef6SBarry Smith 
11187b578ef6SBarry Smith   Level: beginner
11197b578ef6SBarry Smith 
11208434afd1SBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
11217b578ef6SBarry Smith S*/
11228434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeInitialGuessFn)(KSP ksp, Vec x, void *ctx);
11237b578ef6SBarry Smith 
11248434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
11258434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
11268434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
11278434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
11288434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
11298434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
11308434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);
11316c699258SBarry Smith 
113202b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
11339371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
1134953494dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode);
1135557cf195SMatthew G. Knepley 
11362b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
11372b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
11384bf303faSJacob Faibussowitsch 
11394bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
11404bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);
11415d83a8b1SBarry Smith 
11425d83a8b1SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);
1143