xref: /petsc/include/petscksp.h (revision 94a4d4eb42d35d79d8a3878b697c5b4398d9cd44)
1f26ada1bSBarry Smith /*
2f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
3f26ada1bSBarry Smith */
426bd1501SBarry Smith #ifndef PETSCKSP_H
526bd1501SBarry Smith #define PETSCKSP_H
62c8e378dSBarry Smith #include <petscpc.h>
72eac72dbSBarry Smith 
8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
91dbb0a54SBarry Smith 
1028ce4d24SBarry Smith /*S
118f6c3df8SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
128f6c3df8SBarry Smith          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
1328ce4d24SBarry Smith 
1428ce4d24SBarry Smith    Level: beginner
1528ce4d24SBarry Smith 
1695452b02SPatrick Sanan         Notes:
1795452b02SPatrick Sanan     When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
188f6c3df8SBarry Smith        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).
198f6c3df8SBarry Smith 
20bca11509SBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES
2128ce4d24SBarry Smith S*/
22e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
232eac72dbSBarry Smith 
2476bdecfbSBarry Smith /*J
258f6c3df8SBarry Smith     KSPType - String with the name of a PETSc Krylov method.
2628ce4d24SBarry Smith 
2728ce4d24SBarry Smith    Level: beginner
2828ce4d24SBarry Smith 
298f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
3076bdecfbSBarry Smith J*/
3119fd82e9SBarry Smith typedef const char* KSPType;
3282bf6240SBarry Smith #define KSPRICHARDSON "richardson"
336c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3482bf6240SBarry Smith #define KSPCG         "cg"
352c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
362c8fc646SJed Brown #define KSPPIPECG     "pipecg"
37901ccb91SSiegfried Cools #define KSPPIPECGRR   "pipecgrr"
38265663fdSSiegfried Cools #define KSPPIPELCG     "pipelcg"
39b21a8899STyler Chen #define KSPPIPEPRCG    "pipeprcg"
40325e8391SManasi #define KSPPIPECG2     "pipecg2"
41df2a969aSvictorle #define   KSPCGNE       "cgne"
4205de396fSBarry Smith #define   KSPNASH       "nash"
4305de396fSBarry Smith #define   KSPSTCG       "stcg"
4405de396fSBarry Smith #define   KSPGLTR       "gltr"
4505de396fSBarry Smith #define     KSPCGNASH  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"")  "nash"
4605de396fSBarry Smith #define     KSPCGSTCG  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"")  "stcg"
4705de396fSBarry Smith #define     KSPCGGLTR  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
48640f4f53SPatrick Sanan #define KSPFCG        "fcg"
49390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
5082bf6240SBarry Smith #define KSPGMRES      "gmres"
51483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
529dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
539dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
544f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
5561b468f9SJed Brown #define   KSPPGMRES     "pgmres"
5682bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5782bf6240SBarry Smith #define KSPBCGS       "bcgs"
58e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5918ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
60c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
61f0037002Svictorle #define   KSPBCGSL      "bcgsl"
62f154af2dSSiegfried Cools #define   KSPPIPEBCGS   "pipebcgs"
6382bf6240SBarry Smith #define KSPCGS        "cgs"
6482bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
6582bf6240SBarry Smith #define KSPCR         "cr"
662c8fc646SJed Brown #define KSPPIPECR     "pipecr"
6782bf6240SBarry Smith #define KSPLSQR       "lsqr"
6882bf6240SBarry Smith #define KSPPREONLY    "preonly"
6982bf6240SBarry Smith #define KSPQCG        "qcg"
70c9cf9b11SBarry Smith #define KSPBICG       "bicg"
71b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
7201c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
7321356919SSatish Balay #define KSPLCD        "lcd"
741d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
7558ad63f4SBarry Smith #define KSPGCR        "gcr"
76fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
77e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
78e4d80e07Szianekhodja #define KSPCGLS       "cgls"
79329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
8038cfc46eSPierre Jolivet #define KSPHPDDM      "hpddm"
812eac72dbSBarry Smith 
828ba1e511SMatthew Knepley /* Logging support */
83014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
84ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
855358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
868ba1e511SMatthew Knepley 
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
8819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
89c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
945ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat);
955ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPSetMatSolveBlockSize(KSP,PetscInt);
965ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPGetMatSolveBlockSize(KSP,PetscInt*);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
98ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
10023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
1017d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
10225c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
103c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
1042eac72dbSBarry Smith 
105140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
106ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
107bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
10830de9b25SBarry Smith 
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
112c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
1187d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
119c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
121c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
126734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1272a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
12825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
1292eac72dbSBarry Smith 
130d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
131d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
132d4211eb9SBarry Smith 
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
135aabeff55SBarry Smith 
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
140*94a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);
142*94a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*);
143*94a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool);
1444b0e389bSBarry Smith 
145fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
146fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
147fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
148fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
149fa0ddf94SBarry Smith 
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
151cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
153014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
156285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
157b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
158b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
159b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
160b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
161014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
1624a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
163f3b08a26SMatthew G. Knepley /*
164f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
165f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
166f3b08a26SMatthew G. Knepley */
167f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
168f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
169f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
170f3b08a26SMatthew G. Knepley 
1710a71e23eSMatthew Knepley 
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
1742eac72dbSBarry Smith 
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
17858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
179b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
18058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
182d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
183d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1847d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
1854b0e389bSBarry Smith 
186640f4f53SPatrick Sanan /*E
187640f4f53SPatrick Sanan 
18806137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
189640f4f53SPatrick Sanan 
19093f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
19193f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
192640f4f53SPatrick Sanan 
1932b26979fSBarry Smith    Level: intermediate
19406137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
195640f4f53SPatrick Sanan 
196640f4f53SPatrick Sanan E*/
19793f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
19806137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
199640f4f53SPatrick Sanan 
200640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
201640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
202640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
203640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
20406137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
20506137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
206640f4f53SPatrick Sanan 
207390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
208390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
209390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
210390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
211390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
212390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
213390d8e47SPatrick Sanan 
214fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
215fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
216fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
217fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
218fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
219fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
220fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
221fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
22283f0b094SBarry Smith 
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
226e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal);
2279f236934SBarry Smith 
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2331d73ed98SBarry Smith 
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2361d73ed98SBarry Smith 
237483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
238483d6965SPatrick Sanan 
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
24258ad63f4SBarry Smith 
24304d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
24404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
24504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
246568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
247e82af88dSprj- 
248e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
249e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
2505ea7661aSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
251d55ab951SPierre Jolivet /*E
252d55ab951SPierre Jolivet     KSPHPDDMType - Type of Krylov method used by KSPHPDDM
253d55ab951SPierre Jolivet 
254d55ab951SPierre Jolivet     Level: intermediate
255d55ab951SPierre Jolivet 
256d55ab951SPierre Jolivet     Values:
257d55ab951SPierre Jolivet +   KSP_HPDDM_TYPE_GMRES (default)
258d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGMRES
259d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_CG
260d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BCG
261d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_GCRODR
262d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGCRODR
263d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BFBCG
264d55ab951SPierre Jolivet -   KSP_HPDDM_TYPE_PREONLY
265d55ab951SPierre Jolivet 
266d55ab951SPierre Jolivet .seealso: KSPHPDDM, KSPHPDDMSetType()
267d55ab951SPierre Jolivet E*/
268d55ab951SPierre Jolivet typedef enum {
269d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES = 0,
270d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES = 1,
271d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG = 2,
272d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG = 3,
273d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR = 4,
274d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
275d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG = 6,
276d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
277d55ab951SPierre Jolivet } KSPHPDDMType;
278d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
279d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
280d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
281e82af88dSprj- 
282b49fd9e1SBarry Smith /*E
283b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
284b49fd9e1SBarry Smith 
285b49fd9e1SBarry Smith    Level: advanced
286b49fd9e1SBarry Smith 
287bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
288e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
289b49fd9e1SBarry Smith 
290b49fd9e1SBarry Smith E*/
29178d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2926a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2931f7e983dSSatish Balay /*MC
2941957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2958c5b8ba0SBarry Smith 
2968c5b8ba0SBarry Smith    Level: advanced
2978c5b8ba0SBarry Smith 
2988c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2998c5b8ba0SBarry Smith 
300bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
301e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
3028c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3038c5b8ba0SBarry Smith M*/
3048c5b8ba0SBarry Smith 
3051f7e983dSSatish Balay /*MC
3068c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
3078c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
3088c5b8ba0SBarry Smith           poor orthogonality.
3098c5b8ba0SBarry Smith 
3108c5b8ba0SBarry Smith    Level: advanced
3118c5b8ba0SBarry Smith 
3128c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
3138c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
3148c5b8ba0SBarry Smith 
315bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
316e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
3178c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3188c5b8ba0SBarry Smith M*/
3198c5b8ba0SBarry Smith 
3201f7e983dSSatish Balay /*MC
3218c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3228c5b8ba0SBarry Smith 
3238c5b8ba0SBarry Smith    Level: advanced
3248c5b8ba0SBarry Smith 
3258c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
3268c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
3278c5b8ba0SBarry Smith 
3288c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
3298c5b8ba0SBarry Smith 
330bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
331e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
3328c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3338c5b8ba0SBarry Smith M*/
3348c5b8ba0SBarry Smith 
335014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
33708480c60SBarry Smith 
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
340014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
341c38d4ed2SBarry Smith 
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
343014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
345121fd945SKris Buschelman 
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
349e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
350d9492815SBarry Smith 
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
35287d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
3542eac72dbSBarry Smith 
355fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
356fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
357fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
358fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
359390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
360390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
361fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
362fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
363fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
364fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
365e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
366e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
367e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
369fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
37084cb2905SBarry Smith 
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
373c01c455dSBarry Smith 
37423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
37523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3801eb62cbbSBarry Smith 
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
3851f7f0c4fSBarry Smith 
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
38755849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
388fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
38919a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer);
3901b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
391*94a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer);
3921b2b9847SBarry Smith 
39319a666eeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);}
3941b2b9847SBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}
39555849f57SBarry Smith 
39655849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3971eb62cbbSBarry Smith 
398dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
3990e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
401884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
402db9b2ab1SHong Zhang 
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
40568ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
40683ab6a24SBarry Smith 
40728ce4d24SBarry Smith /*E
4088a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
4098a4b9c5cSBarry Smith        test routines.
4108a4b9c5cSBarry Smith 
4118a4b9c5cSBarry Smith    Level: advanced
4128a4b9c5cSBarry Smith 
413a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
4141718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
415a3f661c8SBarry Smith 
41695452b02SPatrick Sanan    Notes:
41795452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4188a4b9c5cSBarry Smith 
41994b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
4201718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
4218a4b9c5cSBarry Smith E*/
4229e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
4239e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
424014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
425a21b2a99SBarry Smith 
4261f7e983dSSatish Balay /*MC
4279793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
4288c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
4298c5b8ba0SBarry Smith           be based on a norm of a residual etc.
4308c5b8ba0SBarry Smith 
4318c5b8ba0SBarry Smith    Level: advanced
4328c5b8ba0SBarry Smith 
4331957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
4348c5b8ba0SBarry Smith 
435ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
4368c5b8ba0SBarry Smith M*/
4378c5b8ba0SBarry Smith 
4381f7e983dSSatish Balay /*MC
4391957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
4408c5b8ba0SBarry Smith        convergence test routine.
4418c5b8ba0SBarry Smith 
4428c5b8ba0SBarry Smith    Level: advanced
4438c5b8ba0SBarry Smith 
4449793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4458c5b8ba0SBarry Smith M*/
4468c5b8ba0SBarry Smith 
4471f7e983dSSatish Balay /*MC
448ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
4498c5b8ba0SBarry Smith        convergence test routine.
4508c5b8ba0SBarry Smith 
4518c5b8ba0SBarry Smith    Level: advanced
4528c5b8ba0SBarry Smith 
4539793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4548c5b8ba0SBarry Smith M*/
4558c5b8ba0SBarry Smith 
4561f7e983dSSatish Balay /*MC
457ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
458390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
4598c5b8ba0SBarry Smith 
4608c5b8ba0SBarry Smith    Level: advanced
4618c5b8ba0SBarry Smith 
4629793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
4638c5b8ba0SBarry Smith M*/
4648c5b8ba0SBarry Smith 
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4708a4b9c5cSBarry Smith 
471f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
4728a4b9c5cSBarry Smith /*E
4731957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
47428ce4d24SBarry Smith 
47528ce4d24SBarry Smith    Level: beginner
47628ce4d24SBarry Smith 
47795452b02SPatrick Sanan    Notes:
47895452b02SPatrick Sanan     See KSPGetConvergedReason() for explanation of each value
47928ce4d24SBarry Smith 
48095452b02SPatrick Sanan    Developer Notes:
48195452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4824d0a8057SBarry Smith 
4834d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4844d0a8057SBarry Smith       any of the values here also change them that array of names.
48586c02ca4SBarry Smith 
4861b2b9847SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView()
48728ce4d24SBarry Smith E*/
488d15094e1SBarry Smith typedef enum {/* converged */
4899ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4909ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
491d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
492d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
493b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4948031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4958031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
496329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
497af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
498d15094e1SBarry Smith               /* diverged */
499b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
500d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
501d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
502d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
503b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
504b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
505b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
5064d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
5076aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
508c0decd05SBarry Smith               KSP_DIVERGED_PC_FAILED           = -11,
509aa4d2078SSatish Balay               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
510d15094e1SBarry Smith 
511d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
512014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
513d15094e1SBarry Smith 
514c838673bSBarry Smith /*MC
515f9fed41fSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
516c838673bSBarry Smith 
517c838673bSBarry Smith    Level: beginner
518c838673bSBarry Smith 
519c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
520c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
521c838673bSBarry Smith        2-norm of the residual for right preconditioning
522c838673bSBarry Smith 
523f9fed41fSBarry Smith    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
524f9fed41fSBarry Smith 
525f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
526c838673bSBarry Smith 
527c838673bSBarry Smith M*/
528c838673bSBarry Smith 
529c838673bSBarry Smith /*MC
530c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
531c838673bSBarry Smith 
532c838673bSBarry Smith    Level: beginner
533c838673bSBarry Smith 
534c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
535c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
536c838673bSBarry Smith        2-norm of the residual for right preconditioning
537c838673bSBarry Smith 
538f9fed41fSBarry Smith    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
539c838673bSBarry Smith 
540f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
541c838673bSBarry Smith 
542c838673bSBarry Smith M*/
543c838673bSBarry Smith 
544c838673bSBarry Smith /*MC
545c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
546c838673bSBarry Smith 
547c838673bSBarry Smith    Level: beginner
548c838673bSBarry Smith 
549c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
550c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
551c838673bSBarry Smith        2-norm of the residual for right preconditioning
552c838673bSBarry Smith 
553c838673bSBarry Smith    Level: beginner
554c838673bSBarry Smith 
555f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
556c838673bSBarry Smith 
557c838673bSBarry Smith M*/
558c838673bSBarry Smith 
559c838673bSBarry Smith /*MC
560c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
561c838673bSBarry Smith       reached
562c838673bSBarry Smith 
563c838673bSBarry Smith    Level: beginner
564c838673bSBarry Smith 
565c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
566c838673bSBarry Smith 
567c838673bSBarry Smith M*/
568c838673bSBarry Smith 
569c838673bSBarry Smith /*MC
5708c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5710059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5724d0a8057SBarry Smith            test routine is set in KSP.
573c838673bSBarry Smith 
574c838673bSBarry Smith    Level: beginner
575c838673bSBarry Smith 
576c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
577c838673bSBarry Smith 
578c838673bSBarry Smith M*/
579c838673bSBarry Smith 
580c838673bSBarry Smith /*MC
581c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5821de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
5831de96524SPierre Jolivet           preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
5841de96524SPierre Jolivet           are colinear.
585c838673bSBarry Smith 
586c838673bSBarry Smith    Level: beginner
587c838673bSBarry Smith 
588c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
589c838673bSBarry Smith 
590c838673bSBarry Smith M*/
591c838673bSBarry Smith 
592c838673bSBarry Smith /*MC
593c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
594c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
595c838673bSBarry Smith 
596c838673bSBarry Smith    Level: beginner
597c838673bSBarry Smith 
598c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
599c838673bSBarry Smith 
600c838673bSBarry Smith M*/
601c838673bSBarry Smith 
602c838673bSBarry Smith /*MC
603c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
604c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
605c838673bSBarry Smith 
606c838673bSBarry Smith    Level: beginner
607c838673bSBarry Smith 
608c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
609c838673bSBarry Smith 
610c838673bSBarry Smith M*/
611c838673bSBarry Smith 
612c838673bSBarry Smith /*MC
613c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
614c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
615c838673bSBarry Smith         be positive definite
616c838673bSBarry Smith 
617c838673bSBarry Smith    Level: beginner
618c838673bSBarry Smith 
61995452b02SPatrick Sanan      Notes:
62095452b02SPatrick Sanan     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
621c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
622c838673bSBarry Smith 
623c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
624c838673bSBarry Smith 
625c838673bSBarry Smith M*/
626c838673bSBarry Smith 
627c838673bSBarry Smith /*MC
628c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
6299fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
6309fc87aa7SBarry Smith      such as PCFIELDSPLIT.
6319fc87aa7SBarry Smith 
6329fc87aa7SBarry Smith    Level: beginner
6339fc87aa7SBarry Smith 
63495452b02SPatrick Sanan     Notes:
63595452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
6369fc87aa7SBarry Smith 
6379fc87aa7SBarry Smith 
6389fc87aa7SBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
6399fc87aa7SBarry Smith 
6409fc87aa7SBarry Smith M*/
6419fc87aa7SBarry Smith 
6429fc87aa7SBarry Smith /*MC
643c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
644c838673bSBarry Smith         while the KSPSolve() is still running.
645c838673bSBarry Smith 
646c838673bSBarry Smith    Level: beginner
647c838673bSBarry Smith 
648c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
649c838673bSBarry Smith 
650c838673bSBarry Smith M*/
651c838673bSBarry Smith 
652014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
653df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
654df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
655014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
6568de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6573eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6588de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
6598de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
6608de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
6618de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
66254b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
6630059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
664014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
665*94a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
666abef13c0SSatish Balay 
66725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
6688ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
66925ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
6708ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
67125ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
6728ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
67325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
6748ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
67525ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
6768ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
67725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
6788ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
6798ea1b3e6SJed Brown 
6800bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
68125ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
682d4fbbf0eSBarry Smith 
68328ce4d24SBarry Smith /*E
68428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
68528ce4d24SBarry Smith 
68628ce4d24SBarry Smith    Level: beginner
68728ce4d24SBarry Smith 
68828ce4d24SBarry Smith .seealso: KSPCGSetType()
68928ce4d24SBarry Smith E*/
6909dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
6916a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
69228ce4d24SBarry Smith 
693014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
694014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);
6958031f4adStmunson 
696dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
697dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
698dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
699fcae7a14Stmunson 
70005de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
70105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
70225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
70325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
7048031f4adStmunson 
705014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
7061d6018f0SLisandro Dalcin 
707014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
708014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
7093369ce9aSBarry Smith 
7109804daf3SBarry Smith #include <petscdrawtypes.h>
711d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
712d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
713d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
714d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
715014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
7162f2e5d10SKris Buschelman 
717014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
718014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
71903919abeSBarry Smith 
720ba36735cSStefano Zampini /*S
721ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
722f8a50e2bSBarry Smith 
723ba36735cSStefano Zampini    Level: beginner
724f8a50e2bSBarry Smith 
725ba36735cSStefano Zampini .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
726ba36735cSStefano Zampini S*/
727ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess;
728ba36735cSStefano Zampini /*J
729ba36735cSStefano Zampini     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
730ba36735cSStefano Zampini 
731ba36735cSStefano Zampini    Level: beginner
732ba36735cSStefano Zampini 
733ba36735cSStefano Zampini .seealso: KSPGuess
734ba36735cSStefano Zampini J*/
735ba36735cSStefano Zampini typedef const char* KSPGuessType;
736ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
737ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
7381d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
739ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
740ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
741ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
742ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
743ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
744ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
745ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
746ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
747ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
748ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
749ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
750ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
751014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
752ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
753ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
754f8a50e2bSBarry Smith 
755470b340bSDmitry Karpeev /*E
756470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
757470b340bSDmitry Karpeev 
758470b340bSDmitry Karpeev     Level: intermediate
759470b340bSDmitry Karpeev 
760470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
761470b340bSDmitry Karpeev E*/
7620c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
763470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
764470b340bSDmitry Karpeev 
765864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
766864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
767864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
768864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
769864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
77092f76d53SAlp Dener 
771014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
772014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
773d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
774bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
775aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
776bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
777470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
778470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
7795bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
7805a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
781470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
782470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
7833f22127dSBarry Smith 
78478e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
78578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
78678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
787864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
788864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
789864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
790864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
791864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
792cd929ea3SAlp Dener 
793cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
794b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
795cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
796cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
797e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
7980ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
799cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
800cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
801cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
802cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
803cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
8042d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
805cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
806cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
807cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
808cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
80992f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
810cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
811cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
812864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
813864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
814cd929ea3SAlp Dener 
815014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
816014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
817014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
818014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
819014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
820fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
82123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
822fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
82323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
82423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
825014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
826014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
827fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
828fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
8296c699258SBarry Smith 
83002b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
8318c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
832e29c0834SMatthew G. Knepley                                            void (**)(PetscInt, PetscInt, PetscInt,
833e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
834e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
835191494d9SMatthew G. Knepley                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
836557cf195SMatthew G. Knepley 
837557cf195SMatthew G. Knepley 
838557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *);
839557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal);
8402eac72dbSBarry Smith #endif
841