xref: /petsc/include/petscksp.h (revision d55ab9515a88c9d141879b24f8f7178823920ea6)
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*);
1600a71e23eSMatthew Knepley 
161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
1632eac72dbSBarry Smith 
164014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
16758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
168b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
16958450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
171d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
172d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1737d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
1744b0e389bSBarry Smith 
175640f4f53SPatrick Sanan /*E
176640f4f53SPatrick Sanan 
17706137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
178640f4f53SPatrick Sanan 
17993f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
18093f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
181640f4f53SPatrick Sanan 
1822b26979fSBarry Smith    Level: intermediate
18306137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
184640f4f53SPatrick Sanan 
185640f4f53SPatrick Sanan E*/
18693f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
18706137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
188640f4f53SPatrick Sanan 
189640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
190640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
191640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
192640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
19306137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
19406137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
195640f4f53SPatrick Sanan 
196390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
197390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
198390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
199390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
200390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
201390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
202390d8e47SPatrick Sanan 
203fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
204fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
205fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
206fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
207fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
208fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
209fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
210fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
21183f0b094SBarry Smith 
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
2159f236934SBarry Smith 
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2211d73ed98SBarry Smith 
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2241d73ed98SBarry Smith 
225483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
226483d6965SPatrick Sanan 
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
23058ad63f4SBarry Smith 
23104d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
23204d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
23304d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
234568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
235e82af88dSprj- 
236e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
237e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
2385ea7661aSPierre 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); }
239*d55ab951SPierre Jolivet /*E
240*d55ab951SPierre Jolivet     KSPHPDDMType - Type of Krylov method used by KSPHPDDM
241*d55ab951SPierre Jolivet 
242*d55ab951SPierre Jolivet     Level: intermediate
243*d55ab951SPierre Jolivet 
244*d55ab951SPierre Jolivet     Values:
245*d55ab951SPierre Jolivet +   KSP_HPDDM_TYPE_GMRES (default)
246*d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGMRES
247*d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_CG
248*d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BCG
249*d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_GCRODR
250*d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BGCRODR
251*d55ab951SPierre Jolivet .   KSP_HPDDM_TYPE_BFBCG
252*d55ab951SPierre Jolivet -   KSP_HPDDM_TYPE_PREONLY
253*d55ab951SPierre Jolivet 
254*d55ab951SPierre Jolivet .seealso: KSPHPDDM, KSPHPDDMSetType()
255*d55ab951SPierre Jolivet E*/
256*d55ab951SPierre Jolivet typedef enum {
257*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES = 0,
258*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES = 1,
259*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG = 2,
260*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG = 3,
261*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR = 4,
262*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
263*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG = 6,
264*d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
265*d55ab951SPierre Jolivet } KSPHPDDMType;
266*d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
267*d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
268*d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
269e82af88dSprj- 
270b49fd9e1SBarry Smith /*E
271b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
272b49fd9e1SBarry Smith 
273b49fd9e1SBarry Smith    Level: advanced
274b49fd9e1SBarry Smith 
275bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
276e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
277b49fd9e1SBarry Smith 
278b49fd9e1SBarry Smith E*/
27978d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2806a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2811f7e983dSSatish Balay /*MC
2821957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2838c5b8ba0SBarry Smith 
2848c5b8ba0SBarry Smith    Level: advanced
2858c5b8ba0SBarry Smith 
2868c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2878c5b8ba0SBarry Smith 
288bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
289e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2908c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2918c5b8ba0SBarry Smith M*/
2928c5b8ba0SBarry Smith 
2931f7e983dSSatish Balay /*MC
2948c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2958c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2968c5b8ba0SBarry Smith           poor orthogonality.
2978c5b8ba0SBarry Smith 
2988c5b8ba0SBarry Smith    Level: advanced
2998c5b8ba0SBarry Smith 
3008c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
3018c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
3028c5b8ba0SBarry Smith 
303bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
304e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
3058c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3068c5b8ba0SBarry Smith M*/
3078c5b8ba0SBarry Smith 
3081f7e983dSSatish Balay /*MC
3098c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3108c5b8ba0SBarry Smith 
3118c5b8ba0SBarry Smith    Level: advanced
3128c5b8ba0SBarry Smith 
3138c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
3148c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
3158c5b8ba0SBarry Smith 
3168c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
3178c5b8ba0SBarry Smith 
318bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
319e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
3208c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
3218c5b8ba0SBarry Smith M*/
3228c5b8ba0SBarry Smith 
323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
32508480c60SBarry Smith 
326014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
327014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
328014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
329c38d4ed2SBarry Smith 
330014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
331014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
332014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
333121fd945SKris Buschelman 
334014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
335014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
337e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
338d9492815SBarry Smith 
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
34087d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
341014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
3422eac72dbSBarry Smith 
343fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
344fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
345fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
346fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
347390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
348390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
349fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
350fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
351fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
352fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
353e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
354e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
355e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
357fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
35884cb2905SBarry Smith 
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
361c01c455dSBarry Smith 
36223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
36323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3681eb62cbbSBarry Smith 
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
3731f7f0c4fSBarry Smith 
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
37555849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
376fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
3772a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3782a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
37955849f57SBarry Smith 
38055849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3811eb62cbbSBarry Smith 
382dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
3830e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
385884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
386db9b2ab1SHong Zhang 
387014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
38968ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
39083ab6a24SBarry Smith 
39128ce4d24SBarry Smith /*E
3928a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3938a4b9c5cSBarry Smith        test routines.
3948a4b9c5cSBarry Smith 
3958a4b9c5cSBarry Smith    Level: advanced
3968a4b9c5cSBarry Smith 
397a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3981718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
399a3f661c8SBarry Smith 
40095452b02SPatrick Sanan    Notes:
40195452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4028a4b9c5cSBarry Smith 
40394b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
4041718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
4058a4b9c5cSBarry Smith E*/
4069e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
4079e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
408014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
409a21b2a99SBarry Smith 
4101f7e983dSSatish Balay /*MC
4119793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
4128c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
4138c5b8ba0SBarry Smith           be based on a norm of a residual etc.
4148c5b8ba0SBarry Smith 
4158c5b8ba0SBarry Smith    Level: advanced
4168c5b8ba0SBarry Smith 
4171957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
4188c5b8ba0SBarry Smith 
419ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
4208c5b8ba0SBarry Smith M*/
4218c5b8ba0SBarry Smith 
4221f7e983dSSatish Balay /*MC
4231957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
4248c5b8ba0SBarry Smith        convergence test routine.
4258c5b8ba0SBarry Smith 
4268c5b8ba0SBarry Smith    Level: advanced
4278c5b8ba0SBarry Smith 
4289793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4298c5b8ba0SBarry Smith M*/
4308c5b8ba0SBarry Smith 
4311f7e983dSSatish Balay /*MC
432ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
4338c5b8ba0SBarry Smith        convergence test routine.
4348c5b8ba0SBarry Smith 
4358c5b8ba0SBarry Smith    Level: advanced
4368c5b8ba0SBarry Smith 
4379793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
4388c5b8ba0SBarry Smith M*/
4398c5b8ba0SBarry Smith 
4401f7e983dSSatish Balay /*MC
441ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
442390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
4438c5b8ba0SBarry Smith 
4448c5b8ba0SBarry Smith    Level: advanced
4458c5b8ba0SBarry Smith 
4469793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
4478c5b8ba0SBarry Smith M*/
4488c5b8ba0SBarry Smith 
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4548a4b9c5cSBarry Smith 
455f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
4568a4b9c5cSBarry Smith /*E
4571957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
45828ce4d24SBarry Smith 
45928ce4d24SBarry Smith    Level: beginner
46028ce4d24SBarry Smith 
46195452b02SPatrick Sanan    Notes:
46295452b02SPatrick Sanan     See KSPGetConvergedReason() for explanation of each value
46328ce4d24SBarry Smith 
46495452b02SPatrick Sanan    Developer Notes:
46595452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
4664d0a8057SBarry Smith 
4674d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4684d0a8057SBarry Smith       any of the values here also change them that array of names.
46986c02ca4SBarry Smith 
470c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
47128ce4d24SBarry Smith E*/
472d15094e1SBarry Smith typedef enum {/* converged */
4739ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4749ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
475d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
476d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
477b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4788031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4798031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
480329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
481af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
482d15094e1SBarry Smith               /* diverged */
483b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
484d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
485d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
486d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
487b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
488b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
489b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4904d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4916aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
492c0decd05SBarry Smith               KSP_DIVERGED_PC_FAILED           = -11,
493aa4d2078SSatish Balay               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
494d15094e1SBarry Smith 
495d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
496014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
497d15094e1SBarry Smith 
498c838673bSBarry Smith /*MC
499f9fed41fSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
500c838673bSBarry Smith 
501c838673bSBarry Smith    Level: beginner
502c838673bSBarry Smith 
503c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
504c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
505c838673bSBarry Smith        2-norm of the residual for right preconditioning
506c838673bSBarry Smith 
507f9fed41fSBarry Smith    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
508f9fed41fSBarry Smith 
509f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
510c838673bSBarry Smith 
511c838673bSBarry Smith M*/
512c838673bSBarry Smith 
513c838673bSBarry Smith /*MC
514c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
515c838673bSBarry Smith 
516c838673bSBarry Smith    Level: beginner
517c838673bSBarry Smith 
518c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
519c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
520c838673bSBarry Smith        2-norm of the residual for right preconditioning
521c838673bSBarry Smith 
522f9fed41fSBarry Smith    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
523c838673bSBarry Smith 
524f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
525c838673bSBarry Smith 
526c838673bSBarry Smith M*/
527c838673bSBarry Smith 
528c838673bSBarry Smith /*MC
529c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
530c838673bSBarry Smith 
531c838673bSBarry Smith    Level: beginner
532c838673bSBarry Smith 
533c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
534c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
535c838673bSBarry Smith        2-norm of the residual for right preconditioning
536c838673bSBarry Smith 
537c838673bSBarry Smith    Level: beginner
538c838673bSBarry Smith 
539f9fed41fSBarry Smith .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
540c838673bSBarry Smith 
541c838673bSBarry Smith M*/
542c838673bSBarry Smith 
543c838673bSBarry Smith /*MC
544c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
545c838673bSBarry Smith       reached
546c838673bSBarry Smith 
547c838673bSBarry Smith    Level: beginner
548c838673bSBarry Smith 
549c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
550c838673bSBarry Smith 
551c838673bSBarry Smith M*/
552c838673bSBarry Smith 
553c838673bSBarry Smith /*MC
5548c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5550059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5564d0a8057SBarry Smith            test routine is set in KSP.
557c838673bSBarry Smith 
558c838673bSBarry Smith    Level: beginner
559c838673bSBarry Smith 
560c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
561c838673bSBarry Smith 
562c838673bSBarry Smith M*/
563c838673bSBarry Smith 
564c838673bSBarry Smith /*MC
565c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5661de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
5671de96524SPierre Jolivet           preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
5681de96524SPierre Jolivet           are colinear.
569c838673bSBarry Smith 
570c838673bSBarry Smith    Level: beginner
571c838673bSBarry Smith 
572c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
573c838673bSBarry Smith 
574c838673bSBarry Smith M*/
575c838673bSBarry Smith 
576c838673bSBarry Smith /*MC
577c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
578c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
579c838673bSBarry Smith 
580c838673bSBarry Smith    Level: beginner
581c838673bSBarry Smith 
582c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
583c838673bSBarry Smith 
584c838673bSBarry Smith M*/
585c838673bSBarry Smith 
586c838673bSBarry Smith /*MC
587c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
588c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
589c838673bSBarry Smith 
590c838673bSBarry Smith    Level: beginner
591c838673bSBarry Smith 
592c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
593c838673bSBarry Smith 
594c838673bSBarry Smith M*/
595c838673bSBarry Smith 
596c838673bSBarry Smith /*MC
597c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
598c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
599c838673bSBarry Smith         be positive definite
600c838673bSBarry Smith 
601c838673bSBarry Smith    Level: beginner
602c838673bSBarry Smith 
60395452b02SPatrick Sanan      Notes:
60495452b02SPatrick Sanan     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
605c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
606c838673bSBarry Smith 
607c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
608c838673bSBarry Smith 
609c838673bSBarry Smith M*/
610c838673bSBarry Smith 
611c838673bSBarry Smith /*MC
612c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
6139fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
6149fc87aa7SBarry Smith      such as PCFIELDSPLIT.
6159fc87aa7SBarry Smith 
6169fc87aa7SBarry Smith    Level: beginner
6179fc87aa7SBarry Smith 
61895452b02SPatrick Sanan     Notes:
61995452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
6209fc87aa7SBarry Smith 
6219fc87aa7SBarry Smith 
6229fc87aa7SBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
6239fc87aa7SBarry Smith 
6249fc87aa7SBarry Smith M*/
6259fc87aa7SBarry Smith 
6269fc87aa7SBarry Smith /*MC
627c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
628c838673bSBarry Smith         while the KSPSolve() is still running.
629c838673bSBarry Smith 
630c838673bSBarry Smith    Level: beginner
631c838673bSBarry Smith 
632c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
633c838673bSBarry Smith 
634c838673bSBarry Smith M*/
635c838673bSBarry Smith 
636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
637df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
638df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
6408de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6413eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
6428de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
6438de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
6448de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
6458de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
64654b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
6470059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
648014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
649abef13c0SSatish Balay 
65025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
6518ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
65225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
6538ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
65425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
6558ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
65625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
6578ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
65825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
6598ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
66025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
6618ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
6628ea1b3e6SJed Brown 
6630bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
66425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
665d4fbbf0eSBarry Smith 
66628ce4d24SBarry Smith /*E
66728ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
66828ce4d24SBarry Smith 
66928ce4d24SBarry Smith    Level: beginner
67028ce4d24SBarry Smith 
67128ce4d24SBarry Smith .seealso: KSPCGSetType()
67228ce4d24SBarry Smith E*/
6739dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
6746a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
67528ce4d24SBarry Smith 
676014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
677014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
6788031f4adStmunson 
679dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
680dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
681dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
682fcae7a14Stmunson 
68305de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
68405de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
68525ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
68625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
6878031f4adStmunson 
688014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6891d6018f0SLisandro Dalcin 
690014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
691014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
6923369ce9aSBarry Smith 
6939804daf3SBarry Smith #include <petscdrawtypes.h>
694d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
695d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
696d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
697d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
698014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6992f2e5d10SKris Buschelman 
700014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
701014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
70203919abeSBarry Smith 
703ba36735cSStefano Zampini /*S
704ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
705f8a50e2bSBarry Smith 
706ba36735cSStefano Zampini    Level: beginner
707f8a50e2bSBarry Smith 
708ba36735cSStefano Zampini .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
709ba36735cSStefano Zampini S*/
710ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess;
711ba36735cSStefano Zampini /*J
712ba36735cSStefano Zampini     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
713ba36735cSStefano Zampini 
714ba36735cSStefano Zampini    Level: beginner
715ba36735cSStefano Zampini 
716ba36735cSStefano Zampini .seealso: KSPGuess
717ba36735cSStefano Zampini J*/
718ba36735cSStefano Zampini typedef const char* KSPGuessType;
719ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
720ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
7211d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
722ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
723ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
724ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
725ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
726ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
727ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
728ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
729ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
730ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
731ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
732ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
733ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
734014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
735ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
736ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
737f8a50e2bSBarry Smith 
738470b340bSDmitry Karpeev /*E
739470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
740470b340bSDmitry Karpeev 
741470b340bSDmitry Karpeev     Level: intermediate
742470b340bSDmitry Karpeev 
743470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
744470b340bSDmitry Karpeev E*/
7450c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
746470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
747470b340bSDmitry Karpeev 
748864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
749864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
750864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
751864588a7SAlp Dener               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
752864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
75392f76d53SAlp Dener 
754014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
755014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
756d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
757bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
758aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
759bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
760470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
761470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
7625bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
7635a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
764470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
765470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
7663f22127dSBarry Smith 
76778e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
76878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
76978e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
770864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
771864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
772864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
773864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
774864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
775cd929ea3SAlp Dener 
776cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
777b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
778cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
779cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
780e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
7810ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
782cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
783cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
784cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
785cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
786cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
7872d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
788cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
789cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
790cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
791cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
79292f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
793cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
794cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
795864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
796864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
797cd929ea3SAlp Dener 
798014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
799014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
800014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
801014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
802014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
803fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
80423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
805fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
80623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
80723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
808014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
809014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
810fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
811fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
8126c699258SBarry Smith 
81302b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
8148c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
815e29c0834SMatthew G. Knepley                                            void (**)(PetscInt, PetscInt, PetscInt,
816e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
817e29c0834SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
818191494d9SMatthew G. Knepley                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
8192eac72dbSBarry Smith #endif
820