Actual source code: petscksp.h

  1: /*
  2:    Defines the interface functions for the Krylov subspace accelerators.
  3: */
  4: #pragma once

  6: #include <petscpc.h>

  8: /* SUBMANSEC = KSP */

 10: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
 11: PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void);

 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).

 17:    Level: beginner

 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.

 23:    Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object

 25:    The `PC` object is used to control preconditioners in PETSc.

 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).

 30: .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
 31: S*/
 32: typedef struct _p_KSP *KSP;

 34: /*J
 35:    KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides.

 37:    Level: beginner

 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"

 94: /* Logging support */
 95: PETSC_EXTERN PetscClassId KSP_CLASSID;
 96: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
 97: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 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);

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));

134: /*S
135:   KSPMonitorRegisterFn - A function prototype for functions provided to `KSPMonitorRegister()`

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

143:   Level: beginner

145:   Note:
146:   This is a `KSPMonitorFn` specialized for a context of `PetscViewerAndFormat`

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);

152: /*S
153:   KSPMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to `KSPMonitorRegister()`

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

161:   Level: beginner

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);

167: /*S
168:   KSPMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to `KSPMonitorRegister()`

170:   Calling Sequence:
171: . vf - a `PetscViewerAndFormat` object to be destroyed, including any context

173:   Level: beginner

175: .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn`
176: S*/
177: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterDestroyFn(PetscViewerAndFormat **result);

179: PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, KSPMonitorRegisterFn *, KSPMonitorRegisterCreateFn *, KSPMonitorRegisterDestroyFn *);

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: }

207: /*S
208:   KSPPSolveFn - A function prototype for functions provided to `KSPSetPreSolve()` and `KSPSetPostSolve()`

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()`

216:   Level: intermediate

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);

222: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, KSPPSolveFn *, PetscCtx);
223: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, KSPPSolveFn *, PetscCtx);

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 *);

230: /*S
231:   KSPMonitorFn - A function prototype for functions provided to `KSPMonitorSet()`

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()`

239:   Level: beginner

241: .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()`
242: S*/
243: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorFn(KSP ksp, PetscInt it, PetscReal rnorm, PetscCtx ctx);

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);

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);

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 *);

274: /*S
275:   PCMGCoarseSpaceConstructorFn - A function prototype for functions registered with `PCMGRegisterCoarseSpaceConstructor()`

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

286:   Level: beginner

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);

292: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
293: PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn *);
294: PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn **);

296: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
297: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);

299: /*E
300:   KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV`

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

307:   Level: intermediate

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;

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[]);

330: /*E
331:   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods

333:   Values:
334: + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to `mmax`) stored directions
335: - `KSP_FCD_TRUNC_TYPE_NOTAY`    - uses the last `max(1,mod(i,mmax))` stored directions at iteration i = 0, 1, ...

337:    Level: intermediate

339:   Note:
340:   Function such as `KSPFCGSetMmax()`, `KSPPIPEGCRSetNMax()`, `KSPPIPEGCRSetNMax()`, and `KSPPIPEFCGSetNMax()` may be
341:   used to provide `nmax` or they may be provided with the option database.

343: .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`,
344:           `KSPPIPEGCRSetTruncationType()`, `KSPPIPEGCRGetTruncationType()`,
345:           `KSPFCGSetMmax()`, `KSPPIPEGCRSetNMax()`, `KSPPIPEGCRGetNMax()`, `KSPPIPEFCGGetNMax()`
346: E*/
347: typedef enum {
348:   KSP_FCD_TRUNC_TYPE_STANDARD,
349:   KSP_FCD_TRUNC_TYPE_NOTAY
350: } KSPFCDTruncationType;
351: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

353: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
354: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
355: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
356: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
357: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
358: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);

360: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
361: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
362: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
363: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
364: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
365: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);

367: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
368: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
369: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
370: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
371: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
372: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
373: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
374: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);

