Actual source code: petscksp.h

  1: /*
  2:    Defines the interface functions for the Krylov subspace accelerators.
  3: */
  4: #ifndef PETSCKSP_H
  5: #define PETSCKSP_H

  7: #include <petscpc.h>

  9: /* SUBMANSEC = KSP */

 11: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);

 13: /*S
 14:      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
 15:          linear solves in PETSc (even those such as direct factorization-based solvers that do no use Krylov accelerators).

 17:    Level: beginner

 19:    Note:
 20:     When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a
 21:     `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver.

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

 27: /*J
 28:     KSPType - String with the name of a PETSc Krylov method.

 30:    Level: beginner

 32: .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
 33: J*/
 34: typedef const char *KSPType;
 35: #define KSPRICHARDSON "richardson"
 36: #define KSPCHEBYSHEV  "chebyshev"
 37: #define KSPCG         "cg"
 38: #define KSPGROPPCG    "groppcg"
 39: #define KSPPIPECG     "pipecg"
 40: #define KSPPIPECGRR   "pipecgrr"
 41: #define KSPPIPELCG    "pipelcg"
 42: #define KSPPIPEPRCG   "pipeprcg"
 43: #define KSPPIPECG2    "pipecg2"
 44: #define KSPCGNE       "cgne"
 45: #define KSPNASH       "nash"
 46: #define KSPSTCG       "stcg"
 47: #define KSPGLTR       "gltr"
 48: #define KSPCGNASH     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash"
 49: #define KSPCGSTCG     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg"
 50: #define KSPCGGLTR     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
 51: #define KSPFCG        "fcg"
 52: #define KSPPIPEFCG    "pipefcg"
 53: #define KSPGMRES      "gmres"
 54: #define KSPPIPEFGMRES "pipefgmres"
 55: #define KSPFGMRES     "fgmres"
 56: #define KSPLGMRES     "lgmres"
 57: #define KSPDGMRES     "dgmres"
 58: #define KSPPGMRES     "pgmres"
 59: #define KSPTCQMR      "tcqmr"
 60: #define KSPBCGS       "bcgs"
 61: #define KSPIBCGS      "ibcgs"
 62: #define KSPQMRCGS     "qmrcgs"
 63: #define KSPFBCGS      "fbcgs"
 64: #define KSPFBCGSR     "fbcgsr"
 65: #define KSPBCGSL      "bcgsl"
 66: #define KSPPIPEBCGS   "pipebcgs"
 67: #define KSPCGS        "cgs"
 68: #define KSPTFQMR      "tfqmr"
 69: #define KSPCR         "cr"
 70: #define KSPPIPECR     "pipecr"
 71: #define KSPLSQR       "lsqr"
 72: #define KSPPREONLY    "preonly"
 73: #define KSPNONE       "none"
 74: #define KSPQCG        "qcg"
 75: #define KSPBICG       "bicg"
 76: #define KSPMINRES     "minres"
 77: #define KSPSYMMLQ     "symmlq"
 78: #define KSPLCD        "lcd"
 79: #define KSPPYTHON     "python"
 80: #define KSPGCR        "gcr"
 81: #define KSPPIPEGCR    "pipegcr"
 82: #define KSPTSIRM      "tsirm"
 83: #define KSPCGLS       "cgls"
 84: #define KSPFETIDP     "fetidp"
 85: #define KSPHPDDM      "hpddm"

 87: /* Logging support */
 88: PETSC_EXTERN PetscClassId KSP_CLASSID;
 89: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
 90: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 92: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
 93: PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
 94: PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
 95: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
 96: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
 97: PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
 98: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
 99: PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
100: PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
101: PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat);
102: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
103: PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
104: {
105:   return KSPSetMatSolveBatchSize(ksp, n);
106: }
107: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
108: PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
109: {
110:   return KSPGetMatSolveBatchSize(ksp, n);
111: }
112: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
113: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
114: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
115: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
116: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
117: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
118: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);

