xref: /petsc/include/petscksp.h (revision 6d5076fa512d0d8e4e12a3685978afdbfee2b1a4)
1 /*
2    Defines the interface functions for the Krylov subspace accelerators.
3 */
4 #ifndef PETSCKSP_H
5 #define PETSCKSP_H
6 
7 #include <petscpc.h>
8 
9 /* SUBMANSEC = KSP */
10 
11 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
12 
13 /*S
14      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
15          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
16 
17    Level: beginner
18 
19    Note:
20     When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a
21     `KSPType` of `KSPPREONLY`, meaning that only application of the preconditioner is used as the linear solver.
22 
23 .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
24 S*/
25 typedef struct _p_KSP *KSP;
26 
27 /*J
28     KSPType - String with the name of a PETSc Krylov method.
29 
30    Level: beginner
31 
32 .seealso: [](chapter_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
33 J*/
34 typedef const char *KSPType;
35 #define KSPRICHARDSON "richardson"
36 #define KSPCHEBYSHEV  "chebyshev"
37 #define KSPCG         "cg"
38 #define KSPGROPPCG    "groppcg"
39 #define KSPPIPECG     "pipecg"
40 #define KSPPIPECGRR   "pipecgrr"
41 #define KSPPIPELCG    "pipelcg"
42 #define KSPPIPEPRCG   "pipeprcg"
43 #define KSPPIPECG2    "pipecg2"
44 #define KSPCGNE       "cgne"
45 #define KSPNASH       "nash"
46 #define KSPSTCG       "stcg"
47 #define KSPGLTR       "gltr"
48 #define KSPCGNASH     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash"
49 #define KSPCGSTCG     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg"
50 #define KSPCGGLTR     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
51 #define KSPFCG        "fcg"
52 #define KSPPIPEFCG    "pipefcg"
53 #define KSPGMRES      "gmres"
54 #define KSPPIPEFGMRES "pipefgmres"
55 #define KSPFGMRES     "fgmres"
56 #define KSPLGMRES     "lgmres"
57 #define KSPDGMRES     "dgmres"
58 #define KSPPGMRES     "pgmres"
59 #define KSPTCQMR      "tcqmr"
60 #define KSPBCGS       "bcgs"
61 #define KSPIBCGS      "ibcgs"
62 #define KSPQMRCGS     "qmrcgs"
63 #define KSPFBCGS      "fbcgs"
64 #define KSPFBCGSR     "fbcgsr"
65 #define KSPBCGSL      "bcgsl"
66 #define KSPPIPEBCGS   "pipebcgs"
67 #define KSPCGS        "cgs"
68 #define KSPTFQMR      "tfqmr"
69 #define KSPCR         "cr"
70 #define KSPPIPECR     "pipecr"
71 #define KSPLSQR       "lsqr"
72 #define KSPPREONLY    "preonly"
73 #define KSPNONE       "none"
74 #define KSPQCG        "qcg"
75 #define KSPBICG       "bicg"
76 #define KSPMINRES     "minres"
77 #define KSPSYMMLQ     "symmlq"
78 #define KSPLCD        "lcd"
79 #define KSPPYTHON     "python"
80 #define KSPGCR        "gcr"
81 #define KSPPIPEGCR    "pipegcr"
82 #define KSPTSIRM      "tsirm"
83 #define KSPCGLS       "cgls"
84 #define KSPFETIDP     "fetidp"
85 #define KSPHPDDM      "hpddm"
86 
87 /* Logging support */
88 PETSC_EXTERN PetscClassId KSP_CLASSID;
89 PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
90 PETSC_EXTERN PetscClassId DMKSP_CLASSID;
91 
92 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
93 PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
94 PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
95 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
96 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
97 PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
98 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
99 PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
100 PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
101 PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
102 PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
103 {
104   return KSPSetMatSolveBatchSize(ksp, n);
105 }
106 PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
107 PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
108 {
109   return KSPGetMatSolveBatchSize(ksp, n);
110 }
111 PETSC_EXTERN PetscErrorCode KSPReset(KSP);
112 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
113 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
114 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
115 PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
116 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
117 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
118 
119 PETSC_EXTERN PetscFunctionList KSPList;
120 PETSC_EXTERN PetscFunctionList KSPGuessList;
121 PETSC_EXTERN PetscFunctionList KSPMonitorList;
122 PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
123 PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
124 PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
125 PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
126 
127 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
128 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
129 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
130 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
131 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
132 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
133 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
134 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
135 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
136 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
137 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
138 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
139 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
140 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
141 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
142 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
143 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
144 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
145 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
146 PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
147 {
148   return KSPCreateVecs(ksp, n, x, m, y);
149 }
150 
151 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
152 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
153 
154 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
155 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
156 
157 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
158 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
159 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
160 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
161 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
162 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool);
163 PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
164 PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool);
165 
166 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
167 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
168 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
169 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
170 
171 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
172 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
173 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
174 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
175 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
176 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
177 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
178 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
179 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
180 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
181 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
182 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
183 PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
184 /*
185   PCMGCoarseList contains the list of coarse space constructor currently registered
186   These are added with PCMGRegisterCoarseSpaceConstructor()
187 */
188 PETSC_EXTERN PetscFunctionList PCMGCoarseList;
189 PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
190 PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
191 
192 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
193 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
194 
195 /*E
196   KSPChebyshevKind - Which Chebyshev polynomial to use
197 
198   Values:
199 + `KSP_CHEBYSHEV_FIRST`      - "classic" first-kind Chebyshev polynomial
200 . `KSP_CHEBYSHEV_FOURTH`     - fourth-kind Chebyshev polynomial
201 - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial
202 
203   Level: intermediate
204 
205 .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind`
206 E*/
207 typedef enum {
208   KSP_CHEBYSHEV_FIRST,
209   KSP_CHEBYSHEV_FOURTH,
210   KSP_CHEBYSHEV_OPT_FOURTH
211 } KSPChebyshevKind;
212 
213 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
214 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
215 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
216 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
217 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
218 PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind);
219 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
220 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
221 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
222 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
223 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
224 
225 /*E
226 
227   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
228 
229   Values:
230 + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
231 - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
232 
233    Level: intermediate
234 
235 .seealso: [](chapter_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
236 E*/
237 typedef enum {
238   KSP_FCD_TRUNC_TYPE_STANDARD,
239   KSP_FCD_TRUNC_TYPE_NOTAY
240 } KSPFCDTruncationType;
241 PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
242 
243 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
244 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
245 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
246 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
247 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
248 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
249 
250 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
251 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
252 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
253 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
254 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
255 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
256 
257 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
258 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
259 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
260 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
261 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
262 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
263 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
264 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
265 
266 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
267 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
268 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
269 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
270 
271 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
272 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
273 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
274 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
275 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
276 
277 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
278 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
279 
280 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
281 
282 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
283 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
284 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
285 
286 PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
287 PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
288 PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);
289 
290 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
291 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
292 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
293 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
294 
295 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
296 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
297 #if PetscDefined(HAVE_HPDDM)
298 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
299 {
300   return KSPHPDDMSetDeflationMat(ksp, U);
301 }
302 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
303 {
304   return KSPHPDDMGetDeflationMat(ksp, U);
305 }
306 #endif
307 PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
308 {
309   return KSPMatSolve(ksp, B, X);
310 }
311 /*E
312     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
313 
314     Level: intermediate
315 
316     Values:
317 +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
318 .   `KSP_HPDDM_TYPE_BGMRES` - block GMRES
319 .   `KSP_HPDDM_TYPE_CG` - Conjugate Gradient
320 .   `KSP_HPDDM_TYPE_BCG` - block CG
321 .   `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
322 .   `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR
323 .   `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG
324 -   `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only
325 
326 .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
327 E*/
328 typedef enum {
329   KSP_HPDDM_TYPE_GMRES   = 0,
330   KSP_HPDDM_TYPE_BGMRES  = 1,
331   KSP_HPDDM_TYPE_CG      = 2,
332   KSP_HPDDM_TYPE_BCG     = 3,
333   KSP_HPDDM_TYPE_GCRODR  = 4,
334   KSP_HPDDM_TYPE_BGCRODR = 5,
335   KSP_HPDDM_TYPE_BFBCG   = 6,
336   KSP_HPDDM_TYPE_PREONLY = 7
337 } KSPHPDDMType;
338 PETSC_EXTERN const char *const KSPHPDDMTypes[];
339 
340 /*E
341     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
342 
343     Level: intermediate
344 
345     Values:
346 +   `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16`
347 .   `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single`
348 .   `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double`
349 -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
350 
351 .seealso: [](chapter_ksp), `KSP`, `KSPHPDDM`
352 E*/
353 typedef enum {
354   KSP_HPDDM_PRECISION_HALF      = 0,
355   KSP_HPDDM_PRECISION_SINGLE    = 1,
356   KSP_HPDDM_PRECISION_DOUBLE    = 2,
357   KSP_HPDDM_PRECISION_QUADRUPLE = 3
358 } KSPHPDDMPrecision;
359 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
360 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
361 
362 /*E
363     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
364 
365    Level: advanced
366 
367     Values:
368 +   `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt
369 .   `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
370 -   `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps
371 
372 .seealso: [](chapter_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
373           `KSPGMRESGetOrthogonalization()`,
374           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
375 E*/
376 typedef enum {
377   KSP_GMRES_CGS_REFINE_NEVER,
378   KSP_GMRES_CGS_REFINE_IFNEEDED,
379   KSP_GMRES_CGS_REFINE_ALWAYS
380 } KSPGMRESCGSRefinementType;
381 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
382 /*MC
383     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
384 
385    Level: advanced
386 
387    Note:
388    Possibly unstable, but the fastest to compute
389 
390 .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
391           `KSP`, `KSPGMRESGetOrthogonalization()`,
392           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
393           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
394 M*/
395 
396 /*MC
397     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
398           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
399           poor orthogonality.
400 
401    Level: advanced
402 
403    Note:
404    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
405    estimate the orthogonality but is more stable.
406 
407 .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
408           `KSP`, `KSPGMRESGetOrthogonalization()`,
409           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
410           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
411 M*/
412 
413 /*MC
414     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
415 
416    Level: advanced
417 
418    Notes:
419    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
420    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
421 
422    You should only use this if you absolutely know that the iterative refinement is needed.
423 
424 .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
425           `KSP`, `KSPGMRESGetOrthogonalization()`,
426           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
427           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
428 M*/
429 
430 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
431 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
432 
433 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
434 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
435 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
436 
437 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
438 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
439 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
440 
441 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
442 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
443 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
444 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
445 
446 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
447 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
448 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
449 
450 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
451 PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
452 PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
453 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
454 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
455 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
456 PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
457 PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
458 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
459 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
460 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
461 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
462 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
463 PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
464 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
465 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
466 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
467 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
468 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
469 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
470 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
471 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
472 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
473 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
474 {
475   return KSPMonitorResidual(ksp, n, rnorm, vf);
476 }
477 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
478 {
479   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
480 }
481 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
482 {
483   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
484 }
485 
486 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
487 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
488 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
489 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
490 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
491 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
492 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
493 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
494 
495 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
496 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
497 
498 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
499 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
500 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
501 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
502 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
503 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
504 
505 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
506 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
507 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
508 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
509 
510 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
511 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
512 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
513 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
514 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
515 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
516 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
517 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
518 
519 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
520 {
521   return KSPConvergedReasonView(ksp, v);
522 }
523 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
524 {
525   return KSPConvergedReasonViewFromOptions(ksp);
526 }
527 
528 #define KSP_FILE_CLASSID 1211223
529 
530 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
531 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
532 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
533 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
534 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
535 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
536 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
537 
538 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
539 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
540 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
541 
542 /*E
543     KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
544        test routines.
545 
546     Values:
547 +   `KSP_NORM_DEFAULT` - use the default for the current `KSPType`
548 .   `KSP_NORM_NONE` - use no norm calculation
549 .   `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm
550 .   `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
551 -   `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator)
552 
553    Level: advanced
554 
555    Note:
556    Each solver only supports a subset of these and some may support different ones
557    depending on left or right preconditioning, see `KSPSetPCSide()`
558 
559    Developer Note:
560    This must match the values in petsc/finclude/petscksp.h
561 
562 .seealso: [](chapter_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
563           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
564 E*/
565 typedef enum {
566   KSP_NORM_DEFAULT          = -1,
567   KSP_NORM_NONE             = 0,
568   KSP_NORM_PRECONDITIONED   = 1,
569   KSP_NORM_UNPRECONDITIONED = 2,
570   KSP_NORM_NATURAL          = 3
571 } KSPNormType;
572 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
573 PETSC_EXTERN const char *const *const KSPNormTypes;
574 
575 /*MC
576     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
577           possibly save some computation but means the convergence test cannot
578           be based on a norm of a residual etc.
579 
580    Level: advanced
581 
582     Note:
583     Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
584 
585 .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
586 M*/
587 
588 /*MC
589     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
590        convergence test routine.
591 
592    Level: advanced
593 
594 .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
595 M*/
596 
597 /*MC
598     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
599        convergence test routine.
600 
601    Level: advanced
602 
603 .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
604 M*/
605 
606 /*MC
607     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
608        convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
609 
610    Level: advanced
611 
612 .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
613 M*/
614 
615 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
616 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
617 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
618 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
619 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
620 
621 #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_NEG_CURVE (since version 3.19)")
622 #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_STEP_LENGTH (since version 3.19)")
623 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
624 /*E
625     KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
626 
627    Level: beginner
628 
629    Values:
630 +  `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations
631 .  `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations
632 .  `KSP_CONVERGED_RTOL` - requested decrease in the residual
633 .  `KSP_CONVERGED_ATOL` - requested absolute value in the residual
634 .  `KSP_CONVERGED_ITS` - requested number of iterations
635 .  `KSP_CONVERGED_NEG_CURVE` - see note below
636 .  `KSP_CONVERGED_STEP_LENGTH` - see note below
637 .  `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred.
638 .  `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES
639 .  `KSP_DIVERGED_ITS` - requested number of iterations
640 .  `KSP_DIVERGED_DTOL` - large increase in the residual norm
641 .  `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method
642 .  `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method
643 .  `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
644 .  `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite
645 .  `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation
646 .  `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite
647 -  `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason
648 
649    Note:
650    The values  `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`,
651    `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
652 
653    Developer Notes:
654    This must match the values in petsc/finclude/petscksp.h
655 
656    The string versions of these are `KSPConvergedReasons`; if you change
657    any of the values here also change them that array of names.
658 
659 .seealso: [](chapter_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
660 E*/
661 typedef enum { /* converged */
662   KSP_CONVERGED_RTOL_NORMAL               = 1,
663   KSP_CONVERGED_ATOL_NORMAL               = 9,
664   KSP_CONVERGED_RTOL                      = 2,
665   KSP_CONVERGED_ATOL                      = 3,
666   KSP_CONVERGED_ITS                       = 4,
667   KSP_CONVERGED_NEG_CURVE                 = 5,
668   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
669   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
670   KSP_CONVERGED_STEP_LENGTH               = 6,
671   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
672   /* diverged */
673   KSP_DIVERGED_NULL                      = -2,
674   KSP_DIVERGED_ITS                       = -3,
675   KSP_DIVERGED_DTOL                      = -4,
676   KSP_DIVERGED_BREAKDOWN                 = -5,
677   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
678   KSP_DIVERGED_NONSYMMETRIC              = -7,
679   KSP_DIVERGED_INDEFINITE_PC             = -8,
680   KSP_DIVERGED_NANORINF                  = -9,
681   KSP_DIVERGED_INDEFINITE_MAT            = -10,
682   KSP_DIVERGED_PC_FAILED                 = -11,
683   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
684 
685   KSP_CONVERGED_ITERATING = 0
686 } KSPConvergedReason;
687 PETSC_EXTERN const char *const *KSPConvergedReasons;
688 
689 /*MC
690      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called
691 
692    Level: beginner
693 
694    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
695    for left preconditioning it is the 2-norm of the preconditioned residual, and the
696    2-norm of the residual for right preconditioning
697 
698    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
699 
700 .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
701 
702 M*/
703 
704 /*MC
705      KSP_CONVERGED_ATOL - norm(r) <= atol
706 
707    Level: beginner
708 
709    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
710        for left preconditioning it is the 2-norm of the preconditioned residual, and the
711        2-norm of the residual for right preconditioning
712 
713    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
714 
715 .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
716 
717 M*/
718 
719 /*MC
720      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
721 
722    Level: beginner
723 
724    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
725        for left preconditioning it is the 2-norm of the preconditioned residual, and the
726        2-norm of the residual for right preconditioning
727 
728    Level: beginner
729 
730 .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
731 
732 M*/
733 
734 /*MC
735      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
736       reached
737 
738    Level: beginner
739 
740 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
741 
742 M*/
743 
744 /*MC
745      KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
746            the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
747            test routine is set in KSP.
748 
749    Level: beginner
750 
751 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
752 
753 M*/
754 
755 /*MC
756      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
757           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
758           preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
759           are colinear.
760 
761    Level: beginner
762 
763 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
764 
765 M*/
766 
767 /*MC
768      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
769           method could not continue to enlarge the Krylov space.
770 
771    Level: beginner
772 
773 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
774 
775 M*/
776 
777 /*MC
778      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
779         symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
780 
781    Level: beginner
782 
783 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
784 
785 M*/
786 
787 /*MC
788      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
789         positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
790         be positive definite
791 
792    Level: beginner
793 
794      Note:
795     This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
796   the `PCICC` preconditioner to generate a positive definite preconditioner
797 
798 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
799 
800 M*/
801 
802 /*MC
803      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
804      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
805      such as `PCFIELDSPLIT`.
806 
807    Level: beginner
808 
809     Note:
810     Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
811 
812 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
813 
814 M*/
815 
816 /*MC
817      KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
818         while `KSPSolve()` is still running.
819 
820    Level: beginner
821 
822 .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
823 
824 M*/
825 
826 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
827 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
828 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
829 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
830 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
831 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
832 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
833 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
834 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
835 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
836 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
837 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
838 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
839 PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
840 PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
841 PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
842 PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
843 
844 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void)
845 { /* never called */
846 }
847 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
848 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void)
849 { /* never called */
850 }
851 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
852 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void)
853 { /* never called */
854 }
855 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
856 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void)
857 { /* never called */
858 }
859 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
860 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void)
861 { /* never called */
862 }
863 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
864 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void)
865 { /* never called */
866 }
867 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
868 
869 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
870 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
871 {
872   return KSPComputeOperator(A, NULL, B);
873 }
874 
875 /*E
876     KSPCGType - Determines what type of `KSPCG` to use
877 
878    Level: beginner
879 
880    Values:
881  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
882  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
883 
884 .seealso: [](chapter_ksp), `KSP`, `KSPCGSetType()`
885 E*/
886 typedef enum {
887   KSP_CG_SYMMETRIC = 0,
888   KSP_CG_HERMITIAN = 1
889 } KSPCGType;
890 PETSC_EXTERN const char *const KSPCGTypes[];
891 
892 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
893 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
894 
895 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
896 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
897 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
898 
899 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
900 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
901 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
902 {
903   return KSPGLTRGetMinEig(ksp, x);
904 }
905 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
906 {
907   return KSPGLTRGetLambda(ksp, x);
908 }
909 
910 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
911 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
912 
913 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
914 PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
915 PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
916 
917 #include <petscdrawtypes.h>
918 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
919 
920 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
921 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
922 
923 /*S
924      KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
925 
926    Level: intermediate
927 
928 .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
929 S*/
930 typedef struct _p_KSPGuess *KSPGuess;
931 /*J
932     KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
933 
934    Level: intermediate
935 
936    Values:
937  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
938  - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition
939 
940 .seealso: [](chapter_ksp), `KSP`, `KSPGuess`
941 J*/
942 typedef const char *KSPGuessType;
943 #define KSPGUESSFISCHER "fischer"
944 #define KSPGUESSPOD     "pod"
945 
946 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
947 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
948 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
949 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
950 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
951 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
952 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
953 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
954 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
955 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
956 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
957 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
958 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
959 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
960 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
961 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
962 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
963 
964 /*E
965     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
966 
967     Level: intermediate
968 
969 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
970           `MatCreateSchurComplementPmat()`
971 E*/
972 typedef enum {
973   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
974   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
975   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
976   MAT_SCHUR_COMPLEMENT_AINV_FULL
977 } MatSchurComplementAinvType;
978 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
979 
980 typedef enum {
981   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
982   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
983   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
984   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
985 } MatLMVMSymBroydenScaleType;
986 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
987 
988 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
989 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
990 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
991 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
992 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
993 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
994 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
995 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
996 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
997 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
998 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
999 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
1000 
1001 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1002 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1003 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1004 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1005 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1006 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1007 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1008 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1009 
1010 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1011 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1012 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1013 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1014 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1015 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1016 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1017 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1018 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1019 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1020 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1021 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1022 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1023 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1024 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1025 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1026 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1027 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1028 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1029 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1030 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1031 
1032 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1033 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1034 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1035 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1036 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1037 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
1038 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1039 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
1040 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1041 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
1042 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1043 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1044 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1045 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1046 
1047 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1048 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);
1049 
1050 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1051 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1052 #endif
1053