xref: /petsc/include/petscksp.h (revision 3220ff8572602716d60bdda8b51773ebceb3c8ea)
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);
KSPSetMatSolveBlockSize(KSP ksp,PetscInt n)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 *);
KSPGetMatSolveBlockSize(KSP ksp,PetscInt * n)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, PetscCtx 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 **);
KSPGetVecs(KSP ksp,PetscInt n,Vec ** x,PetscInt m,Vec ** y)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, PetscCtx ctx);
221 
222 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, KSPPSolveFn *, PetscCtx);
223 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, KSPPSolveFn *, PetscCtx);
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, PetscCtx ctx);
244 
245 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
246 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, KSPMonitorFn *, PetscCtx, PetscCtxDestroyFn *);
247 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
248 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, PetscCtxRt);
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 gradient/directions 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   Note:
341   Function such as `KSPFCGSetMmax()`, `KSPPIPEGCRSetNMax()`, `KSPPIPEGCRSetNMax()`, and `KSPPIPEFCGSetNMax()` may be
342   used to provide `nmax` or they may be provided with the option database.
343 
344 .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`,
345           `KSPPIPEGCRSetTruncationType()`, `KSPPIPEGCRGetTruncationType()`,
346           `KSPFCGSetMmax()`, `KSPPIPEGCRSetNMax()`, `KSPPIPEGCRGetNMax()`, `KSPPIPEFCGGetNMax()`
347 E*/
348 typedef enum {
349   KSP_FCD_TRUNC_TYPE_STANDARD,
350   KSP_FCD_TRUNC_TYPE_NOTAY
351 } KSPFCDTruncationType;
352 PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
353 
354 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
355 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
356 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
357 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
358 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
359 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
360 
361 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
362 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
363 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
364 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
365 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
366 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
367 
368 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
369 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
370 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
371 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
372 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
373 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
374 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
375 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
376 
377 /*S
378   KSPFlexibleModifyPCFn - A prototype of a function used to modify the preconditioner during the use of flexible `KSP` methods, such as `KSPFGMRES`
379 
380   Calling Sequence:
381 + ksp       - the `KSP` context being used.
382 . total_its - the total number of iterations that have occurred.
383 . local_its - the number of iterations since last restart if applicable
384 . res_norm  - the current residual norm
385 - ctx       - optional context variable set with `KSPFlexibleSetModifyPC()`
386 
387   Level: beginner
388 
389 .seealso: [](ch_ksp), `KSP`, `KSPFlexibleSetModifyPC()`
390 S*/
391 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPFlexibleModifyPCFn(KSP ksp, PetscInt total_its, PetscInt local_its, PetscReal res_norm, PetscCtx ctx);
392 
393 PETSC_EXTERN PetscErrorCode KSPFlexibleSetModifyPC(KSP, KSPFlexibleModifyPCFn *, PetscCtx, PetscCtxDestroyFn *);
394 
395 PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleSetModifyPC()", )
KSPPIPEGCRSetModifyPC(KSP ksp,KSPFlexibleModifyPCFn * fun,PetscCtx ctx,PetscCtxDestroyFn * dfun)396 static inline PetscErrorCode KSPPIPEGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fun, PetscCtx ctx, PetscCtxDestroyFn *dfun)
397 {
398   return KSPFlexibleSetModifyPC(ksp, fun, ctx, dfun);
399 }
400 
401 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
402 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
403 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
404 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
405 
406 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
407 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
408 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
409 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
410 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
411 
412 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
413 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
414 
415 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
416 
417 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
418 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
419 
420 PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleSetModifyPC()", )
KSPGCRSetModifyPC(KSP ksp,KSPFlexibleModifyPCFn * fun,PetscCtx ctx,PetscCtxDestroyFn * dfun)421 static inline PetscErrorCode KSPGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fun, PetscCtx ctx, PetscCtxDestroyFn *dfun)
422 {
423   return KSPFlexibleSetModifyPC(ksp, fun, ctx, dfun);
424 }
425 
426 PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
427 PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
428 PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);
429 
430 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
431 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
432 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
433 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
434 
435 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
436 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
437 #if PetscDefined(HAVE_HPDDM)
KSPHPDDMSetDeflationSpace(KSP ksp,Mat U)438 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
439 {
440   return KSPHPDDMSetDeflationMat(ksp, U);
441 }
KSPHPDDMGetDeflationSpace(KSP ksp,Mat * U)442 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
443 {
444   return KSPHPDDMGetDeflationMat(ksp, U);
445 }
446 #endif
KSPHPDDMMatSolve(KSP ksp,Mat B,Mat X)447 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
448 {
449   return KSPMatSolve(ksp, B, X);
450 }
451 /*E
452     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
453 
454     Values:
455 +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
456 .   `KSP_HPDDM_TYPE_BGMRES`          - block GMRES
457 .   `KSP_HPDDM_TYPE_CG`              - Conjugate Gradient
458 .   `KSP_HPDDM_TYPE_BCG`             - block CG
459 .   `KSP_HPDDM_TYPE_GCRODR`          - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
460 .   `KSP_HPDDM_TYPE_BGCRODR`         - block GCRODR
461 .   `KSP_HPDDM_TYPE_BFBCG`           - breakdown-free BCG
462 -   `KSP_HPDDM_TYPE_PREONLY`         - apply the preconditioner only
463 
464     Level: intermediate
465 
466 .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
467 E*/
468 typedef enum {
469   KSP_HPDDM_TYPE_GMRES   = 0,
470   KSP_HPDDM_TYPE_BGMRES  = 1,
471   KSP_HPDDM_TYPE_CG      = 2,
472   KSP_HPDDM_TYPE_BCG     = 3,
473   KSP_HPDDM_TYPE_GCRODR  = 4,
474   KSP_HPDDM_TYPE_BGCRODR = 5,
475   KSP_HPDDM_TYPE_BFBCG   = 6,
476   KSP_HPDDM_TYPE_PREONLY = 7
477 } KSPHPDDMType;
478 PETSC_EXTERN const char *const KSPHPDDMTypes[];
479 
480 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
481 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
482 
483 /*E
484    KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers
485 
486    Values:
487 +  `KSP_GMRES_CGS_REFINE_NEVER`    - one step of classical Gram-Schmidt
488 .  `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
489 -  `KSP_GMRES_CGS_REFINE_ALWAYS`   - always perform two steps
490 
491    Level: advanced
492 
493 .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
494           `KSPGMRESGetOrthogonalization()`,
495           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
496 E*/
497 typedef enum {
498   KSP_GMRES_CGS_REFINE_NEVER,
499   KSP_GMRES_CGS_REFINE_IFNEEDED,
500   KSP_GMRES_CGS_REFINE_ALWAYS
501 } KSPGMRESCGSRefinementType;
502 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
503 
504 /*MC
505    KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
506 
507    Level: advanced
508 
509    Note:
510    Possibly unstable, but the fastest to compute
511 
512 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
513           `KSP`, `KSPGMRESGetOrthogonalization()`,
514           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
515           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
516 M*/
517 
518 /*MC
519     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
520           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
521           poor orthogonality.
522 
523    Level: advanced
524 
525    Note:
526    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
527    estimate the orthogonality but is more stable.
528 
529 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
530           `KSP`, `KSPGMRESGetOrthogonalization()`,
531           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
532           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
533 M*/
534 
535 /*MC
536    KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process.
537 
538    Level: advanced
539 
540    Notes:
541    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
542    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
543 
544    You should only use this if you absolutely know that the iterative refinement is needed.
545 
546 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
547           `KSP`, `KSPGMRESGetOrthogonalization()`,
548           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
549           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
550 M*/
551 
552 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
553 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
554 
555 PETSC_EXTERN KSPFlexibleModifyPCFn KSPFlexibleModifyPCNoChange;
556 PETSC_EXTERN KSPFlexibleModifyPCFn KSPFlexibleModifyPCKSP;
557 
558 PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleModifyPCNoChange()", )
KSPFGMRESModifyPCNoChange(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,PetscCtx ctx)559 static inline PetscErrorCode KSPFGMRESModifyPCNoChange(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx)
560 {
561   return KSPFlexibleModifyPCNoChange(ksp, total_its, loc_its, res_norm, ctx);
562 }
563 
564 PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleModifyPCKSP()", )
KSPFGMRESModifyPCKSP(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,PetscCtx ctx)565 static inline PetscErrorCode KSPFGMRESModifyPCKSP(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx)
566 {
567   return KSPFlexibleModifyPCKSP(ksp, total_its, loc_its, res_norm, ctx);
568 }
569 
570 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
571 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
572 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
573 
574 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
575 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
576 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
577 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
578 
579 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
580 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
581 
582 PETSC_EXTERN PetscErrorCode       KSPMonitorSetFromOptions(KSP, const char[], const char[], PetscCtx);
583 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidual;
584 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualView;
KSPMonitorResidualDraw(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)585 PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorResidualDraw()", ) static inline PetscErrorCode KSPMonitorResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
586 {
587   return KSPMonitorResidualView(ksp, n, rnorm, vf);
588 }
589 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDrawLG;
590 PETSC_EXTERN PetscErrorCode       KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
591 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualShort;
592 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualRange;
593 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidual;
594 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualView;
KSPMonitorTrueResidualDraw(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)595 PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorTrueResidualDraw()", ) static inline PetscErrorCode KSPMonitorTrueResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
596 {
597   return KSPMonitorTrueResidualView(ksp, n, rnorm, vf);
598 }
599 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDrawLG;
600 PETSC_EXTERN PetscErrorCode       KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
601 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualMax;
602 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorError;
603 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDraw;
604 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDrawLG;
605 PETSC_EXTERN PetscErrorCode       KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
606 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolution;
607 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDraw;
608 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDrawLG;
609 PETSC_EXTERN PetscErrorCode       KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
610 PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSingularValue;
611 PETSC_EXTERN PetscErrorCode       KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
KSPMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)612 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
613 {
614   return KSPMonitorResidual(ksp, n, rnorm, vf);
615 }
KSPMonitorTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)616 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
617 {
618   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
619 }
KSPMonitorTrueResidualMaxNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)620 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
621 {
622   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
623 }
624 
625 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
626 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, PetscCtx);
627 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(PetscCtxRt);
628 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
629 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
630 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
631 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
632 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(PetscCtxRt);
633 
634 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
635 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
636 
637 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
638 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
639 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
640 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
641 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
642 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
643 
644 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
645 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
646 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
647 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
648 
649 /*S
650   KSPConvergedReasonViewFn - A prototype of a function used with `KSPConvergedReasonViewSet()`
651 
652   Calling Sequence:
653 + ksp - the `KSP` object whose `KSPConvergedReason` is to be viewed
654 - ctx - context used by the function, set with `KSPConvergedReasonViewSet()`
655 
656   Level: beginner
657 
658 .seealso: [](ch_ksp), `KSP`, `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`, `KSPConvergedReasonViewFromOptions()`, `KSPView()`
659 S*/
660 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergedReasonViewFn(KSP ksp, PetscCtx ctx);
661 
662 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
663 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
664 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
665 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
666 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, KSPConvergedReasonViewFn *, void *, PetscCtxDestroyFn *);
667 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
668 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
669 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
670 
KSPReasonView(KSP ksp,PetscViewer v)671 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
672 {
673   return KSPConvergedReasonView(ksp, v);
674 }
KSPReasonViewFromOptions(KSP ksp)675 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
676 {
677   return KSPConvergedReasonViewFromOptions(ksp);
678 }
679 
680 #define KSP_FILE_CLASSID 1211223
681 
682 PETSC_EXTERN PetscErrorCode       KSPLSQRSetExactMatNorm(KSP, PetscBool);
683 PETSC_EXTERN PetscErrorCode       KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
684 PETSC_EXTERN PetscErrorCode       KSPLSQRGetStandardErrorVec(KSP, Vec *);
685 PETSC_EXTERN PetscErrorCode       KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
686 PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidual;
687 PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidualDrawLG;
688 PETSC_EXTERN PetscErrorCode       KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
689 
690 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
691 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
692 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
693 PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);
694 
695 /*E
696    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
697        test routines.
698 
699    Values:
700 +  `KSP_NORM_DEFAULT`          - use the default for the current `KSPType`
701 .  `KSP_NORM_NONE`             - use no norm calculation
702 .  `KSP_NORM_PRECONDITIONED`   - use the preconditioned residual norm
703 .  `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
704 -  `KSP_NORM_NATURAL`          - use the natural norm (the norm induced by the linear operator)
705 
706    Level: advanced
707 
708    Note:
709    Each solver only supports a subset of these and some may support different ones
710    depending on whether left or right preconditioning is used, see `KSPSetPCSide()`
711 
712 .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
713           `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
714 E*/
715 typedef enum {
716   KSP_NORM_DEFAULT          = -1,
717   KSP_NORM_NONE             = 0,
718   KSP_NORM_PRECONDITIONED   = 1,
719   KSP_NORM_UNPRECONDITIONED = 2,
720   KSP_NORM_NATURAL          = 3
721 } KSPNormType;
722 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
723 PETSC_EXTERN const char *const *const KSPNormTypes;
724 
725 /*MC
726    KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
727    possibly save some computation but means the convergence test cannot
728    be based on a norm of a residual etc.
729 
730    Level: advanced
731 
732    Note:
733    Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
734 
735 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
736 M*/
737 
738 /*MC
739    KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
740    convergence test routine.
741 
742    Level: advanced
743 
744 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
745 M*/
746 
747 /*MC
748    KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
749    convergence test routine.
750 
751    Level: advanced
752 
753 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
754 M*/
755 
756 /*MC
757    KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
758    convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
759 
760    Level: advanced
761 
762 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
763 M*/
764 
765 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
766 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
767 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt);
768 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
769 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
770 
771 #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
772 #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
773 #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", )
774 #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", )
775 /*E
776    KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
777 
778    Values:
779 +  `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR`
780 .  `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR`
781 .  `KSP_CONVERGED_RTOL`                  - requested decrease in the residual
782 .  `KSP_CONVERGED_ATOL`                  - requested absolute value in the residual
783 .  `KSP_CONVERGED_ITS`                   - requested number of iterations
784 .  `KSP_CONVERGED_NEG_CURVE`             - see note below
785 .  `KSP_CONVERGED_STEP_LENGTH`           - see note below
786 .  `KSP_CONVERGED_HAPPY_BREAKDOWN`       - happy breakdown (meaning early convergence of the `KSPType` occurred).
787 .  `KSP_CONVERGED_USER`                  - the user has indicated convergence for an arbitrary reason
788 .  `KSP_DIVERGED_NULL`                   - breakdown when solving the Hessenberg system within `KSPGMRES`
789 .  `KSP_DIVERGED_ITS`                    - requested number of iterations
790 .  `KSP_DIVERGED_DTOL`                   - large increase in the residual norm indicating the solution is diverging
791 .  `KSP_DIVERGED_BREAKDOWN`              - breakdown in the Krylov method
792 .  `KSP_DIVERGED_BREAKDOWN_BICG`         - breakdown in the `KSPBCGS` Krylov method
793 .  `KSP_DIVERGED_NONSYMMETRIC`           - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
794 .  `KSP_DIVERGED_INDEFINITE_PC`          - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
795 .  `KSP_DIVERGED_NANORINF`               - a not a number of infinity was detected in a vector during the computation
796 .  `KSP_DIVERGED_INDEFINITE_MAT`         - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
797 .  `KSP_DIVERGED_PC_FAILED`              - the action of the preconditioner failed for some reason
798 -  `KSP_DIVERGED_USER`                   - the user has indicated divergence for an arbitrary reason
799 
800    Level: beginner
801 
802    Note:
803    The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
804    the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
805 
806    Developer Note:
807    The string versions of these are `KSPConvergedReasons`; if you change
808    any of the values here also change them that array of names.
809 
810 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
811 E*/
812 typedef enum { /* converged */
813   KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    = 1,
814   KSP_CONVERGED_RTOL_NORMAL_EQUATIONS     = 1,
815   KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    = 9,
816   KSP_CONVERGED_ATOL_NORMAL_EQUATIONS     = 9,
817   KSP_CONVERGED_RTOL                      = 2,
818   KSP_CONVERGED_ATOL                      = 3,
819   KSP_CONVERGED_ITS                       = 4,
820   KSP_CONVERGED_NEG_CURVE                 = 5,
821   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
822   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
823   KSP_CONVERGED_STEP_LENGTH               = 6,
824   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
825   KSP_CONVERGED_USER                      = 8,
826   /* diverged */
827   KSP_DIVERGED_NULL                      = -2,
828   KSP_DIVERGED_ITS                       = -3,
829   KSP_DIVERGED_DTOL                      = -4,
830   KSP_DIVERGED_BREAKDOWN                 = -5,
831   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
832   KSP_DIVERGED_NONSYMMETRIC              = -7,
833   KSP_DIVERGED_INDEFINITE_PC             = -8,
834   KSP_DIVERGED_NANORINF                  = -9,
835   KSP_DIVERGED_INDEFINITE_MAT            = -10,
836   KSP_DIVERGED_PC_FAILED                 = -11,
837   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
838   KSP_DIVERGED_USER                      = -12,
839 
840   KSP_CONVERGED_ITERATING = 0
841 } KSPConvergedReason;
842 PETSC_EXTERN const char *const *KSPConvergedReasons;
843 
844 /*MC
845    KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called
846 
847    Level: beginner
848 
849    Notes:
850    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
851    for left preconditioning it is the 2-norm of the preconditioned residual, and the
852    2-norm of the residual for right preconditioning
853 
854    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
855 
856 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
857 M*/
858 
859 /*MC
860    KSP_CONVERGED_ATOL - $||r|| \le atol$
861 
862    Level: beginner
863 
864    Notes:
865    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
866    for left preconditioning it is the 2-norm of the preconditioned residual, and the
867    2-norm of the residual for right preconditioning
868 
869    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
870 
871 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
872 M*/
873 
874 /*MC
875    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$
876 
877    Level: beginner
878 
879    Note:
880    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
881    for left preconditioning it is the 2-norm of the preconditioned residual, and the
882    2-norm of the residual for right preconditioning
883 
884 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
885 M*/
886 
887 /*MC
888    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
889    reached
890 
891    Level: beginner
892 
893 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
894 M*/
895 
896 /*MC
897    KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
898    the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
899    test routine is set in `KSP`.
900 
901    Level: beginner
902 
903 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
904 M*/
905 
906 /*MC
907    KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
908    method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
909    preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
910    are collinear.
911 
912    Level: beginner
913 
914 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
915 M*/
916 
917 /*MC
918    KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
919    method could not continue to enlarge the Krylov space.
920 
921    Level: beginner
922 
923 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
924 M*/
925 
926 /*MC
927    KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
928    symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
929 
930    Level: beginner
931 
932 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
933 M*/
934 
935 /*MC
936    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
937    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
938    be symmetric positive definite (SPD).
939 
940    Level: beginner
941 
942    Note:
943    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
944    the `PCICC` preconditioner to generate a positive definite preconditioner
945 
946 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
947 M*/
948 
949 /*MC
950    KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
951    zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
952    such as `PCFIELDSPLIT`.
953 
954    Level: beginner
955 
956    Note:
957    Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
958 
959 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
960 M*/
961 
962 /*MC
963    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
964    while `KSPSolve()` is still running.
965 
966    Level: beginner
967 
968 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
969 M*/
970 
971 /*S
972   KSPConvergenceTestFn - A prototype of a function used with `KSPSetConvergenceTest()`
973 
974   Calling Sequence:
975 + ksp    - iterative solver obtained from `KSPCreate()`
976 . it     - iteration number
977 . rnorm  - (estimated) 2-norm of (preconditioned) residual
978 . reason - the reason why it has converged or diverged
979 - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`
980 
981   Level: beginner
982 
983 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
984 S*/
985 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergenceTestFn(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, PetscCtx ctx);
986 
987 PETSC_EXTERN PetscErrorCode       KSPSetConvergenceTest(KSP, KSPConvergenceTestFn *, void *, PetscCtxDestroyFn *);
988 PETSC_EXTERN PetscErrorCode       KSPGetConvergenceTest(KSP, KSPConvergenceTestFn **, PetscCtxRt, PetscCtxDestroyFn **);
989 PETSC_EXTERN PetscErrorCode       KSPGetAndClearConvergenceTest(KSP, KSPConvergenceTestFn **, PetscCtxRt, PetscCtxDestroyFn **);
990 PETSC_EXTERN PetscErrorCode       KSPGetConvergenceContext(KSP, PetscCtxRt);
991 PETSC_EXTERN KSPConvergenceTestFn KSPConvergedDefault;
992 PETSC_EXTERN KSPConvergenceTestFn KSPLSQRConvergedDefault;
993 PETSC_EXTERN PetscCtxDestroyFn    KSPConvergedDefaultDestroy;
994 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultCreate(void **);
995 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUIRNorm(KSP);
996 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUMIRNorm(KSP);
997 PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
998 PETSC_EXTERN PetscErrorCode       KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
999 PETSC_EXTERN PetscErrorCode       KSPGetConvergedReason(KSP, KSPConvergedReason *);
1000 PETSC_EXTERN PetscErrorCode       KSPGetConvergedReasonString(KSP, const char *[]);
1001 PETSC_EXTERN PetscErrorCode       KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1002 PETSC_EXTERN PetscErrorCode       KSPSetConvergedNegativeCurvature(KSP, PetscBool);
1003 PETSC_EXTERN PetscErrorCode       KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
1004 
KSPDefaultConverged(void)1005 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
1006 { /* never called */
1007 }
1008 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
KSPDefaultConvergedDestroy(void)1009 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
1010 { /* never called */
1011 }
1012 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
KSPDefaultConvergedCreate(void)1013 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
1014 { /* never called */
1015 }
1016 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
KSPDefaultConvergedSetUIRNorm(void)1017 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
1018 { /* never called */
1019 }
1020 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
KSPDefaultConvergedSetUMIRNorm(void)1021 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
1022 { /* never called */
1023 }
1024 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
KSPSkipConverged(void)1025 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
1026 { /* never called */
1027 }
1028 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
1029 
1030 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
KSPComputeExplicitOperator(KSP A,Mat * B)1031 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
1032 {
1033   return KSPComputeOperator(A, PETSC_NULLPTR, B);
1034 }
1035 
1036 /*E
1037    KSPCGType - Determines what type of `KSPCG` to use
1038 
1039    Values:
1040  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
1041  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
1042 
1043    Level: beginner
1044 
1045 .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
1046 E*/
1047 typedef enum {
1048   KSP_CG_SYMMETRIC = 0,
1049   KSP_CG_HERMITIAN = 1
1050 } KSPCGType;
1051 PETSC_EXTERN const char *const KSPCGTypes[];
1052 
1053 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
1054 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
1055 
1056 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
1057 PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
1058 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
1059 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
1060 
1061 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
1062 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
KSPCGGLTRGetMinEig(KSP ksp,PetscReal * x)1063 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
1064 {
1065   return KSPGLTRGetMinEig(ksp, x);
1066 }
KSPCGGLTRGetLambda(KSP ksp,PetscReal * x)1067 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
1068 {
1069   return KSPGLTRGetLambda(ksp, x);
1070 }
1071 
1072 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
1073 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
1074 
1075 PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
1076 PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
1077 
1078 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
1079 
1080 /*S
1081   PCShellPSolveFn - A function prototype for functions provided to `PCShellSetPreSolve()` and `PCShellSetPostSolve()`
1082 
1083   Calling Sequence:
1084 + pc  - the preconditioner `PC` context
1085 . ksp - the `KSP` context
1086 . xin  - input vector
1087 - xout - output vector
1088 
1089   Level: intermediate
1090 
1091 .seealso: [](ch_snes), `KSPPSolveFn`, `KSP`, `PCShellSetPreSolve()`, `PCShellSetPostSolve()`
1092 S*/
1093 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCShellPSolveFn(PC pc, KSP ksp, Vec xim, Vec xout);
1094 
1095 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PCShellPSolveFn *);
1096 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PCShellPSolveFn *);
1097 
1098 /*S
1099    KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
1100 
1101    Level: intermediate
1102 
1103    Note:
1104    These methods generate initial guesses based on a series of previous, related, linear solves. For example,
1105    in implicit time-stepping with `TS`.
1106 
1107 .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
1108 S*/
1109 typedef struct _p_KSPGuess *KSPGuess;
1110 
1111 /*J
1112    KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
1113 
1114    Values:
1115  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
1116  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition (POD)
1117 
1118    Level: intermediate
1119 
1120 .seealso: [](ch_ksp), `KSP`, `KSPGuess`
1121 J*/
1122 typedef const char *KSPGuessType;
1123 #define KSPGUESSFISCHER "fischer"
1124 #define KSPGUESSPOD     "pod"
1125 
1126 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
1127 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
1128 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
1129 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
1130 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
1131 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
1132 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
1133 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
1134 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
1135 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
1136 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
1137 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
1138 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
1139 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
1140 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
1141 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
1142 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
1143 
1144 /*E
1145     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines
1146 
1147     Level: intermediate
1148 
1149 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
1150           `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
1151 E*/
1152 typedef enum {
1153   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
1154   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
1155   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
1156   MAT_SCHUR_COMPLEMENT_AINV_FULL
1157 } MatSchurComplementAinvType;
1158 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
1159 
1160 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
1161 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
1162 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
1163 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1164 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1165 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
1166 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
1167 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
1168 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
1169 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
1170 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
1171 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
1172 
1173 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1174 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1175 PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1176 PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1177 PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
1178 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1179 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1180 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1181 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1182 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1183 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1184 
1185 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1186 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1187 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1188 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1189 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1190 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1191 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1192 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1193 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1194 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1195 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1196 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1197 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1198 PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *);
1199 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1200 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1201 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1202 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1203 PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1204 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1205 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1206 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1207 
1208 /*E
1209   MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix
1210 
1211   Values:
1212 + `MAT_LMVM_MULT_RECURSIVE`     - Use recursive formulas for products and solves
1213 . `MAT_LMVM_MULT_DENSE`         - Use dense formulas for products and solves when possible
1214 - `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
1215 
1216   Level: advanced
1217 
1218   Options Database Keys:
1219 . -mat_lmvm_mult_algorithm  - the algorithm to use for multiplication (recursive, dense, compact_dense)
1220 
1221 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()`
1222 E*/
1223 typedef enum {
1224   MAT_LMVM_MULT_RECURSIVE,
1225   MAT_LMVM_MULT_DENSE,
1226   MAT_LMVM_MULT_COMPACT_DENSE,
1227 } MatLMVMMultAlgorithm;
1228 
1229 PETSC_EXTERN const char *const MatLMVMMultAlgorithms[];
1230 
1231 PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm);
1232 PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *);
1233 
1234 /*E
1235   MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix.
1236 
1237   Values:
1238 + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - no rescaling
1239 . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar rescaling
1240 . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling
1241 . `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE`
1242 - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE`   - let PETSc decide rescaling
1243 
1244   Level: intermediate
1245 
1246 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1247 E*/
1248 typedef enum {
1249   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1250   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1251   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1252   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3,
1253   MAT_LMVM_SYMBROYDEN_SCALE_DECIDE   = 4
1254 } MatLMVMSymBroydenScaleType;
1255 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
1256 
1257 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1258 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *);
1259 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal);
1260 PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *);
1261 PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal);
1262 
1263 /*E
1264   MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.
1265 
1266   Values:
1267 + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1268 - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement
1269 
1270   Level: intermediate
1271 
1272 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1273 E*/
1274 typedef enum {
1275   MAT_LMVM_DENSE_REORDER,
1276   MAT_LMVM_DENSE_INPLACE
1277 } MatLMVMDenseType;
1278 PETSC_EXTERN const char *const MatLMVMDenseTypes[];
1279 
1280 PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);
1281 
1282 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1283 
1284 /*E
1285   KSPDMActive - Indicates if the `DM` attached to the `KSP` should be used to compute the operator, the right-hand side, or the initial guess
1286 
1287   Values:
1288 + `KSP_DMACTIVE_OPERATOR`      - compute the operator
1289 . `KSP_DMACTIVE_RHS`           - compute the right-hand side
1290 . `KSP_DMACTIVE_INITIAL_GUESS` - compute the initial guess
1291 - `KSP_DMACTIVE_ALL`           - compute all of them
1292 
1293   Level: intermediate
1294 
1295 .seealso: [](ch_ksp), `KSP`, `KSPSetDMActive()`, `KSPSetDM()`
1296 E*/
1297 typedef enum {
1298   KSP_DMACTIVE_OPERATOR      = 1,
1299   KSP_DMACTIVE_RHS           = 2,
1300   KSP_DMACTIVE_INITIAL_GUESS = 4,
1301   KSP_DMACTIVE_ALL           = 1 + 2 + 4
1302 } KSPDMActive;
1303 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, KSPDMActive, PetscBool);
1304 
1305 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1306 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, PetscCtx);
1307 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, PetscCtxRt);
1308 
1309 /*S
1310   KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`
1311 
1312   Calling Sequence:
1313 + ksp  - `ksp` context
1314 . b    - output vector
1315 - ctx - [optional] user-defined function context
1316 
1317   Level: beginner
1318 
1319 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
1320 S*/
1321 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeRHSFn(KSP ksp, Vec b, PetscCtx ctx);
1322 
1323 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);
1324 
1325 /*S
1326   KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`
1327 
1328   Calling Sequence:
1329 + ksp - `KSP` context
1330 . A   - the operator that defines the linear system
1331 . P   - an operator from which to build the preconditioner (often the same as `A`)
1332 - ctx - [optional] user-defined function context
1333 
1334   Level: beginner
1335 
1336 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
1337 S*/
1338 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeOperatorsFn(KSP ksp, Mat A, Mat P, PetscCtx ctx);
1339 
1340 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);
1341 
1342 /*S
1343   KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`
1344 
1345   Calling Sequence:
1346 + ksp  - `ksp` context
1347 . x    - output vector
1348 - ctx - [optional] user-defined function context
1349 
1350   Level: beginner
1351 
1352 .seealso: [](ch_ksp), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
1353 S*/
1354 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeInitialGuessFn(KSP ksp, Vec x, PetscCtx ctx);
1355 
1356 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1357 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1358 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1359 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1360 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1361 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1362 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);
1363 
1364 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1365 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode);
1366 PETSC_EXTERN PetscErrorCode DMSwarmProjectGradientFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode);
1367 
1368 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1369 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1370 
1371 PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1372 PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);
1373 
1374 PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);
1375 
1376 #include <petscdstypes.h>
1377 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, PetscPointFn **, InsertMode, Vec);
1378