xref: /petsc/include/petscksp.h (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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
6ac09b921SBarry Smith 
72c8e378dSBarry Smith #include <petscpc.h>
82eac72dbSBarry Smith 
9ac09b921SBarry Smith /* SUBMANSEC = KSP */
10ac09b921SBarry Smith 
11607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
121dbb0a54SBarry Smith 
1328ce4d24SBarry Smith /*S
148f6c3df8SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
158f6c3df8SBarry Smith          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
1628ce4d24SBarry Smith 
1728ce4d24SBarry Smith    Level: beginner
1828ce4d24SBarry Smith 
1987497f52SBarry Smith    Note:
20db3b2ab5SMatthew Knepley     When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a
2187497f52SBarry Smith     `KSPType` of `KSPPREONLY`, meaning that only application of the preconditioner is used as the linear solver.
228f6c3df8SBarry Smith 
23db781477SPatrick Sanan .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
2428ce4d24SBarry Smith S*/
25e2a1c21fSSatish Balay typedef struct _p_KSP *KSP;
262eac72dbSBarry Smith 
2776bdecfbSBarry Smith /*J
288f6c3df8SBarry Smith     KSPType - String with the name of a PETSc Krylov method.
2928ce4d24SBarry Smith 
3028ce4d24SBarry Smith    Level: beginner
3128ce4d24SBarry Smith 
32db781477SPatrick Sanan .seealso: `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
3376bdecfbSBarry Smith J*/
3419fd82e9SBarry Smith typedef const char *KSPType;
3582bf6240SBarry Smith #define KSPRICHARDSON "richardson"
366c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3782bf6240SBarry Smith #define KSPCG         "cg"
382c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
392c8fc646SJed Brown #define KSPPIPECG     "pipecg"
40901ccb91SSiegfried Cools #define KSPPIPECGRR   "pipecgrr"
41265663fdSSiegfried Cools #define KSPPIPELCG    "pipelcg"
42b21a8899STyler Chen #define KSPPIPEPRCG   "pipeprcg"
43325e8391SManasi #define KSPPIPECG2    "pipecg2"
44df2a969aSvictorle #define KSPCGNE       "cgne"
4505de396fSBarry Smith #define KSPNASH       "nash"
4605de396fSBarry Smith #define KSPSTCG       "stcg"
4705de396fSBarry Smith #define KSPGLTR       "gltr"
4805de396fSBarry Smith #define KSPCGNASH     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash"
4905de396fSBarry Smith #define KSPCGSTCG     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg"
5005de396fSBarry Smith #define KSPCGGLTR     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
51640f4f53SPatrick Sanan #define KSPFCG        "fcg"
52390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
5382bf6240SBarry Smith #define KSPGMRES      "gmres"
54483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
559dcbbd2bSBarry Smith #define KSPFGMRES     "fgmres"
569dcbbd2bSBarry Smith #define KSPLGMRES     "lgmres"
574f8e6cd9SSatish Balay #define KSPDGMRES     "dgmres"
5861b468f9SJed Brown #define KSPPGMRES     "pgmres"
5982bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
6082bf6240SBarry Smith #define KSPBCGS       "bcgs"
61e1c61ce8SBarry Smith #define KSPIBCGS      "ibcgs"
62345ecf0bSXiangmin Jiao #define KSPQMRCGS     "qmrcgs"
6318ac38e6SHong Zhang #define KSPFBCGS      "fbcgs"
64c2b71004SHong Zhang #define KSPFBCGSR     "fbcgsr"
65f0037002Svictorle #define KSPBCGSL      "bcgsl"
66f154af2dSSiegfried Cools #define KSPPIPEBCGS   "pipebcgs"
6782bf6240SBarry Smith #define KSPCGS        "cgs"
6882bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
6982bf6240SBarry Smith #define KSPCR         "cr"
702c8fc646SJed Brown #define KSPPIPECR     "pipecr"
7182bf6240SBarry Smith #define KSPLSQR       "lsqr"
7282bf6240SBarry Smith #define KSPPREONLY    "preonly"
733c2be86cSBarry Smith #define KSPNONE       "none"
7482bf6240SBarry Smith #define KSPQCG        "qcg"
75c9cf9b11SBarry Smith #define KSPBICG       "bicg"
76b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
7701c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
7821356919SSatish Balay #define KSPLCD        "lcd"
791d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
8058ad63f4SBarry Smith #define KSPGCR        "gcr"
81fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
82e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
83e4d80e07Szianekhodja #define KSPCGLS       "cgls"
84329cd39dSStefano Zampini #define KSPFETIDP     "fetidp"
8538cfc46eSPierre Jolivet #define KSPHPDDM      "hpddm"
862eac72dbSBarry Smith 
878ba1e511SMatthew Knepley /* Logging support */
88014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
89ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
905358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
918ba1e511SMatthew Knepley 
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
9319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
94c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
9975f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
1005ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
1013e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
102*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) {
103*9371c9d4SSatish Balay   return KSPSetMatSolveBatchSize(ksp, n);
104*9371c9d4SSatish Balay }
1053e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
106*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) {
107*9371c9d4SSatish Balay   return KSPGetMatSolveBatchSize(ksp, n);
108*9371c9d4SSatish Balay }
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
110ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
11223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
1137d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
11425c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
115c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
1162eac72dbSBarry Smith 
117140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
118ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList;
119798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList;
120798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
121798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
122bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
123798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
12430de9b25SBarry Smith 
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
128c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
1347d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
135c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
137c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
142734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
1432a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
144*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) {
145*9371c9d4SSatish Balay   return KSPCreateVecs(ksp, n, x, m, y);
146*9371c9d4SSatish Balay }
1472eac72dbSBarry Smith 
148d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
149d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
150d4211eb9SBarry Smith 
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
153aabeff55SBarry Smith 
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
1573ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
15894a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool);
16094a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
16194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool);
1624b0e389bSBarry Smith 
163fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
164fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
165fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
166fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
167fa0ddf94SBarry Smith 
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode    PCKSPGetKSP(PC, KSP *);
169cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode    PCKSPSetKSP(PC, KSP);
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode    PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
171014dd563SJed Brown PETSC_EXTERN PetscErrorCode    PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode    PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode    PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
174285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode    PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
175b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode    PCMGGetSmoother(PC, PetscInt, KSP *);
176b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode    PCMGGetSmootherDown(PC, PetscInt, KSP *);
177b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode    PCMGGetSmootherUp(PC, PetscInt, KSP *);
178b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSolve(PC, KSP *);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode    PCGalerkinGetKSP(PC, KSP *);
1804a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode    PCDeflationGetCoarseKSP(PC, KSP *);
181f3b08a26SMatthew G. Knepley /*
182f3b08a26SMatthew G. Knepley   PCMGCoarseList contains the list of coarse space constructor currently registered
183f3b08a26SMatthew G. Knepley   These are added with PCMGRegisterCoarseSpaceConstructor()
184f3b08a26SMatthew G. Knepley */
185f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList;
1862b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
1872b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
188f3b08a26SMatthew G. Knepley 
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
1912eac72dbSBarry Smith 
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
19558450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
196b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
19758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
199d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
200d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
2017d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
2024b0e389bSBarry Smith 
203640f4f53SPatrick Sanan /*E
204640f4f53SPatrick Sanan 
20506137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
206640f4f53SPatrick Sanan 
20767b8a455SSatish Balay   Values:
20887497f52SBarry Smith $ `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
20987497f52SBarry Smith $ `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
210640f4f53SPatrick Sanan 
2112b26979fSBarry Smith    Level: intermediate
212640f4f53SPatrick Sanan 
21387497f52SBarry Smith .seealso `:` `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
214640f4f53SPatrick Sanan E*/
215*9371c9d4SSatish Balay typedef enum {
216*9371c9d4SSatish Balay   KSP_FCD_TRUNC_TYPE_STANDARD,
217*9371c9d4SSatish Balay   KSP_FCD_TRUNC_TYPE_NOTAY
218*9371c9d4SSatish Balay } KSPFCDTruncationType;
21906137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
220640f4f53SPatrick Sanan 
221640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
222640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
223640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
224640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
22506137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
22606137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
227640f4f53SPatrick Sanan 
228390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
229390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
230390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
231390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
232390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
233390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
234390d8e47SPatrick Sanan 
235fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
236fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
237fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
238fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
239fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
240fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
241fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
242fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
24383f0b094SBarry Smith 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
247e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
2489f236934SBarry Smith 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
2541d73ed98SBarry Smith 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2571d73ed98SBarry Smith 
258483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
259483d6965SPatrick Sanan 
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
26358ad63f4SBarry Smith 
26404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
26504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
26604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
267568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
268e82af88dSprj- 
2696cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
2706cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
2716cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM)
272*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) {
273*9371c9d4SSatish Balay   return KSPHPDDMSetDeflationMat(ksp, U);
274*9371c9d4SSatish Balay }
275*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) {
276*9371c9d4SSatish Balay   return KSPHPDDMGetDeflationMat(ksp, U);
277*9371c9d4SSatish Balay }
2786cac28cbSPierre Jolivet #endif
279*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) {
280*9371c9d4SSatish Balay   return KSPMatSolve(ksp, B, X);
281*9371c9d4SSatish Balay }
282d55ab951SPierre Jolivet /*E
28387497f52SBarry Smith     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
284d55ab951SPierre Jolivet 
285d55ab951SPierre Jolivet     Level: intermediate
286d55ab951SPierre Jolivet 
287d55ab951SPierre Jolivet     Values:
28887497f52SBarry Smith $   `KSP_HPDDM_TYPE_GMRES` (default)
28987497f52SBarry Smith $   `KSP_HPDDM_TYPE_BGMRES`
29087497f52SBarry Smith $   `KSP_HPDDM_TYPE_CG`
29187497f52SBarry Smith $   `KSP_HPDDM_TYPE_BCG`
29287497f52SBarry Smith $   `KSP_HPDDM_TYPE_GCRODR`
29387497f52SBarry Smith $   `KSP_HPDDM_TYPE_BGCRODR`
29487497f52SBarry Smith $   `KSP_HPDDM_TYPE_BFBCG`
29587497f52SBarry Smith $   `KSP_HPDDM_TYPE_PREONLY`
296d55ab951SPierre Jolivet 
297db781477SPatrick Sanan .seealso: `KSPHPDDM`, `KSPHPDDMSetType()`
298d55ab951SPierre Jolivet E*/
299d55ab951SPierre Jolivet typedef enum {
300d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GMRES   = 0,
301d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGMRES  = 1,
302d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_CG      = 2,
303d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BCG     = 3,
304d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_GCRODR  = 4,
305d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BGCRODR = 5,
306d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_BFBCG   = 6,
307d55ab951SPierre Jolivet   KSP_HPDDM_TYPE_PREONLY = 7
308d55ab951SPierre Jolivet } KSPHPDDMType;
309d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[];
3102dd49c90SPierre Jolivet /*E
31187497f52SBarry Smith     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
3122dd49c90SPierre Jolivet 
3132dd49c90SPierre Jolivet     Level: intermediate
3142dd49c90SPierre Jolivet 
3152dd49c90SPierre Jolivet     Values:
31687497f52SBarry Smith $   `KSP_HPDDM_PRECISION_HALF` (currently unsupported)
31787497f52SBarry Smith $   `KSP_HPDDM_PRECISION_SINGLE` (default when PETSc is configured --with-precision=single)
31887497f52SBarry Smith $   `KSP_HPDDM_PRECISION_DOUBLE` (default when PETSc is configured --with-precision=double)
31987497f52SBarry Smith $   `KSP_HPDDM_PRECISION_QUADRUPLE` (default when PETSc is configured --with-precision=__float128)
3202dd49c90SPierre Jolivet 
321db781477SPatrick Sanan .seealso: `KSPHPDDM`
3222dd49c90SPierre Jolivet E*/
3232dd49c90SPierre Jolivet typedef enum {
3242dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_HALF      = 0,
3252dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_SINGLE    = 1,
3262dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_DOUBLE    = 2,
3272dd49c90SPierre Jolivet   KSP_HPDDM_PRECISION_QUADRUPLE = 3
3282dd49c90SPierre Jolivet } KSPHPDDMPrecision;
329d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
330d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
331e82af88dSprj- 
332b49fd9e1SBarry Smith /*E
333b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
334b49fd9e1SBarry Smith 
335b49fd9e1SBarry Smith    Level: advanced
336b49fd9e1SBarry Smith 
337db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
338db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
339b49fd9e1SBarry Smith 
340b49fd9e1SBarry Smith E*/
341*9371c9d4SSatish Balay typedef enum {
342*9371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_NEVER,
343*9371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_IFNEEDED,
344*9371c9d4SSatish Balay   KSP_GMRES_CGS_REFINE_ALWAYS
345*9371c9d4SSatish Balay } KSPGMRESCGSRefinementType;
3466a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
3471f7e983dSSatish Balay /*MC
3481957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
3498c5b8ba0SBarry Smith 
3508c5b8ba0SBarry Smith    Level: advanced
3518c5b8ba0SBarry Smith 
35287497f52SBarry Smith    Note:
35387497f52SBarry Smith    Possibly unstable, but the fastest to compute
3548c5b8ba0SBarry Smith 
355db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
356db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
357db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3588c5b8ba0SBarry Smith M*/
3598c5b8ba0SBarry Smith 
3601f7e983dSSatish Balay /*MC
3618c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
3628c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
3638c5b8ba0SBarry Smith           poor orthogonality.
3648c5b8ba0SBarry Smith 
3658c5b8ba0SBarry Smith    Level: advanced
3668c5b8ba0SBarry Smith 
36787497f52SBarry Smith    Note: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
3688c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
3698c5b8ba0SBarry Smith 
370db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
371db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
372db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3738c5b8ba0SBarry Smith M*/
3748c5b8ba0SBarry Smith 
3751f7e983dSSatish Balay /*MC
3768c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
3778c5b8ba0SBarry Smith 
3788c5b8ba0SBarry Smith    Level: advanced
3798c5b8ba0SBarry Smith 
38087497f52SBarry Smith    Notes:
38187497f52SBarry Smith    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
38287497f52SBarry Smith    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
3838c5b8ba0SBarry Smith 
3848c5b8ba0SBarry Smith    You should only use this if you absolutely know that the iterative refinement is needed.
3858c5b8ba0SBarry Smith 
386db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
387db781477SPatrick Sanan           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
388db781477SPatrick Sanan           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
3898c5b8ba0SBarry Smith M*/
3908c5b8ba0SBarry Smith 
391014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
392014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
39308480c60SBarry Smith 
394014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
395014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
397c38d4ed2SBarry Smith 
398014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
401121fd945SKris Buschelman 
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
405e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
406d9492815SBarry Smith 
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
40887d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
4102eac72dbSBarry Smith 
411798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
412798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
413798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
414798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
415798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
416798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
417798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
418798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
419798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
420798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
421798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
422798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
423798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
424798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
425798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
426798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
427798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
428798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
429798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
430798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
431798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
432fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
433798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
434*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {
435*9371c9d4SSatish Balay   return KSPMonitorResidual(ksp, n, rnorm, vf);
436*9371c9d4SSatish Balay }
437*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {
438*9371c9d4SSatish Balay   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
439*9371c9d4SSatish Balay }
440*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {
441*9371c9d4SSatish Balay   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
442*9371c9d4SSatish Balay }
443798534f6SMatthew G. Knepley 
444798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
445390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp, PetscInt its, PetscReal fnorm, void *dummy);
446390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
447e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
448e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
449e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
45084cb2905SBarry Smith 
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
453c01c455dSBarry Smith 
45423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
45523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
4601eb62cbbSBarry Smith 
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
4651f7f0c4fSBarry Smith 
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
46755849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
468fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
46919a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
470c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
4711b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
472c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
47394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
4741b2b9847SBarry Smith 
475*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) {
476*9371c9d4SSatish Balay   return KSPConvergedReasonView(ksp, v);
477*9371c9d4SSatish Balay }
478*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {
479*9371c9d4SSatish Balay   return KSPConvergedReasonViewFromOptions(ksp);
480*9371c9d4SSatish Balay }
48155849f57SBarry Smith 
48255849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
4831eb62cbbSBarry Smith 
484dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
4850e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
487884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
488798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
489798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
490798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
491db9b2ab1SHong Zhang 
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
49468ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
49583ab6a24SBarry Smith 
49628ce4d24SBarry Smith /*E
4978a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
4988a4b9c5cSBarry Smith        test routines.
4998a4b9c5cSBarry Smith 
5008a4b9c5cSBarry Smith    Level: advanced
5018a4b9c5cSBarry Smith 
502a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
50387497f52SBarry Smith    depending on left or right preconditioning, see `KSPSetPCSide()`
504a3f661c8SBarry Smith 
50587497f52SBarry Smith    Developer Note:
50695452b02SPatrick Sanan     this must match petsc/finclude/petscksp.h
5078a4b9c5cSBarry Smith 
508db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
509db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
5108a4b9c5cSBarry Smith E*/
511*9371c9d4SSatish Balay typedef enum {
512*9371c9d4SSatish Balay   KSP_NORM_DEFAULT          = -1,
513*9371c9d4SSatish Balay   KSP_NORM_NONE             = 0,
514*9371c9d4SSatish Balay   KSP_NORM_PRECONDITIONED   = 1,
515*9371c9d4SSatish Balay   KSP_NORM_UNPRECONDITIONED = 2,
516*9371c9d4SSatish Balay   KSP_NORM_NATURAL          = 3
517*9371c9d4SSatish Balay } KSPNormType;
5189e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
519014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes;
520a21b2a99SBarry Smith 
5211f7e983dSSatish Balay /*MC
5229793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
5238c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
5248c5b8ba0SBarry Smith           be based on a norm of a residual etc.
5258c5b8ba0SBarry Smith 
5268c5b8ba0SBarry Smith    Level: advanced
5278c5b8ba0SBarry Smith 
52887497f52SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
5298c5b8ba0SBarry Smith 
530db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
5318c5b8ba0SBarry Smith M*/
5328c5b8ba0SBarry Smith 
5331f7e983dSSatish Balay /*MC
5341957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
5358c5b8ba0SBarry Smith        convergence test routine.
5368c5b8ba0SBarry Smith 
5378c5b8ba0SBarry Smith    Level: advanced
5388c5b8ba0SBarry Smith 
539db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
5408c5b8ba0SBarry Smith M*/
5418c5b8ba0SBarry Smith 
5421f7e983dSSatish Balay /*MC
543ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
5448c5b8ba0SBarry Smith        convergence test routine.
5458c5b8ba0SBarry Smith 
5468c5b8ba0SBarry Smith    Level: advanced
5478c5b8ba0SBarry Smith 
548db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
5498c5b8ba0SBarry Smith M*/
5508c5b8ba0SBarry Smith 
5511f7e983dSSatish Balay /*MC
552ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
55387497f52SBarry Smith        convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
5548c5b8ba0SBarry Smith 
5558c5b8ba0SBarry Smith    Level: advanced
5568c5b8ba0SBarry Smith 
557db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
5588c5b8ba0SBarry Smith M*/
5598c5b8ba0SBarry Smith 
560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
5658a4b9c5cSBarry Smith 
566f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
5678a4b9c5cSBarry Smith /*E
5681957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
56928ce4d24SBarry Smith 
57028ce4d24SBarry Smith    Level: beginner
57128ce4d24SBarry Smith 
57295452b02SPatrick Sanan    Notes:
57387497f52SBarry Smith     See `KSPGetConvergedReason()` for explanation of each value
57428ce4d24SBarry Smith 
57595452b02SPatrick Sanan    Developer Notes:
57687497f52SBarry Smith    This must match petsc/finclude/petscksp.h
5774d0a8057SBarry Smith 
57887497f52SBarry Smith    The string versions of these are `KSPConvergedReasons`; if you change
5794d0a8057SBarry Smith    any of the values here also change them that array of names.
58086c02ca4SBarry Smith 
581db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
58228ce4d24SBarry Smith E*/
583d15094e1SBarry Smith typedef enum { /* converged */
5849ecb1286SBarry Smith   KSP_CONVERGED_RTOL_NORMAL              = 1,
5859ecb1286SBarry Smith   KSP_CONVERGED_ATOL_NORMAL              = 9,
586d15094e1SBarry Smith   KSP_CONVERGED_RTOL                     = 2,
587d15094e1SBarry Smith   KSP_CONVERGED_ATOL                     = 3,
588b335793eSSatish Balay   KSP_CONVERGED_ITS                      = 4,
5898031f4adStmunson   KSP_CONVERGED_CG_NEG_CURVE             = 5,
5908031f4adStmunson   KSP_CONVERGED_CG_CONSTRAINED           = 6,
591329f5518SBarry Smith   KSP_CONVERGED_STEP_LENGTH              = 7,
592af3d8002SBarry Smith   KSP_CONVERGED_HAPPY_BREAKDOWN          = 8,
593d15094e1SBarry Smith   /* diverged */
594b3cc6726SBarry Smith   KSP_DIVERGED_NULL                      = -2,
595d15094e1SBarry Smith   KSP_DIVERGED_ITS                       = -3,
596d15094e1SBarry Smith   KSP_DIVERGED_DTOL                      = -4,
597d15094e1SBarry Smith   KSP_DIVERGED_BREAKDOWN                 = -5,
598b4ac9ba4SBarry Smith   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
599b4ac9ba4SBarry Smith   KSP_DIVERGED_NONSYMMETRIC              = -7,
600b4ac9ba4SBarry Smith   KSP_DIVERGED_INDEFINITE_PC             = -8,
6014d51c080SBarry Smith   KSP_DIVERGED_NANORINF                  = -9,
6026aee1118SBarry Smith   KSP_DIVERGED_INDEFINITE_MAT            = -10,
603c0decd05SBarry Smith   KSP_DIVERGED_PC_FAILED                 = -11,
604aa4d2078SSatish Balay   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
605d15094e1SBarry Smith 
606*9371c9d4SSatish Balay   KSP_CONVERGED_ITERATING = 0
607*9371c9d4SSatish Balay } KSPConvergedReason;
608014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons;
609d15094e1SBarry Smith 
610c838673bSBarry Smith /*MC
61187497f52SBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called
612c838673bSBarry Smith 
613c838673bSBarry Smith    Level: beginner
614c838673bSBarry Smith 
61587497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
616c838673bSBarry Smith    for left preconditioning it is the 2-norm of the preconditioned residual, and the
617c838673bSBarry Smith    2-norm of the residual for right preconditioning
618c838673bSBarry Smith 
61987497f52SBarry Smith    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
620f9fed41fSBarry Smith 
621db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
622c838673bSBarry Smith 
623c838673bSBarry Smith M*/
624c838673bSBarry Smith 
625c838673bSBarry Smith /*MC
626c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
627c838673bSBarry Smith 
628c838673bSBarry Smith    Level: beginner
629c838673bSBarry Smith 
63087497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
631c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
632c838673bSBarry Smith        2-norm of the residual for right preconditioning
633c838673bSBarry Smith 
63487497f52SBarry Smith    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
635c838673bSBarry Smith 
636db781477SPatrick Sanan .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
637c838673bSBarry Smith 
638c838673bSBarry Smith M*/
639c838673bSBarry Smith 
640c838673bSBarry Smith /*MC
641c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
642c838673bSBarry Smith 
643c838673bSBarry Smith    Level: beginner
644c838673bSBarry Smith 
64587497f52SBarry Smith    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
646c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
647c838673bSBarry Smith        2-norm of the residual for right preconditioning
648c838673bSBarry Smith 
649c838673bSBarry Smith    Level: beginner
650c838673bSBarry Smith 
651db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
652c838673bSBarry Smith 
653c838673bSBarry Smith M*/
654c838673bSBarry Smith 
655c838673bSBarry Smith /*MC
656c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
657c838673bSBarry Smith       reached
658c838673bSBarry Smith 
659c838673bSBarry Smith    Level: beginner
660c838673bSBarry Smith 
661db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
662c838673bSBarry Smith 
663c838673bSBarry Smith M*/
664c838673bSBarry Smith 
665c838673bSBarry Smith /*MC
66687497f52SBarry Smith      KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
66787497f52SBarry Smith            the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
6684d0a8057SBarry Smith            test routine is set in KSP.
669c838673bSBarry Smith 
670c838673bSBarry Smith    Level: beginner
671c838673bSBarry Smith 
672db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
673c838673bSBarry Smith 
674c838673bSBarry Smith M*/
675c838673bSBarry Smith 
676c838673bSBarry Smith /*MC
677c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
6781de96524SPierre Jolivet           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
67987497f52SBarry Smith           preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
6801de96524SPierre Jolivet           are colinear.
681c838673bSBarry Smith 
682c838673bSBarry Smith    Level: beginner
683c838673bSBarry Smith 
684db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
685c838673bSBarry Smith 
686c838673bSBarry Smith M*/
687c838673bSBarry Smith 
688c838673bSBarry Smith /*MC
68987497f52SBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
690c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
691c838673bSBarry Smith 
692c838673bSBarry Smith    Level: beginner
693c838673bSBarry Smith 
694db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
695c838673bSBarry Smith 
696c838673bSBarry Smith M*/
697c838673bSBarry Smith 
698c838673bSBarry Smith /*MC
699c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
70087497f52SBarry Smith         symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
701c838673bSBarry Smith 
702c838673bSBarry Smith    Level: beginner
703c838673bSBarry Smith 
704db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
705c838673bSBarry Smith 
706c838673bSBarry Smith M*/
707c838673bSBarry Smith 
708c838673bSBarry Smith /*MC
709c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
71087497f52SBarry Smith         positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
711c838673bSBarry Smith         be positive definite
712c838673bSBarry Smith 
713c838673bSBarry Smith    Level: beginner
714c838673bSBarry Smith 
71587497f52SBarry Smith      Note:
71687497f52SBarry Smith     This can happen with the `PCICC` preconditioner, use -pc_factor_shift_positive_definite to force
71787497f52SBarry Smith   the `PCICC` preconditioner to generate a positive definite preconditioner
718c838673bSBarry Smith 
719db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
720c838673bSBarry Smith 
721c838673bSBarry Smith M*/
722c838673bSBarry Smith 
723c838673bSBarry Smith /*MC
724c0decd05SBarry Smith      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
7259fc87aa7SBarry Smith      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
72687497f52SBarry Smith      such as `PCFIELDSPLIT`.
7279fc87aa7SBarry Smith 
7289fc87aa7SBarry Smith    Level: beginner
7299fc87aa7SBarry Smith 
73095452b02SPatrick Sanan     Notes:
73195452b02SPatrick Sanan     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
7329fc87aa7SBarry Smith 
733db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
7349fc87aa7SBarry Smith 
7359fc87aa7SBarry Smith M*/
7369fc87aa7SBarry Smith 
7379fc87aa7SBarry Smith /*MC
73887497f52SBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call `KSPGetConvergedReason()`
73987497f52SBarry Smith         while the `KSPSolve()` is still running.
740c838673bSBarry Smith 
741c838673bSBarry Smith    Level: beginner
742c838673bSBarry Smith 
743db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
744c838673bSBarry Smith 
745c838673bSBarry Smith M*/
746c838673bSBarry Smith 
747014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
748df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
749df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
7503ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
7518de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
7523eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
7538de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
7548de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
7558de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
7568de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
75754b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
7580059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
759014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
760c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
76194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
762abef13c0SSatish Balay 
763*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) { /* never called */
764*9371c9d4SSatish Balay }
7658ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
766*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) { /* never called */
767*9371c9d4SSatish Balay }
7688ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
769*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) { /* never called */
770*9371c9d4SSatish Balay }
7718ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
772*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) { /* never called */
773*9371c9d4SSatish Balay }
7748ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
775*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */
776*9371c9d4SSatish Balay }
7778ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
778*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) { /* never called */
779*9371c9d4SSatish Balay }
7808ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
7818ea1b3e6SJed Brown 
7820bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
783*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) {
784*9371c9d4SSatish Balay   return KSPComputeOperator(A, NULL, B);
785*9371c9d4SSatish Balay }
786d4fbbf0eSBarry Smith 
78728ce4d24SBarry Smith /*E
78828ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
78928ce4d24SBarry Smith 
79028ce4d24SBarry Smith    Level: beginner
79128ce4d24SBarry Smith 
792db781477SPatrick Sanan .seealso: `KSPCGSetType()`
79328ce4d24SBarry Smith E*/
794*9371c9d4SSatish Balay typedef enum {
795*9371c9d4SSatish Balay   KSP_CG_SYMMETRIC = 0,
796*9371c9d4SSatish Balay   KSP_CG_HERMITIAN = 1
797*9371c9d4SSatish Balay } KSPCGType;
7986a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
79928ce4d24SBarry Smith 
800014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
801014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
8028031f4adStmunson 
803dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
804dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
805dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
806fcae7a14Stmunson 
80705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
80805de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
809*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) {
810*9371c9d4SSatish Balay   return KSPGLTRGetMinEig(ksp, x);
811*9371c9d4SSatish Balay }
812*9371c9d4SSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) {
813*9371c9d4SSatish Balay   return KSPGLTRGetLambda(ksp, x);
814*9371c9d4SSatish Balay }
8158031f4adStmunson 
816014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
817ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
8181d6018f0SLisandro Dalcin 
819f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
820014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
821014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
8223369ce9aSBarry Smith 
8239804daf3SBarry Smith #include <petscdrawtypes.h>
824014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
8252f2e5d10SKris Buschelman 
826014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
827014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
82803919abeSBarry Smith 
829ba36735cSStefano Zampini /*S
830ba36735cSStefano Zampini      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
831f8a50e2bSBarry Smith 
832ba36735cSStefano Zampini    Level: beginner
833f8a50e2bSBarry Smith 
834db781477SPatrick Sanan .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
835ba36735cSStefano Zampini S*/
836ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess;
837ba36735cSStefano Zampini /*J
83887497f52SBarry Smith     KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
839ba36735cSStefano Zampini 
840ba36735cSStefano Zampini    Level: beginner
841ba36735cSStefano Zampini 
842db781477SPatrick Sanan .seealso: `KSPGuess`
843ba36735cSStefano Zampini J*/
844ba36735cSStefano Zampini typedef const char         *KSPGuessType;
845ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer"
846ba36735cSStefano Zampini #define KSPGUESSPOD     "pod"
8471d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
848ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
849ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
850ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
851ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
852ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
853ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
854ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
8558410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
856ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
857ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
858ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
859ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
860ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
861014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
862ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
863ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
864f8a50e2bSBarry Smith 
865470b340bSDmitry Karpeev /*E
866470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
867470b340bSDmitry Karpeev 
868470b340bSDmitry Karpeev     Level: intermediate
869470b340bSDmitry Karpeev 
870db781477SPatrick Sanan .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()`
871470b340bSDmitry Karpeev E*/
872*9371c9d4SSatish Balay typedef enum {
873*9371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
874*9371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
875*9371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
876*9371c9d4SSatish Balay   MAT_SCHUR_COMPLEMENT_AINV_FULL
877*9371c9d4SSatish Balay } MatSchurComplementAinvType;
878470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
879470b340bSDmitry Karpeev 
880*9371c9d4SSatish Balay typedef enum {
881*9371c9d4SSatish Balay   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
882864588a7SAlp Dener   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
883864588a7SAlp Dener   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
884*9371c9d4SSatish Balay   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
885*9371c9d4SSatish Balay } MatLMVMSymBroydenScaleType;
886864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
88792f76d53SAlp Dener 
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
890d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
891bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
892aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
893bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
894470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
895470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
8965bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
8975a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
898470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
899470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
9003f22127dSBarry Smith 
90178e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
90278e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
90378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
904864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
905864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
906864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
907864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
908864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
909cd929ea3SAlp Dener 
910cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
911b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
912cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
913cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
914e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
9150ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
916cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
917cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
918cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
919cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
920cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
9212d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
922cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
923cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
924cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
925cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
92692f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
927cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
928cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
929864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
930864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
931cd929ea3SAlp Dener 
932014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
933014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
934014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
935014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
936014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
937fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
93823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
939fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
94023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
94123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
942014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
943014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
944fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
945fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
9466c699258SBarry Smith 
94702b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
948*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
949557cf195SMatthew G. Knepley 
9502b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
9512b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
9522eac72dbSBarry Smith #endif
953