xref: /petsc/include/petscksp.h (revision 2bd2df7ba2215fb8a722e042edc3791533c0b150)
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
62c8e378dSBarry Smith #include <petscpc.h>
72eac72dbSBarry Smith 
8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
91dbb0a54SBarry Smith 
1028ce4d24SBarry Smith /*S
118f6c3df8SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
128f6c3df8SBarry Smith          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
1328ce4d24SBarry Smith 
1428ce4d24SBarry Smith    Level: beginner
1528ce4d24SBarry Smith 
1695452b02SPatrick Sanan         Notes:
17db3b2ab5SMatthew Knepley     When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a
18db3b2ab5SMatthew Knepley        KSPType of KSPPREONLY, meaning that only application of the preconditioner is used as the linear solver.
198f6c3df8SBarry Smith 
20bca11509SBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES
2128ce4d24SBarry Smith S*/
22e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
232eac72dbSBarry Smith 
2476bdecfbSBarry Smith /*J
258f6c3df8SBarry Smith     KSPType - String with the name of a PETSc Krylov method.
2628ce4d24SBarry Smith 
2728ce4d24SBarry Smith    Level: beginner
2828ce4d24SBarry Smith 
298f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
3076bdecfbSBarry Smith J*/
3119fd82e9SBarry Smith typedef const char* KSPType;
3282bf6240SBarry Smith #define KSPRICHARDSON "richardson"
336c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3482bf6240SBarry Smith #define KSPCG         "cg"
352c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
362c8fc646SJed Brown #define KSPPIPECG     "pipecg"
37901ccb91SSiegfried Cools #define KSPPIPECGRR   "pipecgrr"
38265663fdSSiegfried Cools #define KSPPIPELCG     "pipelcg"
39b21a8899STyler Chen #define KSPPIPEPRCG    "pipeprcg"
40325e8391SManasi #define KSPPIPECG2     "pipecg2"
41df2a969aSvictorle #define   KSPCGNE       "cgne"
4205de396fSBarry Smith #define   KSPNASH       "nash"
4305de396fSBarry Smith #define   KSPSTCG       "stcg"
4405de396fSBarry Smith #define   KSPGLTR       "gltr"
4505de396fSBarry Smith #define     KSPCGNASH  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"")  "nash"
4605de396fSBarry Smith #define     KSPCGSTCG  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"")  "stcg"
4705de396fSBarry Smith #define     KSPCGGLTR  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
48640f4f53SPatrick Sanan #define KSPFCG        "fcg"
49390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
5082bf6240SBarry Smith #define KSPGMRES      "gmres"
51483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
529dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
539dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
544f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
5561b468f9SJed Brown #define   KSPPGMRES     "pgmres"
5682bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5782bf6240SBarry Smith #define KSPBCGS       "bcgs"
58e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5918ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
60c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
61f0037002Svictorle #define   KSPBCGSL      "bcgsl"
62f154af2dSSiegfried Cools #define   KSPPIPEBCGS   "pipebcgs"
6382bf6240SBarry Smith #define KSPCGS        "cgs"
6482bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
6582bf6240SBarry Smith #define KSPCR         "cr"
662c8fc646SJed Brown #define KSPPIPECR     "pipecr"
6782bf6240SBarry Smith #define KSPLSQR       "lsqr"
6882bf6240SBarry Smith #define KSPPREONLY    "preonly"
6982bf6240SBarry Smith #define KSPQCG        "qcg"
70c9cf9b11SBarry Smith #define KSPBICG       "bicg"
71b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
7201c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
7321356919SSatish Balay #define KSPLCD        "lcd"
741d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
7558ad63f4SBarry Smith #define KSPGCR        "gcr"
76fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
77e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
78e4d80e07Szianekhodja #define KSPCGLS       "cgls"
79329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
8038cfc46eSPierre Jolivet #define KSPHPDDM      "hpddm"
812eac72dbSBarry Smith 
828ba1e511SMatthew Knepley /* Logging support */
83014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
84ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
855358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
868ba1e511SMatthew Knepley 
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
8819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
89c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
9475f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP,PetscBool);
955ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat);
965ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPSetMatSolveBlockSize(KSP,PetscInt);
975ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPGetMatSolveBlockSize(KSP,PetscInt*);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
99ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
10123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
1027d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
10325c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
104c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
1052eac72dbSBarry Smith 
106140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
107ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
108798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList;
109798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
110798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
111bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
112798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
11330de9b25SBarry Smith 
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
117c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
1237d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
124c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
126c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
131734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1322a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
13325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
1342eac72dbSBarry Smith 
135d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
136d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
137d4211eb9SBarry Smith 
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
140aabeff55SBarry Smith 
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
14594a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*);
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);
14794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*);
14894a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool);
1494b0e389bSBarry Smith 
150fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
151fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
152fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
153fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
154fa0ddf94SBarry Smith 
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
156cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
157014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
159014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
160014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
161285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
162b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
163b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
164b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
165b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
1674a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
168f3b08a26SMatthew G. Knepley /*
169f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
170f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
171f3b08a26SMatthew G. Knepley */
172f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
173f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
174f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
175f3b08a26SMatthew G. Knepley 
1760a71e23eSMatthew Knepley 
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
1792eac72dbSBarry Smith 
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
18358450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
184b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
18558450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
187d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
188d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1897d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
1904b0e389bSBarry Smith 
191640f4f53SPatrick Sanan /*E
192640f4f53SPatrick Sanan 
19306137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
194640f4f53SPatrick Sanan 
19593f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
19693f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
197640f4f53SPatrick Sanan 
1982b26979fSBarry Smith    Level: intermediate
19906137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
200640f4f53SPatrick Sanan 
201640f4f53SPatrick Sanan E*/
20293f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
20306137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
204640f4f53SPatrick Sanan 
205640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
206640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
207640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
208640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
20906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
21006137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
211640f4f53SPatrick Sanan 
212390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
213390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
214390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
215390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
216390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
217390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
218390d8e47SPatrick Sanan 
219fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
220fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
221fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
222fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
223fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
224fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
225fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
226fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
22783f0b094SBarry Smith 
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
231e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal);
2329f236934SBarry Smith 
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2381d73ed98SBarry Smith 
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2411d73ed98SBarry Smith 
242483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
243483d6965SPatrick Sanan 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
24758ad63f4SBarry Smith 
24804d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
24904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
25004d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
251568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
252e82af88dSprj- 
253e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
254e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
2555ea7661aSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
256d55ab951SPierre Jolivet /*E
257d55ab951SPierre Jolivet     KSPHPDDMType - Type of Krylov method used by KSPHPDDM
258d55ab951SPierre Jolivet 
259d55ab951SPierre Jolivet     Level: intermediate
260d55ab951SPierre Jolivet 
261d55ab951SPierre Jolivet     Values:
262d55ab951SPierre Jolivet +   KSP_HPDDM_TYPE_GMRES (default)
263d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGMRES
264d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_CG
265d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BCG
266d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_GCRODR
267d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGCRODR
268d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BFBCG
269d55ab951SPierre Jolivet -   KSP_HPDDM_TYPE_PREONLY
270d55ab951SPierre Jolivet 
271d55ab951SPierre Jolivet .seealso: KSPHPDDM, KSPHPDDMSetType()
272d55ab951SPierre Jolivet E*/
273d55ab951SPierre Jolivet typedef enum {
274d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES = 0,
275d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES = 1,
276d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG = 2,
277d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG = 3,
278d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR = 4,
279d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
280d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG = 6,
281d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
282d55ab951SPierre Jolivet } KSPHPDDMType;
283d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
284d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
285d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
286e82af88dSprj- 
287b49fd9e1SBarry Smith /*E
288b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
289b49fd9e1SBarry Smith 
290b49fd9e1SBarry Smith    Level: advanced
291b49fd9e1SBarry Smith 
292bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
293e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
294b49fd9e1SBarry Smith 
295b49fd9e1SBarry Smith E*/
29678d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2976a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2981f7e983dSSatish Balay /*MC
2991957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
3008c5b8ba0SBarry Smith 
3018c5b8ba0SBarry Smith    Level: advanced
3028c5b8ba0SBarry Smith 
3038c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
3048c5b8ba0SBarry Smith 
305bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
306e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
3078c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3088c5b8ba0SBarry Smith M*/
3098c5b8ba0SBarry Smith 
3101f7e983dSSatish Balay /*MC
3118c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
3128c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
3138c5b8ba0SBarry Smith           poor orthogonality.
3148c5b8ba0SBarry Smith 
3158c5b8ba0SBarry Smith    Level: advanced
3168c5b8ba0SBarry Smith 
3178c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
3188c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
3198c5b8ba0SBarry Smith 
320bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
321e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
3228c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3238c5b8ba0SBarry Smith M*/
3248c5b8ba0SBarry Smith 
3251f7e983dSSatish Balay /*MC
3268c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3278c5b8ba0SBarry Smith 
3288c5b8ba0SBarry Smith    Level: advanced
3298c5b8ba0SBarry Smith 
3308c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
3318c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
3328c5b8ba0SBarry Smith 
3338c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
3348c5b8ba0SBarry Smith 
335bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
336e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
3378c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3388c5b8ba0SBarry Smith M*/
3398c5b8ba0SBarry Smith 
340014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
341014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
34208480c60SBarry Smith 
343014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
345014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
346c38d4ed2SBarry Smith 
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
350121fd945SKris Buschelman 
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
354e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
355d9492815SBarry Smith 
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
35787d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
3592eac72dbSBarry Smith 
360798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*);
361798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *);
362798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
363798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
364798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
365798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
366798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
367798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
368798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
369798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
370798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
371798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
372798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
373798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
374798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
375798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
376798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
377798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
378798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
379798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
380798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
381fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
382798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
383*2bd2df7bSSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);}
384*2bd2df7bSSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);}
385*2bd2df7bSSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);}
386798534f6SMatthew G. Knepley 
387798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
388390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
389390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
390e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
391e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
392e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
39384cb2905SBarry Smith 
394014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
395014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
396c01c455dSBarry Smith 
39723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
39823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
4031eb62cbbSBarry Smith 
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
4081f7f0c4fSBarry Smith 
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
41055849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
411fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
41219a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer);
4131b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
41494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer);
4151b2b9847SBarry Smith 
41619a666eeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);}
4171b2b9847SBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}
41855849f57SBarry Smith 
41955849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
4201eb62cbbSBarry Smith 
421dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
4220e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
423014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
424884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
425798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
426798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
427798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
428db9b2ab1SHong Zhang 
429014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
430014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
43168ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
43283ab6a24SBarry Smith 
43328ce4d24SBarry Smith /*E
4348a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
4358a4b9c5cSBarry Smith        test routines.
4368a4b9c5cSBarry Smith 
4378a4b9c5cSBarry Smith    Level: advanced
4388a4b9c5cSBarry Smith 
439a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
4401718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
441a3f661c8SBarry Smith 
44295452b02SPatrick Sanan    Notes:
44395452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4448a4b9c5cSBarry Smith 
44594b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
4461718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
4478a4b9c5cSBarry Smith E*/
4489e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
4499e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
450014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
451a21b2a99SBarry Smith 
4521f7e983dSSatish Balay /*MC
4539793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
4548c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
4558c5b8ba0SBarry Smith           be based on a norm of a residual etc.
4568c5b8ba0SBarry Smith 
4578c5b8ba0SBarry Smith    Level: advanced
4588c5b8ba0SBarry Smith 
4591957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
4608c5b8ba0SBarry Smith 
461ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
4628c5b8ba0SBarry Smith M*/
4638c5b8ba0SBarry Smith 
4641f7e983dSSatish Balay /*MC
4651957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
4668c5b8ba0SBarry Smith        convergence test routine.
4678c5b8ba0SBarry Smith 
4688c5b8ba0SBarry Smith    Level: advanced
4698c5b8ba0SBarry Smith 
4709793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4718c5b8ba0SBarry Smith M*/
4728c5b8ba0SBarry Smith 
4731f7e983dSSatish Balay /*MC
474ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
4758c5b8ba0SBarry Smith        convergence test routine.
4768c5b8ba0SBarry Smith 
4778c5b8ba0SBarry Smith    Level: advanced
4788c5b8ba0SBarry Smith 
4799793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4808c5b8ba0SBarry Smith M*/
4818c5b8ba0SBarry Smith 
4821f7e983dSSatish Balay /*MC
483ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
484390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
4858c5b8ba0SBarry Smith 
4868c5b8ba0SBarry Smith    Level: advanced
4878c5b8ba0SBarry Smith 
4889793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
4898c5b8ba0SBarry Smith M*/
4908c5b8ba0SBarry Smith 
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4968a4b9c5cSBarry Smith 
497f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
4988a4b9c5cSBarry Smith /*E
4991957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
50028ce4d24SBarry Smith 
50128ce4d24SBarry Smith    Level: beginner
50228ce4d24SBarry Smith 
50395452b02SPatrick Sanan    Notes:
50495452b02SPatrick Sanan     See KSPGetConvergedReason() for explanation of each value
50528ce4d24SBarry Smith 
50695452b02SPatrick Sanan    Developer Notes:
50795452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
5084d0a8057SBarry Smith 
5094d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
5104d0a8057SBarry Smith       any of the values here also change them that array of names.
51186c02ca4SBarry Smith 
5121b2b9847SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView()
51328ce4d24SBarry Smith E*/
514d15094e1SBarry Smith typedef enum {/* converged */
5159ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
5169ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
517d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
518d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
519b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
5208031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
5218031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
522329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
523af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
524d15094e1SBarry Smith               /* diverged */
525b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
526d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
527d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
528d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
529b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
530b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
531b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
5324d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
5336aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
534c0decd05SBarry Smith               KSP_DIVERGED_PC_FAILED           = -11,
535aa4d2078SSatish Balay               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
536d15094e1SBarry Smith 
537d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
538014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
539d15094e1SBarry Smith 
540c838673bSBarry Smith /*MC
541f9fed41fSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
542c838673bSBarry Smith 
543c838673bSBarry Smith    Level: beginner
544c838673bSBarry Smith 
545c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
546c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
547c838673bSBarry Smith        2-norm of the residual for right preconditioning
548c838673bSBarry Smith 
549f9fed41fSBarry Smith    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
550f9fed41fSBarry Smith 
551f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
552c838673bSBarry Smith 
553c838673bSBarry Smith M*/
554c838673bSBarry Smith 
555c838673bSBarry Smith /*MC
556c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
557c838673bSBarry Smith 
558c838673bSBarry Smith    Level: beginner
559c838673bSBarry Smith 
560c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
561c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
562c838673bSBarry Smith        2-norm of the residual for right preconditioning
563c838673bSBarry Smith 
564f9fed41fSBarry Smith    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
565c838673bSBarry Smith 
566f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
567c838673bSBarry Smith 
568c838673bSBarry Smith M*/
569c838673bSBarry Smith 
570c838673bSBarry Smith /*MC
571c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
572c838673bSBarry Smith 
573c838673bSBarry Smith    Level: beginner
574c838673bSBarry Smith 
575c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
576c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
577c838673bSBarry Smith        2-norm of the residual for right preconditioning
578c838673bSBarry Smith 
579c838673bSBarry Smith    Level: beginner
580c838673bSBarry Smith 
581f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
582c838673bSBarry Smith 
583c838673bSBarry Smith M*/
584c838673bSBarry Smith 
585c838673bSBarry Smith /*MC
586c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
587c838673bSBarry Smith       reached
588c838673bSBarry Smith 
589c838673bSBarry Smith    Level: beginner
590c838673bSBarry Smith 
591c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
592c838673bSBarry Smith 
593c838673bSBarry Smith M*/
594c838673bSBarry Smith 
595c838673bSBarry Smith /*MC
5968c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5970059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5984d0a8057SBarry Smith            test routine is set in KSP.
599c838673bSBarry Smith 
600c838673bSBarry Smith    Level: beginner
601c838673bSBarry Smith 
602c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
603c838673bSBarry Smith 
604c838673bSBarry Smith M*/
605c838673bSBarry Smith 
606c838673bSBarry Smith /*MC
607c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
6081de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
6091de96524SPierre Jolivet           preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
6101de96524SPierre Jolivet           are colinear.
611c838673bSBarry Smith 
612c838673bSBarry Smith    Level: beginner
613c838673bSBarry Smith 
614c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
615c838673bSBarry Smith 
616c838673bSBarry Smith M*/
617c838673bSBarry Smith 
618c838673bSBarry Smith /*MC
619c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
620c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
621c838673bSBarry Smith 
622c838673bSBarry Smith    Level: beginner
623c838673bSBarry Smith 
624c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
625c838673bSBarry Smith 
626c838673bSBarry Smith M*/
627c838673bSBarry Smith 
628c838673bSBarry Smith /*MC
629c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
630c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
631c838673bSBarry Smith 
632c838673bSBarry Smith    Level: beginner
633c838673bSBarry Smith 
634c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
635c838673bSBarry Smith 
636c838673bSBarry Smith M*/
637c838673bSBarry Smith 
638c838673bSBarry Smith /*MC
639c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
640c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
641c838673bSBarry Smith         be positive definite
642c838673bSBarry Smith 
643c838673bSBarry Smith    Level: beginner
644c838673bSBarry Smith 
64595452b02SPatrick Sanan      Notes:
64695452b02SPatrick Sanan     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
647c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
648c838673bSBarry Smith 
649c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
650c838673bSBarry Smith 
651c838673bSBarry Smith M*/
652c838673bSBarry Smith 
653c838673bSBarry Smith /*MC
654c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
6559fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
6569fc87aa7SBarry Smith      such as PCFIELDSPLIT.
6579fc87aa7SBarry Smith 
6589fc87aa7SBarry Smith    Level: beginner
6599fc87aa7SBarry Smith 
66095452b02SPatrick Sanan     Notes:
66195452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
6629fc87aa7SBarry Smith 
6639fc87aa7SBarry Smith 
6649fc87aa7SBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
6659fc87aa7SBarry Smith 
6669fc87aa7SBarry Smith M*/
6679fc87aa7SBarry Smith 
6689fc87aa7SBarry Smith /*MC
669c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
670c838673bSBarry Smith         while the KSPSolve() is still running.
671c838673bSBarry Smith 
672c838673bSBarry Smith    Level: beginner
673c838673bSBarry Smith 
674c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
675c838673bSBarry Smith 
676c838673bSBarry Smith M*/
677c838673bSBarry Smith 
678014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
679df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
680df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
681014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
6828de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6833eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6848de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
6858de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
6868de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
6878de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
68854b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
6890059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
690014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
69194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
692abef13c0SSatish Balay 
69325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
6948ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
69525ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
6968ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
69725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
6988ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
69925ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
7008ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
70125ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
7028ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
70325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
7048ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
7058ea1b3e6SJed Brown 
7060bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
70725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
708d4fbbf0eSBarry Smith 
70928ce4d24SBarry Smith /*E
71028ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
71128ce4d24SBarry Smith 
71228ce4d24SBarry Smith    Level: beginner
71328ce4d24SBarry Smith 
71428ce4d24SBarry Smith .seealso: KSPCGSetType()
71528ce4d24SBarry Smith E*/
7169dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
7176a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
71828ce4d24SBarry Smith 
719014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
720014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);
7218031f4adStmunson 
722dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
723dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
724dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
725fcae7a14Stmunson 
72605de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
72705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
72825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
72925ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
7308031f4adStmunson 
731014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
7321d6018f0SLisandro Dalcin 
733014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
734014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
7353369ce9aSBarry Smith 
7369804daf3SBarry Smith #include <petscdrawtypes.h>
737014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
7382f2e5d10SKris Buschelman 
739014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
740014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
74103919abeSBarry Smith 
742ba36735cSStefano Zampini /*S
743ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
744f8a50e2bSBarry Smith 
745ba36735cSStefano Zampini    Level: beginner
746f8a50e2bSBarry Smith 
747ba36735cSStefano Zampini .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
748ba36735cSStefano Zampini S*/
749ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess;
750ba36735cSStefano Zampini /*J
751ba36735cSStefano Zampini     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
752ba36735cSStefano Zampini 
753ba36735cSStefano Zampini    Level: beginner
754ba36735cSStefano Zampini 
755ba36735cSStefano Zampini .seealso: KSPGuess
756ba36735cSStefano Zampini J*/
757ba36735cSStefano Zampini typedef const char* KSPGuessType;
758ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
759ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
7601d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
761ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
762ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
763ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
764ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
765ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
766ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
767ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
768ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
769ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
770ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
771ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
772ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
773014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
774ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
775ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
776f8a50e2bSBarry Smith 
777470b340bSDmitry Karpeev /*E
778470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
779470b340bSDmitry Karpeev 
780470b340bSDmitry Karpeev     Level: intermediate
781470b340bSDmitry Karpeev 
782470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
783470b340bSDmitry Karpeev E*/
7840c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
785470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
786470b340bSDmitry Karpeev 
787864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
788864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
789864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
790864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
791864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
79292f76d53SAlp Dener 
793014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
794014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
795d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
796bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
797aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
798bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
799470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
800470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
8015bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
8025a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
803470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
804470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
8053f22127dSBarry Smith 
80678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
80778e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
80878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
809864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
810864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
811864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
812864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
813864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
814cd929ea3SAlp Dener 
815cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
816b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
817cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
818cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
819e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
8200ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
821cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
822cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
823cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
824cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
825cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
8262d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
827cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
828cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
829cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
830cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
83192f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
832cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
833cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
834864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
835864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
836cd929ea3SAlp Dener 
837014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
838014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
839014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
840014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
841014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
842fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
84323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
844fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
84523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
84623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
847014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
848014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
849fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
850fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
8516c699258SBarry Smith 
85202b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
8538c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
854e29c0834SMatthew G. Knepley                                            void (**)(PetscInt, PetscInt, PetscInt,
855e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
856e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
857191494d9SMatthew G. Knepley                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
858557cf195SMatthew G. Knepley 
859557cf195SMatthew G. Knepley 
860557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *);
861557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal);
8622eac72dbSBarry Smith #endif
863