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 KSPSetMatSolveBatchSize(KSP, PetscInt);
102: PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
103: {
104: return KSPSetMatSolveBatchSize(ksp, n);
105: }
106: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
107: PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
108: {
109: return KSPGetMatSolveBatchSize(ksp, n);
110: }
111: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
112: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
113: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
114: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
115: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
116: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
117: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
119: PETSC_EXTERN PetscFunctionList KSPList;
120: PETSC_EXTERN PetscFunctionList KSPGuessList;
121: PETSC_EXTERN PetscFunctionList KSPMonitorList;
122: PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
123: PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
124: PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP));
125: PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
127: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
128: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
129: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
130: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
131: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
132: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
133: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
134: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
135: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
136: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
137: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
138: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
139: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
140: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
141: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
142: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
143: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
144: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
145: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
146: PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
147: {
148: return KSPCreateVecs(ksp, n, x, m, y);
149: }
151: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
152: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
154: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
155: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
157: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
158: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
159: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
160: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
161: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
162: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool);
163: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
164: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool);
166: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
167: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
168: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
169: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
171: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
172: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
173: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
174: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
175: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
176: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
177: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
178: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
179: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
180: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
181: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
182: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
183: PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
184: /*
185: PCMGCoarseList contains the list of coarse space constructor currently registered
186: These are added with PCMGRegisterCoarseSpaceConstructor()
187: */
188: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
189: PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
190: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
192: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
193: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
195: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
196: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
197: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
198: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
199: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
200: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
201: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
202: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
203: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
204: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
206: /*E
208: KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
210: Values:
211: + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
212: - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
214: Level: intermediate
216: .seealso: [](chapter_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
217: E*/
218: typedef enum {
219: KSP_FCD_TRUNC_TYPE_STANDARD,
220: KSP_FCD_TRUNC_TYPE_NOTAY
221: } KSPFCDTruncationType;
222: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
224: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
225: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
226: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
227: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
228: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
229: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
231: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
232: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
233: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
234: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
235: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
236: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
238: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
239: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
240: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
241: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
242: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
243: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
244: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
245: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
247: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
248: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
249: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
250: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
252: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
253: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
254: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
255: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
256: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
258: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
259: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
261: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
263: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
264: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
265: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
267: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
268: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
269: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
270: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
272: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
273: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
274: #if PetscDefined(HAVE_HPDDM)
275: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
276: {
277: return KSPHPDDMSetDeflationMat(ksp, U);
278: }
279: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
280: {
281: return KSPHPDDMGetDeflationMat(ksp, U);
282: }
283: #endif
284: PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
285: {
286: return KSPMatSolve(ksp, B, X);
287: }
288: /*E
289: KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
291: Level: intermediate
293: Values:
294: + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
295: . `KSP_HPDDM_TYPE_BGMRES` - block GMRES
296: . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient
297: . `KSP_HPDDM_TYPE_BCG` - block CG
298: . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
299: . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR
300: . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG
301: - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only
303: .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
304: E*/
305: typedef enum {
306: KSP_HPDDM_TYPE_GMRES = 0,
307: KSP_HPDDM_TYPE_BGMRES = 1,
308: KSP_HPDDM_TYPE_CG = 2,
309: KSP_HPDDM_TYPE_BCG = 3,
310: KSP_HPDDM_TYPE_GCRODR = 4,
311: KSP_HPDDM_TYPE_BGCRODR = 5,
312: KSP_HPDDM_TYPE_BFBCG = 6,
313: KSP_HPDDM_TYPE_PREONLY = 7
314: } KSPHPDDMType;
315: PETSC_EXTERN const char *const KSPHPDDMTypes[];
317: /*E
318: KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
320: Level: intermediate
322: Values:
323: + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16`
324: . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single`
325: . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double`
326: - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
328: .seealso: [](chapter_ksp), `KSP`, `KSPHPDDM`
329: E*/
330: typedef enum {
331: KSP_HPDDM_PRECISION_HALF = 0,
332: KSP_HPDDM_PRECISION_SINGLE = 1,
333: KSP_HPDDM_PRECISION_DOUBLE = 2,
334: KSP_HPDDM_PRECISION_QUADRUPLE = 3
335: } KSPHPDDMPrecision;
336: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
337: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
339: /*E
340: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
342: Level: advanced
344: Values:
345: + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt
346: . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
347: - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps
349: .seealso: [](chapter_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
350: `KSPGMRESGetOrthogonalization()`,
351: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
352: E*/
353: typedef enum {
354: KSP_GMRES_CGS_REFINE_NEVER,
355: KSP_GMRES_CGS_REFINE_IFNEEDED,
356: KSP_GMRES_CGS_REFINE_ALWAYS
357: } KSPGMRESCGSRefinementType;
358: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
359: /*MC
360: KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
362: Level: advanced
364: Note:
365: Possibly unstable, but the fastest to compute
367: .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
368: `KSP`, `KSPGMRESGetOrthogonalization()`,
369: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
370: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
371: M*/
373: /*MC
374: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
375: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
376: poor orthogonality.
378: Level: advanced
380: Note:
381: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
382: estimate the orthogonality but is more stable.
384: .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
385: `KSP`, `KSPGMRESGetOrthogonalization()`,
386: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
387: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
388: M*/
390: /*MC
391: KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
393: Level: advanced
395: Notes:
396: This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
397: but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
399: You should only use this if you absolutely know that the iterative refinement is needed.
401: .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
402: `KSP`, `KSPGMRESGetOrthogonalization()`,
403: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
404: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
405: M*/
407: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
408: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
410: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
411: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
412: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
414: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
415: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
416: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
418: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
419: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
420: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
421: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
423: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
424: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
425: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
427: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
428: PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
429: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
430: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
431: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
432: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
433: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
434: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
435: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
436: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
437: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
438: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
439: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
440: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
441: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
442: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
443: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
444: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
445: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
446: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
447: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
448: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
449: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
450: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
451: {
452: return KSPMonitorResidual(ksp, n, rnorm, vf);
453: }
454: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
455: {
456: return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
457: }
458: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
459: {
460: return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
461: }
463: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
464: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
465: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
466: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
467: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
468: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
469: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
470: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
472: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
473: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
475: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
476: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
477: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
478: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
479: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
480: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
482: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
483: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
484: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
485: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
487: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
488: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
489: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
490: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
491: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
492: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
493: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
494: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
496: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
497: {
498: return KSPConvergedReasonView(ksp, v);
499: }
500: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
501: {
502: return KSPConvergedReasonViewFromOptions(ksp);
503: }
505: #define KSP_FILE_CLASSID 1211223
507: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
508: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
509: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
510: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
511: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
512: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
513: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
515: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
516: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
517: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
519: /*E
520: KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
521: test routines.
523: Values:
524: + `KSP_NORM_DEFAULT` - use the default for the current `KSPType`
525: . `KSP_NORM_NONE` - use no norm calculation
526: . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm
527: . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
528: - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator)
530: Level: advanced
532: Note:
533: Each solver only supports a subset of these and some may support different ones
534: depending on left or right preconditioning, see `KSPSetPCSide()`
536: Developer Note:
537: This must match the values in petsc/finclude/petscksp.h
539: .seealso: [](chapter_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
540: `KSPSetConvergenceTest()`, `KSPSetPCSide()`
541: E*/
542: typedef enum {
543: KSP_NORM_DEFAULT = -1,
544: KSP_NORM_NONE = 0,
545: KSP_NORM_PRECONDITIONED = 1,
546: KSP_NORM_UNPRECONDITIONED = 2,
547: KSP_NORM_NATURAL = 3
548: } KSPNormType;
549: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
550: PETSC_EXTERN const char *const *const KSPNormTypes;
552: /*MC
553: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
554: possibly save some computation but means the convergence test cannot
555: be based on a norm of a residual etc.
557: Level: advanced
559: Note:
560: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
562: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
563: M*/
565: /*MC
566: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
567: convergence test routine.
569: Level: advanced
571: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
572: M*/
574: /*MC
575: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
576: convergence test routine.
578: Level: advanced
580: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
581: M*/
583: /*MC
584: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
585: convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
587: Level: advanced
589: .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
590: M*/
592: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
593: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
594: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
595: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
596: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
598: #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
599: /*E
600: KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
602: Level: beginner
604: Values:
605: + `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations
606: . `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations
607: . `KSP_CONVERGED_RTOL` - requested decrease in the residual
608: . `KSP_CONVERGED_ATOL` - requested absolute value in the residual
609: . `KSP_CONVERGED_ITS` - requested number of iterations
610: . `KSP_CONVERGED_CG_NEG_CURVE` - see note below
611: . `KSP_CONVERGED_CG_CONSTRAINED` - see note below
612: . `KSP_CONVERGED_STEP_LENGTH` - see note below
613: . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred.
614: . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES
615: . `KSP_DIVERGED_ITS` - requested number of iterations
616: . `KSP_DIVERGED_DTOL` - large increase in the residual norm
617: . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method
618: . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method
619: . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
620: . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite
621: . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation
622: . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite
623: - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason
625: Note:
626: The values `KSP_CONVERGED_CG_NEG_CURVE`, `KSP_CONVERGED_CG_CONSTRAINED`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`,
627: `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
629: Developer Notes:
630: This must match the values in petsc/finclude/petscksp.h
632: The string versions of these are `KSPConvergedReasons`; if you change
633: any of the values here also change them that array of names.
635: .seealso: [](chapter_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
636: E*/
637: typedef enum { /* converged */
638: KSP_CONVERGED_RTOL_NORMAL = 1,
639: KSP_CONVERGED_ATOL_NORMAL = 9,
640: KSP_CONVERGED_RTOL = 2,
641: KSP_CONVERGED_ATOL = 3,
642: KSP_CONVERGED_ITS = 4,
643: KSP_CONVERGED_CG_NEG_CURVE = 5,
644: KSP_CONVERGED_CG_CONSTRAINED = 6,
645: KSP_CONVERGED_STEP_LENGTH = 7,
646: KSP_CONVERGED_HAPPY_BREAKDOWN = 8,
647: /* diverged */
648: KSP_DIVERGED_NULL = -2,
649: KSP_DIVERGED_ITS = -3,
650: KSP_DIVERGED_DTOL = -4,
651: KSP_DIVERGED_BREAKDOWN = -5,
652: KSP_DIVERGED_BREAKDOWN_BICG = -6,
653: KSP_DIVERGED_NONSYMMETRIC = -7,
654: KSP_DIVERGED_INDEFINITE_PC = -8,
655: KSP_DIVERGED_NANORINF = -9,
656: KSP_DIVERGED_INDEFINITE_MAT = -10,
657: KSP_DIVERGED_PC_FAILED = -11,
658: KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
660: KSP_CONVERGED_ITERATING = 0
661: } KSPConvergedReason;
662: PETSC_EXTERN const char *const *KSPConvergedReasons;
664: /*MC
665: KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called
667: Level: beginner
669: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
670: for left preconditioning it is the 2-norm of the preconditioned residual, and the
671: 2-norm of the residual for right preconditioning
673: See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
675: .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
677: M*/
679: /*MC
680: KSP_CONVERGED_ATOL - norm(r) <= atol
682: Level: beginner
684: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
685: for left preconditioning it is the 2-norm of the preconditioned residual, and the
686: 2-norm of the residual for right preconditioning
688: See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
690: .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
692: M*/
694: /*MC
695: KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
697: Level: beginner
699: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
700: for left preconditioning it is the 2-norm of the preconditioned residual, and the
701: 2-norm of the residual for right preconditioning
703: Level: beginner
705: .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
707: M*/
709: /*MC
710: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
711: reached
713: Level: beginner
715: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
717: M*/
719: /*MC
720: KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
721: the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
722: test routine is set in KSP.
724: Level: beginner
726: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
728: M*/
730: /*MC
731: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
732: method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
733: preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
734: are colinear.
736: Level: beginner
738: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
740: M*/
742: /*MC
743: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
744: method could not continue to enlarge the Krylov space.
746: Level: beginner
748: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
750: M*/
752: /*MC
753: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
754: symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
756: Level: beginner
758: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
760: M*/
762: /*MC
763: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
764: positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
765: be positive definite
767: Level: beginner
769: Note:
770: This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
771: the `PCICC` preconditioner to generate a positive definite preconditioner
773: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
775: M*/
777: /*MC
778: KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
779: zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
780: such as `PCFIELDSPLIT`.
782: Level: beginner
784: Note:
785: Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
787: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
789: M*/
791: /*MC
792: KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
793: while `KSPSolve()` is still running.
795: Level: beginner
797: .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
799: M*/
801: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
802: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
803: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
804: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
805: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
806: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
807: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
808: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
809: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
810: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
811: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
812: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
813: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
814: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
815: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
817: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void)
818: { /* never called */
819: }
820: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
821: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void)
822: { /* never called */
823: }
824: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
825: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void)
826: { /* never called */
827: }
828: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
829: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void)
830: { /* never called */
831: }
832: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
833: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void)
834: { /* never called */
835: }
836: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
837: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void)
838: { /* never called */
839: }
840: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
842: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
843: PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
844: {
845: return KSPComputeOperator(A, NULL, B);
846: }
848: /*E
849: KSPCGType - Determines what type of `KSPCG` to use
851: Level: beginner
853: Values:
854: + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
855: - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
857: .seealso: [](chapter_ksp), `KSP`, `KSPCGSetType()`
858: E*/
859: typedef enum {
860: KSP_CG_SYMMETRIC = 0,
861: KSP_CG_HERMITIAN = 1
862: } KSPCGType;
863: PETSC_EXTERN const char *const KSPCGTypes[];
865: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
866: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
868: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
869: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
870: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
872: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
873: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
874: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
875: {
876: return KSPGLTRGetMinEig(ksp, x);
877: }
878: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
879: {
880: return KSPGLTRGetLambda(ksp, x);
881: }
883: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
884: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
886: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
887: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
888: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
890: #include <petscdrawtypes.h>
891: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
893: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
894: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
896: /*S
897: KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
899: Level: intermediate
901: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
902: S*/
903: typedef struct _p_KSPGuess *KSPGuess;
904: /*J
905: KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
907: Level: intermediate
909: Values:
910: + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
911: - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition
913: .seealso: [](chapter_ksp), `KSP`, `KSPGuess`
914: J*/
915: typedef const char *KSPGuessType;
916: #define KSPGUESSFISCHER "fischer"
917: #define KSPGUESSPOD "pod"
919: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
920: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
921: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
922: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
923: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
924: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
925: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
926: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
927: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
928: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
929: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
930: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
931: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
932: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
933: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
934: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
935: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
937: /*E
938: MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
940: Level: intermediate
942: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
943: `MatCreateSchurComplementPmat()`
944: E*/
945: typedef enum {
946: MAT_SCHUR_COMPLEMENT_AINV_DIAG,
947: MAT_SCHUR_COMPLEMENT_AINV_LUMP,
948: MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
949: MAT_SCHUR_COMPLEMENT_AINV_FULL
950: } MatSchurComplementAinvType;
951: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
953: typedef enum {
954: MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0,
955: MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1,
956: MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
957: MAT_LMVM_SYMBROYDEN_SCALE_USER = 3
958: } MatLMVMSymBroydenScaleType;
959: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
961: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
962: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
963: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
964: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
965: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
966: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
967: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
968: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
969: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
970: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
971: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
972: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
974: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
975: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
976: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
977: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
978: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
979: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
980: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
981: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
983: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
984: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
985: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
986: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
987: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
988: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
989: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
990: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
991: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
992: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
993: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
994: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
995: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
996: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
997: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
998: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
999: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1000: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1001: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1002: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1003: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1005: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1006: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1007: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1008: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1009: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1010: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
1011: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1012: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
1013: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
1014: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
1015: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1016: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1017: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
1018: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
1020: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1021: 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);
1023: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1024: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1025: #endif