376: /*S
377:   KSPFlexibleModifyPCFn - A prototype of a function used to modify the preconditioner during the use of flexible `KSP` methods, such as `KSPFGMRES`

379:   Calling Sequence:
380: + ksp       - the `KSP` context being used.
381: . total_its - the total number of iterations that have occurred.
382: . local_its - the number of iterations since last restart if applicable
383: . res_norm  - the current residual norm
384: - ctx       - optional context variable set with `KSPFlexibleSetModifyPC()`

386:   Level: beginner

388: .seealso: [](ch_ksp), `KSP`, `KSPFlexibleSetModifyPC()`
389: S*/
390: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPFlexibleModifyPCFn(KSP ksp, PetscInt total_its, PetscInt local_its, PetscReal res_norm, PetscCtx ctx);

392: PETSC_EXTERN PetscErrorCode KSPFlexibleSetModifyPC(KSP, KSPFlexibleModifyPCFn *, PetscCtx, PetscCtxDestroyFn *);

394: PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleSetModifyPC()", )
395: static inline PetscErrorCode KSPPIPEGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fun, PetscCtx ctx, PetscCtxDestroyFn *dfun)
396: {
397:   return KSPFlexibleSetModifyPC(ksp, fun, ctx, dfun);
398: }

400: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
401: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
402: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
403: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);

405: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
406: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
407: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
408: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
409: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);

411: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
412: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

414: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);

416: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
417: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);

419: PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleSetModifyPC()", )
420: static inline PetscErrorCode KSPGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fun, PetscCtx ctx, PetscCtxDestroyFn *dfun)
421: {
422:   return KSPFlexibleSetModifyPC(ksp, fun, ctx, dfun);
423: }

425: PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
426: PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
427: PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);

429: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
430: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
431: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
432: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);

434: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
435: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
436: #if PetscDefined(HAVE_HPDDM)
437: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
438: {
439:   return KSPHPDDMSetDeflationMat(ksp, U);
440: }
441: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
442: {
443:   return KSPHPDDMGetDeflationMat(ksp, U);
444: }
445: #endif
446: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
447: {
448:   return KSPMatSolve(ksp, B, X);
449: }
450: /*E
451:     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`

453:     Values:
454: +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
455: .   `KSP_HPDDM_TYPE_BGMRES`          - block GMRES
456: .   `KSP_HPDDM_TYPE_CG`              - Conjugate Gradient
457: .   `KSP_HPDDM_TYPE_BCG`             - block CG
458: .   `KSP_HPDDM_TYPE_GCRODR`          - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
459: .   `KSP_HPDDM_TYPE_BGCRODR`         - block GCRODR
460: .   `KSP_HPDDM_TYPE_BFBCG`           - breakdown-free BCG
461: -   `KSP_HPDDM_TYPE_PREONLY`         - apply the preconditioner only

463:     Level: intermediate

465: .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
466: E*/
467: typedef enum {
468:   KSP_HPDDM_TYPE_GMRES   = 0,
469:   KSP_HPDDM_TYPE_BGMRES  = 1,
470:   KSP_HPDDM_TYPE_CG      = 2,
471:   KSP_HPDDM_TYPE_BCG     = 3,
472:   KSP_HPDDM_TYPE_GCRODR  = 4,
473:   KSP_HPDDM_TYPE_BGCRODR = 5,
474:   KSP_HPDDM_TYPE_BFBCG   = 6,
475:   KSP_HPDDM_TYPE_PREONLY = 7
476: } KSPHPDDMType;
477: PETSC_EXTERN const char *const KSPHPDDMTypes[];

479: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
480: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);