120: PETSC_EXTERN PetscFunctionList KSPList;
121: PETSC_EXTERN PetscFunctionList KSPGuessList;
122: PETSC_EXTERN PetscFunctionList KSPMonitorList;
123: PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
124: PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
125: PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
126: PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));

128: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
129: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
130: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
131: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
132: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
133: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
134: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
135: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
136: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
137: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
138: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
139: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
140: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
141: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
142: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
143: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
144: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
145: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
146: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
147: PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
148: {
149:   return KSPCreateVecs(ksp, n, x, m, y);
150: }

152: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
153: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);

155: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
156: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);

158: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
159: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
160: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
161: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
162: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
163: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool);
164: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
165: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool);

167: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
168: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
169: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
170: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);

172: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
173: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
174: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
175: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
176: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
177: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
178: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
179: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
180: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
181: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
182: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
183: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
184: PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
185: /*
186:   PCMGCoarseList contains the list of coarse space constructor currently registered
187:   These are added with PCMGRegisterCoarseSpaceConstructor()
188: */
189: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
190: PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
191: PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));

193: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
194: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);

196: /*E
197:   KSPChebyshevKind - Which kind of Chebyshev polynomial to use

199:   Values:
200: + `KSP_CHEBYSHEV_FIRST`      - "classic" first-kind Chebyshev polynomial
201: . `KSP_CHEBYSHEV_FOURTH`     - fourth-kind Chebyshev polynomial
202: - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial

204:   Level: intermediate

206: .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()`
207: E*/
208: typedef enum {
209:   KSP_CHEBYSHEV_FIRST,
210:   KSP_CHEBYSHEV_FOURTH,
211:   KSP_CHEBYSHEV_OPT_FOURTH
212: } KSPChebyshevKind;

214: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
215: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
216: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
217: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
218: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
219: PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind);
220: PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *);
221: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
222: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
223: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
224: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
225: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);

227: /*E

229:   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods

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

235:    Level: intermediate

237: .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
238: E*/
239: typedef enum {
240:   KSP_FCD_TRUNC_TYPE_STANDARD,
241:   KSP_FCD_TRUNC_TYPE_NOTAY
242: } KSPFCDTruncationType;
243: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

245: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
246: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
247: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
248: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
249: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
250: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);

252: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
253: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
254: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
255: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
256: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
257: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);

259: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
260: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
261: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
262: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
263: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
264: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
265: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
266: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);

268: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
269: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
270: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
271: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);

273: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
274: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
275: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
276: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
277: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);

279: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
280: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

282: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);

284: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
285: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
286: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

288: PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
289: PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
290: PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);

292: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
293: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
294: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
295: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);

297: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
298: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
299: #if PetscDefined(HAVE_HPDDM)
300: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
301: {
302:   return KSPHPDDMSetDeflationMat(ksp, U);
303: }
304: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
305: {
306:   return KSPHPDDMGetDeflationMat(ksp, U);
307: }
308: #endif
309: PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
310: {
311:   return KSPMatSolve(ksp, B, X);
312: }
313: /*E
314:     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`

316:     Values:
317: +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
318: .   `KSP_HPDDM_TYPE_BGMRES` - block GMRES
319: .   `KSP_HPDDM_TYPE_CG` - Conjugate Gradient
320: .   `KSP_HPDDM_TYPE_BCG` - block CG
321: .   `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
322: .   `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR
323: .   `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG
324: -   `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only

326:     Level: intermediate

328: .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
329: E*/
330: typedef enum {
331:   KSP_HPDDM_TYPE_GMRES   = 0,
332:   KSP_HPDDM_TYPE_BGMRES  = 1,
333:   KSP_HPDDM_TYPE_CG      = 2,
334:   KSP_HPDDM_TYPE_BCG     = 3,
335:   KSP_HPDDM_TYPE_GCRODR  = 4,
336:   KSP_HPDDM_TYPE_BGCRODR = 5,
337:   KSP_HPDDM_TYPE_BFBCG   = 6,
338:   KSP_HPDDM_TYPE_PREONLY = 7
339: } KSPHPDDMType;
340: PETSC_EXTERN const char *const KSPHPDDMTypes[];

