xref: /petsc/include/petscksp.h (revision 901ccb91967f7f6c9e6e8ad0b9d61032d623db9c)
1f26ada1bSBarry Smith /*
2f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
3f26ada1bSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCKSP_H
50a835dfdSSatish Balay #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 
1628ce4d24SBarry Smith   Concepts: Krylov methods
1728ce4d24SBarry Smith 
188f6c3df8SBarry Smith         Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
198f6c3df8SBarry Smith        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).
208f6c3df8SBarry Smith 
218f6c3df8SBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy()
2228ce4d24SBarry Smith S*/
23e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
242eac72dbSBarry Smith 
2576bdecfbSBarry Smith /*J
268f6c3df8SBarry Smith     KSPType - String with the name of a PETSc Krylov method.
2728ce4d24SBarry Smith 
2828ce4d24SBarry Smith    Level: beginner
2928ce4d24SBarry Smith 
308f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
3176bdecfbSBarry Smith J*/
3219fd82e9SBarry Smith typedef const char* KSPType;
3382bf6240SBarry Smith #define KSPRICHARDSON "richardson"
346c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3582bf6240SBarry Smith #define KSPCG         "cg"
362c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
372c8fc646SJed Brown #define KSPPIPECG     "pipecg"
38*901ccb91SSiegfried Cools #define KSPPIPECGRR   "pipecgrr"
39df2a969aSvictorle #define   KSPCGNE       "cgne"
40fcae7a14Stmunson #define   KSPNASH       "nash"
4180e17431SMatthew Knepley #define   KSPSTCG       "stcg"
428031f4adStmunson #define   KSPGLTR       "gltr"
43640f4f53SPatrick Sanan #define KSPFCG        "fcg"
4482bf6240SBarry Smith #define KSPGMRES      "gmres"
459dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
469dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
474f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
4861b468f9SJed Brown #define   KSPPGMRES     "pgmres"
4982bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5082bf6240SBarry Smith #define KSPBCGS       "bcgs"
51e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5218ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
53c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
54f0037002Svictorle #define   KSPBCGSL      "bcgsl"
5582bf6240SBarry Smith #define KSPCGS        "cgs"
5682bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
5782bf6240SBarry Smith #define KSPCR         "cr"
582c8fc646SJed Brown #define KSPPIPECR     "pipecr"
5982bf6240SBarry Smith #define KSPLSQR       "lsqr"
6082bf6240SBarry Smith #define KSPPREONLY    "preonly"
6182bf6240SBarry Smith #define KSPQCG        "qcg"
62c9cf9b11SBarry Smith #define KSPBICG       "bicg"
63b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
6401c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
6521356919SSatish Balay #define KSPLCD        "lcd"
661d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
6758ad63f4SBarry Smith #define KSPGCR        "gcr"
68e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
69e4d80e07Szianekhodja #define KSPCGLS       "cgls"
702eac72dbSBarry Smith 
718ba1e511SMatthew Knepley /* Logging support */
72014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
735358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
748ba1e511SMatthew Knepley 
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
7619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
77c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
8423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
8525c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
862eac72dbSBarry Smith 
87140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
88bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
8930de9b25SBarry Smith 
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
93c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
1017d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool );
102c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
104c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
109734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1102a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1112a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
1122eac72dbSBarry Smith 
113d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
114d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
115d4211eb9SBarry Smith 
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
118aabeff55SBarry Smith 
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1254b0e389bSBarry Smith 
126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
129fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
130fa0ddf94SBarry Smith 
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
139b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1410a71e23eSMatthew Knepley 
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1442eac72dbSBarry Smith 
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
14858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
14982fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool);
150ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
15158450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
153d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
154d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1557d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]);
1564b0e389bSBarry Smith 
157640f4f53SPatrick Sanan /*E
158640f4f53SPatrick Sanan 
159bd242c34SBarry Smith   KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG
160640f4f53SPatrick Sanan 
1619662b6c4SBarry Smith   KSP_FCG_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
1629662b6c4SBarry Smith   KSP_FCG_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
163640f4f53SPatrick Sanan 
1642b26979fSBarry Smith    Level: intermediate
1652b26979fSBarry Smith 
166bd242c34SBarry Smith .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
167640f4f53SPatrick Sanan 
168640f4f53SPatrick Sanan E*/
1699662b6c4SBarry Smith typedef enum {KSP_FCG_TRUNC_TYPE_STANDARD,KSP_FCG_TRUNC_TYPE_NOTAY} KSPFCGTruncationType;
170bd242c34SBarry Smith PETSC_EXTERN const char *const KSPFCGTruncationTypes[];
171640f4f53SPatrick Sanan 
172640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
173640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
174640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
175640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
176bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType);
177bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*);
178640f4f53SPatrick Sanan 
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1829f236934SBarry Smith 
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1881d73ed98SBarry Smith 
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
1911d73ed98SBarry Smith 
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
19558ad63f4SBarry Smith 
196b49fd9e1SBarry Smith /*E
197b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
198b49fd9e1SBarry Smith 
199b49fd9e1SBarry Smith    Level: advanced
200b49fd9e1SBarry Smith 
201bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
202e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
203b49fd9e1SBarry Smith 
204b49fd9e1SBarry Smith E*/
20578d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2066a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2071f7e983dSSatish Balay /*MC
2081957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2098c5b8ba0SBarry Smith 
2108c5b8ba0SBarry Smith    Level: advanced
2118c5b8ba0SBarry Smith 
2128c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2138c5b8ba0SBarry Smith 
214bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
215e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2168c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2178c5b8ba0SBarry Smith M*/
2188c5b8ba0SBarry Smith 
2191f7e983dSSatish Balay /*MC
2208c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2218c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2228c5b8ba0SBarry Smith           poor orthogonality.
2238c5b8ba0SBarry Smith 
2248c5b8ba0SBarry Smith    Level: advanced
2258c5b8ba0SBarry Smith 
2268c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2278c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2288c5b8ba0SBarry Smith 
229bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
230e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2318c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2328c5b8ba0SBarry Smith M*/
2338c5b8ba0SBarry Smith 
2341f7e983dSSatish Balay /*MC
2358c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2368c5b8ba0SBarry Smith 
2378c5b8ba0SBarry Smith    Level: advanced
2388c5b8ba0SBarry Smith 
2398c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2408c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2418c5b8ba0SBarry Smith 
2428c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2438c5b8ba0SBarry Smith 
244bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
245e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2468c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2478c5b8ba0SBarry Smith M*/
2488c5b8ba0SBarry Smith 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
25108480c60SBarry Smith 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
255c38d4ed2SBarry Smith 
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
259121fd945SKris Buschelman 
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
263e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
264d9492815SBarry Smith 
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2672eac72dbSBarry Smith 
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
272390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
273390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
275816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
278e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
279e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
280e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
28284cb2905SBarry Smith 
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
285c01c455dSBarry Smith 
28623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
28723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
2927287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
2937287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
2941eb62cbbSBarry Smith 
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
2991f7f0c4fSBarry Smith 
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
30155849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
302685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
3032a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3042a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
30555849f57SBarry Smith 
30655849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3071eb62cbbSBarry Smith 
308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
310db9b2ab1SHong Zhang 
311014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
312014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
31368ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
31483ab6a24SBarry Smith 
31528ce4d24SBarry Smith /*E
3168a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3178a4b9c5cSBarry Smith        test routines.
3188a4b9c5cSBarry Smith 
3198a4b9c5cSBarry Smith    Level: advanced
3208a4b9c5cSBarry Smith 
321a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3221718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
323a3f661c8SBarry Smith 
324af0996ceSBarry Smith    Notes: this must match petsc/finclude/petscksp.h
3258a4b9c5cSBarry Smith 
32694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3271718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3288a4b9c5cSBarry Smith E*/
3299e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3309e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
331014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
332a21b2a99SBarry Smith 
3331f7e983dSSatish Balay /*MC
3349793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3358c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3368c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3378c5b8ba0SBarry Smith 
3388c5b8ba0SBarry Smith    Level: advanced
3398c5b8ba0SBarry Smith 
3401957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3418c5b8ba0SBarry Smith 
342ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3438c5b8ba0SBarry Smith M*/
3448c5b8ba0SBarry Smith 
3451f7e983dSSatish Balay /*MC
3461957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3478c5b8ba0SBarry Smith        convergence test routine.
3488c5b8ba0SBarry Smith 
3498c5b8ba0SBarry Smith    Level: advanced
3508c5b8ba0SBarry Smith 
3519793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3528c5b8ba0SBarry Smith M*/
3538c5b8ba0SBarry Smith 
3541f7e983dSSatish Balay /*MC
355ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3568c5b8ba0SBarry Smith        convergence test routine.
3578c5b8ba0SBarry Smith 
3588c5b8ba0SBarry Smith    Level: advanced
3598c5b8ba0SBarry Smith 
3609793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3618c5b8ba0SBarry Smith M*/
3628c5b8ba0SBarry Smith 
3631f7e983dSSatish Balay /*MC
364ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
365640f4f53SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG
3668c5b8ba0SBarry Smith 
3678c5b8ba0SBarry Smith    Level: advanced
3688c5b8ba0SBarry Smith 
3699793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3708c5b8ba0SBarry Smith M*/
3718c5b8ba0SBarry Smith 
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
3778a4b9c5cSBarry Smith 
3788a4b9c5cSBarry Smith /*E
3791957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
38028ce4d24SBarry Smith 
38128ce4d24SBarry Smith    Level: beginner
38228ce4d24SBarry Smith 
3834d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
38428ce4d24SBarry Smith 
385af0996ceSBarry Smith    Developer notes: this must match petsc/finclude/petscksp.h
3864d0a8057SBarry Smith 
3874d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3884d0a8057SBarry Smith       any of the values here also change them that array of names.
38986c02ca4SBarry Smith 
390c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
39128ce4d24SBarry Smith E*/
392d15094e1SBarry Smith typedef enum {/* converged */
3939ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
3949ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
395d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
396d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
397b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3988031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3998031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
400329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
401af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
402d15094e1SBarry Smith               /* diverged */
403b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
404d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
405d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
406d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
407b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
408b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
409b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4104d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4116aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
412422a814eSBarry Smith               KSP_DIVERGED_PCSETUP_FAILED      = -11,
413d15094e1SBarry Smith 
414d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
415014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
416d15094e1SBarry Smith 
417c838673bSBarry Smith /*MC
418c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
419c838673bSBarry Smith 
420c838673bSBarry Smith    Level: beginner
421c838673bSBarry Smith 
422c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
423c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
424c838673bSBarry Smith        2-norm of the residual for right preconditioning
425c838673bSBarry Smith 
426c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
427c838673bSBarry Smith 
428c838673bSBarry Smith M*/
429c838673bSBarry Smith 
430c838673bSBarry Smith /*MC
431c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
432c838673bSBarry Smith 
433c838673bSBarry Smith    Level: beginner
434c838673bSBarry Smith 
435c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
436c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
437c838673bSBarry Smith        2-norm of the residual for right preconditioning
438c838673bSBarry Smith 
439c838673bSBarry Smith    Level: beginner
440c838673bSBarry Smith 
441c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
442c838673bSBarry Smith 
443c838673bSBarry Smith M*/
444c838673bSBarry Smith 
445c838673bSBarry Smith /*MC
446c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
447c838673bSBarry Smith 
448c838673bSBarry Smith    Level: beginner
449c838673bSBarry Smith 
450c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
451c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
452c838673bSBarry Smith        2-norm of the residual for right preconditioning
453c838673bSBarry Smith 
454c838673bSBarry Smith    Level: beginner
455c838673bSBarry Smith 
456c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
457c838673bSBarry Smith 
458c838673bSBarry Smith M*/
459c838673bSBarry Smith 
460c838673bSBarry Smith /*MC
461c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
462c838673bSBarry Smith       reached
463c838673bSBarry Smith 
464c838673bSBarry Smith    Level: beginner
465c838673bSBarry Smith 
466c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
467c838673bSBarry Smith 
468c838673bSBarry Smith M*/
469c838673bSBarry Smith 
470c838673bSBarry Smith /*MC
4718c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4720059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
4734d0a8057SBarry Smith            test routine is set in KSP.
474c838673bSBarry Smith 
475c838673bSBarry Smith    Level: beginner
476c838673bSBarry Smith 
477c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
478c838673bSBarry Smith 
479c838673bSBarry Smith M*/
480c838673bSBarry Smith 
481c838673bSBarry Smith /*MC
482c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
4833014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
4843014e516SBarry Smith           preconditioner.
485c838673bSBarry Smith 
486c838673bSBarry Smith    Level: beginner
487c838673bSBarry Smith 
488c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
489c838673bSBarry Smith 
490c838673bSBarry Smith M*/
491c838673bSBarry Smith 
492c838673bSBarry Smith /*MC
493c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
494c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
495c838673bSBarry Smith 
496c838673bSBarry Smith    Level: beginner
497c838673bSBarry Smith 
498c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
499c838673bSBarry Smith 
500c838673bSBarry Smith M*/
501c838673bSBarry Smith 
502c838673bSBarry Smith /*MC
503c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
504c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
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
513c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
514c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
515c838673bSBarry Smith         be positive definite
516c838673bSBarry Smith 
517c838673bSBarry Smith    Level: beginner
518c838673bSBarry Smith 
5192401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
520c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
521c838673bSBarry Smith 
522c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
523c838673bSBarry Smith 
524c838673bSBarry Smith M*/
525c838673bSBarry Smith 
526c838673bSBarry Smith /*MC
527c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
528c838673bSBarry Smith         while the KSPSolve() is still running.
529c838673bSBarry Smith 
530c838673bSBarry Smith    Level: beginner
531c838673bSBarry Smith 
532c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
533c838673bSBarry Smith 
534c838673bSBarry Smith M*/
535c838673bSBarry Smith 
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
5388de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
5408de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
5418de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
5428de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
5438de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
5440059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
546abef13c0SSatish Balay 
5478ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
5488ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
5498ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
5508ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
5518ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
5528ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
5538ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
5548ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
5558ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
5568ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
5578ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
5588ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
5598ea1b3e6SJed Brown 
560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
561d4fbbf0eSBarry Smith 
56228ce4d24SBarry Smith /*E
56328ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
56428ce4d24SBarry Smith 
56528ce4d24SBarry Smith    Level: beginner
56628ce4d24SBarry Smith 
56728ce4d24SBarry Smith .seealso: KSPCGSetType()
56828ce4d24SBarry Smith E*/
5699dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5706a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
57128ce4d24SBarry Smith 
572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
5748031f4adStmunson 
575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
578fcae7a14Stmunson 
579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
582e559a7feSLois Curfman McInnes 
583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
5888031f4adStmunson 
589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
5901d6018f0SLisandro Dalcin 
591014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
592014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
5933369ce9aSBarry Smith 
5949804daf3SBarry Smith #include <petscdrawtypes.h>
595d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
596d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
597d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
598d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6002f2e5d10SKris Buschelman 
601014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
60303919abeSBarry Smith 
604f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
605ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
606f8a50e2bSBarry Smith 
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
613f8a50e2bSBarry Smith 
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
617f8a50e2bSBarry Smith 
618470b340bSDmitry Karpeev /*E
619470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
620470b340bSDmitry Karpeev 
621470b340bSDmitry Karpeev     Level: intermediate
622470b340bSDmitry Karpeev 
623470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
624470b340bSDmitry Karpeev E*/
625470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
626470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
627470b340bSDmitry Karpeev 
628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
629014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
630d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
631bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
632aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
633bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
634470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
635470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
6365bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
6375a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
638470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
639470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
6403f22127dSBarry Smith 
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
646fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
64723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
648fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
64923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
65023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
651014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
652014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
653fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
654fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6556c699258SBarry Smith 
65602b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
65702b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec);
6582eac72dbSBarry Smith #endif
659