482: /*E
483:    KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers

485:    Values:
486: +  `KSP_GMRES_CGS_REFINE_NEVER`    - one step of classical Gram-Schmidt
487: .  `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
488: -  `KSP_GMRES_CGS_REFINE_ALWAYS`   - always perform two steps

490:    Level: advanced

492: .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
493:           `KSPGMRESGetOrthogonalization()`,
494:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
495: E*/
496: typedef enum {
497:   KSP_GMRES_CGS_REFINE_NEVER,
498:   KSP_GMRES_CGS_REFINE_IFNEEDED,
499:   KSP_GMRES_CGS_REFINE_ALWAYS
500: } KSPGMRESCGSRefinementType;
501: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];

503: /*MC
504:    KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

506:    Level: advanced

508:    Note:
509:    Possibly unstable, but the fastest to compute

511: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
512:           `KSP`, `KSPGMRESGetOrthogonalization()`,
513:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
514:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
515: M*/

517: /*MC
518:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
519:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
520:           poor orthogonality.

522:    Level: advanced

524:    Note:
525:    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
526:    estimate the orthogonality but is more stable.

528: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
529:           `KSP`, `KSPGMRESGetOrthogonalization()`,
530:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
531:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
532: M*/

534: /*MC
535:    KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process.

537:    Level: advanced

539:    Notes:
540:    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
541:    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.

543:    You should only use this if you absolutely know that the iterative refinement is needed.

545: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
546:           `KSP`, `KSPGMRESGetOrthogonalization()`,
547:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
548:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
549: M*/

551: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
552: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);

554: PETSC_EXTERN KSPFlexibleModifyPCFn KSPFlexibleModifyPCNoChange;
555: PETSC_EXTERN KSPFlexibleModifyPCFn KSPFlexibleModifyPCKSP;

557: PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleModifyPCNoChange()", )
558: static inline PetscErrorCode KSPFGMRESModifyPCNoChange(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx)
559: {
560:   return KSPFlexibleModifyPCNoChange(ksp, total_its, loc_its, res_norm, ctx);
561: }

563: PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleModifyPCKSP()", )
564: static inline PetscErrorCode KSPFGMRESModifyPCKSP(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx)
565: {
566:   return KSPFlexibleModifyPCKSP(ksp, total_its, loc_its, res_norm, ctx);
567: }

569: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
570: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
571: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);

573: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
574: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
575: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
576: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);

578: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
579: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);

581: PETSC_EXTERN PetscErrorCode       KSPMonitorSetFromOptions(KSP, const char[], const char[], PetscCtx);
582: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidual;
583: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualView;
584: PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorResidualDraw()", ) static inline PetscErrorCode KSPMonitorResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
585: {
586:   return KSPMonitorResidualView(ksp, n, rnorm, vf);
587: }
588: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDrawLG;
589: PETSC_EXTERN PetscErrorCode       KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
590: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualShort;
591: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualRange;
592: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidual;
593: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualView;
594: PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorTrueResidualDraw()", ) static inline PetscErrorCode KSPMonitorTrueResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
595: {
596:   return KSPMonitorTrueResidualView(ksp, n, rnorm, vf);
597: }
598: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDrawLG;
599: PETSC_EXTERN PetscErrorCode       KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
600: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualMax;
601: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorError;
602: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDraw;
603: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDrawLG;
604: PETSC_EXTERN PetscErrorCode       KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
605: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolution;
606: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDraw;
607: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDrawLG;
608: PETSC_EXTERN PetscErrorCode       KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
609: PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSingularValue;
610: PETSC_EXTERN PetscErrorCode       KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **);
611: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
612: {
613:   return KSPMonitorResidual(ksp, n, rnorm, vf);
614: }
615: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
616: {
617:   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
618: }
619: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
620: {
621:   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
622: }

624: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
625: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, PetscCtx);
626: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(PetscCtxRt);
627: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
628: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
629: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
630: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
631: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(PetscCtxRt);

633: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
634: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);

636: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
637: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
638: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
639: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
640: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
641: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);

643: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
644: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
645: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
646: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);

