xref: /petsc/include/petscksp.h (revision 4d4d2bdc375ba73538fdff3fed6ff2932e6875ae)
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 KSPMonitorResidualDraw;
577 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDrawLG;
578 PETSC_EXTERN PetscErrorCode       KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
579 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualShort;
580 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualRange;
581 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidual;
582 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDraw;
583 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDrawLG;
584 PETSC_EXTERN PetscErrorCode       KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
585 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualMax;
586 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorError;
587 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDraw;
588 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDrawLG;
589 PETSC_EXTERN PetscErrorCode       KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
590 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolution;
591 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDraw;
592 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDrawLG;
593 PETSC_EXTERN PetscErrorCode       KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
594 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSingularValue;
595 PETSC_EXTERN PetscErrorCode       KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
596 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
597 {
598   return KSPMonitorResidual(ksp, n, rnorm, vf);
599 }
600 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
601 {
602   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
603 }
604 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
605 {
606   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
607 }
608 
609 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
610 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
611 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
612 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
613 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
614 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
615 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
616 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
617 
618 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
619 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
620 
621 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
622 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
623 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
624 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
625 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
626 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
627 
628 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
629 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
630 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
631 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
632 
633 /*S
634   KSPConvergedReasonViewFn - A prototype of a function used with `KSPConvergedReasonViewSet()`
635 
636   Calling Sequence:
637 + ksp - the `KSP` object whose `KSPConvergedReason` is to be viewed
638 - ctx - context used by the function, set with `KSPConvergedReasonViewSet()`
639 
640   Level: beginner
641 
642 .seealso: [](ch_ksp), `KSP`, `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`, `KSPConvergedReasonViewFromOptions()`, `KSPView()`
643 S*/
644 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergedReasonViewFn(KSP ksp, void *ctx);
645 
646 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
647 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
648 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
649 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
650 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, KSPConvergedReasonViewFn *, void *, PetscCtxDestroyFn *);
651 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
652 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
653 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
654 
655 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
656 {
657   return KSPConvergedReasonView(ksp, v);
658 }
659 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
660 {
661   return KSPConvergedReasonViewFromOptions(ksp);
662 }
663 
664 #define KSP_FILE_CLASSID 1211223
665 
666 PETSC_EXTERN PetscErrorCode       KSPLSQRSetExactMatNorm(KSP, PetscBool);
667 PETSC_EXTERN PetscErrorCode       KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
668 PETSC_EXTERN PetscErrorCode       KSPLSQRGetStandardErrorVec(KSP, Vec *);
669 PETSC_EXTERN PetscErrorCode       KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
670 PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidual;
671 PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidualDrawLG;
672 PETSC_EXTERN PetscErrorCode       KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
673 
674 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
675 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
676 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
677 PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);
678 
679 /*E
680    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
681        test routines.
682 
683    Values:
684 +  `KSP_NORM_DEFAULT`          - use the default for the current `KSPType`
685 .  `KSP_NORM_NONE`             - use no norm calculation
686 .  `KSP_NORM_PRECONDITIONED`   - use the preconditioned residual norm
687 .  `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
688 -  `KSP_NORM_NATURAL`          - use the natural norm (the norm induced by the linear operator)
689 
690    Level: advanced
691 
692    Note:
693    Each solver only supports a subset of these and some may support different ones
694    depending on whether left or right preconditioning is used, see `KSPSetPCSide()`
695 
696 .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
697           `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
698 E*/
699 typedef enum {
700   KSP_NORM_DEFAULT          = -1,
701   KSP_NORM_NONE             = 0,
702   KSP_NORM_PRECONDITIONED   = 1,
703   KSP_NORM_UNPRECONDITIONED = 2,
704   KSP_NORM_NATURAL          = 3
705 } KSPNormType;
706 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
707 PETSC_EXTERN const char *const *const KSPNormTypes;
708 
709 /*MC
710    KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
711    possibly save some computation but means the convergence test cannot
712    be based on a norm of a residual etc.
713 
714    Level: advanced
715 
716    Note:
717    Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
718 
719 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
720 M*/
721 
722 /*MC
723    KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
724    convergence test routine.
725 
726    Level: advanced
727 
728 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
729 M*/
730 
731 /*MC
732    KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
733    convergence test routine.
734 
735    Level: advanced
736 
737 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
738 M*/
739 
740 /*MC
741    KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
742    convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
743 
744    Level: advanced
745 
746 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
747 M*/
748 
749 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
750 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
751 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt);
752 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
753 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
754 
755 #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
756 #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
757 #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", )
758 #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", )
759 /*E
760    KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
761 
762    Values:
763 +  `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR`
764 .  `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR`
765 .  `KSP_CONVERGED_RTOL`                  - requested decrease in the residual
766 .  `KSP_CONVERGED_ATOL`                  - requested absolute value in the residual
767 .  `KSP_CONVERGED_ITS`                   - requested number of iterations
768 .  `KSP_CONVERGED_NEG_CURVE`             - see note below
769 .  `KSP_CONVERGED_STEP_LENGTH`           - see note below
770 .  `KSP_CONVERGED_HAPPY_BREAKDOWN`       - happy breakdown (meaning early convergence of the `KSPType` occurred).
771 .  `KSP_DIVERGED_NULL`                   - breakdown when solving the Hessenberg system within `KSPGMRES`
772 .  `KSP_DIVERGED_ITS`                    - requested number of iterations
773 .  `KSP_DIVERGED_DTOL`                   - large increase in the residual norm indicating the solution is diverging
774 .  `KSP_DIVERGED_BREAKDOWN`              - breakdown in the Krylov method
775 .  `KSP_DIVERGED_BREAKDOWN_BICG`         - breakdown in the `KSPBCGS` Krylov method
776 .  `KSP_DIVERGED_NONSYMMETRIC`           - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
777 .  `KSP_DIVERGED_INDEFINITE_PC`          - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
778 .  `KSP_DIVERGED_NANORINF`               - a not a number of infinity was detected in a vector during the computation
779 .  `KSP_DIVERGED_INDEFINITE_MAT`         - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
780 -  `KSP_DIVERGED_PC_FAILED`              - the action of the preconditioner failed for some reason
781 
782    Level: beginner
783 
784    Note:
785    The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
786    the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
787 
788    Developer Note:
789    The string versions of these are `KSPConvergedReasons`; if you change
790    any of the values here also change them that array of names.
791 
792 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
793 E*/
794 typedef enum { /* converged */
795   KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    = 1,
796   KSP_CONVERGED_RTOL_NORMAL_EQUATIONS     = 1,
797   KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    = 9,
798   KSP_CONVERGED_ATOL_NORMAL_EQUATIONS     = 9,
799   KSP_CONVERGED_RTOL                      = 2,
800   KSP_CONVERGED_ATOL                      = 3,
801   KSP_CONVERGED_ITS                       = 4,
802   KSP_CONVERGED_NEG_CURVE                 = 5,
803   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
804   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
805   KSP_CONVERGED_STEP_LENGTH               = 6,
806   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
807   /* diverged */
808   KSP_DIVERGED_NULL                      = -2,
809   KSP_DIVERGED_ITS                       = -3,
810   KSP_DIVERGED_DTOL                      = -4,
811   KSP_DIVERGED_BREAKDOWN                 = -5,
812   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
813   KSP_DIVERGED_NONSYMMETRIC              = -7,
814   KSP_DIVERGED_INDEFINITE_PC             = -8,
815   KSP_DIVERGED_NANORINF                  = -9,
816   KSP_DIVERGED_INDEFINITE_MAT            = -10,
817   KSP_DIVERGED_PC_FAILED                 = -11,
818   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
819 
820   KSP_CONVERGED_ITERATING = 0
821 } KSPConvergedReason;
822 PETSC_EXTERN const char *const *KSPConvergedReasons;
823 
824 /*MC
825    KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called
826 
827    Level: beginner
828 
829    Notes:
830    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
831    for left preconditioning it is the 2-norm of the preconditioned residual, and the
832    2-norm of the residual for right preconditioning
833 
834    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
835 
836 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
837 M*/
838 
839 /*MC
840    KSP_CONVERGED_ATOL - $||r|| \le atol$
841 
842    Level: beginner
843 
844    Notes:
845    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
846    for left preconditioning it is the 2-norm of the preconditioned residual, and the
847    2-norm of the residual for right preconditioning
848 
849    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
850 
851 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
852 M*/
853 
854 /*MC
855    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$
856 
857    Level: beginner
858 
859    Note:
860    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
861    for left preconditioning it is the 2-norm of the preconditioned residual, and the
862    2-norm of the residual for right preconditioning
863 
864 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
865 M*/
866 
867 /*MC
868    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
869    reached
870 
871    Level: beginner
872 
873 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
874 M*/
875 
876 /*MC
877    KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
878    the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
879    test routine is set in `KSP`.
880 
881    Level: beginner
882 
883 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
884 M*/
885 
886 /*MC
887    KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
888    method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
889    preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
890    are collinear.
891 
892    Level: beginner
893 
894 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
895 M*/
896 
897 /*MC
898    KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
899    method could not continue to enlarge the Krylov space.
900 
901    Level: beginner
902 
903 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
904 M*/
905 
906 /*MC
907    KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
908    symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
909 
910    Level: beginner
911 
912 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
913 M*/
914 
915 /*MC
916    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
917    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
918    be symmetric positive definite (SPD).
919 
920    Level: beginner
921 
922    Note:
923    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
924    the `PCICC` preconditioner to generate a positive definite preconditioner
925 
926 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
927 M*/
928 
929 /*MC
930    KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
931    zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
932    such as `PCFIELDSPLIT`.
933 
934    Level: beginner
935 
936    Note:
937    Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
938 
939 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
940 M*/
941 
942 /*MC
943    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
944    while `KSPSolve()` is still running.
945 
946    Level: beginner
947 
948 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
949 M*/
950 
951 /*S
952   KSPConvergenceTestFn - A prototype of a function used with `KSPSetConvergenceTest()`
953 
954   Calling Sequence:
955 + ksp    - iterative solver obtained from `KSPCreate()`
956 . it     - iteration number
957 . rnorm  - (estimated) 2-norm of (preconditioned) residual
958 . reason - the reason why it has converged or diverged
959 - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`
960 
961   Level: beginner
962 
963 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
964 S*/
965 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergenceTestFn(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx);
966 
967 PETSC_EXTERN PetscErrorCode       KSPSetConvergenceTest(KSP, KSPConvergenceTestFn *, void *, PetscCtxDestroyFn *);
968 PETSC_EXTERN PetscErrorCode       KSPGetConvergenceTest(KSP, KSPConvergenceTestFn **, void **, PetscCtxDestroyFn **);
969 PETSC_EXTERN PetscErrorCode       KSPGetAndClearConvergenceTest(KSP, KSPConvergenceTestFn **, void **, PetscCtxDestroyFn **);
970 PETSC_EXTERN PetscErrorCode       KSPGetConvergenceContext(KSP, void *);
971 PETSC_EXTERN KSPConvergenceTestFn KSPConvergedDefault;
972 PETSC_EXTERN KSPConvergenceTestFn KSPLSQRConvergedDefault;
973 PETSC_EXTERN PetscCtxDestroyFn    KSPConvergedDefaultDestroy;
974 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultCreate(void **);
975 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUIRNorm(KSP);
976 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUMIRNorm(KSP);
977 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
978 PETSC_EXTERN PetscErrorCode       KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
979 PETSC_EXTERN PetscErrorCode       KSPGetConvergedReason(KSP, KSPConvergedReason *);
980 PETSC_EXTERN PetscErrorCode       KSPGetConvergedReasonString(KSP, const char *[]);
981 PETSC_EXTERN PetscErrorCode       KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
982 PETSC_EXTERN PetscErrorCode       KSPSetConvergedNegativeCurvature(KSP, PetscBool);
983 PETSC_EXTERN PetscErrorCode       KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
984 
985 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
986 { /* never called */
987 }
988 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
989 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
990 { /* never called */
991 }
992 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
993 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
994 { /* never called */
995 }
996 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
997 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
998 { /* never called */
999 }
1000 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
1001 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
1002 { /* never called */
1003 }
1004 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
1005 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
1006 { /* never called */
1007 }
1008 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
1009 
1010 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
1011 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
1012 {
1013   return KSPComputeOperator(A, PETSC_NULLPTR, B);
1014 }
1015 
1016 /*E
1017    KSPCGType - Determines what type of `KSPCG` to use
1018 
1019    Values:
1020  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
1021  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
1022 
1023    Level: beginner
1024 
1025 .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
1026 E*/
1027 typedef enum {
1028   KSP_CG_SYMMETRIC = 0,
1029   KSP_CG_HERMITIAN = 1
1030 } KSPCGType;
1031 PETSC_EXTERN const char *const KSPCGTypes[];
1032 
1033 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
1034 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
1035 
1036 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
1037 PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
1038 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
1039 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
1040 
1041 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
1042 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
1043 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
1044 {
1045   return KSPGLTRGetMinEig(ksp, x);
1046 }
1047 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
1048 {
1049   return KSPGLTRGetLambda(ksp, x);
1050 }
1051 
1052 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
1053 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
1054 
1055 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
1056 PETSC_EXTERN PetscErrorCode PCSetPostSolve(PC, PetscErrorCode (*)(PC, KSP));
1057 PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
1058 PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
1059 
1060 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
1061 
1062 /*S
1063   PCShellPSolveFn - A function prototype for functions provided to `PCShellSetPreSolve()` and `PCShellSetPostSolve()`
1064 
1065   Calling Sequence:
1066 + pc  - the preconditioner `PC` context
1067 . ksp - the `KSP` context
1068 . xin  - input vector
1069 - xout - output vector
1070 
1071   Level: intermediate
1072 
1073 .seealso: [](ch_snes), `KSPPSolveFn`, `KSP`, `PCShellSetPreSolve()`, `PCShellSetPostSolve()`
1074 S*/
1075 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCShellPSolveFn(PC pc, KSP ksp, Vec xim, Vec xout);
1076 
1077 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PCShellPSolveFn *);
1078 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PCShellPSolveFn *);
1079 
1080 /*S
1081    KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
1082 
1083    Level: intermediate
1084 
1085    Note:
1086    These methods generate initial guesses based on a series of previous, related, linear solves. For example,
1087    in implicit time-stepping with `TS`.
1088 
1089 .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
1090 S*/
1091 typedef struct _p_KSPGuess *KSPGuess;
1092 
1093 /*J
1094    KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
1095 
1096    Values:
1097  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
1098  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition (POD)
1099 
1100    Level: intermediate
1101 
1102 .seealso: [](ch_ksp), `KSP`, `KSPGuess`
1103 J*/
1104 typedef const char *KSPGuessType;
1105 #define KSPGUESSFISCHER "fischer"
1106 #define KSPGUESSPOD     "pod"
1107 
1108 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
1109 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
1110 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
1111 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
1112 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
1113 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
1114 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
1115 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
1116 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
1117 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
1118 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
1119 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
1120 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
1121 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
1122 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
1123 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
1124 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
1125 
1126 /*E
1127     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines
1128 
1129     Level: intermediate
1130 
1131 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
1132           `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
1133 E*/
1134 typedef enum {
1135   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
1136   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
1137   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
1138   MAT_SCHUR_COMPLEMENT_AINV_FULL
1139 } MatSchurComplementAinvType;
1140 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
1141 
1142 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
1143 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
1144 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
1145 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1146 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1147 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
1148 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
1149 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
1150 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
1151 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
1152 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
1153 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
1154 
1155 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1156 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1157 PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1158 PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1159 PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
1160 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1161 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1162 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1163 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1164 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1165 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1166 
1167 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1168 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1169 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1170 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1171 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1172 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1173 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1174 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1175 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1176 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1177 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1178 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1179 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1180 PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *);
1181 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1182 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1183 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1184 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1185 PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1186 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1187 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1188 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1189 
1190 /*E
1191   MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix
1192 
1193   Values:
1194 + `MAT_LMVM_MULT_RECURSIVE`     - Use recursive formulas for products and solves
1195 . `MAT_LMVM_MULT_DENSE`         - Use dense formulas for products and solves when possible
1196 - `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
1197 
1198   Level: advanced
1199 
1200   Options Database Keys:
1201 . -mat_lmvm_mult_algorithm  - the algorithm to use for multiplication (recursive, dense, compact_dense)
1202 
1203 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()`
1204 E*/
1205 typedef enum {
1206   MAT_LMVM_MULT_RECURSIVE,
1207   MAT_LMVM_MULT_DENSE,
1208   MAT_LMVM_MULT_COMPACT_DENSE,
1209 } MatLMVMMultAlgorithm;
1210 
1211 PETSC_EXTERN const char *const MatLMVMMultAlgorithms[];
1212 
1213 PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm);
1214 PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *);
1215 
1216 /*E
1217   MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix.
1218 
1219   Values:
1220 + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - no rescaling
1221 . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar rescaling
1222 . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling
1223 . `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE`
1224 - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE`   - let PETSc decide rescaling
1225 
1226   Level: intermediate
1227 
1228 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1229 E*/
1230 typedef enum {
1231   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1232   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1233   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1234   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3,
1235   MAT_LMVM_SYMBROYDEN_SCALE_DECIDE   = 4
1236 } MatLMVMSymBroydenScaleType;
1237 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
1238 
1239 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1240 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *);
1241 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal);
1242 PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *);
1243 PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal);
1244 
1245 /*E
1246   MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.
1247 
1248   Values:
1249 + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1250 - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement
1251 
1252   Level: intermediate
1253 
1254 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1255 E*/
1256 typedef enum {
1257   MAT_LMVM_DENSE_REORDER,
1258   MAT_LMVM_DENSE_INPLACE
1259 } MatLMVMDenseType;
1260 PETSC_EXTERN const char *const MatLMVMDenseTypes[];
1261 
1262 PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);
1263 
1264 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1265 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1266 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1267 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1268 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1269 
1270 /*S
1271   KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`
1272 
1273   Calling Sequence:
1274 + ksp  - `ksp` context
1275 . b    - output vector
1276 - ctx - [optional] user-defined function context
1277 
1278   Level: beginner
1279 
1280 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
1281 S*/
1282 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeRHSFn(KSP ksp, Vec b, void *ctx);
1283 
1284 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);
1285 
1286 /*S
1287   KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`
1288 
1289   Calling Sequence:
1290 + ksp - `KSP` context
1291 . A   - the operator that defines the linear system
1292 . P   - an operator from which to build the preconditioner (often the same as `A`)
1293 - ctx - [optional] user-defined function context
1294 
1295   Level: beginner
1296 
1297 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
1298 S*/
1299 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeOperatorsFn(KSP ksp, Mat A, Mat P, void *ctx);
1300 
1301 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);
1302 
1303 /*S
1304   KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`
1305 
1306   Calling Sequence:
1307 + ksp  - `ksp` context
1308 . x    - output vector
1309 - ctx - [optional] user-defined function context
1310 
1311   Level: beginner
1312 
1313 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
1314 S*/
1315 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeInitialGuessFn(KSP ksp, Vec x, void *ctx);
1316 
1317 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1318 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1319 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1320 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1321 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1322 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1323 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);
1324 
1325 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1326 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode);
1327 
1328 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1329 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1330 
1331 PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1332 PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);
1333 
1334 PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);
1335 
1336 #include <petscdstypes.h>
1337 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, PetscPointFn **, InsertMode, Vec);
1338