xref: /petsc/include/petscksp.h (revision 2a28d964fe6fe359fc6e4247643cfa1e8e64ec79)
1f26ada1bSBarry Smith /*
2f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
3f26ada1bSBarry Smith */
426bd1501SBarry Smith #ifndef PETSCKSP_H
526bd1501SBarry Smith #define PETSCKSP_H
6ac09b921SBarry Smith 
72c8e378dSBarry Smith #include <petscpc.h>
82eac72dbSBarry Smith 
9ac09b921SBarry Smith /* SUBMANSEC = KSP */
10ac09b921SBarry Smith 
11607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
121dbb0a54SBarry Smith 
1328ce4d24SBarry Smith /*S
148f6c3df8SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
158f6c3df8SBarry Smith          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
1628ce4d24SBarry Smith 
1728ce4d24SBarry Smith    Level: beginner
1828ce4d24SBarry Smith 
1987497f52SBarry Smith    Note:
20a4d1885cSBarry Smith     When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a
2187497f52SBarry Smith     `KSPType` of `KSPPREONLY`, meaning that only application of the preconditioner is used as the linear solver.
228f6c3df8SBarry Smith 
23a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
2428ce4d24SBarry Smith S*/
25e2a1c21fSSatish Balay typedef struct _p_KSP *KSP;
262eac72dbSBarry Smith 
2776bdecfbSBarry Smith /*J
288f6c3df8SBarry Smith     KSPType - String with the name of a PETSc Krylov method.
2928ce4d24SBarry Smith 
3028ce4d24SBarry Smith    Level: beginner
3128ce4d24SBarry Smith 
32a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
3376bdecfbSBarry Smith J*/
3419fd82e9SBarry Smith typedef const char *KSPType;
3582bf6240SBarry Smith #define KSPRICHARDSON "richardson"
366c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3782bf6240SBarry Smith #define KSPCG         "cg"
382c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
392c8fc646SJed Brown #define KSPPIPECG     "pipecg"
40901ccb91SSiegfried Cools #define KSPPIPECGRR   "pipecgrr"
41265663fdSSiegfried Cools #define KSPPIPELCG    "pipelcg"
42b21a8899STyler Chen #define KSPPIPEPRCG   "pipeprcg"
43325e8391SManasi #define KSPPIPECG2    "pipecg2"
44df2a969aSvictorle #define KSPCGNE       "cgne"
4505de396fSBarry Smith #define KSPNASH       "nash"
4605de396fSBarry Smith #define KSPSTCG       "stcg"
4705de396fSBarry Smith #define KSPGLTR       "gltr"
4805de396fSBarry Smith #define KSPCGNASH     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash"
4905de396fSBarry Smith #define KSPCGSTCG     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg"
5005de396fSBarry Smith #define KSPCGGLTR     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
51640f4f53SPatrick Sanan #define KSPFCG        "fcg"
52390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
5382bf6240SBarry Smith #define KSPGMRES      "gmres"
54483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
559dcbbd2bSBarry Smith #define KSPFGMRES     "fgmres"
569dcbbd2bSBarry Smith #define KSPLGMRES     "lgmres"
574f8e6cd9SSatish Balay #define KSPDGMRES     "dgmres"
5861b468f9SJed Brown #define KSPPGMRES     "pgmres"
5982bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
6082bf6240SBarry Smith #define KSPBCGS       "bcgs"
61e1c61ce8SBarry Smith #define KSPIBCGS      "ibcgs"
62345ecf0bSXiangmin Jiao #define KSPQMRCGS     "qmrcgs"
6318ac38e6SHong Zhang #define KSPFBCGS      "fbcgs"
64c2b71004SHong Zhang #define KSPFBCGSR     "fbcgsr"
65f0037002Svictorle #define KSPBCGSL      "bcgsl"
66f154af2dSSiegfried Cools #define KSPPIPEBCGS   "pipebcgs"
6782bf6240SBarry Smith #define KSPCGS        "cgs"
6882bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
6982bf6240SBarry Smith #define KSPCR         "cr"
702c8fc646SJed Brown #define KSPPIPECR     "pipecr"
7182bf6240SBarry Smith #define KSPLSQR       "lsqr"
7282bf6240SBarry Smith #define KSPPREONLY    "preonly"
733c2be86cSBarry Smith #define KSPNONE       "none"
7482bf6240SBarry Smith #define KSPQCG        "qcg"
75c9cf9b11SBarry Smith #define KSPBICG       "bicg"
76b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
7701c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
7821356919SSatish Balay #define KSPLCD        "lcd"
791d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
8058ad63f4SBarry Smith #define KSPGCR        "gcr"
81fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
82e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
83e4d80e07Szianekhodja #define KSPCGLS       "cgls"
84329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
8538cfc46eSPierre Jolivet #define KSPHPDDM      "hpddm"
862eac72dbSBarry Smith 
878ba1e511SMatthew Knepley /* Logging support */
88014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
89ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
905358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
918ba1e511SMatthew Knepley 
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
9319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
94c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
9975f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
1005ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
1013e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
102d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
103d71ae5a4SJacob Faibussowitsch {
1049371c9d4SSatish Balay   return KSPSetMatSolveBatchSize(ksp, n);
1059371c9d4SSatish Balay }
1063e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
107d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
108d71ae5a4SJacob Faibussowitsch {
1099371c9d4SSatish Balay   return KSPGetMatSolveBatchSize(ksp, n);
1109371c9d4SSatish Balay }
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
112ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
11423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
1157d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
11625c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
117c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
1182eac72dbSBarry Smith 
119140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
120ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
121798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList;
122798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
123798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
124bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
125798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
12630de9b25SBarry Smith 
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
130c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
1367d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
137c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
139c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
144734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
1452a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
146d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
147d71ae5a4SJacob Faibussowitsch {
1489371c9d4SSatish Balay   return KSPCreateVecs(ksp, n, x, m, y);
1499371c9d4SSatish Balay }
1502eac72dbSBarry Smith 
151d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
152d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
153d4211eb9SBarry Smith 
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
156aabeff55SBarry Smith 
157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
1603ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
16194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool);
16394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
16494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool);
1654b0e389bSBarry Smith 
166fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
167fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
168fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
169fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
170fa0ddf94SBarry Smith 
171014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
172cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
177285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
178b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
179b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
180b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
181b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
1834a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
184f3b08a26SMatthew G. Knepley /*
185f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
186f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
187f3b08a26SMatthew G. Knepley */
188f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
1892b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
1902b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
191f3b08a26SMatthew G. Knepley 
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
1942eac72dbSBarry Smith 
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
19858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
199b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
20058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
202d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
203d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
2047d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
2054b0e389bSBarry Smith 
206640f4f53SPatrick Sanan /*E
207640f4f53SPatrick Sanan 
20806137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
209640f4f53SPatrick Sanan 
21067b8a455SSatish Balay   Values:
211a1cb98faSBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
212a1cb98faSBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
213640f4f53SPatrick Sanan 
2142b26979fSBarry Smith    Level: intermediate
215640f4f53SPatrick Sanan 
216a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
217640f4f53SPatrick Sanan E*/
2189371c9d4SSatish Balay typedef enum {
2199371c9d4SSatish Balay   KSP_FCD_TRUNC_TYPE_STANDARD,
2209371c9d4SSatish Balay   KSP_FCD_TRUNC_TYPE_NOTAY
2219371c9d4SSatish Balay } KSPFCDTruncationType;
22206137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
223640f4f53SPatrick Sanan 
224640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
225640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
226640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
227640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
22806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
22906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
230640f4f53SPatrick Sanan 
231390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
232390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
233390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
234390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
235390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
236390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
237390d8e47SPatrick Sanan 
238fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
239fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
240fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
241fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
242fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
243fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
244fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
245fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
24683f0b094SBarry Smith 
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
250e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
2519f236934SBarry Smith 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
2571d73ed98SBarry Smith 
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2601d73ed98SBarry Smith 
261483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
262483d6965SPatrick Sanan 
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
26658ad63f4SBarry Smith 
26704d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
26804d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
26904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
270568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
271e82af88dSprj- 
2726cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
2736cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
2746cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM)
275d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
276d71ae5a4SJacob Faibussowitsch {
2779371c9d4SSatish Balay   return KSPHPDDMSetDeflationMat(ksp, U);
2789371c9d4SSatish Balay }
279d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
280d71ae5a4SJacob Faibussowitsch {
2819371c9d4SSatish Balay   return KSPHPDDMGetDeflationMat(ksp, U);
2829371c9d4SSatish Balay }
2836cac28cbSPierre Jolivet #endif
284d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
285d71ae5a4SJacob Faibussowitsch {
2869371c9d4SSatish Balay   return KSPMatSolve(ksp, B, X);
2879371c9d4SSatish Balay }
288d55ab951SPierre Jolivet /*E
28987497f52SBarry Smith     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
290d55ab951SPierre Jolivet 
291d55ab951SPierre Jolivet     Level: intermediate
292d55ab951SPierre Jolivet 
293d55ab951SPierre Jolivet     Values:
294a4d1885cSBarry Smith +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
295a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BGMRES` - block GMRES
296a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_CG` - Conjugate Gradient
297a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BCG` - block CG
298a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
299a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR
300a4d1885cSBarry Smith .   `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG
301a4d1885cSBarry Smith -   `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only
302d55ab951SPierre Jolivet 
303a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
304d55ab951SPierre Jolivet E*/
305d55ab951SPierre Jolivet typedef enum {
306d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES   = 0,
307d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES  = 1,
308d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG      = 2,
309d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG     = 3,
310d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR  = 4,
311d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
312d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG   = 6,
313d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
314d55ab951SPierre Jolivet } KSPHPDDMType;
315d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
316a4d1885cSBarry Smith 
3172dd49c90SPierre Jolivet /*E
31887497f52SBarry Smith     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
3192dd49c90SPierre Jolivet 
3202dd49c90SPierre Jolivet     Level: intermediate
3212dd49c90SPierre Jolivet 
3222dd49c90SPierre Jolivet     Values:
323a4d1885cSBarry Smith +   `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16`
324a4d1885cSBarry Smith .   `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single`
325a4d1885cSBarry Smith .   `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double`
326a4d1885cSBarry Smith -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
3272dd49c90SPierre Jolivet 
328a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPHPDDM`
3292dd49c90SPierre Jolivet E*/
3302dd49c90SPierre Jolivet typedef enum {
3312dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_HALF      = 0,
3322dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_SINGLE    = 1,
3332dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_DOUBLE    = 2,
3342dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_QUADRUPLE = 3
3352dd49c90SPierre Jolivet } KSPHPDDMPrecision;
336d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
337d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
338e82af88dSprj- 
339b49fd9e1SBarry Smith /*E
340b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
341b49fd9e1SBarry Smith 
342b49fd9e1SBarry Smith    Level: advanced
343b49fd9e1SBarry Smith 
344a4d1885cSBarry Smith     Values:
345a4d1885cSBarry Smith +   `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt
346a4d1885cSBarry Smith .   `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
347a4d1885cSBarry Smith -   `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps
348b49fd9e1SBarry Smith 
349a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
350a4d1885cSBarry Smith           `KSPGMRESGetOrthogonalization()`,
351a4d1885cSBarry Smith           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
352b49fd9e1SBarry Smith E*/
3539371c9d4SSatish Balay typedef enum {
3549371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_NEVER,
3559371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_IFNEEDED,
3569371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_ALWAYS
3579371c9d4SSatish Balay } KSPGMRESCGSRefinementType;
3586a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
3591f7e983dSSatish Balay /*MC
3601957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
3618c5b8ba0SBarry Smith 
3628c5b8ba0SBarry Smith    Level: advanced
3638c5b8ba0SBarry Smith 
36487497f52SBarry Smith    Note:
36587497f52SBarry Smith    Possibly unstable, but the fastest to compute
3668c5b8ba0SBarry Smith 
367a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
368a4d1885cSBarry Smith           `KSP`, `KSPGMRESGetOrthogonalization()`,
369db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
370db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3718c5b8ba0SBarry Smith M*/
3728c5b8ba0SBarry Smith 
3731f7e983dSSatish Balay /*MC
3748c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
3758c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
3768c5b8ba0SBarry Smith           poor orthogonality.
3778c5b8ba0SBarry Smith 
3788c5b8ba0SBarry Smith    Level: advanced
3798c5b8ba0SBarry Smith 
380a4d1885cSBarry Smith    Note:
381a4d1885cSBarry Smith    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
3828c5b8ba0SBarry Smith    estimate the orthogonality but is more stable.
3838c5b8ba0SBarry Smith 
384a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
385a4d1885cSBarry Smith           `KSP`, `KSPGMRESGetOrthogonalization()`,
386db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
387db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3888c5b8ba0SBarry Smith M*/
3898c5b8ba0SBarry Smith 
3901f7e983dSSatish Balay /*MC
3918c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3928c5b8ba0SBarry Smith 
3938c5b8ba0SBarry Smith    Level: advanced
3948c5b8ba0SBarry Smith 
39587497f52SBarry Smith    Notes:
39687497f52SBarry Smith    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
39787497f52SBarry Smith    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
3988c5b8ba0SBarry Smith 
3998c5b8ba0SBarry Smith    You should only use this if you absolutely know that the iterative refinement is needed.
4008c5b8ba0SBarry Smith 
401a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
402a4d1885cSBarry Smith           `KSP`, `KSPGMRESGetOrthogonalization()`,
403db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
404db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
4058c5b8ba0SBarry Smith M*/
4068c5b8ba0SBarry Smith 
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
40908480c60SBarry Smith 
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
413c38d4ed2SBarry Smith 
414014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
415014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
416014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
417121fd945SKris Buschelman 
418014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
419014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
420014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
421e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
422d9492815SBarry Smith 
423014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
42487d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
425014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
4262eac72dbSBarry Smith 
427798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
428798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
429798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
430798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
431798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
432798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
433798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
434798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
435798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
436798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
437798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
438798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
439798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
440798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
441798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
442798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
443798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
444798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
445798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
446798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
447798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
448fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
449798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
450d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
451d71ae5a4SJacob Faibussowitsch {
4529371c9d4SSatish Balay   return KSPMonitorResidual(ksp, n, rnorm, vf);
4539371c9d4SSatish Balay }
454d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
455d71ae5a4SJacob Faibussowitsch {
4569371c9d4SSatish Balay   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
4579371c9d4SSatish Balay }
458d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
459d71ae5a4SJacob Faibussowitsch {
4609371c9d4SSatish Balay   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
4619371c9d4SSatish Balay }
462798534f6SMatthew G. Knepley 
463798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
464341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
465341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
466341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
467341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
468e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
469e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
470e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
47184cb2905SBarry Smith 
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
474c01c455dSBarry Smith 
47523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
47623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
4811eb62cbbSBarry Smith 
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
4861f7f0c4fSBarry Smith 
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
48855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
489fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
49019a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
491c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
4921b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
493c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
49494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
4951b2b9847SBarry Smith 
496d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
497d71ae5a4SJacob Faibussowitsch {
4989371c9d4SSatish Balay   return KSPConvergedReasonView(ksp, v);
4999371c9d4SSatish Balay }
500d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
501d71ae5a4SJacob Faibussowitsch {
5029371c9d4SSatish Balay   return KSPConvergedReasonViewFromOptions(ksp);
5039371c9d4SSatish Balay }
50455849f57SBarry Smith 
50555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
5061eb62cbbSBarry Smith 
507dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
5080e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
509014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
510884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
511798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
512798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
513798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
514db9b2ab1SHong Zhang 
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
51768ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
51883ab6a24SBarry Smith 
51928ce4d24SBarry Smith /*E
520a4d1885cSBarry Smith     KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
5218a4b9c5cSBarry Smith        test routines.
5228a4b9c5cSBarry Smith 
523a4d1885cSBarry Smith     Values:
524a4d1885cSBarry Smith +   `KSP_NORM_DEFAULT` - use the default for the current `KSPType`
525a4d1885cSBarry Smith .   `KSP_NORM_NONE` - use no norm calculation
526a4d1885cSBarry Smith .   `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm
527a4d1885cSBarry Smith .   `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
528a4d1885cSBarry Smith -   `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator)
529a4d1885cSBarry Smith 
5308a4b9c5cSBarry Smith    Level: advanced
5318a4b9c5cSBarry Smith 
532a4d1885cSBarry Smith    Note:
533a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
53487497f52SBarry Smith    depending on left or right preconditioning, see `KSPSetPCSide()`
535a3f661c8SBarry Smith 
53687497f52SBarry Smith    Developer Note:
537a4d1885cSBarry Smith    This must match the values in petsc/finclude/petscksp.h
5388a4b9c5cSBarry Smith 
539a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
540db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
5418a4b9c5cSBarry Smith E*/
5429371c9d4SSatish Balay typedef enum {
5439371c9d4SSatish Balay   KSP_NORM_DEFAULT          = -1,
5449371c9d4SSatish Balay   KSP_NORM_NONE             = 0,
5459371c9d4SSatish Balay   KSP_NORM_PRECONDITIONED   = 1,
5469371c9d4SSatish Balay   KSP_NORM_UNPRECONDITIONED = 2,
5479371c9d4SSatish Balay   KSP_NORM_NATURAL          = 3
5489371c9d4SSatish Balay } KSPNormType;
5499e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
550014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes;
551a21b2a99SBarry Smith 
5521f7e983dSSatish Balay /*MC
5539793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
5548c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
5558c5b8ba0SBarry Smith           be based on a norm of a residual etc.
5568c5b8ba0SBarry Smith 
5578c5b8ba0SBarry Smith    Level: advanced
5588c5b8ba0SBarry Smith 
559a4d1885cSBarry Smith     Note:
560a4d1885cSBarry Smith     Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
5618c5b8ba0SBarry Smith 
562a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
5638c5b8ba0SBarry Smith M*/
5648c5b8ba0SBarry Smith 
5651f7e983dSSatish Balay /*MC
5661957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
5678c5b8ba0SBarry Smith        convergence test routine.
5688c5b8ba0SBarry Smith 
5698c5b8ba0SBarry Smith    Level: advanced
5708c5b8ba0SBarry Smith 
571a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
5728c5b8ba0SBarry Smith M*/
5738c5b8ba0SBarry Smith 
5741f7e983dSSatish Balay /*MC
575ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
5768c5b8ba0SBarry Smith        convergence test routine.
5778c5b8ba0SBarry Smith 
5788c5b8ba0SBarry Smith    Level: advanced
5798c5b8ba0SBarry Smith 
580a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
5818c5b8ba0SBarry Smith M*/
5828c5b8ba0SBarry Smith 
5831f7e983dSSatish Balay /*MC
584ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
58587497f52SBarry Smith        convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
5868c5b8ba0SBarry Smith 
5878c5b8ba0SBarry Smith    Level: advanced
5888c5b8ba0SBarry Smith 
589a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
5908c5b8ba0SBarry Smith M*/
5918c5b8ba0SBarry Smith 
592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
593014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
594014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
5978a4b9c5cSBarry Smith 
5984a221d59SStefano Zampini #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_NEG_CURVE (since version 3.19)")
5994a221d59SStefano Zampini #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_STEP_LENGTH (since version 3.19)")
600f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
6018a4b9c5cSBarry Smith /*E
602a4d1885cSBarry Smith     KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
60328ce4d24SBarry Smith 
60428ce4d24SBarry Smith    Level: beginner
60528ce4d24SBarry Smith 
606a4d1885cSBarry Smith    Values:
607a4d1885cSBarry Smith +  `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations
608a4d1885cSBarry Smith .  `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations
609a4d1885cSBarry Smith .  `KSP_CONVERGED_RTOL` - requested decrease in the residual
610a4d1885cSBarry Smith .  `KSP_CONVERGED_ATOL` - requested absolute value in the residual
611a4d1885cSBarry Smith .  `KSP_CONVERGED_ITS` - requested number of iterations
6124a221d59SStefano Zampini .  `KSP_CONVERGED_NEG_CURVE` - see note below
613a4d1885cSBarry Smith .  `KSP_CONVERGED_STEP_LENGTH` - see note below
614a4d1885cSBarry Smith .  `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred.
615a4d1885cSBarry Smith .  `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES
616a4d1885cSBarry Smith .  `KSP_DIVERGED_ITS` - requested number of iterations
617a4d1885cSBarry Smith .  `KSP_DIVERGED_DTOL` - large increase in the residual norm
618a4d1885cSBarry Smith .  `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method
619a4d1885cSBarry Smith .  `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method
620a4d1885cSBarry Smith .  `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
621a4d1885cSBarry Smith .  `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite
622a4d1885cSBarry Smith .  `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation
623a4d1885cSBarry Smith .  `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite
624a4d1885cSBarry Smith -  `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason
625a4d1885cSBarry Smith 
626a4d1885cSBarry Smith    Note:
6274a221d59SStefano Zampini    The values  `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`,
628a4d1885cSBarry Smith    `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
62928ce4d24SBarry Smith 
63095452b02SPatrick Sanan    Developer Notes:
631a4d1885cSBarry Smith    This must match the values in petsc/finclude/petscksp.h
6324d0a8057SBarry Smith 
63387497f52SBarry Smith    The string versions of these are `KSPConvergedReasons`; if you change
6344d0a8057SBarry Smith    any of the values here also change them that array of names.
63586c02ca4SBarry Smith 
636a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
63728ce4d24SBarry Smith E*/
638d15094e1SBarry Smith typedef enum { /* converged */
6399ecb1286SBarry Smith   KSP_CONVERGED_RTOL_NORMAL               = 1,
6409ecb1286SBarry Smith   KSP_CONVERGED_ATOL_NORMAL               = 9,
641d15094e1SBarry Smith   KSP_CONVERGED_RTOL                      = 2,
642d15094e1SBarry Smith   KSP_CONVERGED_ATOL                      = 3,
643b335793eSSatish Balay   KSP_CONVERGED_ITS                       = 4,
6444a221d59SStefano Zampini   KSP_CONVERGED_NEG_CURVE                 = 5,
6454a221d59SStefano Zampini   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
6464a221d59SStefano Zampini   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
6474a221d59SStefano Zampini   KSP_CONVERGED_STEP_LENGTH               = 6,
6484a221d59SStefano Zampini   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
649d15094e1SBarry Smith   /* diverged */
650b3cc6726SBarry Smith   KSP_DIVERGED_NULL                      = -2,
651d15094e1SBarry Smith   KSP_DIVERGED_ITS                       = -3,
652d15094e1SBarry Smith   KSP_DIVERGED_DTOL                      = -4,
653d15094e1SBarry Smith   KSP_DIVERGED_BREAKDOWN                 = -5,
654b4ac9ba4SBarry Smith   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
655b4ac9ba4SBarry Smith   KSP_DIVERGED_NONSYMMETRIC              = -7,
656b4ac9ba4SBarry Smith   KSP_DIVERGED_INDEFINITE_PC             = -8,
6574d51c080SBarry Smith   KSP_DIVERGED_NANORINF                  = -9,
6586aee1118SBarry Smith   KSP_DIVERGED_INDEFINITE_MAT            = -10,
659c0decd05SBarry Smith   KSP_DIVERGED_PC_FAILED                 = -11,
660aa4d2078SSatish Balay   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
661d15094e1SBarry Smith 
6629371c9d4SSatish Balay   KSP_CONVERGED_ITERATING = 0
6639371c9d4SSatish Balay } KSPConvergedReason;
664014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons;
665d15094e1SBarry Smith 
666c838673bSBarry Smith /*MC
66787497f52SBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called
668c838673bSBarry Smith 
669c838673bSBarry Smith    Level: beginner
670c838673bSBarry Smith 
67187497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
672c838673bSBarry Smith    for left preconditioning it is the 2-norm of the preconditioned residual, and the
673c838673bSBarry Smith    2-norm of the residual for right preconditioning
674c838673bSBarry Smith 
67587497f52SBarry Smith    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
676f9fed41fSBarry Smith 
677a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
678c838673bSBarry Smith 
679c838673bSBarry Smith M*/
680c838673bSBarry Smith 
681c838673bSBarry Smith /*MC
682c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
683c838673bSBarry Smith 
684c838673bSBarry Smith    Level: beginner
685c838673bSBarry Smith 
68687497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
687c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
688c838673bSBarry Smith        2-norm of the residual for right preconditioning
689c838673bSBarry Smith 
69087497f52SBarry Smith    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
691c838673bSBarry Smith 
692a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
693c838673bSBarry Smith 
694c838673bSBarry Smith M*/
695c838673bSBarry Smith 
696c838673bSBarry Smith /*MC
697c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
698c838673bSBarry Smith 
699c838673bSBarry Smith    Level: beginner
700c838673bSBarry Smith 
70187497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
702c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
703c838673bSBarry Smith        2-norm of the residual for right preconditioning
704c838673bSBarry Smith 
705c838673bSBarry Smith    Level: beginner
706c838673bSBarry Smith 
707a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
708c838673bSBarry Smith 
709c838673bSBarry Smith M*/
710c838673bSBarry Smith 
711c838673bSBarry Smith /*MC
712c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
713c838673bSBarry Smith       reached
714c838673bSBarry Smith 
715c838673bSBarry Smith    Level: beginner
716c838673bSBarry Smith 
717a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
718c838673bSBarry Smith 
719c838673bSBarry Smith M*/
720c838673bSBarry Smith 
721c838673bSBarry Smith /*MC
72287497f52SBarry Smith      KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
72387497f52SBarry Smith            the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
7244d0a8057SBarry Smith            test routine is set in KSP.
725c838673bSBarry Smith 
726c838673bSBarry Smith    Level: beginner
727c838673bSBarry Smith 
728a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
729c838673bSBarry Smith 
730c838673bSBarry Smith M*/
731c838673bSBarry Smith 
732c838673bSBarry Smith /*MC
733c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
7341de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
73587497f52SBarry Smith           preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
7361de96524SPierre Jolivet           are colinear.
737c838673bSBarry Smith 
738c838673bSBarry Smith    Level: beginner
739c838673bSBarry Smith 
740a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
741c838673bSBarry Smith 
742c838673bSBarry Smith M*/
743c838673bSBarry Smith 
744c838673bSBarry Smith /*MC
74587497f52SBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
746c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
747c838673bSBarry Smith 
748c838673bSBarry Smith    Level: beginner
749c838673bSBarry Smith 
750a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
751c838673bSBarry Smith 
752c838673bSBarry Smith M*/
753c838673bSBarry Smith 
754c838673bSBarry Smith /*MC
755c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
75687497f52SBarry Smith         symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
757c838673bSBarry Smith 
758c838673bSBarry Smith    Level: beginner
759c838673bSBarry Smith 
760a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
761c838673bSBarry Smith 
762c838673bSBarry Smith M*/
763c838673bSBarry Smith 
764c838673bSBarry Smith /*MC
765c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
76687497f52SBarry Smith         positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
767c838673bSBarry Smith         be positive definite
768c838673bSBarry Smith 
769c838673bSBarry Smith    Level: beginner
770c838673bSBarry Smith 
77187497f52SBarry Smith      Note:
772a4d1885cSBarry Smith     This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
77387497f52SBarry Smith   the `PCICC` preconditioner to generate a positive definite preconditioner
774c838673bSBarry Smith 
775a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
776c838673bSBarry Smith 
777c838673bSBarry Smith M*/
778c838673bSBarry Smith 
779c838673bSBarry Smith /*MC
780c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
7819fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
78287497f52SBarry Smith      such as `PCFIELDSPLIT`.
7839fc87aa7SBarry Smith 
7849fc87aa7SBarry Smith    Level: beginner
7859fc87aa7SBarry Smith 
786a4d1885cSBarry Smith     Note:
787a4d1885cSBarry Smith     Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
7889fc87aa7SBarry Smith 
789a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
7909fc87aa7SBarry Smith 
7919fc87aa7SBarry Smith M*/
7929fc87aa7SBarry Smith 
7939fc87aa7SBarry Smith /*MC
794a4d1885cSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
795a4d1885cSBarry Smith         while `KSPSolve()` is still running.
796c838673bSBarry Smith 
797c838673bSBarry Smith    Level: beginner
798c838673bSBarry Smith 
799a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
800c838673bSBarry Smith 
801c838673bSBarry Smith M*/
802c838673bSBarry Smith 
803014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
804df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
805df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
8063ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
8078de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
8083eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
8098de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
8108de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
8118de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
8128de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
81354b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
8140059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
815014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
816c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
81794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
818*2a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
819*2a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
820abef13c0SSatish Balay 
821d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void)
822d71ae5a4SJacob Faibussowitsch { /* never called */
8239371c9d4SSatish Balay }
8248ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
825d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void)
826d71ae5a4SJacob Faibussowitsch { /* never called */
8279371c9d4SSatish Balay }
8288ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
829d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void)
830d71ae5a4SJacob Faibussowitsch { /* never called */
8319371c9d4SSatish Balay }
8328ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
833d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void)
834d71ae5a4SJacob Faibussowitsch { /* never called */
8359371c9d4SSatish Balay }
8368ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
837d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void)
838d71ae5a4SJacob Faibussowitsch { /* never called */
8399371c9d4SSatish Balay }
8408ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
841d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void)
842d71ae5a4SJacob Faibussowitsch { /* never called */
8439371c9d4SSatish Balay }
8448ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
8458ea1b3e6SJed Brown 
8460bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
847d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
848d71ae5a4SJacob Faibussowitsch {
8499371c9d4SSatish Balay   return KSPComputeOperator(A, NULL, B);
8509371c9d4SSatish Balay }
851d4fbbf0eSBarry Smith 
85228ce4d24SBarry Smith /*E
853a4d1885cSBarry Smith     KSPCGType - Determines what type of `KSPCG` to use
85428ce4d24SBarry Smith 
85528ce4d24SBarry Smith    Level: beginner
85628ce4d24SBarry Smith 
857a4d1885cSBarry Smith    Values:
858a4d1885cSBarry Smith  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
859a4d1885cSBarry Smith  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
860a4d1885cSBarry Smith 
861a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPCGSetType()`
86228ce4d24SBarry Smith E*/
8639371c9d4SSatish Balay typedef enum {
8649371c9d4SSatish Balay   KSP_CG_SYMMETRIC = 0,
8659371c9d4SSatish Balay   KSP_CG_HERMITIAN = 1
8669371c9d4SSatish Balay } KSPCGType;
8676a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
86828ce4d24SBarry Smith 
869014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
870014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
8718031f4adStmunson 
872dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
873dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
874dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
875fcae7a14Stmunson 
87605de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
87705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
878d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
879d71ae5a4SJacob Faibussowitsch {
8809371c9d4SSatish Balay   return KSPGLTRGetMinEig(ksp, x);
8819371c9d4SSatish Balay }
882d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
883d71ae5a4SJacob Faibussowitsch {
8849371c9d4SSatish Balay   return KSPGLTRGetLambda(ksp, x);
8859371c9d4SSatish Balay }
8868031f4adStmunson 
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
888ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
8891d6018f0SLisandro Dalcin 
890f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
891014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
892014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
8933369ce9aSBarry Smith 
8949804daf3SBarry Smith #include <petscdrawtypes.h>
895014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
8962f2e5d10SKris Buschelman 
897014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
898014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
89903919abeSBarry Smith 
900ba36735cSStefano Zampini /*S
901a4d1885cSBarry Smith      KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
902f8a50e2bSBarry Smith 
903a4d1885cSBarry Smith    Level: intermediate
904f8a50e2bSBarry Smith 
905a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
906ba36735cSStefano Zampini S*/
907ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess;
908ba36735cSStefano Zampini /*J
90987497f52SBarry Smith     KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
910ba36735cSStefano Zampini 
911a4d1885cSBarry Smith    Level: intermediate
912ba36735cSStefano Zampini 
913a4d1885cSBarry Smith    Values:
914a4d1885cSBarry Smith  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
915a4d1885cSBarry Smith  - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition
916a4d1885cSBarry Smith 
917a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPGuess`
918ba36735cSStefano Zampini J*/
919ba36735cSStefano Zampini typedef const char *KSPGuessType;
920ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
921ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
922a4d1885cSBarry Smith 
9231d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
924ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
925ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
926ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
927ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
928ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
929ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
930ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
9318410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
932ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
933ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
934ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
935ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
936ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
937014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
938ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
939ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
940f8a50e2bSBarry Smith 
941470b340bSDmitry Karpeev /*E
942470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
943470b340bSDmitry Karpeev 
944470b340bSDmitry Karpeev     Level: intermediate
945470b340bSDmitry Karpeev 
946a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
947a4d1885cSBarry Smith           `MatCreateSchurComplementPmat()`
948470b340bSDmitry Karpeev E*/
9499371c9d4SSatish Balay typedef enum {
9509371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
9519371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
9529371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
9539371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_FULL
9549371c9d4SSatish Balay } MatSchurComplementAinvType;
955470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
956470b340bSDmitry Karpeev 
9579371c9d4SSatish Balay typedef enum {
9589371c9d4SSatish Balay   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
959864588a7SAlp Dener   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
960864588a7SAlp Dener   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
9619371c9d4SSatish Balay   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
9629371c9d4SSatish Balay } MatLMVMSymBroydenScaleType;
963864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
96492f76d53SAlp Dener 
965014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
966014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
967d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
968bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
969aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
970bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
971470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
972470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
9735bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
9745a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
975470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
976470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
9773f22127dSBarry Smith 
97878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
97978e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
98078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
981864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
982864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
983864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
984864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
985864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
986cd929ea3SAlp Dener 
987cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
988b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
989cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
990cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
991e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
9920ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
993cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
994cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
995cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
996cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
997cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
9982d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
999cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1000cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1001cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1002cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
100392f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1004cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1005cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1006864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1007864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1008cd929ea3SAlp Dener 
1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1010014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1012014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1013014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1014fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
101523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1016fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
101723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
101823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1021fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1022fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
10236c699258SBarry Smith 
102402b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
10259371c9d4SSatish 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);
1026557cf195SMatthew G. Knepley 
10272b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
10282b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
10292eac72dbSBarry Smith #endif
1030