648: /*S
649:   KSPConvergedReasonViewFn - A prototype of a function used with `KSPConvergedReasonViewSet()`

651:   Calling Sequence:
652: + ksp - the `KSP` object whose `KSPConvergedReason` is to be viewed
653: - ctx - context used by the function, set with `KSPConvergedReasonViewSet()`

655:   Level: beginner

657: .seealso: [](ch_ksp), `KSP`, `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`, `KSPConvergedReasonViewFromOptions()`, `KSPView()`
658: S*/
659: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergedReasonViewFn(KSP ksp, PetscCtx ctx);

661: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
662: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
663: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
664: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
665: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, KSPConvergedReasonViewFn *, void *, PetscCtxDestroyFn *);
666: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
667: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
668: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);

670: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
671: {
672:   return KSPConvergedReasonView(ksp, v);
673: }
674: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
675: {
676:   return KSPConvergedReasonViewFromOptions(ksp);
677: }

679: #define KSP_FILE_CLASSID 1211223

681: PETSC_EXTERN PetscErrorCode       KSPLSQRSetExactMatNorm(KSP, PetscBool);
682: PETSC_EXTERN PetscErrorCode       KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
683: PETSC_EXTERN PetscErrorCode       KSPLSQRGetStandardErrorVec(KSP, Vec *);
684: PETSC_EXTERN PetscErrorCode       KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
685: PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidual;
686: PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidualDrawLG;
687: PETSC_EXTERN PetscErrorCode       KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);

689: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
690: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
691: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
692: PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);

694: /*E
695:    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
696:        test routines.

698:    Values:
699: +  `KSP_NORM_DEFAULT`          - use the default for the current `KSPType`
700: .  `KSP_NORM_NONE`             - use no norm calculation
701: .  `KSP_NORM_PRECONDITIONED`   - use the preconditioned residual norm
702: .  `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
703: -  `KSP_NORM_NATURAL`          - use the natural norm (the norm induced by the linear operator)

705:    Level: advanced

707:    Note:
708:    Each solver only supports a subset of these and some may support different ones
709:    depending on whether left or right preconditioning is used, see `KSPSetPCSide()`

711: .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
712:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
713: E*/
714: typedef enum {
715:   KSP_NORM_DEFAULT          = -1,
716:   KSP_NORM_NONE             = 0,
717:   KSP_NORM_PRECONDITIONED   = 1,
718:   KSP_NORM_UNPRECONDITIONED = 2,
719:   KSP_NORM_NATURAL          = 3
720: } KSPNormType;
721: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
722: PETSC_EXTERN const char *const *const KSPNormTypes;

724: /*MC
725:    KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
726:    possibly save some computation but means the convergence test cannot
727:    be based on a norm of a residual etc.

729:    Level: advanced

731:    Note:
732:    Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored

734: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
735: M*/

737: /*MC
738:    KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
739:    convergence test routine.

741:    Level: advanced

743: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
744: M*/

746: /*MC
747:    KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
748:    convergence test routine.

750:    Level: advanced

752: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
753: M*/

755: /*MC
756:    KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
757:    convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`

759:    Level: advanced

761: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
762: M*/

764: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
765: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
766: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt);
767: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
768: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);

