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 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`, meaning that only application of the preconditioner is used as the linear solver.
23: .seealso: [](chapter_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: [](chapter_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 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: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind`
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 KSPChebyshevEstEigGetKSP(KSP, KSP *);
221: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
222: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
223: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
224: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
226: /*E
228: KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
230: Values:
231: + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
232: - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
234: Level: intermediate
236: .seealso: [](chapter_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
237: E*/
238: typedef enum {
239: KSP_FCD_TRUNC_TYPE_STANDARD,
240: KSP_FCD_TRUNC_TYPE_NOTAY
241: } KSPFCDTruncationType;
242: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
244: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
245: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
246: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
247: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
248: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
249: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
251: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
252: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
253: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
254: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
255: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
256: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
258: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
259: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
260: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
261: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
262: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
263: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
264: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
265: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
267: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
268: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
269: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
270: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
272: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
273: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
274: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
275: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
276: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
278: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
279: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
281: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
283: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
284: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
285: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
287: PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
288: PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
289: PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);
291: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
292: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
293: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
294: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
296: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
297: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
298: #if PetscDefined(HAVE_HPDDM)
299: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
300: {
301: return KSPHPDDMSetDeflationMat(ksp, U);
302: }
303: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
304: {
305: return KSPHPDDMGetDeflationMat(ksp, U);
306: }
307: #endif
308: PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
309: {
310: return KSPMatSolve(ksp, B, X);
311: }
312: /*E
313: KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
315: Level: intermediate
317: Values:
318: + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
319: . `KSP_HPDDM_TYPE_BGMRES` - block GMRES
320: . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient
321: . `KSP_HPDDM_TYPE_BCG` - block CG
322: . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
323: . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR
324: . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG
325: - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only
327: .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
328: E*/
329: typedef enum {
330: KSP_HPDDM_TYPE_GMRES = 0,
331: KSP_HPDDM_TYPE_BGMRES = 1,
332: KSP_HPDDM_TYPE_CG = 2,
333: KSP_HPDDM_TYPE_BCG = 3,
334: KSP_HPDDM_TYPE_GCRODR = 4,
335: KSP_HPDDM_TYPE_BGCRODR = 5,
336: KSP_HPDDM_TYPE_BFBCG = 6,
337: KSP_HPDDM_TYPE_PREONLY = 7
338: } KSPHPDDMType;
339: PETSC_EXTERN const char *const KSPHPDDMTypes[];
341: /*E
342: KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
344: Level: intermediate
346: Values:
347: + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16`
348: . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single`
349: . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double`
350: - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
352: .seealso: [](chapter_ksp), `KSP`, `KSPHPDDM`
353: E*/
354: typedef enum {
355: KSP_HPDDM_PRECISION_HALF = 0,
356: KSP_HPDDM_PRECISION_SINGLE = 1,
357: KSP_HPDDM_PRECISION_DOUBLE = 2,
358: KSP_HPDDM_PRECISION_QUADRUPLE = 3
359: } KSPHPDDMPrecision;
360: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
361: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
363: /*E
364: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
366: Level: advanced
368: Values:
369: + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt
370: . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
371: - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps
373: .seealso: [](chapter_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
374: `KSPGMRESGetOrthogonalization()`,
375: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
376: E*/
377: typedef enum {
378: KSP_GMRES_CGS_REFINE_NEVER,
379: KSP_GMRES_CGS_REFINE_IFNEEDED,
380: KSP_GMRES_CGS_REFINE_ALWAYS
381: } KSPGMRESCGSRefinementType;
382: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
383: /*MC
384: KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
386: Level: advanced
388: Note:
389: Possibly unstable, but the fastest to compute
391: .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
392: `KSP`, `KSPGMRESGetOrthogonalization()`,
393: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
394: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
395: M*/
397: /*MC
398: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
399: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
400: poor orthogonality.
402: Level: advanced
404: Note:
405: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
406: estimate the orthogonality but is more stable.
408: .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
409: `KSP`, `KSPGMRESGetOrthogonalization()`,
410: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
411: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
412: M*/
414: /*MC
415: KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
417: Level: advanced
419: Notes:
420: This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
421: but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
423: You should only use this if you absolutely know that the iterative refinement is needed.
425: .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
426: `KSP`, `KSPGMRESGetOrthogonalization()`,
427: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
428: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
429: M*/
431: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
432: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
434: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
435: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
436: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
438: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
439: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
440: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
442: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
443: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
444: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
445: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
447: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
448: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
449: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
451: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
452: PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
453: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
454: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
455: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
456: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
457: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
458: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
459: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
460: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
461: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
462: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
463: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
464: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
465: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
466: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
467: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
468: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
469: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
470: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
471: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
472: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
473: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
474: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
475: {
476: return KSPMonitorResidual(ksp, n, rnorm, vf);
477: }
478: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
479: {
480: return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
481: }
482: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
483: {
484: return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
485: }
487: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
488: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
489: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
490: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
491: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
492: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
493: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
494: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
496: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
497: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
499: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
500: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
501: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
502: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
503: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
504: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
506: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
507: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
508: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
509: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
511: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
512: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
513: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
514: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
515: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
516: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
517: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
518: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
520: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
521: {
522: return KSPConvergedReasonView(ksp, v);
523: }
524: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
525: {
526: return KSPConvergedReasonViewFromOptions(ksp);
527: }
529: #define KSP_FILE_CLASSID 1211223
531: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
532: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
533: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
534: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
535: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
536: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
537: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
539: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
540: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
541: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
543: /*E
544: KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
545: test routines.
547: Values:
548: + `KSP_NORM_DEFAULT` - use the default for the current `KSPType`
549: . `KSP_NORM_NONE` - use no norm calculation
550: . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm
551: . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
552: - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator)
554: Level: advanced
556: Note:
557: Each solver only supports a subset of these and some may support different ones
558: depending on left or right preconditioning, see `KSPSetPCSide()`
560: Developer Note:
561: This must match the values in petsc/finclude/petscksp.h
563: .seealso: [](chapter_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
564: `KSPSetConvergenceTest()`, `KSPSetPCSide()`
565: E*/
566: typedef enum {
567: KSP_NORM_DEFAULT = -1,
568: KSP_NORM_NONE = 0,
569: KSP_NORM_PRECONDITIONED = 1,
570: KSP_NORM_UNPRECONDITIONED = 2,
571: KSP_NORM_NATURAL = 3
572: } KSPNormType;
573: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
574: PETSC_EXTERN const char *const *const KSPNormTypes;
576: /*MC
577: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
578: possibly save some computation but means the convergence test cannot
579: be based on a norm of a residual etc.
581: Level: advanced
583: Note:
584: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
586: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
587: M*/
589: /*MC
590: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
591: convergence test routine.
593: Level: advanced
595: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
596: M*/
598: /*MC
599: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
600: convergence test routine.
602: Level: advanced
604: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
605: M*/
607: /*MC
608: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
609: convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
611: Level: advanced
613: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
614: M*/
616: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
617: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
618: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
619: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
620: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
622: #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_NEG_CURVE (since version 3.19)")
623: #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_STEP_LENGTH (since version 3.19)")
624: #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
625: /*E
626: KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
628: Level: beginner
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: Note:
651: The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`,
652: `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
654: Developer Notes:
655: This must match the values in petsc/finclude/petscksp.h
657: The string versions of these are `KSPConvergedReasons`; if you change
658: any of the values here also change them that array of names.
660: .seealso: [](chapter_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
661: E*/
662: typedef enum { /* converged */
663: KSP_CONVERGED_RTOL_NORMAL = 1,
664: KSP_CONVERGED_ATOL_NORMAL = 9,
665: KSP_CONVERGED_RTOL = 2,
666: KSP_CONVERGED_ATOL = 3,
667: KSP_CONVERGED_ITS = 4,
668: KSP_CONVERGED_NEG_CURVE = 5,
669: KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5,
670: KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
671: KSP_CONVERGED_STEP_LENGTH = 6,
672: KSP_CONVERGED_HAPPY_BREAKDOWN = 7,
673: /* diverged */
674: KSP_DIVERGED_NULL = -2,
675: KSP_DIVERGED_ITS = -3,
676: KSP_DIVERGED_DTOL = -4,
677: KSP_DIVERGED_BREAKDOWN = -5,
678: KSP_DIVERGED_BREAKDOWN_BICG = -6,
679: KSP_DIVERGED_NONSYMMETRIC = -7,
680: KSP_DIVERGED_INDEFINITE_PC = -8,
681: KSP_DIVERGED_NANORINF = -9,
682: KSP_DIVERGED_INDEFINITE_MAT = -10,
683: KSP_DIVERGED_PC_FAILED = -11,
684: KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
686: KSP_CONVERGED_ITERATING = 0
687: } KSPConvergedReason;
688: PETSC_EXTERN const char *const *KSPConvergedReasons;
690: /*MC
691: KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called
693: Level: beginner
695: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
696: for left preconditioning it is the 2-norm of the preconditioned residual, and the
697: 2-norm of the residual for right preconditioning
699: See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
701: .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
703: M*/
705: /*MC
706: KSP_CONVERGED_ATOL - norm(r) <= atol
708: Level: beginner
710: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
711: for left preconditioning it is the 2-norm of the preconditioned residual, and the
712: 2-norm of the residual for right preconditioning
714: See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
716: .seealso: [](chapter_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: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
733: M*/
735: /*MC
736: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
737: reached
739: Level: beginner
741: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
743: M*/
745: /*MC
746: KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
747: the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
748: test routine is set in KSP.
750: Level: beginner
752: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
754: M*/
756: /*MC
757: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
758: method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
759: preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
760: are colinear.
762: Level: beginner
764: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
766: M*/
768: /*MC
769: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
770: method could not continue to enlarge the Krylov space.
772: Level: beginner
774: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
776: M*/
778: /*MC
779: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
780: symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
782: Level: beginner
784: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
786: M*/
788: /*MC
789: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
790: positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
791: be positive definite
793: Level: beginner
795: Note:
796: This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
797: the `PCICC` preconditioner to generate a positive definite preconditioner
799: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
801: M*/
803: /*MC
804: KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
805: zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
806: such as `PCFIELDSPLIT`.
808: Level: beginner
810: Note:
811: Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
813: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
815: M*/
817: /*MC
818: KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
819: while `KSPSolve()` is still running.
821: Level: beginner
823: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
825: M*/
827: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
828: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
829: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
830: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
831: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
832: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
833: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
834: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
835: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
836: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
837: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
838: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
839: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
840: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
841: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
842: PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
843: PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
845: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void)
846: { /* never called */
847: }
848: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
849: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void)
850: { /* never called */
851: }
852: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
853: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void)
854: { /* never called */
855: }
856: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
857: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void)
858: { /* never called */
859: }
860: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
861: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void)
862: { /* never called */
863: }
864: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
865: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void)
866: { /* never called */
867: }
868: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
870: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
871: PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
872: {
873: return KSPComputeOperator(A, NULL, B);
874: }
876: /*E
877: KSPCGType - Determines what type of `KSPCG` to use
879: Level: beginner
881: Values:
882: + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
883: - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
885: .seealso: [](chapter_ksp), `KSP`, `KSPCGSetType()`
886: E*/
887: typedef enum {
888: KSP_CG_SYMMETRIC = 0,
889: KSP_CG_HERMITIAN = 1
890: } KSPCGType;
891: PETSC_EXTERN const char *const KSPCGTypes[];
893: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
894: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
896: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
897: PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
898: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
899: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
901: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
902: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
903: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
904: {
905: return KSPGLTRGetMinEig(ksp, x);
906: }
907: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
908: {
909: return KSPGLTRGetLambda(ksp, x);
910: }
912: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
913: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
915: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
916: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
917: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
919: #include <petscdrawtypes.h>
920: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
922: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
923: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
925: /*S
926: KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
928: Level: intermediate
930: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
931: S*/
932: typedef struct _p_KSPGuess *KSPGuess;
933: /*J
934: KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
936: Level: intermediate
938: Values:
939: + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
940: - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition
942: .seealso: [](chapter_ksp), `KSP`, `KSPGuess`
943: J*/
944: typedef const char *KSPGuessType;
945: #define KSPGUESSFISCHER "fischer"
946: #define KSPGUESSPOD "pod"
948: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
949: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
950: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
951: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
952: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
953: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
954: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
955: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
956: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
957: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
958: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
959: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
960: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
961: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
962: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
963: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
964: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
966: /*E
967: MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
969: Level: intermediate
971: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
972: `MatCreateSchurComplementPmat()`
973: E*/
974: typedef enum {
975: MAT_SCHUR_COMPLEMENT_AINV_DIAG,
976: MAT_SCHUR_COMPLEMENT_AINV_LUMP,
977: MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
978: MAT_SCHUR_COMPLEMENT_AINV_FULL
979: } MatSchurComplementAinvType;
980: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
982: typedef enum {
983: MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0,
984: MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1,
985: MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
986: MAT_LMVM_SYMBROYDEN_SCALE_USER = 3
987: } MatLMVMSymBroydenScaleType;
988: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
990: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
991: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
992: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
993: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
994: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
995: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
996: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
997: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
998: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
999: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
1000: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
1001: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
1003: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1004: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1005: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1006: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1007: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1008: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1009: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1010: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1012: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1013: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1014: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1015: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1016: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1017: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1018: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1019: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1020: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1021: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1022: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1023: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1024: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1025: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1026: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1027: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1028: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1029: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1030: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1031: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1032: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1034: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1035: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1036: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1037: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1038: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1039: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
1040: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1041: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
1042: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1043: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
1044: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1045: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1046: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1047: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1049: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1050: 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);
1052: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1053: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1054: #endif