xref: /petsc/include/petscksp.h (revision 734794cf6305c19bac7ec37970b3aba1264ea7c2)
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"
38df2a969aSvictorle #define   KSPCGNE       "cgne"
39fcae7a14Stmunson #define   KSPNASH       "nash"
4080e17431SMatthew Knepley #define   KSPSTCG       "stcg"
418031f4adStmunson #define   KSPGLTR       "gltr"
42640f4f53SPatrick Sanan #define KSPFCG        "fcg"
4382bf6240SBarry Smith #define KSPGMRES      "gmres"
449dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
459dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
464f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
4761b468f9SJed Brown #define   KSPPGMRES     "pgmres"
4882bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
4982bf6240SBarry Smith #define KSPBCGS       "bcgs"
50e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5118ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
52c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
53f0037002Svictorle #define   KSPBCGSL      "bcgsl"
5482bf6240SBarry Smith #define KSPCGS        "cgs"
5582bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
5682bf6240SBarry Smith #define KSPCR         "cr"
572c8fc646SJed Brown #define KSPPIPECR     "pipecr"
5882bf6240SBarry Smith #define KSPLSQR       "lsqr"
5982bf6240SBarry Smith #define KSPPREONLY    "preonly"
6082bf6240SBarry Smith #define KSPQCG        "qcg"
61c9cf9b11SBarry Smith #define KSPBICG       "bicg"
62b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
6301c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
6421356919SSatish Balay #define KSPLCD        "lcd"
651d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
6658ad63f4SBarry Smith #define KSPGCR        "gcr"
672eac72dbSBarry Smith 
688ba1e511SMatthew Knepley /* Logging support */
69014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
705358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
718ba1e511SMatthew Knepley 
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
7319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
74c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
8123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
8225c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
832eac72dbSBarry Smith 
84140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
85bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
8630de9b25SBarry Smith 
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
88014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
90c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
98c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
100c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
105*734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
1082a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1092a7a6963SBarry 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);}
1102eac72dbSBarry Smith 
111d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
112d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
113d4211eb9SBarry Smith 
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
116aabeff55SBarry Smith 
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1234b0e389bSBarry Smith 
124fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
125fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
128fa0ddf94SBarry Smith 
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
134b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
135b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1390a71e23eSMatthew Knepley 
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1422eac72dbSBarry Smith 
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
14658450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
147ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
14858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
150d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
151d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1524b0e389bSBarry Smith 
153640f4f53SPatrick Sanan /*E
154640f4f53SPatrick Sanan 
155bd242c34SBarry Smith   KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG
156640f4f53SPatrick Sanan 
1579662b6c4SBarry Smith   KSP_FCG_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
1589662b6c4SBarry Smith   KSP_FCG_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
159640f4f53SPatrick Sanan 
1602b26979fSBarry Smith    Level: intermediate
1612b26979fSBarry Smith 
162bd242c34SBarry Smith .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
163640f4f53SPatrick Sanan 
164640f4f53SPatrick Sanan E*/
1659662b6c4SBarry Smith typedef enum {KSP_FCG_TRUNC_TYPE_STANDARD,KSP_FCG_TRUNC_TYPE_NOTAY} KSPFCGTruncationType;
166bd242c34SBarry Smith PETSC_EXTERN const char *const KSPFCGTruncationTypes[];
167640f4f53SPatrick Sanan 
168640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
169640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
170640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
171640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
172bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType);
173bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*);
174640f4f53SPatrick Sanan 
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1789f236934SBarry Smith 
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1841d73ed98SBarry Smith 
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
1871d73ed98SBarry Smith 
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
19158ad63f4SBarry Smith 
192b49fd9e1SBarry Smith /*E
193b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
194b49fd9e1SBarry Smith 
195b49fd9e1SBarry Smith    Level: advanced
196b49fd9e1SBarry Smith 
197bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
198e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
199b49fd9e1SBarry Smith 
200b49fd9e1SBarry Smith E*/
20178d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2026a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2031f7e983dSSatish Balay /*MC
2041957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2058c5b8ba0SBarry Smith 
2068c5b8ba0SBarry Smith    Level: advanced
2078c5b8ba0SBarry Smith 
2088c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2098c5b8ba0SBarry Smith 
210bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
211e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2128c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2138c5b8ba0SBarry Smith M*/
2148c5b8ba0SBarry Smith 
2151f7e983dSSatish Balay /*MC
2168c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2178c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2188c5b8ba0SBarry Smith           poor orthogonality.
2198c5b8ba0SBarry Smith 
2208c5b8ba0SBarry Smith    Level: advanced
2218c5b8ba0SBarry Smith 
2228c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2238c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2248c5b8ba0SBarry Smith 
225bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
226e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2278c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2288c5b8ba0SBarry Smith M*/
2298c5b8ba0SBarry Smith 
2301f7e983dSSatish Balay /*MC
2318c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2328c5b8ba0SBarry Smith 
2338c5b8ba0SBarry Smith    Level: advanced
2348c5b8ba0SBarry Smith 
2358c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2368c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2378c5b8ba0SBarry Smith 
2388c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2398c5b8ba0SBarry Smith 
240bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
241e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2428c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2438c5b8ba0SBarry Smith M*/
2448c5b8ba0SBarry Smith 
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
24708480c60SBarry Smith 
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
251c38d4ed2SBarry Smith 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
255121fd945SKris Buschelman 
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
259e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
260d9492815SBarry Smith 
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2632eac72dbSBarry Smith 
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
268390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
269390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
271816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
274e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
275e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
276e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
27884cb2905SBarry Smith 
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
281c01c455dSBarry Smith 
28223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
28323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
2887287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
2897287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
2901eb62cbbSBarry Smith 
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
2951f7f0c4fSBarry Smith 
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
29755849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
298ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
2992a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3002a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
30155849f57SBarry Smith 
30255849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3031eb62cbbSBarry Smith 
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
306db9b2ab1SHong Zhang 
307014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
308014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
30983ab6a24SBarry Smith 
31028ce4d24SBarry Smith /*E
3118a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3128a4b9c5cSBarry Smith        test routines.
3138a4b9c5cSBarry Smith 
3148a4b9c5cSBarry Smith    Level: advanced
3158a4b9c5cSBarry Smith 
316a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3171718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
318a3f661c8SBarry Smith 
3198fb42493SBarry Smith    Notes: this must match petsc-finclude/petscksp.h
3208a4b9c5cSBarry Smith 
32194b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3221718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3238a4b9c5cSBarry Smith E*/
3249e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3259e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
326014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
327a21b2a99SBarry Smith 
3281f7e983dSSatish Balay /*MC
3299793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3308c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3318c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3328c5b8ba0SBarry Smith 
3338c5b8ba0SBarry Smith    Level: advanced
3348c5b8ba0SBarry Smith 
3351957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3368c5b8ba0SBarry Smith 
337ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3388c5b8ba0SBarry Smith M*/
3398c5b8ba0SBarry Smith 
3401f7e983dSSatish Balay /*MC
3411957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3428c5b8ba0SBarry Smith        convergence test routine.
3438c5b8ba0SBarry Smith 
3448c5b8ba0SBarry Smith    Level: advanced
3458c5b8ba0SBarry Smith 
3469793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3478c5b8ba0SBarry Smith M*/
3488c5b8ba0SBarry Smith 
3491f7e983dSSatish Balay /*MC
350ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3518c5b8ba0SBarry Smith        convergence test routine.
3528c5b8ba0SBarry Smith 
3538c5b8ba0SBarry Smith    Level: advanced
3548c5b8ba0SBarry Smith 
3559793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3568c5b8ba0SBarry Smith M*/
3578c5b8ba0SBarry Smith 
3581f7e983dSSatish Balay /*MC
359ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
360640f4f53SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG
3618c5b8ba0SBarry Smith 
3628c5b8ba0SBarry Smith    Level: advanced
3638c5b8ba0SBarry Smith 
3649793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3658c5b8ba0SBarry Smith M*/
3668c5b8ba0SBarry Smith 
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
3728a4b9c5cSBarry Smith 
3738a4b9c5cSBarry Smith /*E
3741957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
37528ce4d24SBarry Smith 
37628ce4d24SBarry Smith    Level: beginner
37728ce4d24SBarry Smith 
3784d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
37928ce4d24SBarry Smith 
3808fb42493SBarry Smith    Developer notes: this must match petsc-finclude/petscksp.h
3814d0a8057SBarry Smith 
3824d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3834d0a8057SBarry Smith       any of the values here also change them that array of names.
38486c02ca4SBarry Smith 
385c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
38628ce4d24SBarry Smith E*/
387d15094e1SBarry Smith typedef enum {/* converged */
3889ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
3899ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
390d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
391d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
392b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3938031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3948031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
395329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
396af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
397d15094e1SBarry Smith               /* diverged */
398b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
399d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
400d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
401d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
402b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
403b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
404b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4054d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4066aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
407d15094e1SBarry Smith 
408d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
409014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
410d15094e1SBarry Smith 
411c838673bSBarry Smith /*MC
412c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
413c838673bSBarry Smith 
414c838673bSBarry Smith    Level: beginner
415c838673bSBarry Smith 
416c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
417c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
418c838673bSBarry Smith        2-norm of the residual for right preconditioning
419c838673bSBarry Smith 
420c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
421c838673bSBarry Smith 
422c838673bSBarry Smith M*/
423c838673bSBarry Smith 
424c838673bSBarry Smith /*MC
425c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
426c838673bSBarry Smith 
427c838673bSBarry Smith    Level: beginner
428c838673bSBarry Smith 
429c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
430c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
431c838673bSBarry Smith        2-norm of the residual for right preconditioning
432c838673bSBarry Smith 
433c838673bSBarry Smith    Level: beginner
434c838673bSBarry Smith 
435c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
436c838673bSBarry Smith 
437c838673bSBarry Smith M*/
438c838673bSBarry Smith 
439c838673bSBarry Smith /*MC
440c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
441c838673bSBarry Smith 
442c838673bSBarry Smith    Level: beginner
443c838673bSBarry Smith 
444c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
445c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
446c838673bSBarry Smith        2-norm of the residual for right preconditioning
447c838673bSBarry Smith 
448c838673bSBarry Smith    Level: beginner
449c838673bSBarry Smith 
450c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
451c838673bSBarry Smith 
452c838673bSBarry Smith M*/
453c838673bSBarry Smith 
454c838673bSBarry Smith /*MC
455c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
456c838673bSBarry Smith       reached
457c838673bSBarry Smith 
458c838673bSBarry Smith    Level: beginner
459c838673bSBarry Smith 
460c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
461c838673bSBarry Smith 
462c838673bSBarry Smith M*/
463c838673bSBarry Smith 
464c838673bSBarry Smith /*MC
4658c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4660059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
4674d0a8057SBarry Smith            test routine is set in KSP.
468c838673bSBarry Smith 
469c838673bSBarry Smith    Level: beginner
470c838673bSBarry Smith 
471c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
472c838673bSBarry Smith 
473c838673bSBarry Smith M*/
474c838673bSBarry Smith 
475c838673bSBarry Smith /*MC
476c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
4773014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
4783014e516SBarry Smith           preconditioner.
479c838673bSBarry Smith 
480c838673bSBarry Smith    Level: beginner
481c838673bSBarry Smith 
482c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
483c838673bSBarry Smith 
484c838673bSBarry Smith M*/
485c838673bSBarry Smith 
486c838673bSBarry Smith /*MC
487c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
488c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
489c838673bSBarry Smith 
490c838673bSBarry Smith    Level: beginner
491c838673bSBarry Smith 
492c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
493c838673bSBarry Smith 
494c838673bSBarry Smith M*/
495c838673bSBarry Smith 
496c838673bSBarry Smith /*MC
497c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
498c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
499c838673bSBarry Smith 
500c838673bSBarry Smith    Level: beginner
501c838673bSBarry Smith 
502c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
503c838673bSBarry Smith 
504c838673bSBarry Smith M*/
505c838673bSBarry Smith 
506c838673bSBarry Smith /*MC
507c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
508c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
509c838673bSBarry Smith         be positive definite
510c838673bSBarry Smith 
511c838673bSBarry Smith    Level: beginner
512c838673bSBarry Smith 
5132401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
514c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
515c838673bSBarry Smith 
516c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
517c838673bSBarry Smith 
518c838673bSBarry Smith M*/
519c838673bSBarry Smith 
520c838673bSBarry Smith /*MC
521c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
522c838673bSBarry Smith         while the KSPSolve() is still running.
523c838673bSBarry Smith 
524c838673bSBarry Smith    Level: beginner
525c838673bSBarry Smith 
526c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
527c838673bSBarry Smith 
528c838673bSBarry Smith M*/
529c838673bSBarry Smith 
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
5328de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
5348de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
5358de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
5368de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
5378de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
5380059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
540abef13c0SSatish Balay 
5418ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
5428ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
5438ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
5448ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
5458ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
5468ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
5478ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
5488ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
5498ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
5508ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
5518ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
5528ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
5538ea1b3e6SJed Brown 
554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
555d4fbbf0eSBarry Smith 
55628ce4d24SBarry Smith /*E
55728ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
55828ce4d24SBarry Smith 
55928ce4d24SBarry Smith    Level: beginner
56028ce4d24SBarry Smith 
56128ce4d24SBarry Smith .seealso: KSPCGSetType()
56228ce4d24SBarry Smith E*/
5639dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5646a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
56528ce4d24SBarry Smith 
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
5688031f4adStmunson 
569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
570014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
572fcae7a14Stmunson 
573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
576e559a7feSLois Curfman McInnes 
577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
5828031f4adStmunson 
583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
5841d6018f0SLisandro Dalcin 
585014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
586014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
5873369ce9aSBarry Smith 
5889804daf3SBarry Smith #include <petscdrawtypes.h>
589435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscObject**);
590435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,PetscObject*);
591435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscObject**);
592bbe26c46SLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(const char[],const char[],int,int,int,int,PetscObject**);
593435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,PetscObject*);
594435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscObject**);
595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
5962f2e5d10SKris Buschelman 
597014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
598014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
59903919abeSBarry Smith 
600f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
601ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
602f8a50e2bSBarry Smith 
603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
609f8a50e2bSBarry Smith 
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
613f8a50e2bSBarry Smith 
614470b340bSDmitry Karpeev /*E
615470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
616470b340bSDmitry Karpeev 
617470b340bSDmitry Karpeev     Level: intermediate
618470b340bSDmitry Karpeev 
619470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
620470b340bSDmitry Karpeev E*/
621470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
622470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
623470b340bSDmitry Karpeev 
624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
625014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
626d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
627bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
628aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
629bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
630470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
631470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
6325bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
6335a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
634470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
635470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
6363f22127dSBarry Smith 
637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
642fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
64323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
644fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
64523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
64623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
647014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
648014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
649fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
650fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6516c699258SBarry Smith 
65202b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
65302b02e71SToby 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);
6542eac72dbSBarry Smith #endif
655