770: #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
771: #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
772: #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", )
773: #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", )
774: /*E
775:    KSPConvergedReason - reason a Krylov method was determined to have converged or diverged

777:    Values:
778: +  `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR`
779: .  `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR`
780: .  `KSP_CONVERGED_RTOL`                  - requested decrease in the residual
781: .  `KSP_CONVERGED_ATOL`                  - requested absolute value in the residual
782: .  `KSP_CONVERGED_ITS`                   - requested number of iterations
783: .  `KSP_CONVERGED_NEG_CURVE`             - see note below
784: .  `KSP_CONVERGED_STEP_LENGTH`           - see note below
785: .  `KSP_CONVERGED_HAPPY_BREAKDOWN`       - happy breakdown (meaning early convergence of the `KSPType` occurred).
786: .  `KSP_CONVERGED_USER`                  - the user has indicated convergence for an arbitrary reason
787: .  `KSP_DIVERGED_NULL`                   - breakdown when solving the Hessenberg system within `KSPGMRES`
788: .  `KSP_DIVERGED_ITS`                    - requested number of iterations
789: .  `KSP_DIVERGED_DTOL`                   - large increase in the residual norm indicating the solution is diverging
790: .  `KSP_DIVERGED_BREAKDOWN`              - breakdown in the Krylov method
791: .  `KSP_DIVERGED_BREAKDOWN_BICG`         - breakdown in the `KSPBCGS` Krylov method
792: .  `KSP_DIVERGED_NONSYMMETRIC`           - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
793: .  `KSP_DIVERGED_INDEFINITE_PC`          - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
794: .  `KSP_DIVERGED_NANORINF`               - a not a number of infinity was detected in a vector during the computation
795: .  `KSP_DIVERGED_INDEFINITE_MAT`         - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
796: .  `KSP_DIVERGED_PC_FAILED`              - the action of the preconditioner failed for some reason
797: -  `KSP_DIVERGED_USER`                   - the user has indicated divergence for an arbitrary reason

799:    Level: beginner

801:    Note:
802:    The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
803:    the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.

805:    Developer Note:
806:    The string versions of these are `KSPConvergedReasons`; if you change
807:    any of the values here also change them that array of names.

809: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
810: E*/
811: typedef enum { /* converged */
812:   KSP_CONVERGED_RTOL_NORMAL_DEPRECATED    = 1,
813:   KSP_CONVERGED_RTOL_NORMAL_EQUATIONS     = 1,
814:   KSP_CONVERGED_ATOL_NORMAL_DEPRECATED    = 9,
815:   KSP_CONVERGED_ATOL_NORMAL_EQUATIONS     = 9,
816:   KSP_CONVERGED_RTOL                      = 2,
817:   KSP_CONVERGED_ATOL                      = 3,
818:   KSP_CONVERGED_ITS                       = 4,
819:   KSP_CONVERGED_NEG_CURVE                 = 5,
820:   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
821:   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
822:   KSP_CONVERGED_STEP_LENGTH               = 6,
823:   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
824:   KSP_CONVERGED_USER                      = 8,
825:   /* diverged */
826:   KSP_DIVERGED_NULL                      = -2,
827:   KSP_DIVERGED_ITS                       = -3,
828:   KSP_DIVERGED_DTOL                      = -4,
829:   KSP_DIVERGED_BREAKDOWN                 = -5,
830:   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
831:   KSP_DIVERGED_NONSYMMETRIC              = -7,
832:   KSP_DIVERGED_INDEFINITE_PC             = -8,
833:   KSP_DIVERGED_NANORINF                  = -9,
834:   KSP_DIVERGED_INDEFINITE_MAT            = -10,
835:   KSP_DIVERGED_PC_FAILED                 = -11,
836:   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
837:   KSP_DIVERGED_USER                      = -12,

839:   KSP_CONVERGED_ITERATING = 0
840: } KSPConvergedReason;
841: PETSC_EXTERN const char *const *KSPConvergedReasons;

843: /*MC
844:    KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called

846:    Level: beginner

848:    Notes:
849:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
850:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
851:    2-norm of the residual for right preconditioning

853:    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.

855: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
856: M*/

858: /*MC
859:    KSP_CONVERGED_ATOL - $||r|| \le atol$

861:    Level: beginner

863:    Notes:
864:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
865:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
866:    2-norm of the residual for right preconditioning

868:    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.

870: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
871: M*/

873: /*MC
874:    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$

876:    Level: beginner

878:    Note:
879:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
880:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
881:    2-norm of the residual for right preconditioning

883: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
884: M*/

886: /*MC
887:    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
888:    reached

890:    Level: beginner

892: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
893: M*/

895: /*MC
896:    KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
897:    the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
898:    test routine is set in `KSP`.

900:    Level: beginner

902: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
903: M*/