342: /*E
343:     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`

345:     Values:
346: +   `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16`
347: .   `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single`
348: .   `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double`
349: -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`

351:     Level: intermediate

353: .seealso: [](ch_ksp), `KSP`, `KSPHPDDM`
354: E*/
355: typedef enum {
356:   KSP_HPDDM_PRECISION_HALF      = 0,
357:   KSP_HPDDM_PRECISION_SINGLE    = 1,
358:   KSP_HPDDM_PRECISION_DOUBLE    = 2,
359:   KSP_HPDDM_PRECISION_QUADRUPLE = 3
360: } KSPHPDDMPrecision;
361: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
362: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);

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

367:     Values:
368: +   `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt
369: .   `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
370: -   `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps

372:    Level: advanced

374: .seealso: [](ch_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
375:           `KSPGMRESGetOrthogonalization()`,
376:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
377: E*/
378: typedef enum {
379:   KSP_GMRES_CGS_REFINE_NEVER,
380:   KSP_GMRES_CGS_REFINE_IFNEEDED,
381:   KSP_GMRES_CGS_REFINE_ALWAYS
382: } KSPGMRESCGSRefinementType;
383: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];

385: /*MC
386:     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

388:    Level: advanced

390:    Note:
391:    Possibly unstable, but the fastest to compute

393: .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
394:           `KSP`, `KSPGMRESGetOrthogonalization()`,
395:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
396:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
397: M*/

399: /*MC
400:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
401:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
402:           poor orthogonality.

404:    Level: advanced

406:    Note:
407:    This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
408:    estimate the orthogonality but is more stable.

410: .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
411:           `KSP`, `KSPGMRESGetOrthogonalization()`,
412:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
413:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
414: M*/

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

419:    Level: advanced

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

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

427: .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
428:           `KSP`, `KSPGMRESGetOrthogonalization()`,
429:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
430:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
431: M*/

433: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
434: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);

436: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
437: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
438: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

440: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
441: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
442: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);

444: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
445: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
446: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
447: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);

449: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
450: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
451: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

453: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
454: PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
455: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
456: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
457: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
458: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
459: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
460: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
461: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
462: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
463: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
464: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
465: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
466: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
467: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
468: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
469: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
470: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
471: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
472: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
473: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
474: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
475: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
476: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
477: {
478:   return KSPMonitorResidual(ksp, n, rnorm, vf);
479: }
480: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
481: {
482:   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
483: }
484: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
485: {
486:   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
487: }

489: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
490: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
491: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
492: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
493: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
494: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
495: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
496: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);

498: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
499: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);

501: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
502: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
503: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
504: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
505: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
506: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);

508: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
509: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
510: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
511: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);

513: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
514: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
515: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
516: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
517: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
518: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
519: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
520: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);

522: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
523: {
524:   return KSPConvergedReasonView(ksp, v);
525: }
526: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
527: {
528:   return KSPConvergedReasonViewFromOptions(ksp);
529: }

531: #define KSP_FILE_CLASSID 1211223

533: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
534: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
535: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
536: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
537: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
538: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
539: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);

541: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
542: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
543: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);

545: /*E
546:     KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
547:        test routines.

549:     Values:
550: +   `KSP_NORM_DEFAULT` - use the default for the current `KSPType`
551: .   `KSP_NORM_NONE` - use no norm calculation
552: .   `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm
553: .   `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
554: -   `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator)

556:    Level: advanced

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

562:    Developer Note:
563:    This must match the values in petsc/finclude/petscksp.h

