xref: /petsc/include/petscksp.h (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1 /*
2    Defines the interface functions for the Krylov subspace accelerators.
3 */
4 #pragma once
5 
6 #include <petscpc.h>
7 
8 /* SUBMANSEC = KSP */
9 
10 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
11 PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void);
12 
13 /*S
14    KSP - Abstract PETSc object that manages the linear solves in PETSc (even those such as direct factorization-based solvers that
15          do not use Krylov accelerators).
16 
17    Level: beginner
18 
19    Notes:
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` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver.
22 
23    Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object
24 
25    The `PC` object is used to control preconditioners in PETSc.
26 
27   `KSP` can also be used to solve some least squares problems (over or under-determined linear systems), using, for example, `KSPLSQR`, see `PETSCREGRESSORLINEAR`
28   for additional methods that can be used to solve least squares problems and other linear regressions).
29 
30 .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
31 S*/
32 typedef struct _p_KSP *KSP;
33 
34 /*J
35    KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides.
36 
37    Level: beginner
38 
39 .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
40 J*/
41 typedef const char *KSPType;
42 #define KSPRICHARDSON "richardson"
43 #define KSPCHEBYSHEV  "chebyshev"
44 #define KSPCG         "cg"
45 #define KSPGROPPCG    "groppcg"
46 #define KSPPIPECG     "pipecg"
47 #define KSPPIPECGRR   "pipecgrr"
48 #define KSPPIPELCG    "pipelcg"
49 #define KSPPIPEPRCG   "pipeprcg"
50 #define KSPPIPECG2    "pipecg2"
51 #define KSPCGNE       "cgne"
52 #define KSPNASH       "nash"
53 #define KSPSTCG       "stcg"
54 #define KSPGLTR       "gltr"
55 #define KSPCGNASH     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash"
56 #define KSPCGSTCG     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg"
57 #define KSPCGGLTR     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr"
58 #define KSPFCG        "fcg"
59 #define KSPPIPEFCG    "pipefcg"
60 #define KSPGMRES      "gmres"
61 #define KSPPIPEFGMRES "pipefgmres"
62 #define KSPFGMRES     "fgmres"
63 #define KSPLGMRES     "lgmres"
64 #define KSPDGMRES     "dgmres"
65 #define KSPPGMRES     "pgmres"
66 #define KSPTCQMR      "tcqmr"
67 #define KSPBCGS       "bcgs"
68 #define KSPIBCGS      "ibcgs"
69 #define KSPQMRCGS     "qmrcgs"
70 #define KSPFBCGS      "fbcgs"
71 #define KSPFBCGSR     "fbcgsr"
72 #define KSPBCGSL      "bcgsl"
73 #define KSPPIPEBCGS   "pipebcgs"
74 #define KSPCGS        "cgs"
75 #define KSPTFQMR      "tfqmr"
76 #define KSPCR         "cr"
77 #define KSPPIPECR     "pipecr"
78 #define KSPLSQR       "lsqr"
79 #define KSPPREONLY    "preonly"
80 #define KSPNONE       "none"
81 #define KSPQCG        "qcg"
82 #define KSPBICG       "bicg"
83 #define KSPMINRES     "minres"
84 #define KSPSYMMLQ     "symmlq"
85 #define KSPLCD        "lcd"
86 #define KSPPYTHON     "python"
87 #define KSPGCR        "gcr"
88 #define KSPPIPEGCR    "pipegcr"
89 #define KSPTSIRM      "tsirm"
90 #define KSPCGLS       "cgls"
91 #define KSPFETIDP     "fetidp"
92 #define KSPHPDDM      "hpddm"
93 
94 /* Logging support */
95 PETSC_EXTERN PetscClassId KSP_CLASSID;
96 PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
97 PETSC_EXTERN PetscClassId DMKSP_CLASSID;
98 
99 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
100 PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
101 PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
102 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
103 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
104 PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
105 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
106 PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
107 PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
108 PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat);
109 PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
110 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
111 {
112   return KSPSetMatSolveBatchSize(ksp, n);
113 }
114 PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
115 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
116 {
117   return KSPGetMatSolveBatchSize(ksp, n);
118 }
119 PETSC_EXTERN PetscErrorCode KSPReset(KSP);
120 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
121 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
122 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
123 PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
124 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
125 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
126 
127 PETSC_EXTERN PetscFunctionList KSPList;
128 PETSC_EXTERN PetscFunctionList KSPGuessList;
129 PETSC_EXTERN PetscFunctionList KSPMonitorList;
130 PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
131 PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
132 PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
133 
134 /*S
135   KSPMonitorRegisterFn - A function prototype for functions provided to `KSPMonitorRegister()`
136 
137   Calling Sequence:
138 + ksp   - iterative solver obtained from `KSPCreate()`
139 . it    - iteration number
140 . rnorm - (estimated) 2-norm of (preconditioned) residual
141 - ctx   - `PetscViewerAndFormat` object
142 
143   Level: beginner
144 
145   Note:
146   This is a `KSPMonitorFn` specialized for a context of `PetscViewerAndFormat`
147 
148 .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn`, `KSPMonitorRegisterDestroyFn`
149 S*/
150 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterFn(KSP ksp, PetscInt it, PetscReal rnorm, PetscViewerAndFormat *ctx);
151 
152 /*S
153   KSPMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to `KSPMonitorRegister()`
154 
155   Calling Sequence:
156 + viewer - the viewer to be used with the `KSPMonitorRegisterFn`
157 . format - the format of the viewer
158 . ctx    - a context for the monitor
159 - result - a `PetscViewerAndFormat` object
160 
161   Level: beginner
162 
163 .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterDestroyFn`
164 S*/
165 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterCreateFn(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **result);
166 
167 /*S
168   KSPMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to `KSPMonitorRegister()`
169 
170   Calling Sequence:
171 . vf - a `PetscViewerAndFormat` object to be destroyed, including any context
172 
173   Level: beginner
174 
175 .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn`
176 S*/
177 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterDestroyFn(PetscViewerAndFormat **result);
178 
179 PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, KSPMonitorRegisterFn *, KSPMonitorRegisterCreateFn *, KSPMonitorRegisterDestroyFn *);
180 
181 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
182 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
183 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
184 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
185 PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt);
186 PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *);
187 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
188 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
189 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
190 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
191 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
192 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
193 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
194 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
195 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
196 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
197 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
198 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
199 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
200 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
201 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
202 PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
203 {
204   return KSPCreateVecs(ksp, n, x, m, y);
205 }
206 
207 /*S
208   KSPPSolveFn - A function prototype for functions provided to `KSPSetPreSolve()` and `KSPSetPostSolve()`
209 
210   Calling Sequence:
211 + ksp - the `KSP` context
212 . rhs - the right-hand side vector
213 . x   - the solution vector
214 - ctx - optional context that was provided with `KSPSetPreSolve()` or `KSPSetPostSolve()`
215 
216   Level: intermediate
217 
218 .seealso: [](ch_snes), `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`, `PCShellPSolveFn`
219 S*/
220 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPPSolveFn(KSP ksp, Vec rhs, Vec x, void *ctx);
221 
222 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, KSPPSolveFn *, void *);
223 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, KSPPSolveFn *, void *);
224 
225 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
226 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
227 PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt);
228 PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *);
229 
230 /*S
231   KSPMonitorFn - A function prototype for functions provided to `KSPMonitorSet()`
232 
233   Calling Sequence:
234 + ksp   - iterative solver obtained from `KSPCreate()`
235 . it    - iteration number
236 . rnorm - (estimated) 2-norm of (preconditioned) residual
237 - ctx   - optional monitoring context, as provided with `KSPMonitorSet()`
238 
239   Level: beginner
240 
241 .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()`
242 S*/
243 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorFn(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx);
244 
245 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
246 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, KSPMonitorFn *, void *, PetscCtxDestroyFn *);
247 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
248 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
249 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
250 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool);
251 PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
252 PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool);
253 
254 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
255 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
256 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
257 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
258 
259 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
260 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
261 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
262 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
263 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
264 PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]);
265 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
266 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
267 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
268 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
269 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
270 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
271 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
272 PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
273 
274 /*S
275   PCMGCoarseSpaceConstructorFn - A function prototype for functions registered with `PCMGRegisterCoarseSpaceConstructor()`
276 
277   Calling Sequence:
278 + pc        - The `PC` object
279 . l         - The multigrid level, 0 is the coarse level
280 . dm        - The `DM` for this level
281 . smooth    - The level smoother
282 . Nc        - The size of the coarse space
283 . initGuess - Basis for an initial guess for the space
284 - coarseSp  - A basis for the computed coarse space
285 
286   Level: beginner
287 
288 .seealso: [](ch_ksp), `PCMGRegisterCoarseSpaceConstructor()`, `PCMGGetCoarseSpaceConstructor()`
289 S*/
290 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCMGCoarseSpaceConstructorFn(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp);
291 
292 PETSC_EXTERN PetscFunctionList PCMGCoarseList;
293 PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn *);
294 PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn **);
295 
296 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
297 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
298 
299 /*E
300   KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV`
301 
302   Values:
303 + `KSP_CHEBYSHEV_FIRST`      - "classic" first-kind Chebyshev polynomial
304 . `KSP_CHEBYSHEV_FOURTH`     - fourth-kind Chebyshev polynomial
305 - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial
306 
307   Level: intermediate
308 
309 .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()`
310 E*/
311 typedef enum {
312   KSP_CHEBYSHEV_FIRST,
313   KSP_CHEBYSHEV_FOURTH,
314   KSP_CHEBYSHEV_OPT_FOURTH
315 } KSPChebyshevKind;
316 
317 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
318 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
319 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
320 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
321 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
322 PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind);
323 PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *);
324 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
325 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
326 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
327 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
328 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
329 
330 /*E
331 
332   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
333 
334   Values:
335 + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
336 - `KSP_FCD_TRUNC_TYPE_NOTAY`    - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
337 
338    Level: intermediate
339 
340 .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
341 E*/
342 typedef enum {
343   KSP_FCD_TRUNC_TYPE_STANDARD,
344   KSP_FCD_TRUNC_TYPE_NOTAY
345 } KSPFCDTruncationType;
346 PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
347 
348 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
349 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
350 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
351 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
352 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
353 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
354 
355 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
356 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
357 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
358 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
359 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
360 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
361 
362 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
363 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
364 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
365 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
366 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
367 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
368 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
369 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
370 
371 /*S
372   KSPFlexibleModifyPCFn - A prototype of a function used to modify the preconditioner during the use of flexible `KSP` methods, such as `KSPFGMRES`
373 
374   Calling Sequence:
375 + ksp       - the `KSP` context being used.
376 . total_its - the total number of iterations that have occurred.
377 . local_its - the number of iterations since last restart if applicable
378 . res_norm  - the current residual norm
379 - ctx       - optional context variable set with `KSPFlexibleSetModifyPC()`, `KSPPIPEGCRSetModifyPC()`, `KSPGCRSetModifyPC()`, `KSPFGMRESSetModifyPC()`
380 
381   Level: beginner
382 
383 .seealso: [](ch_ksp), `KSP`, `KSPFlexibleSetModifyPC()`, `KSPPIPEGCRSetModifyPC()`, `KSPGCRSetModifyPC()`, `KSPFGMRESSetModifyPC()`
384 S*/
385 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPFlexibleModifyPCFn(KSP ksp, PetscInt total_its, PetscInt local_its, PetscReal res_norm, void *ctx);
386 
387 PETSC_EXTERN PetscErrorCode KSPFlexibleSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *);
388 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *);
389 
390 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
391 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
392 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
393 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
394 
395 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
396 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
397 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
398 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
399 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
400 
401 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
402 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
403 
404 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
405 
406 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
407 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
408 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *);
409 
410 PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
411 PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
412 PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);
413 
414 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
415 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
416 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
417 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
418 
419 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
420 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
421 #if PetscDefined(HAVE_HPDDM)
422 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
423 {
424   return KSPHPDDMSetDeflationMat(ksp, U);
425 }
426 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
427 {
428   return KSPHPDDMGetDeflationMat(ksp, U);
429 }
430 #endif
431 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
432 {
433   return KSPMatSolve(ksp, B, X);
434 }
435 /*E
436     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
437 
438     Values:
439 +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
440 .   `KSP_HPDDM_TYPE_BGMRES`          - block GMRES
441 .   `KSP_HPDDM_TYPE_CG`              - Conjugate Gradient
442 .   `KSP_HPDDM_TYPE_BCG`             - block CG
443 .   `KSP_HPDDM_TYPE_GCRODR`          - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
444 .   `KSP_HPDDM_TYPE_BGCRODR`         - block GCRODR
445 .   `KSP_HPDDM_TYPE_BFBCG`           - breakdown-free BCG
446 -   `KSP_HPDDM_TYPE_PREONLY`         - apply the preconditioner only
447 
448     Level: intermediate
449 
450 .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
451 E*/
452 typedef enum {
453   KSP_HPDDM_TYPE_GMRES   = 0,
454   KSP_HPDDM_TYPE_BGMRES  = 1,
455   KSP_HPDDM_TYPE_CG      = 2,
456   KSP_HPDDM_TYPE_BCG     = 3,
457   KSP_HPDDM_TYPE_GCRODR  = 4,
458   KSP_HPDDM_TYPE_BGCRODR = 5,
459   KSP_HPDDM_TYPE_BFBCG   = 6,
460   KSP_HPDDM_TYPE_PREONLY = 7
461 } KSPHPDDMType;
462 PETSC_EXTERN const char *const KSPHPDDMTypes[];
463 
464 /*E
465     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
466 
467     Values:
468 +   `KSP_HPDDM_PRECISION_HALF`      - default when PETSc is configured `--with-precision=__fp16`
469 .   `KSP_HPDDM_PRECISION_SINGLE`    - default when PETSc is configured `--with-precision=single`
470 .   `KSP_HPDDM_PRECISION_DOUBLE`    - default when PETSc is configured `--with-precision=double`
471 -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
472 
473     Level: intermediate
474 
475 .seealso: [](ch_ksp), `KSP`, `KSPHPDDM`
476 E*/
477 typedef enum {
478   KSP_HPDDM_PRECISION_HALF      = 0,
479   KSP_HPDDM_PRECISION_SINGLE    = 1,
480   KSP_HPDDM_PRECISION_DOUBLE    = 2,
481   KSP_HPDDM_PRECISION_QUADRUPLE = 3
482 } KSPHPDDMPrecision;
483 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
484 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
485 
486 /*E
487    KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers
488 
489    Values:
490 +  `KSP_GMRES_CGS_REFINE_NEVER`    - one step of classical Gram-Schmidt
491 .  `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
492 -  `KSP_GMRES_CGS_REFINE_ALWAYS`   - always perform two steps
493 
494    Level: advanced
495 
496 .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
497           `KSPGMRESGetOrthogonalization()`,
498           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
499 E*/
500 typedef enum {
501   KSP_GMRES_CGS_REFINE_NEVER,
502   KSP_GMRES_CGS_REFINE_IFNEEDED,
503   KSP_GMRES_CGS_REFINE_ALWAYS
504 } KSPGMRESCGSRefinementType;
505 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
506 
507 /*MC
508    KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
509 
510    Level: advanced
511 
512    Note:
513    Possibly unstable, but the fastest to compute
514 
515 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
516           `KSP`, `KSPGMRESGetOrthogonalization()`,
517           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
518           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
519 M*/
520 
521 /*MC
522     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
523           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
524           poor orthogonality.
525 
526    Level: advanced
527 
528    Note:
529    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
530    estimate the orthogonality but is more stable.
531 
532 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
533           `KSP`, `KSPGMRESGetOrthogonalization()`,
534           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
535           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
536 M*/
537 
538 /*MC
539    KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process.
540 
541    Level: advanced
542 
543    Notes:
544    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
545    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
546 
547    You should only use this if you absolutely know that the iterative refinement is needed.
548 
549 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
550           `KSP`, `KSPGMRESGetOrthogonalization()`,
551           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
552           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
553 M*/
554 
555 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
556 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
557 
558 PETSC_EXTERN KSPFlexibleModifyPCFn KSPFGMRESModifyPCNoChange;
559 PETSC_EXTERN KSPFlexibleModifyPCFn KSPFGMRESModifyPCKSP;
560 PETSC_EXTERN PetscErrorCode        KSPFGMRESSetModifyPC(KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *);
561 
562 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
563 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
564 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
565 
566 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
567 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
568 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
569 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
570 
571 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
572 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
573 
574 PETSC_EXTERN PetscErrorCode       KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
575 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidual;
576 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualView;
577 PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorResidualDraw()", ) static inline PetscErrorCode KSPMonitorResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
578 {
579   return KSPMonitorResidualView(ksp, n, rnorm, vf);
580 }
581 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDrawLG;
582 PETSC_EXTERN PetscErrorCode       KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
583 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualShort;
584 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualRange;
585 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidual;
586 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualView;
587 PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorTrueResidualDraw()", ) static inline PetscErrorCode KSPMonitorTrueResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
588 {
589   return KSPMonitorTrueResidualView(ksp, n, rnorm, vf);
590 }
591 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDrawLG;
592 PETSC_EXTERN PetscErrorCode       KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
593 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualMax;
594 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorError;
595 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDraw;
596 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDrawLG;
597 PETSC_EXTERN PetscErrorCode       KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
598 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolution;
599 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDraw;
600 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDrawLG;
601 PETSC_EXTERN PetscErrorCode       KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
602 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSingularValue;
603 PETSC_EXTERN PetscErrorCode       KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
604 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
605 {
606   return KSPMonitorResidual(ksp, n, rnorm, vf);
607 }
608 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
609 {
610   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
611 }
612 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
613 {
614   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
615 }
616 
617 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
618 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
619 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
620 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
621 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
622 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
623 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
624 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
625 
626 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
627 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
628 
629 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
630 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
631 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
632 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
633 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
634 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
635 
636 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
637 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
638 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
639 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
640 
641 /*S
642   KSPConvergedReasonViewFn - A prototype of a function used with `KSPConvergedReasonViewSet()`
643 
644   Calling Sequence:
645 + ksp - the `KSP` object whose `KSPConvergedReason` is to be viewed
646 - ctx - context used by the function, set with `KSPConvergedReasonViewSet()`
647 
648   Level: beginner
649 
650 .seealso: [](ch_ksp), `KSP`, `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`, `KSPConvergedReasonViewFromOptions()`, `KSPView()`
651 S*/
652 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergedReasonViewFn(KSP ksp, void *ctx);
653 
654 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
655 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
656 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
657 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
658 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, KSPConvergedReasonViewFn *, void *, PetscCtxDestroyFn *);
659 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
660 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
661 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
662 
663 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
664 {
665   return KSPConvergedReasonView(ksp, v);
666 }
667 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
668 {
669   return KSPConvergedReasonViewFromOptions(ksp);
670 }
671 
672 #define KSP_FILE_CLASSID 1211223
673 
674 PETSC_EXTERN PetscErrorCode       KSPLSQRSetExactMatNorm(KSP, PetscBool);
675 PETSC_EXTERN PetscErrorCode       KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
676 PETSC_EXTERN PetscErrorCode       KSPLSQRGetStandardErrorVec(KSP, Vec *);
677 PETSC_EXTERN PetscErrorCode       KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
678 PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidual;
679 PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidualDrawLG;
680 PETSC_EXTERN PetscErrorCode       KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
681 
682 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
683 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
684 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
685 PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);
686 
687 /*E
688    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
689        test routines.
690 
691    Values:
692 +  `KSP_NORM_DEFAULT`          - use the default for the current `KSPType`
693 .  `KSP_NORM_NONE`             - use no norm calculation
694 .  `KSP_NORM_PRECONDITIONED`   - use the preconditioned residual norm
695 .  `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
696 -  `KSP_NORM_NATURAL`          - use the natural norm (the norm induced by the linear operator)
697 
698    Level: advanced
699 
700    Note:
701    Each solver only supports a subset of these and some may support different ones
702    depending on whether left or right preconditioning is used, see `KSPSetPCSide()`
703 
704 .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
705           `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
706 E*/
707 typedef enum {
708   KSP_NORM_DEFAULT          = -1,
709   KSP_NORM_NONE             = 0,
710   KSP_NORM_PRECONDITIONED   = 1,
711   KSP_NORM_UNPRECONDITIONED = 2,
712   KSP_NORM_NATURAL          = 3
713 } KSPNormType;
714 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
715 PETSC_EXTERN const char *const *const KSPNormTypes;
716 
717 /*MC
718    KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
719    possibly save some computation but means the convergence test cannot
720    be based on a norm of a residual etc.
721 
722    Level: advanced
723 
724    Note:
725    Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
726 
727 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
728 M*/
729 
730 /*MC
731    KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
732    convergence test routine.
733 
734    Level: advanced
735 
736 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
737 M*/
738 
739 /*MC
740    KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
741    convergence test routine.
742 
743    Level: advanced
744 
745 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
746 M*/
747 
748 /*MC
749    KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
750    convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
751 
752    Level: advanced
753 
754 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
755 M*/
756 
757 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
758 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
759 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt);
760 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
761 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
762 
763 #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
764 #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
765 #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", )
766 #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", )
767 /*E
768    KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
769 
770    Values:
771 +  `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR`
772 .  `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR`
773 .  `KSP_CONVERGED_RTOL`                  - requested decrease in the residual
774 .  `KSP_CONVERGED_ATOL`                  - requested absolute value in the residual
775 .  `KSP_CONVERGED_ITS`                   - requested number of iterations
776 .  `KSP_CONVERGED_NEG_CURVE`             - see note below
777 .  `KSP_CONVERGED_STEP_LENGTH`           - see note below
778 .  `KSP_CONVERGED_HAPPY_BREAKDOWN`       - happy breakdown (meaning early convergence of the `KSPType` occurred).
779 .  `KSP_CONVERGED_USER`                  - the user has indicated convergence for an arbitrary reason
780 .  `KSP_DIVERGED_NULL`                   - breakdown when solving the Hessenberg system within `KSPGMRES`
781 .  `KSP_DIVERGED_ITS`                    - requested number of iterations
782 .  `KSP_DIVERGED_DTOL`                   - large increase in the residual norm indicating the solution is diverging
783 .  `KSP_DIVERGED_BREAKDOWN`              - breakdown in the Krylov method
784 .  `KSP_DIVERGED_BREAKDOWN_BICG`         - breakdown in the `KSPBCGS` Krylov method
785 .  `KSP_DIVERGED_NONSYMMETRIC`           - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
786 .  `KSP_DIVERGED_INDEFINITE_PC`          - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
787 .  `KSP_DIVERGED_NANORINF`               - a not a number of infinity was detected in a vector during the computation
788 .  `KSP_DIVERGED_INDEFINITE_MAT`         - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
789 .  `KSP_DIVERGED_PC_FAILED`              - the action of the preconditioner failed for some reason
790 -  `KSP_DIVERGED_USER`                   - the user has indicated divergence for an arbitrary reason
791 
792    Level: beginner
793 
794    Note:
795    The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
796    the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
797 
798    Developer Note:
799    The string versions of these are `KSPConvergedReasons`; if you change
800    any of the values here also change them that array of names.
801 
802 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
803 E*/
804 typedef enum { /* converged */
805   KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    = 1,
806   KSP_CONVERGED_RTOL_NORMAL_EQUATIONS     = 1,
807   KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    = 9,
808   KSP_CONVERGED_ATOL_NORMAL_EQUATIONS     = 9,
809   KSP_CONVERGED_RTOL                      = 2,
810   KSP_CONVERGED_ATOL                      = 3,
811   KSP_CONVERGED_ITS                       = 4,
812   KSP_CONVERGED_NEG_CURVE                 = 5,
813   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
814   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
815   KSP_CONVERGED_STEP_LENGTH               = 6,
816   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
817   KSP_CONVERGED_USER                      = 8,
818   /* diverged */
819   KSP_DIVERGED_NULL                      = -2,
820   KSP_DIVERGED_ITS                       = -3,
821   KSP_DIVERGED_DTOL                      = -4,
822   KSP_DIVERGED_BREAKDOWN                 = -5,
823   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
824   KSP_DIVERGED_NONSYMMETRIC              = -7,
825   KSP_DIVERGED_INDEFINITE_PC             = -8,
826   KSP_DIVERGED_NANORINF                  = -9,
827   KSP_DIVERGED_INDEFINITE_MAT            = -10,
828   KSP_DIVERGED_PC_FAILED                 = -11,
829   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
830   KSP_DIVERGED_USER                      = -12,
831 
832   KSP_CONVERGED_ITERATING = 0
833 } KSPConvergedReason;
834 PETSC_EXTERN const char *const *KSPConvergedReasons;
835 
836 /*MC
837    KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called
838 
839    Level: beginner
840 
841    Notes:
842    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
843    for left preconditioning it is the 2-norm of the preconditioned residual, and the
844    2-norm of the residual for right preconditioning
845 
846    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
847 
848 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
849 M*/
850 
851 /*MC
852    KSP_CONVERGED_ATOL - $||r|| \le atol$
853 
854    Level: beginner
855 
856    Notes:
857    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
858    for left preconditioning it is the 2-norm of the preconditioned residual, and the
859    2-norm of the residual for right preconditioning
860 
861    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
862 
863 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
864 M*/
865 
866 /*MC
867    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$
868 
869    Level: beginner
870 
871    Note:
872    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
873    for left preconditioning it is the 2-norm of the preconditioned residual, and the
874    2-norm of the residual for right preconditioning
875 
876 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
877 M*/
878 
879 /*MC
880    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
881    reached
882 
883    Level: beginner
884 
885 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
886 M*/
887 
888 /*MC
889    KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
890    the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
891    test routine is set in `KSP`.
892 
893    Level: beginner
894 
895 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
896 M*/
897 
898 /*MC
899    KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
900    method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
901    preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
902    are collinear.
903 
904    Level: beginner
905 
906 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
907 M*/
908 
909 /*MC
910    KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
911    method could not continue to enlarge the Krylov space.
912 
913    Level: beginner
914 
915 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
916 M*/
917 
918 /*MC
919    KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
920    symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
921 
922    Level: beginner
923 
924 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
925 M*/
926 
927 /*MC
928    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
929    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
930    be symmetric positive definite (SPD).
931 
932    Level: beginner
933 
934    Note:
935    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
936    the `PCICC` preconditioner to generate a positive definite preconditioner
937 
938 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
939 M*/
940 
941 /*MC
942    KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
943    zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
944    such as `PCFIELDSPLIT`.
945 
946    Level: beginner
947 
948    Note:
949    Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
950 
951 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
952 M*/
953 
954 /*MC
955    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
956    while `KSPSolve()` is still running.
957 
958    Level: beginner
959 
960 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
961 M*/
962 
963 /*S
964   KSPConvergenceTestFn - A prototype of a function used with `KSPSetConvergenceTest()`
965 
966   Calling Sequence:
967 + ksp    - iterative solver obtained from `KSPCreate()`
968 . it     - iteration number
969 . rnorm  - (estimated) 2-norm of (preconditioned) residual
970 . reason - the reason why it has converged or diverged
971 - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`
972 
973   Level: beginner
974 
975 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
976 S*/
977 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergenceTestFn(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx);
978 
979 PETSC_EXTERN PetscErrorCode       KSPSetConvergenceTest(KSP, KSPConvergenceTestFn *, void *, PetscCtxDestroyFn *);
980 PETSC_EXTERN PetscErrorCode       KSPGetConvergenceTest(KSP, KSPConvergenceTestFn **, void **, PetscCtxDestroyFn **);
981 PETSC_EXTERN PetscErrorCode       KSPGetAndClearConvergenceTest(KSP, KSPConvergenceTestFn **, void **, PetscCtxDestroyFn **);
982 PETSC_EXTERN PetscErrorCode       KSPGetConvergenceContext(KSP, void *);
983 PETSC_EXTERN KSPConvergenceTestFn KSPConvergedDefault;
984 PETSC_EXTERN KSPConvergenceTestFn KSPLSQRConvergedDefault;
985 PETSC_EXTERN PetscCtxDestroyFn    KSPConvergedDefaultDestroy;
986 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultCreate(void **);
987 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUIRNorm(KSP);
988 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUMIRNorm(KSP);
989 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
990 PETSC_EXTERN PetscErrorCode       KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
991 PETSC_EXTERN PetscErrorCode       KSPGetConvergedReason(KSP, KSPConvergedReason *);
992 PETSC_EXTERN PetscErrorCode       KSPGetConvergedReasonString(KSP, const char *[]);
993 PETSC_EXTERN PetscErrorCode       KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
994 PETSC_EXTERN PetscErrorCode       KSPSetConvergedNegativeCurvature(KSP, PetscBool);
995 PETSC_EXTERN PetscErrorCode       KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
996 
997 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
998 { /* never called */
999 }
1000 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
1001 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
1002 { /* never called */
1003 }
1004 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
1005 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
1006 { /* never called */
1007 }
1008 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
1009 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
1010 { /* never called */
1011 }
1012 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
1013 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
1014 { /* never called */
1015 }
1016 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
1017 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
1018 { /* never called */
1019 }
1020 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
1021 
1022 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
1023 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
1024 {
1025   return KSPComputeOperator(A, PETSC_NULLPTR, B);
1026 }
1027 
1028 /*E
1029    KSPCGType - Determines what type of `KSPCG` to use
1030 
1031    Values:
1032  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
1033  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
1034 
1035    Level: beginner
1036 
1037 .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
1038 E*/
1039 typedef enum {
1040   KSP_CG_SYMMETRIC = 0,
1041   KSP_CG_HERMITIAN = 1
1042 } KSPCGType;
1043 PETSC_EXTERN const char *const KSPCGTypes[];
1044 
1045 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
1046 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
1047 
1048 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
1049 PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
1050 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
1051 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
1052 
1053 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
1054 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
1055 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
1056 {
1057   return KSPGLTRGetMinEig(ksp, x);
1058 }
1059 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
1060 {
1061   return KSPGLTRGetLambda(ksp, x);
1062 }
1063 
1064 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
1065 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
1066 
1067 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
1068 PETSC_EXTERN PetscErrorCode PCSetPostSolve(PC, PetscErrorCode (*)(PC, KSP));
1069 PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
1070 PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
1071 
1072 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
1073 
1074 /*S
1075   PCShellPSolveFn - A function prototype for functions provided to `PCShellSetPreSolve()` and `PCShellSetPostSolve()`
1076 
1077   Calling Sequence:
1078 + pc  - the preconditioner `PC` context
1079 . ksp - the `KSP` context
1080 . xin  - input vector
1081 - xout - output vector
1082 
1083   Level: intermediate
1084 
1085 .seealso: [](ch_snes), `KSPPSolveFn`, `KSP`, `PCShellSetPreSolve()`, `PCShellSetPostSolve()`
1086 S*/
1087 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCShellPSolveFn(PC pc, KSP ksp, Vec xim, Vec xout);
1088 
1089 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PCShellPSolveFn *);
1090 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PCShellPSolveFn *);
1091 
1092 /*S
1093    KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
1094 
1095    Level: intermediate
1096 
1097    Note:
1098    These methods generate initial guesses based on a series of previous, related, linear solves. For example,
1099    in implicit time-stepping with `TS`.
1100 
1101 .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
1102 S*/
1103 typedef struct _p_KSPGuess *KSPGuess;
1104 
1105 /*J
1106    KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
1107 
1108    Values:
1109  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
1110  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition (POD)
1111 
1112    Level: intermediate
1113 
1114 .seealso: [](ch_ksp), `KSP`, `KSPGuess`
1115 J*/
1116 typedef const char *KSPGuessType;
1117 #define KSPGUESSFISCHER "fischer"
1118 #define KSPGUESSPOD     "pod"
1119 
1120 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
1121 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
1122 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
1123 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
1124 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
1125 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
1126 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
1127 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
1128 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
1129 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
1130 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
1131 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
1132 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
1133 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
1134 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
1135 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
1136 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
1137 
1138 /*E
1139     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines
1140 
1141     Level: intermediate
1142 
1143 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
1144           `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
1145 E*/
1146 typedef enum {
1147   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
1148   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
1149   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
1150   MAT_SCHUR_COMPLEMENT_AINV_FULL
1151 } MatSchurComplementAinvType;
1152 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
1153 
1154 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
1155 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
1156 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
1157 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1158 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1159 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
1160 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
1161 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
1162 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
1163 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
1164 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
1165 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
1166 
1167 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1168 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1169 PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1170 PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1171 PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
1172 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1173 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1174 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1175 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1176 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1177 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1178 
1179 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1180 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1181 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1182 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1183 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1184 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1185 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1186 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1187 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1188 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1189 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1190 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1191 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1192 PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *);
1193 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1194 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1195 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1196 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1197 PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1198 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1199 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1200 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1201 
1202 /*E
1203   MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix
1204 
1205   Values:
1206 + `MAT_LMVM_MULT_RECURSIVE`     - Use recursive formulas for products and solves
1207 . `MAT_LMVM_MULT_DENSE`         - Use dense formulas for products and solves when possible
1208 - `MAT_LMVM_MULT_COMPACT_DENSE` - The same as `MATLMVM_MULT_DENSE`, but go further and ensure products and solves are computed in compact low-rank update form
1209 
1210   Level: advanced
1211 
1212   Options Database Keys:
1213 . -mat_lmvm_mult_algorithm  - the algorithm to use for multiplication (recursive, dense, compact_dense)
1214 
1215 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()`
1216 E*/
1217 typedef enum {
1218   MAT_LMVM_MULT_RECURSIVE,
1219   MAT_LMVM_MULT_DENSE,
1220   MAT_LMVM_MULT_COMPACT_DENSE,
1221 } MatLMVMMultAlgorithm;
1222 
1223 PETSC_EXTERN const char *const MatLMVMMultAlgorithms[];
1224 
1225 PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm);
1226 PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *);
1227 
1228 /*E
1229   MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix.
1230 
1231   Values:
1232 + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - no rescaling
1233 . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar rescaling
1234 . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling
1235 . `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE`
1236 - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE`   - let PETSc decide rescaling
1237 
1238   Level: intermediate
1239 
1240 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1241 E*/
1242 typedef enum {
1243   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1244   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1245   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1246   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3,
1247   MAT_LMVM_SYMBROYDEN_SCALE_DECIDE   = 4
1248 } MatLMVMSymBroydenScaleType;
1249 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
1250 
1251 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1252 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *);
1253 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal);
1254 PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *);
1255 PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal);
1256 
1257 /*E
1258   MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.
1259 
1260   Values:
1261 + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1262 - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement
1263 
1264   Level: intermediate
1265 
1266 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1267 E*/
1268 typedef enum {
1269   MAT_LMVM_DENSE_REORDER,
1270   MAT_LMVM_DENSE_INPLACE
1271 } MatLMVMDenseType;
1272 PETSC_EXTERN const char *const MatLMVMDenseTypes[];
1273 
1274 PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);
1275 
1276 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1277 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1278 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1279 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1280 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1281 
1282 /*S
1283   KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`
1284 
1285   Calling Sequence:
1286 + ksp  - `ksp` context
1287 . b    - output vector
1288 - ctx - [optional] user-defined function context
1289 
1290   Level: beginner
1291 
1292 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
1293 S*/
1294 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeRHSFn(KSP ksp, Vec b, void *ctx);
1295 
1296 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);
1297 
1298 /*S
1299   KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`
1300 
1301   Calling Sequence:
1302 + ksp - `KSP` context
1303 . A   - the operator that defines the linear system
1304 . P   - an operator from which to build the preconditioner (often the same as `A`)
1305 - ctx - [optional] user-defined function context
1306 
1307   Level: beginner
1308 
1309 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
1310 S*/
1311 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeOperatorsFn(KSP ksp, Mat A, Mat P, void *ctx);
1312 
1313 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);
1314 
1315 /*S
1316   KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`
1317 
1318   Calling Sequence:
1319 + ksp  - `ksp` context
1320 . x    - output vector
1321 - ctx - [optional] user-defined function context
1322 
1323   Level: beginner
1324 
1325 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
1326 S*/
1327 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeInitialGuessFn(KSP ksp, Vec x, void *ctx);
1328 
1329 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1330 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1331 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1332 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1333 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1334 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1335 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);
1336 
1337 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1338 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode);
1339 PETSC_EXTERN PetscErrorCode DMSwarmProjectGradientFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode);
1340 
1341 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1342 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1343 
1344 PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1345 PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);
1346 
1347 PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);
1348 
1349 #include <petscdstypes.h>
1350 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, PetscPointFn **, InsertMode, Vec);
1351