905: /*MC
906:    KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
907:    method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
908:    preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
909:    are collinear.

911:    Level: beginner

913: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
914: M*/

916: /*MC
917:    KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
918:    method could not continue to enlarge the Krylov space.

920:    Level: beginner

922: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
923: M*/

925: /*MC
926:    KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
927:    symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry

929:    Level: beginner

931: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
932: M*/

934: /*MC
935:    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
936:    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
937:    be symmetric positive definite (SPD).

939:    Level: beginner

941:    Note:
942:    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
943:    the `PCICC` preconditioner to generate a positive definite preconditioner

945: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
946: M*/

948: /*MC
949:    KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
950:    zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
951:    such as `PCFIELDSPLIT`.

953:    Level: beginner

955:    Note:
956:    Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.

958: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
959: M*/

961: /*MC
962:    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
963:    while `KSPSolve()` is still running.

965:    Level: beginner

967: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
968: M*/

970: /*S
971:   KSPConvergenceTestFn - A prototype of a function used with `KSPSetConvergenceTest()`

973:   Calling Sequence:
974: + ksp    - iterative solver obtained from `KSPCreate()`
975: . it     - iteration number
976: . rnorm  - (estimated) 2-norm of (preconditioned) residual
977: . reason - the reason why it has converged or diverged
978: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

980:   Level: beginner

982: .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
983: S*/
984: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergenceTestFn(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, PetscCtx ctx);

986: PETSC_EXTERN PetscErrorCode       KSPSetConvergenceTest(KSP, KSPConvergenceTestFn *, void *, PetscCtxDestroyFn *);
987: PETSC_EXTERN PetscErrorCode       KSPGetConvergenceTest(KSP, KSPConvergenceTestFn **, PetscCtxRt, PetscCtxDestroyFn **);
988: PETSC_EXTERN PetscErrorCode       KSPGetAndClearConvergenceTest(KSP, KSPConvergenceTestFn **, PetscCtxRt, PetscCtxDestroyFn **);
989: PETSC_EXTERN PetscErrorCode       KSPGetConvergenceContext(KSP, PetscCtxRt);
990: PETSC_EXTERN KSPConvergenceTestFn KSPConvergedDefault;
991: PETSC_EXTERN KSPConvergenceTestFn KSPLSQRConvergedDefault;
992: PETSC_EXTERN PetscCtxDestroyFn    KSPConvergedDefaultDestroy;
993: PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultCreate(void **);
994: PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUIRNorm(KSP);
995: PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetUMIRNorm(KSP);
996: PETSC_EXTERN PetscErrorCode       KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
997: PETSC_EXTERN PetscErrorCode       KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
998: PETSC_EXTERN PetscErrorCode       KSPGetConvergedReason(KSP, KSPConvergedReason *);
999: PETSC_EXTERN PetscErrorCode       KSPGetConvergedReasonString(KSP, const char *[]);
1000: PETSC_EXTERN PetscErrorCode       KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1001: PETSC_EXTERN PetscErrorCode       KSPSetConvergedNegativeCurvature(KSP, PetscBool);
1002: PETSC_EXTERN PetscErrorCode       KSPGetConvergedNegativeCurvature(KSP, PetscBool *);

1004: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
1005: { /* never called */
1006: }
1007: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
1008: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
1009: { /* never called */
1010: }
1011: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
1012: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
1013: { /* never called */
1014: }
1015: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
1016: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
1017: { /* never called */
1018: }
1019: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
1020: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
1021: { /* never called */
1022: }
1023: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
1024: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
1025: { /* never called */
1026: }
1027: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

1029: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
1030: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
1031: {
1032:   return KSPComputeOperator(A, PETSC_NULLPTR, B);
1033: }

