xref: /petsc/include/petscksp.h (revision 557cf1953a2187a2a87a3adbc9420a6cc217e91d)
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"
40df2a969aSvictorle #define   KSPCGNE       "cgne"
4105de396fSBarry Smith #define   KSPNASH       "nash"
4205de396fSBarry Smith #define   KSPSTCG       "stcg"
4305de396fSBarry Smith #define   KSPGLTR       "gltr"
4405de396fSBarry Smith #define     KSPCGNASH  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"")  "nash"
4505de396fSBarry Smith #define     KSPCGSTCG  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"")  "stcg"
4605de396fSBarry Smith #define     KSPCGGLTR  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
47640f4f53SPatrick Sanan #define KSPFCG        "fcg"
48390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
4982bf6240SBarry Smith #define KSPGMRES      "gmres"
50483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
519dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
529dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
534f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
5461b468f9SJed Brown #define   KSPPGMRES     "pgmres"
5582bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5682bf6240SBarry Smith #define KSPBCGS       "bcgs"
57e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5818ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
59c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
60f0037002Svictorle #define   KSPBCGSL      "bcgsl"
61f154af2dSSiegfried Cools #define   KSPPIPEBCGS   "pipebcgs"
6282bf6240SBarry Smith #define KSPCGS        "cgs"
6382bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
6482bf6240SBarry Smith #define KSPCR         "cr"
652c8fc646SJed Brown #define KSPPIPECR     "pipecr"
6682bf6240SBarry Smith #define KSPLSQR       "lsqr"
6782bf6240SBarry Smith #define KSPPREONLY    "preonly"
6882bf6240SBarry Smith #define KSPQCG        "qcg"
69c9cf9b11SBarry Smith #define KSPBICG       "bicg"
70b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
7101c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
7221356919SSatish Balay #define KSPLCD        "lcd"
731d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
7458ad63f4SBarry Smith #define KSPGCR        "gcr"
75fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
76e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
77e4d80e07Szianekhodja #define KSPCGLS       "cgls"
78329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
7938cfc46eSPierre Jolivet #define KSPHPDDM      "hpddm"
802eac72dbSBarry Smith 
818ba1e511SMatthew Knepley /* Logging support */
82014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
83ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
845358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
858ba1e511SMatthew Knepley 
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
8719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
88c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
935ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat);
945ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPSetMatSolveBlockSize(KSP,PetscInt);
955ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPGetMatSolveBlockSize(KSP,PetscInt*);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
97ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
9923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
1007d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
10125c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
102c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
1032eac72dbSBarry Smith 
104140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
105ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
106bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
10730de9b25SBarry Smith 
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
111c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
112014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
1177d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
118c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
120c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
125734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1262a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
12725ef9dfeSBarry 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);}
1282eac72dbSBarry Smith 
129d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
130d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
131d4211eb9SBarry Smith 
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
134aabeff55SBarry Smith 
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1414b0e389bSBarry Smith 
142fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
143fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
144fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
145fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
146fa0ddf94SBarry Smith 
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
148cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
153285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
154b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
155b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
156b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
157b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
1594a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
160f3b08a26SMatthew G. Knepley /*
161f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
162f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
163f3b08a26SMatthew G. Knepley */
164f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
165f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
166f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
167f3b08a26SMatthew G. Knepley 
1680a71e23eSMatthew Knepley 
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
1712eac72dbSBarry Smith 
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
17558450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
176b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
17758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
179d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
180d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1817d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
1824b0e389bSBarry Smith 
183640f4f53SPatrick Sanan /*E
184640f4f53SPatrick Sanan 
18506137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
186640f4f53SPatrick Sanan 
18793f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
18893f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
189640f4f53SPatrick Sanan 
1902b26979fSBarry Smith    Level: intermediate
19106137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
192640f4f53SPatrick Sanan 
193640f4f53SPatrick Sanan E*/
19493f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
19506137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
196640f4f53SPatrick Sanan 
197640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
198640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
199640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
200640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
20106137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
20206137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
203640f4f53SPatrick Sanan 
204390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
205390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
206390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
207390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
208390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
209390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
210390d8e47SPatrick Sanan 
211fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
212fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
213fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
214fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
215fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
216fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
217fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
218fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
21983f0b094SBarry Smith 
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
2239f236934SBarry Smith 
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2291d73ed98SBarry Smith 
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2321d73ed98SBarry Smith 
233483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
234483d6965SPatrick Sanan 
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
23858ad63f4SBarry Smith 
23904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
24004d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
24104d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
242568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
243e82af88dSprj- 
244e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
245e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
2465ea7661aSPierre 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); }
247d55ab951SPierre Jolivet /*E
248d55ab951SPierre Jolivet     KSPHPDDMType - Type of Krylov method used by KSPHPDDM
249d55ab951SPierre Jolivet 
250d55ab951SPierre Jolivet     Level: intermediate
251d55ab951SPierre Jolivet 
252d55ab951SPierre Jolivet     Values:
253d55ab951SPierre Jolivet +   KSP_HPDDM_TYPE_GMRES (default)
254d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGMRES
255d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_CG
256d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BCG
257d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_GCRODR
258d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGCRODR
259d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BFBCG
260d55ab951SPierre Jolivet -   KSP_HPDDM_TYPE_PREONLY
261d55ab951SPierre Jolivet 
262d55ab951SPierre Jolivet .seealso: KSPHPDDM, KSPHPDDMSetType()
263d55ab951SPierre Jolivet E*/
264d55ab951SPierre Jolivet typedef enum {
265d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES = 0,
266d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES = 1,
267d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG = 2,
268d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG = 3,
269d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR = 4,
270d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
271d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG = 6,
272d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
273d55ab951SPierre Jolivet } KSPHPDDMType;
274d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
275d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
276d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
277e82af88dSprj- 
278b49fd9e1SBarry Smith /*E
279b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
280b49fd9e1SBarry Smith 
281b49fd9e1SBarry Smith    Level: advanced
282b49fd9e1SBarry Smith 
283bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
284e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
285b49fd9e1SBarry Smith 
286b49fd9e1SBarry Smith E*/
28778d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2886a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2891f7e983dSSatish Balay /*MC
2901957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2918c5b8ba0SBarry Smith 
2928c5b8ba0SBarry Smith    Level: advanced
2938c5b8ba0SBarry Smith 
2948c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2958c5b8ba0SBarry Smith 
296bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
297e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2988c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2998c5b8ba0SBarry Smith M*/
3008c5b8ba0SBarry Smith 
3011f7e983dSSatish Balay /*MC
3028c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
3038c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
3048c5b8ba0SBarry Smith           poor orthogonality.
3058c5b8ba0SBarry Smith 
3068c5b8ba0SBarry Smith    Level: advanced
3078c5b8ba0SBarry Smith 
3088c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
3098c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
3108c5b8ba0SBarry Smith 
311bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
312e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
3138c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3148c5b8ba0SBarry Smith M*/
3158c5b8ba0SBarry Smith 
3161f7e983dSSatish Balay /*MC
3178c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3188c5b8ba0SBarry Smith 
3198c5b8ba0SBarry Smith    Level: advanced
3208c5b8ba0SBarry Smith 
3218c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
3228c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
3238c5b8ba0SBarry Smith 
3248c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
3258c5b8ba0SBarry Smith 
326bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
327e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
3288c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3298c5b8ba0SBarry Smith M*/
3308c5b8ba0SBarry Smith 
331014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
332014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
33308480c60SBarry Smith 
334014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
335014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
337c38d4ed2SBarry Smith 
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
340014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
341121fd945SKris Buschelman 
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
343014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
345e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
346d9492815SBarry Smith 
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
34887d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
3502eac72dbSBarry Smith 
351fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
352fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
353fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
354fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
355390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
356390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
357fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
358fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
359fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
360fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
361e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
362e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
363e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
365fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
36684cb2905SBarry Smith 
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
369c01c455dSBarry Smith 
37023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
37123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3761eb62cbbSBarry Smith 
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
3811f7f0c4fSBarry Smith 
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
38355849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
384fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
3852a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3862a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
38755849f57SBarry Smith 
38855849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3891eb62cbbSBarry Smith 
390dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
3910e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
392014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
393884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
394db9b2ab1SHong Zhang 
395014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
39768ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
39883ab6a24SBarry Smith 
39928ce4d24SBarry Smith /*E
4008a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
4018a4b9c5cSBarry Smith        test routines.
4028a4b9c5cSBarry Smith 
4038a4b9c5cSBarry Smith    Level: advanced
4048a4b9c5cSBarry Smith 
405a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
4061718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
407a3f661c8SBarry Smith 
40895452b02SPatrick Sanan    Notes:
40995452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4108a4b9c5cSBarry Smith 
41194b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
4121718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
4138a4b9c5cSBarry Smith E*/
4149e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
4159e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
416014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
417a21b2a99SBarry Smith 
4181f7e983dSSatish Balay /*MC
4199793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
4208c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
4218c5b8ba0SBarry Smith           be based on a norm of a residual etc.
4228c5b8ba0SBarry Smith 
4238c5b8ba0SBarry Smith    Level: advanced
4248c5b8ba0SBarry Smith 
4251957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
4268c5b8ba0SBarry Smith 
427ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
4288c5b8ba0SBarry Smith M*/
4298c5b8ba0SBarry Smith 
4301f7e983dSSatish Balay /*MC
4311957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
4328c5b8ba0SBarry Smith        convergence test routine.
4338c5b8ba0SBarry Smith 
4348c5b8ba0SBarry Smith    Level: advanced
4358c5b8ba0SBarry Smith 
4369793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4378c5b8ba0SBarry Smith M*/
4388c5b8ba0SBarry Smith 
4391f7e983dSSatish Balay /*MC
440ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
4418c5b8ba0SBarry Smith        convergence test routine.
4428c5b8ba0SBarry Smith 
4438c5b8ba0SBarry Smith    Level: advanced
4448c5b8ba0SBarry Smith 
4459793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4468c5b8ba0SBarry Smith M*/
4478c5b8ba0SBarry Smith 
4481f7e983dSSatish Balay /*MC
449ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
450390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
4518c5b8ba0SBarry Smith 
4528c5b8ba0SBarry Smith    Level: advanced
4538c5b8ba0SBarry Smith 
4549793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
4558c5b8ba0SBarry Smith M*/
4568c5b8ba0SBarry Smith 
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4628a4b9c5cSBarry Smith 
463f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
4648a4b9c5cSBarry Smith /*E
4651957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
46628ce4d24SBarry Smith 
46728ce4d24SBarry Smith    Level: beginner
46828ce4d24SBarry Smith 
46995452b02SPatrick Sanan    Notes:
47095452b02SPatrick Sanan     See KSPGetConvergedReason() for explanation of each value
47128ce4d24SBarry Smith 
47295452b02SPatrick Sanan    Developer Notes:
47395452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4744d0a8057SBarry Smith 
4754d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4764d0a8057SBarry Smith       any of the values here also change them that array of names.
47786c02ca4SBarry Smith 
478c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
47928ce4d24SBarry Smith E*/
480d15094e1SBarry Smith typedef enum {/* converged */
4819ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4829ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
483d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
484d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
485b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4868031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4878031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
488329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
489af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
490d15094e1SBarry Smith               /* diverged */
491b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
492d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
493d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
494d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
495b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
496b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
497b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4984d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4996aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
500c0decd05SBarry Smith               KSP_DIVERGED_PC_FAILED           = -11,
501aa4d2078SSatish Balay               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
502d15094e1SBarry Smith 
503d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
504014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
505d15094e1SBarry Smith 
506c838673bSBarry Smith /*MC
507f9fed41fSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
508c838673bSBarry Smith 
509c838673bSBarry Smith    Level: beginner
510c838673bSBarry Smith 
511c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
512c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
513c838673bSBarry Smith        2-norm of the residual for right preconditioning
514c838673bSBarry Smith 
515f9fed41fSBarry Smith    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
516f9fed41fSBarry Smith 
517f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
518c838673bSBarry Smith 
519c838673bSBarry Smith M*/
520c838673bSBarry Smith 
521c838673bSBarry Smith /*MC
522c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
523c838673bSBarry Smith 
524c838673bSBarry Smith    Level: beginner
525c838673bSBarry Smith 
526c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
527c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
528c838673bSBarry Smith        2-norm of the residual for right preconditioning
529c838673bSBarry Smith 
530f9fed41fSBarry Smith    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
531c838673bSBarry Smith 
532f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
533c838673bSBarry Smith 
534c838673bSBarry Smith M*/
535c838673bSBarry Smith 
536c838673bSBarry Smith /*MC
537c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
538c838673bSBarry Smith 
539c838673bSBarry Smith    Level: beginner
540c838673bSBarry Smith 
541c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
542c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
543c838673bSBarry Smith        2-norm of the residual for right preconditioning
544c838673bSBarry Smith 
545c838673bSBarry Smith    Level: beginner
546c838673bSBarry Smith 
547f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
548c838673bSBarry Smith 
549c838673bSBarry Smith M*/
550c838673bSBarry Smith 
551c838673bSBarry Smith /*MC
552c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
553c838673bSBarry Smith       reached
554c838673bSBarry Smith 
555c838673bSBarry Smith    Level: beginner
556c838673bSBarry Smith 
557c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
558c838673bSBarry Smith 
559c838673bSBarry Smith M*/
560c838673bSBarry Smith 
561c838673bSBarry Smith /*MC
5628c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5630059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5644d0a8057SBarry Smith            test routine is set in KSP.
565c838673bSBarry Smith 
566c838673bSBarry Smith    Level: beginner
567c838673bSBarry Smith 
568c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
569c838673bSBarry Smith 
570c838673bSBarry Smith M*/
571c838673bSBarry Smith 
572c838673bSBarry Smith /*MC
573c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5741de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
5751de96524SPierre Jolivet           preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
5761de96524SPierre Jolivet           are colinear.
577c838673bSBarry Smith 
578c838673bSBarry Smith    Level: beginner
579c838673bSBarry Smith 
580c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
581c838673bSBarry Smith 
582c838673bSBarry Smith M*/
583c838673bSBarry Smith 
584c838673bSBarry Smith /*MC
585c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
586c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
587c838673bSBarry Smith 
588c838673bSBarry Smith    Level: beginner
589c838673bSBarry Smith 
590c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
591c838673bSBarry Smith 
592c838673bSBarry Smith M*/
593c838673bSBarry Smith 
594c838673bSBarry Smith /*MC
595c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
596c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
597c838673bSBarry Smith 
598c838673bSBarry Smith    Level: beginner
599c838673bSBarry Smith 
600c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
601c838673bSBarry Smith 
602c838673bSBarry Smith M*/
603c838673bSBarry Smith 
604c838673bSBarry Smith /*MC
605c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
606c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
607c838673bSBarry Smith         be positive definite
608c838673bSBarry Smith 
609c838673bSBarry Smith    Level: beginner
610c838673bSBarry Smith 
61195452b02SPatrick Sanan      Notes:
61295452b02SPatrick Sanan     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
613c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
614c838673bSBarry Smith 
615c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
616c838673bSBarry Smith 
617c838673bSBarry Smith M*/
618c838673bSBarry Smith 
619c838673bSBarry Smith /*MC
620c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
6219fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
6229fc87aa7SBarry Smith      such as PCFIELDSPLIT.
6239fc87aa7SBarry Smith 
6249fc87aa7SBarry Smith    Level: beginner
6259fc87aa7SBarry Smith 
62695452b02SPatrick Sanan     Notes:
62795452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
6289fc87aa7SBarry Smith 
6299fc87aa7SBarry Smith 
6309fc87aa7SBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
6319fc87aa7SBarry Smith 
6329fc87aa7SBarry Smith M*/
6339fc87aa7SBarry Smith 
6349fc87aa7SBarry Smith /*MC
635c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
636c838673bSBarry Smith         while the KSPSolve() is still running.
637c838673bSBarry Smith 
638c838673bSBarry Smith    Level: beginner
639c838673bSBarry Smith 
640c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
641c838673bSBarry Smith 
642c838673bSBarry Smith M*/
643c838673bSBarry Smith 
644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
645df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
646df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
647014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
6488de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6493eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6508de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
6518de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
6528de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
6538de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
65454b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
6550059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
656014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
657abef13c0SSatish Balay 
65825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
6598ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
66025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
6618ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
66225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
6638ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
66425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
6658ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
66625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
6678ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
66825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
6698ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
6708ea1b3e6SJed Brown 
6710bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
67225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
673d4fbbf0eSBarry Smith 
67428ce4d24SBarry Smith /*E
67528ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
67628ce4d24SBarry Smith 
67728ce4d24SBarry Smith    Level: beginner
67828ce4d24SBarry Smith 
67928ce4d24SBarry Smith .seealso: KSPCGSetType()
68028ce4d24SBarry Smith E*/
6819dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
6826a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
68328ce4d24SBarry Smith 
684014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
685014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
6868031f4adStmunson 
687dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
688dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
689dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
690fcae7a14Stmunson 
69105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
69205de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
69325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
69425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
6958031f4adStmunson 
696014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6971d6018f0SLisandro Dalcin 
698014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
699014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
7003369ce9aSBarry Smith 
7019804daf3SBarry Smith #include <petscdrawtypes.h>
702d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
703d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
704d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
705d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
706014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
7072f2e5d10SKris Buschelman 
708014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
709014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
71003919abeSBarry Smith 
711ba36735cSStefano Zampini /*S
712ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
713f8a50e2bSBarry Smith 
714ba36735cSStefano Zampini    Level: beginner
715f8a50e2bSBarry Smith 
716ba36735cSStefano Zampini .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
717ba36735cSStefano Zampini S*/
718ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess;
719ba36735cSStefano Zampini /*J
720ba36735cSStefano Zampini     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
721ba36735cSStefano Zampini 
722ba36735cSStefano Zampini    Level: beginner
723ba36735cSStefano Zampini 
724ba36735cSStefano Zampini .seealso: KSPGuess
725ba36735cSStefano Zampini J*/
726ba36735cSStefano Zampini typedef const char* KSPGuessType;
727ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
728ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
7291d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
730ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
731ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
732ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
733ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
734ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
735ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
736ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
737ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
738ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
739ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
740ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
741ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
742014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
743ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
744ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
745f8a50e2bSBarry Smith 
746470b340bSDmitry Karpeev /*E
747470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
748470b340bSDmitry Karpeev 
749470b340bSDmitry Karpeev     Level: intermediate
750470b340bSDmitry Karpeev 
751470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
752470b340bSDmitry Karpeev E*/
7530c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
754470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
755470b340bSDmitry Karpeev 
756864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
757864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
758864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
759864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
760864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
76192f76d53SAlp Dener 
762014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
763014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
764d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
765bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
766aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
767bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
768470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
769470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
7705bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
7715a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
772470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
773470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
7743f22127dSBarry Smith 
77578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
77678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
77778e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
778864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
779864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
780864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
781864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
782864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
783cd929ea3SAlp Dener 
784cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
785b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
786cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
787cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
788e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
7890ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
790cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
791cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
792cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
793cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
794cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
7952d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
796cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
797cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
798cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
799cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
80092f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
801cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
802cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
803864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
804864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
805cd929ea3SAlp Dener 
806014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
807014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
808014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
809014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
810014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
811fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
81223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
813fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
81423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
81523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
816014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
817014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
818fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
819fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
8206c699258SBarry Smith 
82102b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
8228c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
823e29c0834SMatthew G. Knepley                                            void (**)(PetscInt, PetscInt, PetscInt,
824e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
825e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
826191494d9SMatthew G. Knepley                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
827*557cf195SMatthew G. Knepley 
828*557cf195SMatthew G. Knepley 
829*557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *);
830*557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal);
8312eac72dbSBarry Smith #endif
832