565: .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
566:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
567: E*/
568: typedef enum {
569:   KSP_NORM_DEFAULT          = -1,
570:   KSP_NORM_NONE             = 0,
571:   KSP_NORM_PRECONDITIONED   = 1,
572:   KSP_NORM_UNPRECONDITIONED = 2,
573:   KSP_NORM_NATURAL          = 3
574: } KSPNormType;
575: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
576: PETSC_EXTERN const char *const *const KSPNormTypes;

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

583:    Level: advanced

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

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

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

595:    Level: advanced

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

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

604:    Level: advanced

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

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

613:    Level: advanced

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

618: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
619: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
620: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
621: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
622: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);

624: #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_NEG_CURVE (since version 3.19)")
625: #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_STEP_LENGTH (since version 3.19)")
626: #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
627: /*E
628:     KSPConvergedReason - reason a Krylov method was determined to have converged or diverged

630:    Values:
631: +  `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations
632: .  `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations
633: .  `KSP_CONVERGED_RTOL` - requested decrease in the residual
634: .  `KSP_CONVERGED_ATOL` - requested absolute value in the residual
635: .  `KSP_CONVERGED_ITS` - requested number of iterations
636: .  `KSP_CONVERGED_NEG_CURVE` - see note below
637: .  `KSP_CONVERGED_STEP_LENGTH` - see note below
638: .  `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred.
639: .  `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES
640: .  `KSP_DIVERGED_ITS` - requested number of iterations
641: .  `KSP_DIVERGED_DTOL` - large increase in the residual norm
642: .  `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method
643: .  `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method
644: .  `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
645: .  `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite
646: .  `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation
647: .  `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite
648: -  `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason

650:    Level: beginner

652:    Note:
653:    The values  `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`,
654:    `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.

656:    Developer Notes:
657:    This must match the values in petsc/finclude/petscksp.h

659:    The string versions of these are `KSPConvergedReasons`; if you change
660:    any of the values here also change them that array of names.

662: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
663: E*/
664: typedef enum { /* converged */
665:   KSP_CONVERGED_RTOL_NORMAL               = 1,
666:   KSP_CONVERGED_ATOL_NORMAL               = 9,
667:   KSP_CONVERGED_RTOL                      = 2,
668:   KSP_CONVERGED_ATOL                      = 3,
669:   KSP_CONVERGED_ITS                       = 4,
670:   KSP_CONVERGED_NEG_CURVE                 = 5,
671:   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
672:   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
673:   KSP_CONVERGED_STEP_LENGTH               = 6,
674:   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
675:   /* diverged */
676:   KSP_DIVERGED_NULL                      = -2,
677:   KSP_DIVERGED_ITS                       = -3,
678:   KSP_DIVERGED_DTOL                      = -4,
679:   KSP_DIVERGED_BREAKDOWN                 = -5,
680:   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
681:   KSP_DIVERGED_NONSYMMETRIC              = -7,
682:   KSP_DIVERGED_INDEFINITE_PC             = -8,
683:   KSP_DIVERGED_NANORINF                  = -9,
684:   KSP_DIVERGED_INDEFINITE_MAT            = -10,
685:   KSP_DIVERGED_PC_FAILED                 = -11,
686:   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,

688:   KSP_CONVERGED_ITERATING = 0
689: } KSPConvergedReason;
690: PETSC_EXTERN const char *const *KSPConvergedReasons;

692: /*MC
693:      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called

695:    Level: beginner

697:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
698:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
699:    2-norm of the residual for right preconditioning

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

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

706: /*MC
707:      KSP_CONVERGED_ATOL - norm(r) <= atol

709:    Level: beginner

711:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
712:        for left preconditioning it is the 2-norm of the preconditioned residual, and the
713:        2-norm of the residual for right preconditioning

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

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

720: /*MC
721:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

723:    Level: beginner

725:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
726:        for left preconditioning it is the 2-norm of the preconditioned residual, and the
727:        2-norm of the residual for right preconditioning

729:    Level: beginner

731: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
732: M*/

