xref: /petsc/include/petscksp.h (revision 2b26979f2976cb8840d2aed9f03519841c1cab77)
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"
67bbce1389SJed Brown #define KSPSPECEST    "specest"
682eac72dbSBarry Smith 
698ba1e511SMatthew Knepley /* Logging support */
70014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
715358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
728ba1e511SMatthew Knepley 
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
7419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
75c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
8223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
832eac72dbSBarry Smith 
84140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
85014dd563SJed Brown PETSC_EXTERN PetscBool         KSPRegisterAllCalled;
86607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegisterAll(void);
87bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
88607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void);
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 );
101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
103c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
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);
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
149ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
152d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
153d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1544b0e389bSBarry Smith 
155640f4f53SPatrick Sanan /*E
156640f4f53SPatrick Sanan 
157bd242c34SBarry Smith   KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG
158640f4f53SPatrick Sanan 
1599662b6c4SBarry Smith   KSP_FCG_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
1609662b6c4SBarry Smith   KSP_FCG_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
161640f4f53SPatrick Sanan 
162*2b26979fSBarry Smith    Level: intermediate
163*2b26979fSBarry Smith 
164bd242c34SBarry Smith .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
165640f4f53SPatrick Sanan 
166640f4f53SPatrick Sanan E*/
1679662b6c4SBarry Smith typedef enum {KSP_FCG_TRUNC_TYPE_STANDARD,KSP_FCG_TRUNC_TYPE_NOTAY} KSPFCGTruncationType;
168bd242c34SBarry Smith PETSC_EXTERN const char *const KSPFCGTruncationTypes[];
169640f4f53SPatrick Sanan 
170640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
171640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
172640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
173640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
174bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType);
175bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*);
176640f4f53SPatrick Sanan 
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1809f236934SBarry Smith 
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1861d73ed98SBarry Smith 
187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
1891d73ed98SBarry Smith 
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
19358ad63f4SBarry Smith 
194b49fd9e1SBarry Smith /*E
195b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
196b49fd9e1SBarry Smith 
197b49fd9e1SBarry Smith    Level: advanced
198b49fd9e1SBarry Smith 
199bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
200e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
201b49fd9e1SBarry Smith 
202b49fd9e1SBarry Smith E*/
20378d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2046a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2051f7e983dSSatish Balay /*MC
2061957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2078c5b8ba0SBarry Smith 
2088c5b8ba0SBarry Smith    Level: advanced
2098c5b8ba0SBarry Smith 
2108c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2118c5b8ba0SBarry Smith 
212bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
213e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2148c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2158c5b8ba0SBarry Smith M*/
2168c5b8ba0SBarry Smith 
2171f7e983dSSatish Balay /*MC
2188c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2198c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2208c5b8ba0SBarry Smith           poor orthogonality.
2218c5b8ba0SBarry Smith 
2228c5b8ba0SBarry Smith    Level: advanced
2238c5b8ba0SBarry Smith 
2248c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2258c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2268c5b8ba0SBarry Smith 
227bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
228e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2298c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2308c5b8ba0SBarry Smith M*/
2318c5b8ba0SBarry Smith 
2321f7e983dSSatish Balay /*MC
2338c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2348c5b8ba0SBarry Smith 
2358c5b8ba0SBarry Smith    Level: advanced
2368c5b8ba0SBarry Smith 
2378c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2388c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2398c5b8ba0SBarry Smith 
2408c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2418c5b8ba0SBarry Smith 
242bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
243e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2448c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2458c5b8ba0SBarry Smith M*/
2468c5b8ba0SBarry Smith 
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
24908480c60SBarry Smith 
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
253c38d4ed2SBarry Smith 
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
257121fd945SKris Buschelman 
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
261e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
262d9492815SBarry Smith 
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2652eac72dbSBarry Smith 
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
270390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
271390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
273816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
276e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
277e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
278e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
28084cb2905SBarry Smith 
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
283c01c455dSBarry Smith 
28423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
28523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
2907287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
2917287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
2921eb62cbbSBarry Smith 
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
2971f7f0c4fSBarry Smith 
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
29955849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
300ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
3012a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3022a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
30355849f57SBarry Smith 
30455849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3051eb62cbbSBarry Smith 
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
308db9b2ab1SHong Zhang 
309014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
31183ab6a24SBarry Smith 
31228ce4d24SBarry Smith /*E
3138a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3148a4b9c5cSBarry Smith        test routines.
3158a4b9c5cSBarry Smith 
3168a4b9c5cSBarry Smith    Level: advanced
3178a4b9c5cSBarry Smith 
318a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3191718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
320a3f661c8SBarry Smith 
3218fb42493SBarry Smith    Notes: this must match petsc-finclude/petscksp.h
3228a4b9c5cSBarry Smith 
32394b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3241718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3258a4b9c5cSBarry Smith E*/
3269e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3279e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
328014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
329a21b2a99SBarry Smith 
3301f7e983dSSatish Balay /*MC
3319793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3328c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3338c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3348c5b8ba0SBarry Smith 
3358c5b8ba0SBarry Smith    Level: advanced
3368c5b8ba0SBarry Smith 
3371957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3388c5b8ba0SBarry Smith 
339ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3408c5b8ba0SBarry Smith M*/
3418c5b8ba0SBarry Smith 
3421f7e983dSSatish Balay /*MC
3431957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3448c5b8ba0SBarry Smith        convergence test routine.
3458c5b8ba0SBarry Smith 
3468c5b8ba0SBarry Smith    Level: advanced
3478c5b8ba0SBarry Smith 
3489793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3498c5b8ba0SBarry Smith M*/
3508c5b8ba0SBarry Smith 
3511f7e983dSSatish Balay /*MC
352ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3538c5b8ba0SBarry Smith        convergence test routine.
3548c5b8ba0SBarry Smith 
3558c5b8ba0SBarry Smith    Level: advanced
3568c5b8ba0SBarry Smith 
3579793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3588c5b8ba0SBarry Smith M*/
3598c5b8ba0SBarry Smith 
3601f7e983dSSatish Balay /*MC
361ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
362640f4f53SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG
3638c5b8ba0SBarry Smith 
3648c5b8ba0SBarry Smith    Level: advanced
3658c5b8ba0SBarry Smith 
3669793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3678c5b8ba0SBarry Smith M*/
3688c5b8ba0SBarry Smith 
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
3748a4b9c5cSBarry Smith 
3758a4b9c5cSBarry Smith /*E
3761957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
37728ce4d24SBarry Smith 
37828ce4d24SBarry Smith    Level: beginner
37928ce4d24SBarry Smith 
3804d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
38128ce4d24SBarry Smith 
3828fb42493SBarry Smith    Developer notes: this must match petsc-finclude/petscksp.h
3834d0a8057SBarry Smith 
3844d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3854d0a8057SBarry Smith       any of the values here also change them that array of names.
38686c02ca4SBarry Smith 
387c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
38828ce4d24SBarry Smith E*/
389d15094e1SBarry Smith typedef enum {/* converged */
3909ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
3919ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
392d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
393d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
394b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3958031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3968031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
397329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
398af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
399d15094e1SBarry Smith               /* diverged */
400b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
401d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
402d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
403d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
404b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
405b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
406b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4074d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4086aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
409d15094e1SBarry Smith 
410d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
411014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
412d15094e1SBarry Smith 
413c838673bSBarry Smith /*MC
414c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
415c838673bSBarry Smith 
416c838673bSBarry Smith    Level: beginner
417c838673bSBarry Smith 
418c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
419c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
420c838673bSBarry Smith        2-norm of the residual for right preconditioning
421c838673bSBarry Smith 
422c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
423c838673bSBarry Smith 
424c838673bSBarry Smith M*/
425c838673bSBarry Smith 
426c838673bSBarry Smith /*MC
427c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
428c838673bSBarry Smith 
429c838673bSBarry Smith    Level: beginner
430c838673bSBarry Smith 
431c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
432c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
433c838673bSBarry Smith        2-norm of the residual for right preconditioning
434c838673bSBarry Smith 
435c838673bSBarry Smith    Level: beginner
436c838673bSBarry Smith 
437c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
438c838673bSBarry Smith 
439c838673bSBarry Smith M*/
440c838673bSBarry Smith 
441c838673bSBarry Smith /*MC
442c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
443c838673bSBarry Smith 
444c838673bSBarry Smith    Level: beginner
445c838673bSBarry Smith 
446c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
447c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
448c838673bSBarry Smith        2-norm of the residual for right preconditioning
449c838673bSBarry Smith 
450c838673bSBarry Smith    Level: beginner
451c838673bSBarry Smith 
452c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
453c838673bSBarry Smith 
454c838673bSBarry Smith M*/
455c838673bSBarry Smith 
456c838673bSBarry Smith /*MC
457c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
458c838673bSBarry Smith       reached
459c838673bSBarry Smith 
460c838673bSBarry Smith    Level: beginner
461c838673bSBarry Smith 
462c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
463c838673bSBarry Smith 
464c838673bSBarry Smith M*/
465c838673bSBarry Smith 
466c838673bSBarry Smith /*MC
4678c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4680059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
4694d0a8057SBarry Smith            test routine is set in KSP.
470c838673bSBarry Smith 
471c838673bSBarry Smith    Level: beginner
472c838673bSBarry Smith 
473c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
474c838673bSBarry Smith 
475c838673bSBarry Smith M*/
476c838673bSBarry Smith 
477c838673bSBarry Smith /*MC
478c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
4793014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
4803014e516SBarry Smith           preconditioner.
481c838673bSBarry Smith 
482c838673bSBarry Smith    Level: beginner
483c838673bSBarry Smith 
484c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
485c838673bSBarry Smith 
486c838673bSBarry Smith M*/
487c838673bSBarry Smith 
488c838673bSBarry Smith /*MC
489c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
490c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
491c838673bSBarry Smith 
492c838673bSBarry Smith    Level: beginner
493c838673bSBarry Smith 
494c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
495c838673bSBarry Smith 
496c838673bSBarry Smith M*/
497c838673bSBarry Smith 
498c838673bSBarry Smith /*MC
499c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
500c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
501c838673bSBarry Smith 
502c838673bSBarry Smith    Level: beginner
503c838673bSBarry Smith 
504c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
505c838673bSBarry Smith 
506c838673bSBarry Smith M*/
507c838673bSBarry Smith 
508c838673bSBarry Smith /*MC
509c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
510c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
511c838673bSBarry Smith         be positive definite
512c838673bSBarry Smith 
513c838673bSBarry Smith    Level: beginner
514c838673bSBarry Smith 
5152401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
516c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
517c838673bSBarry Smith 
518c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
519c838673bSBarry Smith 
520c838673bSBarry Smith M*/
521c838673bSBarry Smith 
522c838673bSBarry Smith /*MC
523c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
524c838673bSBarry Smith         while the KSPSolve() is still running.
525c838673bSBarry Smith 
526c838673bSBarry Smith    Level: beginner
527c838673bSBarry Smith 
528c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
529c838673bSBarry Smith 
530c838673bSBarry Smith M*/
531c838673bSBarry Smith 
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
5348de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
5368de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
5378de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
5388de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
5398de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
5400059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
542abef13c0SSatish Balay 
5438ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
5448ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
5458ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
5468ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
5478ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
5488ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
5498ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
5508ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
5518ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
5528ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
5538ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
5548ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
5558ea1b3e6SJed Brown 
556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
557d4fbbf0eSBarry Smith 
55828ce4d24SBarry Smith /*E
55928ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
56028ce4d24SBarry Smith 
56128ce4d24SBarry Smith    Level: beginner
56228ce4d24SBarry Smith 
56328ce4d24SBarry Smith .seealso: KSPCGSetType()
56428ce4d24SBarry Smith E*/
5659dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5666a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
56728ce4d24SBarry Smith 
568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
5708031f4adStmunson 
571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
574fcae7a14Stmunson 
575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
578e559a7feSLois Curfman McInnes 
579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
5848031f4adStmunson 
585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
5861d6018f0SLisandro Dalcin 
587014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
588014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
5893369ce9aSBarry Smith 
5909804daf3SBarry Smith #include <petscdrawtypes.h>
591435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscObject**);
592435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,PetscObject*);
593435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscObject**);
594435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscObject**);
595435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,PetscObject*);
596435c5a64SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscObject**);
597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
5982f2e5d10SKris Buschelman 
599014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
600014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
60103919abeSBarry Smith 
602f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
603ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
604f8a50e2bSBarry Smith 
605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
611f8a50e2bSBarry Smith 
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
615f8a50e2bSBarry Smith 
616470b340bSDmitry Karpeev /*E
617470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
618470b340bSDmitry Karpeev 
619470b340bSDmitry Karpeev     Level: intermediate
620470b340bSDmitry Karpeev 
621470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
622470b340bSDmitry Karpeev E*/
623470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
624470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
625470b340bSDmitry Karpeev 
626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
627014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
628d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
629bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
630aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
631bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
632470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
633470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
6345bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
6355a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
636470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
637470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
6383f22127dSBarry Smith 
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
644fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
64523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
646fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
64723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
64823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
649014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
650014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
651fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
652fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6536c699258SBarry Smith 
65402b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
65502b02e71SToby 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);
6562eac72dbSBarry Smith #endif
657