1035: /*E
1036:    KSPCGType - Determines what type of `KSPCG` to use

1038:    Values:
1039:  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
1040:  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian

1042:    Level: beginner

1044: .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
1045: E*/
1046: typedef enum {
1047:   KSP_CG_SYMMETRIC = 0,
1048:   KSP_CG_HERMITIAN = 1
1049: } KSPCGType;
1050: PETSC_EXTERN const char *const KSPCGTypes[];

1052: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
1053: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);

1055: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
1056: PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
1057: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
1058: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);

1060: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
1061: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
1062: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
1063: {
1064:   return KSPGLTRGetMinEig(ksp, x);
1065: }
1066: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
1067: {
1068:   return KSPGLTRGetLambda(ksp, x);
1069: }

1071: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
1072: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);

1074: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
1075: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);

1077: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);

1079: /*S
1080:   PCShellPSolveFn - A function prototype for functions provided to `PCShellSetPreSolve()` and `PCShellSetPostSolve()`

1082:   Calling Sequence:
1083: + pc  - the preconditioner `PC` context
1084: . ksp - the `KSP` context
1085: . xin  - input vector
1086: - xout - output vector

1088:   Level: intermediate

1090: .seealso: [](ch_snes), `KSPPSolveFn`, `KSP`, `PCShellSetPreSolve()`, `PCShellSetPostSolve()`
1091: S*/
1092: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCShellPSolveFn(PC pc, KSP ksp, Vec xim, Vec xout);

1094: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PCShellPSolveFn *);
1095: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PCShellPSolveFn *);

1097: /*S
1098:    KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.

1100:    Level: intermediate

1102:    Note:
1103:    These methods generate initial guesses based on a series of previous, related, linear solves. For example,
1104:    in implicit time-stepping with `TS`.

1106: .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
1107: S*/
1108: typedef struct _p_KSPGuess *KSPGuess;

1110: /*J
1111:    KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.

1113:    Values:
1114:  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
1115:  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition (POD)

1117:    Level: intermediate

1119: .seealso: [](ch_ksp), `KSP`, `KSPGuess`
1120: J*/
1121: typedef const char *KSPGuessType;
1122: #define KSPGUESSFISCHER "fischer"
1123: #define KSPGUESSPOD     "pod"

1125: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
1126: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
1127: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
1128: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
1129: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
1130: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
1131: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
1132: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
1133: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
1134: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
1135: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
1136: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
1137: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
1138: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
1139: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
1140: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
1141: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);

1143: /*E
1144:     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines

1146:     Level: intermediate

1148: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
1149:           `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
1150: E*/
1151: typedef enum {
1152:   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
1153:   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
1154:   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
1155:   MAT_SCHUR_COMPLEMENT_AINV_FULL
1156: } MatSchurComplementAinvType;
1157: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

1159: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
1160: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
1161: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
1162: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1163: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
1164: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
1165: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
1166: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
1167: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
1168: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
1169: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
1170: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);

1172: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1173: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1174: PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1175: PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1176: PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
1177: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1178: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1179: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1180: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1181: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1182: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);

1184: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1185: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1186: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1187: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1188: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1189: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1190: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1191: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1192: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1193: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1194: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1195: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1196: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1197: PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *);
1198: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1199: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1200: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1201: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1202: PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1203: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1204: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1205: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);

1207: /*E
1208:   MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix

1210:   Values:
1211: + `MAT_LMVM_MULT_RECURSIVE`     - Use recursive formulas for products and solves
1212: . `MAT_LMVM_MULT_DENSE`         - Use dense formulas for products and solves when possible
1213: - `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:   Level: advanced

1217:   Options Database Keys:
1218: . -mat_lmvm_mult_algorithm (recursive|dense|compact_dense) - the algorithm to use for multiplication

1220: .seealso: [](ch_matrices), `MATLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()`
1221: E*/
1222: typedef enum {
1223:   MAT_LMVM_MULT_RECURSIVE,
1224:   MAT_LMVM_MULT_DENSE,
1225:   MAT_LMVM_MULT_COMPACT_DENSE,
1226: } MatLMVMMultAlgorithm;