734: /*MC
735:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
736:       reached

738:    Level: beginner

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

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

748:    Level: beginner

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

753: /*MC
754:      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
755:           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
756:           preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
757:           are colinear.

759:    Level: beginner

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

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

768:    Level: beginner

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

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

777:    Level: beginner

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

782: /*MC
783:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
784:         positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
785:         be positive definite

787:    Level: beginner

789:      Note:
790:     This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
791:   the `PCICC` preconditioner to generate a positive definite preconditioner

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

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

801:    Level: beginner

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

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

809: /*MC
810:      KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
811:         while `KSPSolve()` is still running.

813:    Level: beginner

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

818: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
819: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
820: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
821: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
822: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
823: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
824: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
825: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
826: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
827: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
828: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
829: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
830: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
831: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
832: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
833: PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
834: PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);

836: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void)
837: { /* never called */
838: }
839: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
840: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void)
841: { /* never called */
842: }
843: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
844: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void)
845: { /* never called */
846: }
847: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
848: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void)
849: { /* never called */
850: }
851: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
852: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void)
853: { /* never called */
854: }
855: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
856: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void)
857: { /* never called */
858: }
859: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

861: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
862: PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
863: {
864:   return KSPComputeOperator(A, PETSC_NULLPTR, B);
865: }

867: /*E
868:     KSPCGType - Determines what type of `KSPCG` to use

870:    Values:
871:  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
872:  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian

874:    Level: beginner

876: .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
877: E*/
878: typedef enum {
879:   KSP_CG_SYMMETRIC = 0,
880:   KSP_CG_HERMITIAN = 1
881: } KSPCGType;
882: PETSC_EXTERN const char *const KSPCGTypes[];

884: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
885: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);

887: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
888: PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
889: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
890: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);

892: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
893: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
894: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
895: {
896:   return KSPGLTRGetMinEig(ksp, x);
897: }
898: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
899: {
900:   return KSPGLTRGetLambda(ksp, x);
901: }

903: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
904: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);

906: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
907: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
908: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);

910: #include <petscdrawtypes.h>
911: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);

913: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
914: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));

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

919:    Level: intermediate

921: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
922: S*/
923: typedef struct _p_KSPGuess *KSPGuess;

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

928:    Values:
929:  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
930:  - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition

932:    Level: intermediate

934: .seealso: [](ch_ksp), `KSP`, `KSPGuess`
935: J*/
936: typedef const char *KSPGuessType;
937: #define KSPGUESSFISCHER "fischer"
938: #define KSPGUESSPOD     "pod"

940: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
941: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
942: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
943: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
944: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
945: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
946: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
947: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
948: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
949: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
950: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
951: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
952: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
953: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
954: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
955: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
956: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);

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

961:     Level: intermediate

963: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
964:           `MatCreateSchurComplementPmat()`
965: E*/
966: typedef enum {
967:   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
968:   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
969:   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
970:   MAT_SCHUR_COMPLEMENT_AINV_FULL
971: } MatSchurComplementAinvType;
972: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

974: typedef enum {
975:   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
976:   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
977:   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
978:   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
979: } MatLMVMSymBroydenScaleType;
980: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

982: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
983: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
984: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
985: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
986: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
987: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
988: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
989: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
990: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
991: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
992: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
993: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);

995: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
996: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
997: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
998: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
999: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1000: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1001: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1002: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);

1004: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1005: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1006: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1007: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1008: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1009: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1010: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1011: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1012: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1013: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1014: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1015: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1016: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1017: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1018: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1019: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1020: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1021: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1022: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1023: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1024: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);

1026: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1027: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1028: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1029: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1030: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1031: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
1032: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1033: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
1034: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1035: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
1036: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1037: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1038: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1039: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);

1041: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1042: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);

1044: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1045: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1046: #endif