xref: /petsc/include/petscksp.h (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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 
1995452b02SPatrick Sanan         Notes:
20db3b2ab5SMatthew Knepley     When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a
21db3b2ab5SMatthew Knepley        KSPType of KSPPREONLY, meaning that only application of the preconditioner is used as the linear solver.
228f6c3df8SBarry Smith 
23db781477SPatrick Sanan .seealso: `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 
32db781477SPatrick Sanan .seealso: `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);
1029fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp,PetscInt n) {return KSPSetMatSolveBatchSize(ksp,n);}
1033e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP,PetscInt*);
1049fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp,PetscInt *n) {return KSPGetMatSolveBatchSize(ksp,n);}
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
106ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
10823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
1097d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
11025c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
111c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
1122eac72dbSBarry Smith 
113140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
114ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
115798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList;
116798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
117798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
118bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
119798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
12030de9b25SBarry Smith 
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
124c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
1307d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
131c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
133c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
138734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1392a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1409fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
1412eac72dbSBarry Smith 
142d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
143d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
144d4211eb9SBarry Smith 
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
147aabeff55SBarry Smith 
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
1513ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void*);
15294a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*);
153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);
15494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*);
15594a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool);
1564b0e389bSBarry Smith 
157fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
158fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
159fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
160fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
161fa0ddf94SBarry Smith 
162014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
163cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
164014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
165014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
167014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
168285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
169b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
170b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
171b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
172b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
1744a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
175f3b08a26SMatthew G. Knepley /*
176f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
177f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
178f3b08a26SMatthew G. Knepley */
179f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
1802b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*));
1812b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*));
182f3b08a26SMatthew G. Knepley 
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
1852eac72dbSBarry Smith 
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
18958450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
190b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
19158450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
193d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
194d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1957d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
1964b0e389bSBarry Smith 
197640f4f53SPatrick Sanan /*E
198640f4f53SPatrick Sanan 
19906137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
200640f4f53SPatrick Sanan 
20167b8a455SSatish Balay   Values:
20267b8a455SSatish Balay $ KSP_FCD_TRUNC_TYPE_STANDARD - uses all (up to mmax) stored directions
20367b8a455SSatish Balay $ KSP_FCD_TRUNC_TYPE_NOTAY - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
204640f4f53SPatrick Sanan 
2052b26979fSBarry Smith    Level: intermediate
206c2e3fba1SPatrick Sanan .seealso `:` `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
207640f4f53SPatrick Sanan 
208640f4f53SPatrick Sanan E*/
20993f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
21006137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
211640f4f53SPatrick Sanan 
212640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
213640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
214640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
215640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
21606137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
21706137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
218640f4f53SPatrick Sanan 
219390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
220390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
221390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
222390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
223390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
224390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
225390d8e47SPatrick Sanan 
226fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
227fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
228fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
229fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
230fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
231fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
232fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
233fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
23483f0b094SBarry Smith 
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
238e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal);
2399f236934SBarry Smith 
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2451d73ed98SBarry Smith 
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2481d73ed98SBarry Smith 
249483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
250483d6965SPatrick Sanan 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
25458ad63f4SBarry Smith 
25504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
25604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
25704d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
258568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
259e82af88dSprj- 
2606cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP,Mat);
2616cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP,Mat*);
2626cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM)
2636cac28cbSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) { return KSPHPDDMSetDeflationMat(ksp, U); }
2646cac28cbSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) { return KSPHPDDMGetDeflationMat(ksp, U); }
2656cac28cbSPierre Jolivet #endif
2669fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
267d55ab951SPierre Jolivet /*E
268d55ab951SPierre Jolivet     KSPHPDDMType - Type of Krylov method used by KSPHPDDM
269d55ab951SPierre Jolivet 
270d55ab951SPierre Jolivet     Level: intermediate
271d55ab951SPierre Jolivet 
272d55ab951SPierre Jolivet     Values:
27367b8a455SSatish Balay $   KSP_HPDDM_TYPE_GMRES (default)
27467b8a455SSatish Balay $   KSP_HPDDM_TYPE_BGMRES
27567b8a455SSatish Balay $   KSP_HPDDM_TYPE_CG
27667b8a455SSatish Balay $   KSP_HPDDM_TYPE_BCG
27767b8a455SSatish Balay $   KSP_HPDDM_TYPE_GCRODR
27867b8a455SSatish Balay $   KSP_HPDDM_TYPE_BGCRODR
27967b8a455SSatish Balay $   KSP_HPDDM_TYPE_BFBCG
28067b8a455SSatish Balay $   KSP_HPDDM_TYPE_PREONLY
281d55ab951SPierre Jolivet 
282db781477SPatrick Sanan .seealso: `KSPHPDDM`, `KSPHPDDMSetType()`
283d55ab951SPierre Jolivet E*/
284d55ab951SPierre Jolivet typedef enum {
285d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES = 0,
286d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES = 1,
287d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG = 2,
288d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG = 3,
289d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR = 4,
290d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
291d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG = 6,
292d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
293d55ab951SPierre Jolivet } KSPHPDDMType;
294d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
2952dd49c90SPierre Jolivet /*E
2962dd49c90SPierre Jolivet     KSPHPDDMPrecision - Precision of Krylov bases used by KSPHPDDM
2972dd49c90SPierre Jolivet 
2982dd49c90SPierre Jolivet     Level: intermediate
2992dd49c90SPierre Jolivet 
3002dd49c90SPierre Jolivet     Values:
30167b8a455SSatish Balay $   KSP_HPDDM_PRECISION_HALF (currently unsupported)
30267b8a455SSatish Balay $   KSP_HPDDM_PRECISION_SINGLE (default when PETSc is configured --with-precision=single)
30367b8a455SSatish Balay $   KSP_HPDDM_PRECISION_DOUBLE (default when PETSc is configured --with-precision=double)
30467b8a455SSatish Balay $   KSP_HPDDM_PRECISION_QUADRUPLE (currently unsupported)
3052dd49c90SPierre Jolivet 
306db781477SPatrick Sanan .seealso: `KSPHPDDM`
3072dd49c90SPierre Jolivet E*/
3082dd49c90SPierre Jolivet typedef enum {
3092dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_HALF = 0,
3102dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_SINGLE = 1,
3112dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_DOUBLE = 2,
3122dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_QUADRUPLE = 3
3132dd49c90SPierre Jolivet } KSPHPDDMPrecision;
314d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
315d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
316e82af88dSprj- 
317b49fd9e1SBarry Smith /*E
318b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
319b49fd9e1SBarry Smith 
320b49fd9e1SBarry Smith    Level: advanced
321b49fd9e1SBarry Smith 
322db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
323db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
324b49fd9e1SBarry Smith 
325b49fd9e1SBarry Smith E*/
32678d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
3276a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
3281f7e983dSSatish Balay /*MC
3291957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
3308c5b8ba0SBarry Smith 
3318c5b8ba0SBarry Smith    Level: advanced
3328c5b8ba0SBarry Smith 
3338c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
3348c5b8ba0SBarry Smith 
335db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
336db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
337db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3388c5b8ba0SBarry Smith M*/
3398c5b8ba0SBarry Smith 
3401f7e983dSSatish Balay /*MC
3418c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
3428c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
3438c5b8ba0SBarry Smith           poor orthogonality.
3448c5b8ba0SBarry Smith 
3458c5b8ba0SBarry Smith    Level: advanced
3468c5b8ba0SBarry Smith 
3478c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
3488c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
3498c5b8ba0SBarry Smith 
350db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
351db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
352db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3538c5b8ba0SBarry Smith M*/
3548c5b8ba0SBarry Smith 
3551f7e983dSSatish Balay /*MC
3568c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3578c5b8ba0SBarry Smith 
3588c5b8ba0SBarry Smith    Level: advanced
3598c5b8ba0SBarry Smith 
3608c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
3618c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
3628c5b8ba0SBarry Smith 
3638c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
3648c5b8ba0SBarry Smith 
365db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
366db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
367db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3688c5b8ba0SBarry Smith M*/
3698c5b8ba0SBarry Smith 
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
37208480c60SBarry Smith 
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
376c38d4ed2SBarry Smith 
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
380121fd945SKris Buschelman 
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
384e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
385d9492815SBarry Smith 
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
38787d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
3892eac72dbSBarry Smith 
390798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*);
391798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *);
392798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
393798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
394798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
395798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
396798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
397798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
398798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
399798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
400798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
401798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
402798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
403798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
404798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
405798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
406798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
407798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
408798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
409798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
410798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
411fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
412798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
4139fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);}
4149fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);}
4159fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);}
416798534f6SMatthew G. Knepley 
417798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
418390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
419390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
420e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
421e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
422e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
42384cb2905SBarry Smith 
424014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
425014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
426c01c455dSBarry Smith 
42723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
42823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
429014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
430014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
431014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
432014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
4331eb62cbbSBarry Smith 
434014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
435014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
436014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
437014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
4381f7f0c4fSBarry Smith 
439014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
44055849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
441fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
44219a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer);
443c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP,PetscErrorCode (*)(KSP,void*),void *vctx,PetscErrorCode (*)(void**));
4441b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
445c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
44694a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer);
4471b2b9847SBarry Smith 
4489fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);}
4499fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}
45055849f57SBarry Smith 
45155849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
4521eb62cbbSBarry Smith 
453dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
4540e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
456884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
457798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
458798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
459798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
460db9b2ab1SHong Zhang 
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
46368ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
46483ab6a24SBarry Smith 
46528ce4d24SBarry Smith /*E
4668a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
4678a4b9c5cSBarry Smith        test routines.
4688a4b9c5cSBarry Smith 
4698a4b9c5cSBarry Smith    Level: advanced
4708a4b9c5cSBarry Smith 
471a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
4721718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
473a3f661c8SBarry Smith 
47495452b02SPatrick Sanan    Notes:
47595452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4768a4b9c5cSBarry Smith 
477db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
478db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
4798a4b9c5cSBarry Smith E*/
4809e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
4819e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
482014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
483a21b2a99SBarry Smith 
4841f7e983dSSatish Balay /*MC
4859793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
4868c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
4878c5b8ba0SBarry Smith           be based on a norm of a residual etc.
4888c5b8ba0SBarry Smith 
4898c5b8ba0SBarry Smith    Level: advanced
4908c5b8ba0SBarry Smith 
4911957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
4928c5b8ba0SBarry Smith 
493db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
4948c5b8ba0SBarry Smith M*/
4958c5b8ba0SBarry Smith 
4961f7e983dSSatish Balay /*MC
4971957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
4988c5b8ba0SBarry Smith        convergence test routine.
4998c5b8ba0SBarry Smith 
5008c5b8ba0SBarry Smith    Level: advanced
5018c5b8ba0SBarry Smith 
502db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
5038c5b8ba0SBarry Smith M*/
5048c5b8ba0SBarry Smith 
5051f7e983dSSatish Balay /*MC
506ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
5078c5b8ba0SBarry Smith        convergence test routine.
5088c5b8ba0SBarry Smith 
5098c5b8ba0SBarry Smith    Level: advanced
5108c5b8ba0SBarry Smith 
511db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
5128c5b8ba0SBarry Smith M*/
5138c5b8ba0SBarry Smith 
5141f7e983dSSatish Balay /*MC
515ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
516390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
5178c5b8ba0SBarry Smith 
5188c5b8ba0SBarry Smith    Level: advanced
5198c5b8ba0SBarry Smith 
520db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
5218c5b8ba0SBarry Smith M*/
5228c5b8ba0SBarry Smith 
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
5288a4b9c5cSBarry Smith 
529f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
5308a4b9c5cSBarry Smith /*E
5311957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
53228ce4d24SBarry Smith 
53328ce4d24SBarry Smith    Level: beginner
53428ce4d24SBarry Smith 
53595452b02SPatrick Sanan    Notes:
53695452b02SPatrick Sanan     See KSPGetConvergedReason() for explanation of each value
53728ce4d24SBarry Smith 
53895452b02SPatrick Sanan    Developer Notes:
53995452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
5404d0a8057SBarry Smith 
5414d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
5424d0a8057SBarry Smith       any of the values here also change them that array of names.
54386c02ca4SBarry Smith 
544db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
54528ce4d24SBarry Smith E*/
546d15094e1SBarry Smith typedef enum {/* converged */
5479ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
5489ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
549d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
550d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
551b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
5528031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
5538031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
554329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
555af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
556d15094e1SBarry Smith               /* diverged */
557b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
558d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
559d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
560d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
561b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
562b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
563b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
5644d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
5656aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
566c0decd05SBarry Smith               KSP_DIVERGED_PC_FAILED           = -11,
567aa4d2078SSatish Balay               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
568d15094e1SBarry Smith 
569d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
570014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
571d15094e1SBarry Smith 
572c838673bSBarry Smith /*MC
573f9fed41fSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
574c838673bSBarry Smith 
575c838673bSBarry Smith    Level: beginner
576c838673bSBarry Smith 
577c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
578c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
579c838673bSBarry Smith        2-norm of the residual for right preconditioning
580c838673bSBarry Smith 
581f9fed41fSBarry Smith    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
582f9fed41fSBarry Smith 
583db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
584c838673bSBarry Smith 
585c838673bSBarry Smith M*/
586c838673bSBarry Smith 
587c838673bSBarry Smith /*MC
588c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
589c838673bSBarry Smith 
590c838673bSBarry Smith    Level: beginner
591c838673bSBarry Smith 
592c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
593c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
594c838673bSBarry Smith        2-norm of the residual for right preconditioning
595c838673bSBarry Smith 
596f9fed41fSBarry Smith    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
597c838673bSBarry Smith 
598db781477SPatrick Sanan .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
599c838673bSBarry Smith 
600c838673bSBarry Smith M*/
601c838673bSBarry Smith 
602c838673bSBarry Smith /*MC
603c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
604c838673bSBarry Smith 
605c838673bSBarry Smith    Level: beginner
606c838673bSBarry Smith 
607c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
608c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
609c838673bSBarry Smith        2-norm of the residual for right preconditioning
610c838673bSBarry Smith 
611c838673bSBarry Smith    Level: beginner
612c838673bSBarry Smith 
613db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
614c838673bSBarry Smith 
615c838673bSBarry Smith M*/
616c838673bSBarry Smith 
617c838673bSBarry Smith /*MC
618c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
619c838673bSBarry Smith       reached
620c838673bSBarry Smith 
621c838673bSBarry Smith    Level: beginner
622c838673bSBarry Smith 
623db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
624c838673bSBarry Smith 
625c838673bSBarry Smith M*/
626c838673bSBarry Smith 
627c838673bSBarry Smith /*MC
6288c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
6290059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
6304d0a8057SBarry Smith            test routine is set in KSP.
631c838673bSBarry Smith 
632c838673bSBarry Smith    Level: beginner
633c838673bSBarry Smith 
634db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
635c838673bSBarry Smith 
636c838673bSBarry Smith M*/
637c838673bSBarry Smith 
638c838673bSBarry Smith /*MC
639c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
6401de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
6411de96524SPierre Jolivet           preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
6421de96524SPierre Jolivet           are colinear.
643c838673bSBarry Smith 
644c838673bSBarry Smith    Level: beginner
645c838673bSBarry Smith 
646db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
647c838673bSBarry Smith 
648c838673bSBarry Smith M*/
649c838673bSBarry Smith 
650c838673bSBarry Smith /*MC
651c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
652c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
653c838673bSBarry Smith 
654c838673bSBarry Smith    Level: beginner
655c838673bSBarry Smith 
656db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
657c838673bSBarry Smith 
658c838673bSBarry Smith M*/
659c838673bSBarry Smith 
660c838673bSBarry Smith /*MC
661c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
662c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
663c838673bSBarry Smith 
664c838673bSBarry Smith    Level: beginner
665c838673bSBarry Smith 
666db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
667c838673bSBarry Smith 
668c838673bSBarry Smith M*/
669c838673bSBarry Smith 
670c838673bSBarry Smith /*MC
671c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
672c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
673c838673bSBarry Smith         be positive definite
674c838673bSBarry Smith 
675c838673bSBarry Smith    Level: beginner
676c838673bSBarry Smith 
67795452b02SPatrick Sanan      Notes:
67895452b02SPatrick Sanan     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
679c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
680c838673bSBarry Smith 
681db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
682c838673bSBarry Smith 
683c838673bSBarry Smith M*/
684c838673bSBarry Smith 
685c838673bSBarry Smith /*MC
686c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
6879fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
6889fc87aa7SBarry Smith      such as PCFIELDSPLIT.
6899fc87aa7SBarry Smith 
6909fc87aa7SBarry Smith    Level: beginner
6919fc87aa7SBarry Smith 
69295452b02SPatrick Sanan     Notes:
69395452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
6949fc87aa7SBarry Smith 
695db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
6969fc87aa7SBarry Smith 
6979fc87aa7SBarry Smith M*/
6989fc87aa7SBarry Smith 
6999fc87aa7SBarry Smith /*MC
700c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
701c838673bSBarry Smith         while the KSPSolve() is still running.
702c838673bSBarry Smith 
703c838673bSBarry Smith    Level: beginner
704c838673bSBarry Smith 
705db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
706c838673bSBarry Smith 
707c838673bSBarry Smith M*/
708c838673bSBarry Smith 
709014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
710df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
711df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
7123ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void*);
7138de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
7143eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
7158de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
7168de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
7178de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
7188de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
71954b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
7200059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
721014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
722c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP,const char**);
72394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
724abef13c0SSatish Balay 
7259fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) { /* never called */ }
7268ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
7279fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) { /* never called */ }
7288ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
7299fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) { /* never called */ }
7308ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
7319fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
7328ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
7339fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
7348ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
7359fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) { /* never called */ }
7368ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
7378ea1b3e6SJed Brown 
7380bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
7399fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
740d4fbbf0eSBarry Smith 
74128ce4d24SBarry Smith /*E
74228ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
74328ce4d24SBarry Smith 
74428ce4d24SBarry Smith    Level: beginner
74528ce4d24SBarry Smith 
746db781477SPatrick Sanan .seealso: `KSPCGSetType()`
74728ce4d24SBarry Smith E*/
7489dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
7496a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
75028ce4d24SBarry Smith 
751014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
752014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);
7538031f4adStmunson 
754dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
755dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
756dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
757fcae7a14Stmunson 
75805de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
75905de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
7609fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
7619fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
7628031f4adStmunson 
763014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
764*ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP,const char*[]);
7651d6018f0SLisandro Dalcin 
766f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC,PetscErrorCode (*)(PC,KSP));
767014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
768014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
7693369ce9aSBarry Smith 
7709804daf3SBarry Smith #include <petscdrawtypes.h>
771014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
7722f2e5d10SKris Buschelman 
773014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
774014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
77503919abeSBarry Smith 
776ba36735cSStefano Zampini /*S
777ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
778f8a50e2bSBarry Smith 
779ba36735cSStefano Zampini    Level: beginner
780f8a50e2bSBarry Smith 
781db781477SPatrick Sanan .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
782ba36735cSStefano Zampini S*/
783ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess;
784ba36735cSStefano Zampini /*J
785ba36735cSStefano Zampini     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
786ba36735cSStefano Zampini 
787ba36735cSStefano Zampini    Level: beginner
788ba36735cSStefano Zampini 
789db781477SPatrick Sanan .seealso: `KSPGuess`
790ba36735cSStefano Zampini J*/
791ba36735cSStefano Zampini typedef const char* KSPGuessType;
792ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
793ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
7941d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
795ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
796ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
797ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
798ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
799ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
800ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
801ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
8028410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess,PetscReal);
803ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
804ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
805ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
806ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
807ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
808014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
809ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
810ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
811f8a50e2bSBarry Smith 
812470b340bSDmitry Karpeev /*E
813470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
814470b340bSDmitry Karpeev 
815470b340bSDmitry Karpeev     Level: intermediate
816470b340bSDmitry Karpeev 
817db781477SPatrick Sanan .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()`
818470b340bSDmitry Karpeev E*/
8190c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
820470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
821470b340bSDmitry Karpeev 
822864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
823864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
824864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
825864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
826864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
82792f76d53SAlp Dener 
828014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
829014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
830d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
831bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
832aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
833bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
834470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
835470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
8365bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
8375a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
838470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
839470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
8403f22127dSBarry Smith 
84178e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
84278e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
84378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
844864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
845864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
846864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
847864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
848864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
849cd929ea3SAlp Dener 
850cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
851b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
852cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
853cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
854e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
8550ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
856cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
857cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
858cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
859cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
860cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
8612d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
862cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
863cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
864cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
865cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
86692f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
867cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
868cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
869864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
870864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
871cd929ea3SAlp Dener 
872014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
873014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
874014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
875014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
876014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
877fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
87823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
879fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
88023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
88123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
882014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
884fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
885fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
8866c699258SBarry Smith 
88702b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
8888c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
889e29c0834SMatthew G. Knepley                                            void (**)(PetscInt, PetscInt, PetscInt,
890e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
891e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
892191494d9SMatthew G. Knepley                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
893557cf195SMatthew G. Knepley 
8942b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
8952b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
8962eac72dbSBarry Smith #endif
897