1228: PETSC_EXTERN const char *const MatLMVMMultAlgorithms[];

1230: PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm);
1231: PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *);

1233: /*E
1234:   MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix.

1236:   Values:
1237: + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - no rescaling
1238: . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar rescaling
1239: . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling
1240: . `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE`
1241: - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE`   - let PETSc decide rescaling

1243:   Level: intermediate

1245: .seealso: [](ch_matrices), `MATLMVM`, `MatLMVMSymBroydenSetScaleType()`
1246: E*/
1247: typedef enum {
1248:   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1249:   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1250:   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1251:   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3,
1252:   MAT_LMVM_SYMBROYDEN_SCALE_DECIDE   = 4
1253: } MatLMVMSymBroydenScaleType;
1254: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

1256: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1257: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *);
1258: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal);
1259: PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *);
1260: PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal);

1262: /*E
1263:   MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.

1265:   Values:
1266: + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1267: - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement

1269:   Level: intermediate

1271: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1272: E*/
1273: typedef enum {
1274:   MAT_LMVM_DENSE_REORDER,
1275:   MAT_LMVM_DENSE_INPLACE
1276: } MatLMVMDenseType;
1277: PETSC_EXTERN const char *const MatLMVMDenseTypes[];

1279: PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);

1281: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);

1283: /*E
1284:   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:   Values:
1287: + `KSP_DMACTIVE_OPERATOR`      - compute the operator
1288: . `KSP_DMACTIVE_RHS`           - compute the right-hand side
1289: . `KSP_DMACTIVE_INITIAL_GUESS` - compute the initial guess
1290: - `KSP_DMACTIVE_ALL`           - compute all of them

1292:   Level: intermediate

1294: .seealso: [](ch_ksp), `KSP`, `KSPSetDMActive()`, `KSPSetDM()`
1295: E*/
1296: typedef enum {
1297:   KSP_DMACTIVE_OPERATOR      = 1,
1298:   KSP_DMACTIVE_RHS           = 2,
1299:   KSP_DMACTIVE_INITIAL_GUESS = 4,
1300:   KSP_DMACTIVE_ALL           = 1 + 2 + 4
1301: } KSPDMActive;
1302: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, KSPDMActive, PetscBool);

1304: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1305: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, PetscCtx);
1306: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, PetscCtxRt);

1308: /*S
1309:   KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`

1311:   Calling Sequence:
1312: + ksp  - `ksp` context
1313: . b    - output vector
1314: - ctx - [optional] user-defined function context

1316:   Level: beginner

1318: .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
1319: S*/
1320: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeRHSFn(KSP ksp, Vec b, PetscCtx ctx);

1322: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);

1324: /*S
1325:   KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`

1327:   Calling Sequence:
1328: + ksp - `KSP` context
1329: . A   - the operator that defines the linear system
1330: . P   - an operator from which to build the preconditioner (often the same as `A`)
1331: - ctx - [optional] user-defined function context

1333:   Level: beginner

1335: .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
1336: S*/
1337: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeOperatorsFn(KSP ksp, Mat A, Mat P, PetscCtx ctx);

1339: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);

1341: /*S
1342:   KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`

1344:   Calling Sequence:
1345: + ksp  - `ksp` context
1346: . x    - output vector
1347: - ctx - [optional] user-defined function context

1349:   Level: beginner

1351: .seealso: [](ch_ksp), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
1352: S*/
1353: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeInitialGuessFn(KSP ksp, Vec x, PetscCtx ctx);

1355: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1356: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1357: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1358: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1359: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1360: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1361: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);

1363: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1364: PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode);
1365: PETSC_EXTERN PetscErrorCode DMSwarmProjectGradientFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode);

1367: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1368: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);

1370: PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1371: PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);

1373: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);

1375: #include <petscdstypes.h>
1376: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, PetscPointFn **, InsertMode, Vec);