xref: /petsc/include/petscksp.h (revision 4a99276e2822d1294806a2082917fdc729d00cac)
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:
1795452b02SPatrick Sanan     When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
188f6c3df8SBarry Smith        KSPType of KSPPREONLY (meaning application of the preconditioner is only 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"
39df2a969aSvictorle #define   KSPCGNE       "cgne"
4005de396fSBarry Smith #define   KSPNASH       "nash"
4105de396fSBarry Smith #define   KSPSTCG       "stcg"
4205de396fSBarry Smith #define   KSPGLTR       "gltr"
4305de396fSBarry Smith #define     KSPCGNASH  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"")  "nash"
4405de396fSBarry Smith #define     KSPCGSTCG  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"")  "stcg"
4505de396fSBarry Smith #define     KSPCGGLTR  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
46640f4f53SPatrick Sanan #define KSPFCG        "fcg"
47390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
4882bf6240SBarry Smith #define KSPGMRES      "gmres"
49483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
509dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
519dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
524f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
5361b468f9SJed Brown #define   KSPPGMRES     "pgmres"
5482bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5582bf6240SBarry Smith #define KSPBCGS       "bcgs"
56e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5718ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
58c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
59f0037002Svictorle #define   KSPBCGSL      "bcgsl"
60f154af2dSSiegfried Cools #define   KSPPIPEBCGS   "pipebcgs"
6182bf6240SBarry Smith #define KSPCGS        "cgs"
6282bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
6382bf6240SBarry Smith #define KSPCR         "cr"
642c8fc646SJed Brown #define KSPPIPECR     "pipecr"
6582bf6240SBarry Smith #define KSPLSQR       "lsqr"
6682bf6240SBarry Smith #define KSPPREONLY    "preonly"
6782bf6240SBarry Smith #define KSPQCG        "qcg"
68c9cf9b11SBarry Smith #define KSPBICG       "bicg"
69b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
7001c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
7121356919SSatish Balay #define KSPLCD        "lcd"
721d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
7358ad63f4SBarry Smith #define KSPGCR        "gcr"
74fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
75e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
76e4d80e07Szianekhodja #define KSPCGLS       "cgls"
77329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
782eac72dbSBarry Smith 
798ba1e511SMatthew Knepley /* Logging support */
80014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
81ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
825358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
838ba1e511SMatthew Knepley 
84014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
8519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
86c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
88014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
92ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
9423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
9525c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
96c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
972eac72dbSBarry Smith 
98140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
99ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
100bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
10130de9b25SBarry Smith 
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
105c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
1117d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
112c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
114c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
119734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1202a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
12125ef9dfeSBarry 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);}
1222eac72dbSBarry Smith 
123d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
124d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
125d4211eb9SBarry Smith 
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
128aabeff55SBarry Smith 
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1354b0e389bSBarry Smith 
136fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
137fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
138fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
139fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
140fa0ddf94SBarry Smith 
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
142cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
147285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
148b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
149b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
150b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
151b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
153*4a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
1540a71e23eSMatthew Knepley 
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
1572eac72dbSBarry Smith 
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
16158450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
162b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
16358450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
164014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
165d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
166d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1677d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
1684b0e389bSBarry Smith 
169640f4f53SPatrick Sanan /*E
170640f4f53SPatrick Sanan 
17106137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
172640f4f53SPatrick Sanan 
17393f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
17493f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
175640f4f53SPatrick Sanan 
1762b26979fSBarry Smith    Level: intermediate
17706137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
178640f4f53SPatrick Sanan 
179640f4f53SPatrick Sanan E*/
18093f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
18106137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
182640f4f53SPatrick Sanan 
183640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
184640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
185640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
186640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
18706137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
18806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
189640f4f53SPatrick Sanan 
190390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
191390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
192390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
193390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
194390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
195390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
196390d8e47SPatrick Sanan 
197fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
198fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
199fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
200fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
201fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
202fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
203fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
204fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
20583f0b094SBarry Smith 
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
2099f236934SBarry Smith 
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2151d73ed98SBarry Smith 
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2181d73ed98SBarry Smith 
219483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
220483d6965SPatrick Sanan 
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
22458ad63f4SBarry Smith 
22504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
22604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
22704d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
228568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
229b49fd9e1SBarry Smith /*E
230b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
231b49fd9e1SBarry Smith 
232b49fd9e1SBarry Smith    Level: advanced
233b49fd9e1SBarry Smith 
234bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
235e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
236b49fd9e1SBarry Smith 
237b49fd9e1SBarry Smith E*/
23878d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2396a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2401f7e983dSSatish Balay /*MC
2411957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2428c5b8ba0SBarry Smith 
2438c5b8ba0SBarry Smith    Level: advanced
2448c5b8ba0SBarry Smith 
2458c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2468c5b8ba0SBarry Smith 
247bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
248e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2498c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2508c5b8ba0SBarry Smith M*/
2518c5b8ba0SBarry Smith 
2521f7e983dSSatish Balay /*MC
2538c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2548c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2558c5b8ba0SBarry Smith           poor orthogonality.
2568c5b8ba0SBarry Smith 
2578c5b8ba0SBarry Smith    Level: advanced
2588c5b8ba0SBarry Smith 
2598c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2608c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2618c5b8ba0SBarry Smith 
262bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
263e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2648c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2658c5b8ba0SBarry Smith M*/
2668c5b8ba0SBarry Smith 
2671f7e983dSSatish Balay /*MC
2688c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2698c5b8ba0SBarry Smith 
2708c5b8ba0SBarry Smith    Level: advanced
2718c5b8ba0SBarry Smith 
2728c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2738c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2748c5b8ba0SBarry Smith 
2758c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2768c5b8ba0SBarry Smith 
277bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
278e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2798c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2808c5b8ba0SBarry Smith M*/
2818c5b8ba0SBarry Smith 
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
28408480c60SBarry Smith 
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
288c38d4ed2SBarry Smith 
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
292121fd945SKris Buschelman 
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
296e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
297d9492815SBarry Smith 
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
29987d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
3012eac72dbSBarry Smith 
302fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
303fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
304fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
305fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
306390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
307390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
308fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
309fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
310fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
311fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
312e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
313e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
314e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
315014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
316fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
31784cb2905SBarry Smith 
318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
320c01c455dSBarry Smith 
32123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
32223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
326014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3271eb62cbbSBarry Smith 
328014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
329014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
330014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
331014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
3321f7f0c4fSBarry Smith 
333014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
33455849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
335685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
3362a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3372a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
33855849f57SBarry Smith 
33955849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3401eb62cbbSBarry Smith 
341dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
3420e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
343014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
344884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
345db9b2ab1SHong Zhang 
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
34868ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
34983ab6a24SBarry Smith 
35028ce4d24SBarry Smith /*E
3518a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3528a4b9c5cSBarry Smith        test routines.
3538a4b9c5cSBarry Smith 
3548a4b9c5cSBarry Smith    Level: advanced
3558a4b9c5cSBarry Smith 
356a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3571718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
358a3f661c8SBarry Smith 
35995452b02SPatrick Sanan    Notes:
36095452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
3618a4b9c5cSBarry Smith 
36294b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3631718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3648a4b9c5cSBarry Smith E*/
3659e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3669e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
367014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
368a21b2a99SBarry Smith 
3691f7e983dSSatish Balay /*MC
3709793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3718c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3728c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3738c5b8ba0SBarry Smith 
3748c5b8ba0SBarry Smith    Level: advanced
3758c5b8ba0SBarry Smith 
3761957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3778c5b8ba0SBarry Smith 
378ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3798c5b8ba0SBarry Smith M*/
3808c5b8ba0SBarry Smith 
3811f7e983dSSatish Balay /*MC
3821957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3838c5b8ba0SBarry Smith        convergence test routine.
3848c5b8ba0SBarry Smith 
3858c5b8ba0SBarry Smith    Level: advanced
3868c5b8ba0SBarry Smith 
3879793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3888c5b8ba0SBarry Smith M*/
3898c5b8ba0SBarry Smith 
3901f7e983dSSatish Balay /*MC
391ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3928c5b8ba0SBarry Smith        convergence test routine.
3938c5b8ba0SBarry Smith 
3948c5b8ba0SBarry Smith    Level: advanced
3958c5b8ba0SBarry Smith 
3969793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3978c5b8ba0SBarry Smith M*/
3988c5b8ba0SBarry Smith 
3991f7e983dSSatish Balay /*MC
400ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
401390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
4028c5b8ba0SBarry Smith 
4038c5b8ba0SBarry Smith    Level: advanced
4048c5b8ba0SBarry Smith 
4059793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
4068c5b8ba0SBarry Smith M*/
4078c5b8ba0SBarry Smith 
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4138a4b9c5cSBarry Smith 
4148a4b9c5cSBarry Smith /*E
4151957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
41628ce4d24SBarry Smith 
41728ce4d24SBarry Smith    Level: beginner
41828ce4d24SBarry Smith 
41995452b02SPatrick Sanan    Notes:
42095452b02SPatrick Sanan     See KSPGetConvergedReason() for explanation of each value
42128ce4d24SBarry Smith 
42295452b02SPatrick Sanan    Developer Notes:
42395452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4244d0a8057SBarry Smith 
4254d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4264d0a8057SBarry Smith       any of the values here also change them that array of names.
42786c02ca4SBarry Smith 
428c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
42928ce4d24SBarry Smith E*/
43025ef9dfeSBarry Smith #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
431d15094e1SBarry Smith typedef enum {/* converged */
4329ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4339ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
434d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
435d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
436b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4378031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4388031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
439329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
440af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
441d15094e1SBarry Smith               /* diverged */
442b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
443d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
444d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
445d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
446b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
447b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
448b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4494d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4506aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
451c0decd05SBarry Smith               KSP_DIVERGED_PC_FAILED           = -11,
452aa4d2078SSatish Balay               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
453d15094e1SBarry Smith 
454d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
455014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
456d15094e1SBarry Smith 
457c838673bSBarry Smith /*MC
458f9fed41fSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
459c838673bSBarry Smith 
460c838673bSBarry Smith    Level: beginner
461c838673bSBarry Smith 
462c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
463c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
464c838673bSBarry Smith        2-norm of the residual for right preconditioning
465c838673bSBarry Smith 
466f9fed41fSBarry Smith    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
467f9fed41fSBarry Smith 
468f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
469c838673bSBarry Smith 
470c838673bSBarry Smith M*/
471c838673bSBarry Smith 
472c838673bSBarry Smith /*MC
473c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
474c838673bSBarry Smith 
475c838673bSBarry Smith    Level: beginner
476c838673bSBarry Smith 
477c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
478c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
479c838673bSBarry Smith        2-norm of the residual for right preconditioning
480c838673bSBarry Smith 
481f9fed41fSBarry Smith    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
482c838673bSBarry Smith 
483f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
484c838673bSBarry Smith 
485c838673bSBarry Smith M*/
486c838673bSBarry Smith 
487c838673bSBarry Smith /*MC
488c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
489c838673bSBarry Smith 
490c838673bSBarry Smith    Level: beginner
491c838673bSBarry Smith 
492c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
493c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
494c838673bSBarry Smith        2-norm of the residual for right preconditioning
495c838673bSBarry Smith 
496c838673bSBarry Smith    Level: beginner
497c838673bSBarry Smith 
498f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
499c838673bSBarry Smith 
500c838673bSBarry Smith M*/
501c838673bSBarry Smith 
502c838673bSBarry Smith /*MC
503c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
504c838673bSBarry Smith       reached
505c838673bSBarry Smith 
506c838673bSBarry Smith    Level: beginner
507c838673bSBarry Smith 
508c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
509c838673bSBarry Smith 
510c838673bSBarry Smith M*/
511c838673bSBarry Smith 
512c838673bSBarry Smith /*MC
5138c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5140059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5154d0a8057SBarry Smith            test routine is set in KSP.
516c838673bSBarry Smith 
517c838673bSBarry Smith    Level: beginner
518c838673bSBarry Smith 
519c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
520c838673bSBarry Smith 
521c838673bSBarry Smith M*/
522c838673bSBarry Smith 
523c838673bSBarry Smith /*MC
524c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5253014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
5263014e516SBarry Smith           preconditioner.
527c838673bSBarry Smith 
528c838673bSBarry Smith    Level: beginner
529c838673bSBarry Smith 
530c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
531c838673bSBarry Smith 
532c838673bSBarry Smith M*/
533c838673bSBarry Smith 
534c838673bSBarry Smith /*MC
535c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
536c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
537c838673bSBarry Smith 
538c838673bSBarry Smith    Level: beginner
539c838673bSBarry Smith 
540c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
541c838673bSBarry Smith 
542c838673bSBarry Smith M*/
543c838673bSBarry Smith 
544c838673bSBarry Smith /*MC
545c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
546c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
547c838673bSBarry Smith 
548c838673bSBarry Smith    Level: beginner
549c838673bSBarry Smith 
550c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
551c838673bSBarry Smith 
552c838673bSBarry Smith M*/
553c838673bSBarry Smith 
554c838673bSBarry Smith /*MC
555c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
556c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
557c838673bSBarry Smith         be positive definite
558c838673bSBarry Smith 
559c838673bSBarry Smith    Level: beginner
560c838673bSBarry Smith 
56195452b02SPatrick Sanan      Notes:
56295452b02SPatrick Sanan     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
563c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
564c838673bSBarry Smith 
565c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
566c838673bSBarry Smith 
567c838673bSBarry Smith M*/
568c838673bSBarry Smith 
569c838673bSBarry Smith /*MC
570c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
5719fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
5729fc87aa7SBarry Smith      such as PCFIELDSPLIT.
5739fc87aa7SBarry Smith 
5749fc87aa7SBarry Smith    Level: beginner
5759fc87aa7SBarry Smith 
57695452b02SPatrick Sanan     Notes:
57795452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
5789fc87aa7SBarry Smith 
5799fc87aa7SBarry Smith 
5809fc87aa7SBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
5819fc87aa7SBarry Smith 
5829fc87aa7SBarry Smith M*/
5839fc87aa7SBarry Smith 
5849fc87aa7SBarry Smith /*MC
585c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
586c838673bSBarry Smith         while the KSPSolve() is still running.
587c838673bSBarry Smith 
588c838673bSBarry Smith    Level: beginner
589c838673bSBarry Smith 
590c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
591c838673bSBarry Smith 
592c838673bSBarry Smith M*/
593c838673bSBarry Smith 
594014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
595df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
596df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
5988de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
5993eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6008de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
6018de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
6028de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
6038de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
6040059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
606abef13c0SSatish Balay 
60725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
6088ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
60925ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
6108ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
61125ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
6128ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
61325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
6148ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
61525ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
6168ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
61725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
6188ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
6198ea1b3e6SJed Brown 
6200bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
62125ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
622d4fbbf0eSBarry Smith 
62328ce4d24SBarry Smith /*E
62428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
62528ce4d24SBarry Smith 
62628ce4d24SBarry Smith    Level: beginner
62728ce4d24SBarry Smith 
62828ce4d24SBarry Smith .seealso: KSPCGSetType()
62928ce4d24SBarry Smith E*/
6309dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
6316a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
63228ce4d24SBarry Smith 
633014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
6358031f4adStmunson 
636dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
637dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
638dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
639fcae7a14Stmunson 
64005de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
64105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
64225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
64325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
6448031f4adStmunson 
645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6461d6018f0SLisandro Dalcin 
647014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
648014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
6493369ce9aSBarry Smith 
6509804daf3SBarry Smith #include <petscdrawtypes.h>
651d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
652d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
653d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
654d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
655014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6562f2e5d10SKris Buschelman 
657014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
658014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
65903919abeSBarry Smith 
660ba36735cSStefano Zampini /*S
661ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
662f8a50e2bSBarry Smith 
663ba36735cSStefano Zampini    Level: beginner
664f8a50e2bSBarry Smith 
665ba36735cSStefano Zampini .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
666ba36735cSStefano Zampini S*/
667ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess;
668ba36735cSStefano Zampini /*J
669ba36735cSStefano Zampini     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
670ba36735cSStefano Zampini 
671ba36735cSStefano Zampini    Level: beginner
672ba36735cSStefano Zampini 
673ba36735cSStefano Zampini .seealso: KSPGuess
674ba36735cSStefano Zampini J*/
675ba36735cSStefano Zampini typedef const char* KSPGuessType;
676ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
677ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
6781d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
679ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
680ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
681ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
682ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
683ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
684ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
685ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
686ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
687ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
688ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
689ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
690ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
691014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
692ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
693ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
694f8a50e2bSBarry Smith 
695470b340bSDmitry Karpeev /*E
696470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
697470b340bSDmitry Karpeev 
698470b340bSDmitry Karpeev     Level: intermediate
699470b340bSDmitry Karpeev 
700470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
701470b340bSDmitry Karpeev E*/
7020c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
703470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
704470b340bSDmitry Karpeev 
705014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
706014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
707d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
708bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
709aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
710bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
711470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
712470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
7135bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
7145a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
715470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
716470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
7173f22127dSBarry Smith 
71878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
71978e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
72078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
72178e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
72278e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
72378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
72418ad3407SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
72550b47da0SAdam Denchfield PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
726cd929ea3SAlp Dener 
727cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
728b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
729cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
730cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
731e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
7320ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
733cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
734cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
735cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
736cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
737cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
7382d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
739cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
740cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
741cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
742cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
743cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
744cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
745d5ae2380SAlp Dener PETSC_EXTERN PetscErrorCode MatSymBrdnSetDelta(Mat, PetscScalar);
746cd929ea3SAlp Dener 
747014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
748014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
749014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
750014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
751014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
752fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
75323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
754fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
75523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
75623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
757014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
758014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
759fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
760fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
7616c699258SBarry Smith 
76202b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
7638c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
764e29c0834SMatthew G. Knepley                                            void (**)(PetscInt, PetscInt, PetscInt,
765e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
766e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
767191494d9SMatthew G. Knepley                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
7682eac72dbSBarry Smith #endif
769