Actual source code: petscksp.h
1: /*
2: Defines the interface functions for the Krylov subspace accelerators.
3: */
4: #pragma once
6: #include <petscpc.h>
8: /* SUBMANSEC = KSP */
10: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
11: PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void);
13: /*S
14: KSP - Abstract PETSc object that manages the linear solves in PETSc (even those such as direct factorization-based solvers that
15: do not use Krylov accelerators).
17: Level: beginner
19: Notes:
20: When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a
21: `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver.
23: Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object
25: The `PC` object is used to control preconditioners in PETSc.
27: `KSP` can also be used to solve some least squares problems (over or under-determined linear systems), using, for example, `KSPLSQR`, see `PETSCREGRESSORLINEAR`
28: for additional methods that can be used to solve least squares problems and other linear regressions).
30: .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
31: S*/
32: typedef struct _p_KSP *KSP;
34: /*J
35: KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides.
37: Level: beginner
39: .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
40: J*/
41: typedef const char *KSPType;
42: #define KSPRICHARDSON "richardson"
43: #define KSPCHEBYSHEV "chebyshev"
44: #define KSPCG "cg"
45: #define KSPGROPPCG "groppcg"
46: #define KSPPIPECG "pipecg"
47: #define KSPPIPECGRR "pipecgrr"
48: #define KSPPIPELCG "pipelcg"
49: #define KSPPIPEPRCG "pipeprcg"
50: #define KSPPIPECG2 "pipecg2"
51: #define KSPCGNE "cgne"
52: #define KSPNASH "nash"
53: #define KSPSTCG "stcg"
54: #define KSPGLTR "gltr"
55: #define KSPCGNASH PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash"
56: #define KSPCGSTCG PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg"
57: #define KSPCGGLTR PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr"
58: #define KSPFCG "fcg"
59: #define KSPPIPEFCG "pipefcg"
60: #define KSPGMRES "gmres"
61: #define KSPPIPEFGMRES "pipefgmres"
62: #define KSPFGMRES "fgmres"
63: #define KSPLGMRES "lgmres"
64: #define KSPDGMRES "dgmres"
65: #define KSPPGMRES "pgmres"
66: #define KSPTCQMR "tcqmr"
67: #define KSPBCGS "bcgs"
68: #define KSPIBCGS "ibcgs"
69: #define KSPQMRCGS "qmrcgs"
70: #define KSPFBCGS "fbcgs"
71: #define KSPFBCGSR "fbcgsr"
72: #define KSPBCGSL "bcgsl"
73: #define KSPPIPEBCGS "pipebcgs"
74: #define KSPCGS "cgs"
75: #define KSPTFQMR "tfqmr"
76: #define KSPCR "cr"
77: #define KSPPIPECR "pipecr"
78: #define KSPLSQR "lsqr"
79: #define KSPPREONLY "preonly"
80: #define KSPNONE "none"
81: #define KSPQCG "qcg"
82: #define KSPBICG "bicg"
83: #define KSPMINRES "minres"
84: #define KSPSYMMLQ "symmlq"
85: #define KSPLCD "lcd"
86: #define KSPPYTHON "python"
87: #define KSPGCR "gcr"
88: #define KSPPIPEGCR "pipegcr"
89: #define KSPTSIRM "tsirm"
90: #define KSPCGLS "cgls"
91: #define KSPFETIDP "fetidp"
92: #define KSPHPDDM "hpddm"
94: /* Logging support */
95: PETSC_EXTERN PetscClassId KSP_CLASSID;
96: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
97: PETSC_EXTERN PetscClassId DMKSP_CLASSID;
99: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
100: PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
101: PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
102: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
103: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
104: PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
105: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
106: PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
107: PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
108: PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat);
109: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
110: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
111: {
112: return KSPSetMatSolveBatchSize(ksp, n);
113: }
114: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
115: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
116: {
117: return KSPGetMatSolveBatchSize(ksp, n);
118: }
119: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
120: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
121: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
122: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
123: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
124: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
125: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);
127: PETSC_EXTERN PetscFunctionList KSPList;
128: PETSC_EXTERN PetscFunctionList KSPGuessList;
129: PETSC_EXTERN PetscFunctionList KSPMonitorList;
130: PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
131: PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
132: PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP));
133: PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
135: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
136: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
137: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
138: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
139: PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt);
140: PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *);
141: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
142: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
143: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
144: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
145: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
146: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
147: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
148: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
149: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
150: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
151: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
152: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
153: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
154: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
155: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
156: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
157: {
158: return KSPCreateVecs(ksp, n, x, m, y);
159: }
161: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
162: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
164: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
165: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
166: PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt);
167: PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *);
169: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
170: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscCtxDestroyFn *);
171: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
172: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
173: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
174: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool);
175: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
176: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool);
178: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
179: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
180: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
181: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);
183: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
184: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
185: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
186: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
187: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
188: PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]);
189: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
190: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
191: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
192: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
193: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
194: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
195: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
196: PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
197: /*
198: PCMGCoarseList contains the list of coarse space constructor currently registered
199: These are added with PCMGRegisterCoarseSpaceConstructor()
200: */
201: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
202: PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
203: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
205: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
206: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);
208: /*E
209: KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV`
211: Values:
212: + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial
213: . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial
214: - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial
216: Level: intermediate
218: .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()`
219: E*/
220: typedef enum {
221: KSP_CHEBYSHEV_FIRST,
222: KSP_CHEBYSHEV_FOURTH,
223: KSP_CHEBYSHEV_OPT_FOURTH
224: } KSPChebyshevKind;
226: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
227: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
228: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
229: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
230: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
231: PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind);
232: PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *);
233: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
234: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
235: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
236: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
237: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);
239: /*E
241: KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
243: Values:
244: + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
245: - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
247: Level: intermediate
249: .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
250: E*/
251: typedef enum {
252: KSP_FCD_TRUNC_TYPE_STANDARD,
253: KSP_FCD_TRUNC_TYPE_NOTAY
254: } KSPFCDTruncationType;
255: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
257: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
258: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
259: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
260: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
261: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
262: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);
264: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
265: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
266: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
267: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
268: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
269: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);
271: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
272: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
273: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
274: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
275: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
276: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
277: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
278: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
279: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
281: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
282: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
283: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
284: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);
286: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
287: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
288: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
289: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
290: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);
292: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
293: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
295: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);
297: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
298: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
299: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
301: PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
302: PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
303: PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);
305: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
306: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
307: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
308: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);
310: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
311: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
312: #if PetscDefined(HAVE_HPDDM)
313: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
314: {
315: return KSPHPDDMSetDeflationMat(ksp, U);
316: }
317: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
318: {
319: return KSPHPDDMGetDeflationMat(ksp, U);
320: }
321: #endif
322: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
323: {
324: return KSPMatSolve(ksp, B, X);
325: }
326: /*E
327: KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
329: Values:
330: + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
331: . `KSP_HPDDM_TYPE_BGMRES` - block GMRES
332: . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient
333: . `KSP_HPDDM_TYPE_BCG` - block CG
334: . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
335: . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR
336: . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG
337: - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only
339: Level: intermediate
341: .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
342: E*/
343: typedef enum {
344: KSP_HPDDM_TYPE_GMRES = 0,
345: KSP_HPDDM_TYPE_BGMRES = 1,
346: KSP_HPDDM_TYPE_CG = 2,
347: KSP_HPDDM_TYPE_BCG = 3,
348: KSP_HPDDM_TYPE_GCRODR = 4,
349: KSP_HPDDM_TYPE_BGCRODR = 5,
350: KSP_HPDDM_TYPE_BFBCG = 6,
351: KSP_HPDDM_TYPE_PREONLY = 7
352: } KSPHPDDMType;
353: PETSC_EXTERN const char *const KSPHPDDMTypes[];
355: /*E
356: KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
358: Values:
359: + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16`
360: . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single`
361: . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double`
362: - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`
364: Level: intermediate
366: .seealso: [](ch_ksp), `KSP`, `KSPHPDDM`
367: E*/
368: typedef enum {
369: KSP_HPDDM_PRECISION_HALF = 0,
370: KSP_HPDDM_PRECISION_SINGLE = 1,
371: KSP_HPDDM_PRECISION_DOUBLE = 2,
372: KSP_HPDDM_PRECISION_QUADRUPLE = 3
373: } KSPHPDDMPrecision;
374: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
375: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);
377: /*E
378: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers
380: Values:
381: + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt
382: . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
383: - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps
385: Level: advanced
387: .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
388: `KSPGMRESGetOrthogonalization()`,
389: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
390: E*/
391: typedef enum {
392: KSP_GMRES_CGS_REFINE_NEVER,
393: KSP_GMRES_CGS_REFINE_IFNEEDED,
394: KSP_GMRES_CGS_REFINE_ALWAYS
395: } KSPGMRESCGSRefinementType;
396: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
398: /*MC
399: KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
401: Level: advanced
403: Note:
404: Possibly unstable, but the fastest to compute
406: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
407: `KSP`, `KSPGMRESGetOrthogonalization()`,
408: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
409: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
410: M*/
412: /*MC
413: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
414: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
415: poor orthogonality.
417: Level: advanced
419: Note:
420: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
421: estimate the orthogonality but is more stable.
423: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
424: `KSP`, `KSPGMRESGetOrthogonalization()`,
425: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
426: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
427: M*/
429: /*MC
430: KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process.
432: Level: advanced
434: Notes:
435: This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
436: but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
438: You should only use this if you absolutely know that the iterative refinement is needed.
440: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
441: `KSP`, `KSPGMRESGetOrthogonalization()`,
442: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
443: `KSPGMRESModifiedGramSchmidtOrthogonalization()`
444: M*/
446: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
447: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);
449: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
450: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
451: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));
453: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
454: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
455: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);
457: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
458: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
459: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
460: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);
462: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
463: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
465: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
466: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
467: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
468: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
469: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
470: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
471: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
472: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
473: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
474: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
475: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
476: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
477: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
478: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
479: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
480: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
481: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
482: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
483: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
484: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
485: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
486: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
487: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
488: {
489: return KSPMonitorResidual(ksp, n, rnorm, vf);
490: }
491: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
492: {
493: return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
494: }
495: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
496: {
497: return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
498: }
500: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
501: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
502: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
503: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
504: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
505: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
506: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
507: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);
509: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
510: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);
512: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
513: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
514: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
515: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
516: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
517: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);
519: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
520: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
521: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
522: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);
524: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
525: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
526: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
527: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
528: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *, PetscCtxDestroyFn *);
529: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
530: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
531: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);
533: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
534: {
535: return KSPConvergedReasonView(ksp, v);
536: }
537: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
538: {
539: return KSPConvergedReasonViewFromOptions(ksp);
540: }
542: #define KSP_FILE_CLASSID 1211223
544: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
545: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
546: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
547: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
548: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
549: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
550: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
552: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
553: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
554: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
555: PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);
557: /*E
558: KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
559: test routines.
561: Values:
562: + `KSP_NORM_DEFAULT` - use the default for the current `KSPType`
563: . `KSP_NORM_NONE` - use no norm calculation
564: . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm
565: . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
566: - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator)
568: Level: advanced
570: Note:
571: Each solver only supports a subset of these and some may support different ones
572: depending on whether left or right preconditioning is used, see `KSPSetPCSide()`
574: .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
575: `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
576: E*/
577: typedef enum {
578: KSP_NORM_DEFAULT = -1,
579: KSP_NORM_NONE = 0,
580: KSP_NORM_PRECONDITIONED = 1,
581: KSP_NORM_UNPRECONDITIONED = 2,
582: KSP_NORM_NATURAL = 3
583: } KSPNormType;
584: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
585: PETSC_EXTERN const char *const *const KSPNormTypes;
587: /*MC
588: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
589: possibly save some computation but means the convergence test cannot
590: be based on a norm of a residual etc.
592: Level: advanced
594: Note:
595: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
597: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
598: M*/
600: /*MC
601: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
602: convergence test routine.
604: Level: advanced
606: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
607: M*/
609: /*MC
610: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
611: convergence test routine.
613: Level: advanced
615: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
616: M*/
618: /*MC
619: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
620: convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
622: Level: advanced
624: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
625: M*/
627: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
628: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
629: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt);
630: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
631: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);
633: #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
634: #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
635: #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", )
636: #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", )
637: /*E
638: KSPConvergedReason - reason a Krylov method was determined to have converged or diverged
640: Values:
641: + `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR`
642: . `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR`
643: . `KSP_CONVERGED_RTOL` - requested decrease in the residual
644: . `KSP_CONVERGED_ATOL` - requested absolute value in the residual
645: . `KSP_CONVERGED_ITS` - requested number of iterations
646: . `KSP_CONVERGED_NEG_CURVE` - see note below
647: . `KSP_CONVERGED_STEP_LENGTH` - see note below
648: . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred).
649: . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within `KSPGMRES`
650: . `KSP_DIVERGED_ITS` - requested number of iterations
651: . `KSP_DIVERGED_DTOL` - large increase in the residual norm indicating the solution is diverging
652: . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method
653: . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBCGS` Krylov method
654: . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
655: . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
656: . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation
657: . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG`
658: - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason
660: Level: beginner
662: Note:
663: The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
664: the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.
666: Developer Note:
667: The string versions of these are `KSPConvergedReasons`; if you change
668: any of the values here also change them that array of names.
670: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
671: E*/
672: typedef enum { /* converged */
673: KSP_CONVERGED_RTOL_NORMAL_DEPRECATED = 1,
674: KSP_CONVERGED_RTOL_NORMAL_EQUATIONS = 1,
675: KSP_CONVERGED_ATOL_NORMAL_DEPRECATED = 9,
676: KSP_CONVERGED_ATOL_NORMAL_EQUATIONS = 9,
677: KSP_CONVERGED_RTOL = 2,
678: KSP_CONVERGED_ATOL = 3,
679: KSP_CONVERGED_ITS = 4,
680: KSP_CONVERGED_NEG_CURVE = 5,
681: KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5,
682: KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
683: KSP_CONVERGED_STEP_LENGTH = 6,
684: KSP_CONVERGED_HAPPY_BREAKDOWN = 7,
685: /* diverged */
686: KSP_DIVERGED_NULL = -2,
687: KSP_DIVERGED_ITS = -3,
688: KSP_DIVERGED_DTOL = -4,
689: KSP_DIVERGED_BREAKDOWN = -5,
690: KSP_DIVERGED_BREAKDOWN_BICG = -6,
691: KSP_DIVERGED_NONSYMMETRIC = -7,
692: KSP_DIVERGED_INDEFINITE_PC = -8,
693: KSP_DIVERGED_NANORINF = -9,
694: KSP_DIVERGED_INDEFINITE_MAT = -10,
695: KSP_DIVERGED_PC_FAILED = -11,
696: KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
698: KSP_CONVERGED_ITERATING = 0
699: } KSPConvergedReason;
700: PETSC_EXTERN const char *const *KSPConvergedReasons;
702: /*MC
703: KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called
705: Level: beginner
707: Notes:
708: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
709: for left preconditioning it is the 2-norm of the preconditioned residual, and the
710: 2-norm of the residual for right preconditioning
712: See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
714: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
715: M*/
717: /*MC
718: KSP_CONVERGED_ATOL - $||r|| \le atol$
720: Level: beginner
722: Notes:
723: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
724: for left preconditioning it is the 2-norm of the preconditioned residual, and the
725: 2-norm of the residual for right preconditioning
727: See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
729: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
730: M*/
732: /*MC
733: KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$
735: Level: beginner
737: Note:
738: See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
739: for left preconditioning it is the 2-norm of the preconditioned residual, and the
740: 2-norm of the residual for right preconditioning
742: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
743: M*/
745: /*MC
746: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
747: reached
749: Level: beginner
751: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
752: M*/
754: /*MC
755: KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
756: the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
757: test routine is set in `KSP`.
759: Level: beginner
761: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
762: M*/
764: /*MC
765: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
766: method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
767: preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
768: are collinear.
770: Level: beginner
772: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
773: M*/
775: /*MC
776: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
777: method could not continue to enlarge the Krylov space.
779: Level: beginner
781: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
782: M*/
784: /*MC
785: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
786: symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
788: Level: beginner
790: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
791: M*/
793: /*MC
794: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
795: positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
796: be symmetric positive definite (SPD).
798: Level: beginner
800: Note:
801: This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
802: the `PCICC` preconditioner to generate a positive definite preconditioner
804: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
805: M*/
807: /*MC
808: KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
809: zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
810: such as `PCFIELDSPLIT`.
812: Level: beginner
814: Note:
815: Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.
817: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
818: M*/
820: /*MC
821: KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
822: while `KSPSolve()` is still running.
824: Level: beginner
826: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
827: M*/
829: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
830: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
831: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
832: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
833: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
834: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
835: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
836: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
837: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
838: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
839: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
840: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
841: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
842: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
843: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
844: PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
845: PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);
847: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
848: { /* never called */
849: }
850: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
851: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
852: { /* never called */
853: }
854: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
855: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
856: { /* never called */
857: }
858: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
859: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
860: { /* never called */
861: }
862: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
863: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
864: { /* never called */
865: }
866: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
867: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
868: { /* never called */
869: }
870: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
872: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
873: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
874: {
875: return KSPComputeOperator(A, PETSC_NULLPTR, B);
876: }
878: /*E
879: KSPCGType - Determines what type of `KSPCG` to use
881: Values:
882: + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
883: - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian
885: Level: beginner
887: .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
888: E*/
889: typedef enum {
890: KSP_CG_SYMMETRIC = 0,
891: KSP_CG_HERMITIAN = 1
892: } KSPCGType;
893: PETSC_EXTERN const char *const KSPCGTypes[];
895: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
896: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);
898: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
899: PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
900: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
901: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);
903: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
904: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
905: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
906: {
907: return KSPGLTRGetMinEig(ksp, x);
908: }
909: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
910: {
911: return KSPGLTRGetLambda(ksp, x);
912: }
914: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
915: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);
917: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
918: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
919: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);
921: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);
923: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
924: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
926: /*S
927: KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.
929: Level: intermediate
931: Note:
932: These methods generate initial guesses based on a series of previous, related, linear solves. For example,
933: in implicit time-stepping with `TS`.
935: .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
936: S*/
937: typedef struct _p_KSPGuess *KSPGuess;
939: /*J
940: KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
942: Values:
943: + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
944: - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition (POD)
946: Level: intermediate
948: .seealso: [](ch_ksp), `KSP`, `KSPGuess`
949: J*/
950: typedef const char *KSPGuessType;
951: #define KSPGUESSFISCHER "fischer"
952: #define KSPGUESSPOD "pod"
954: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
955: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
956: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
957: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
958: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
959: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
960: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
961: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
962: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
963: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
964: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
965: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
966: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
967: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
968: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
969: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
970: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);
972: /*E
973: MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
975: Level: intermediate
977: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
978: `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
979: E*/
980: typedef enum {
981: MAT_SCHUR_COMPLEMENT_AINV_DIAG,
982: MAT_SCHUR_COMPLEMENT_AINV_LUMP,
983: MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
984: MAT_SCHUR_COMPLEMENT_AINV_FULL
985: } MatSchurComplementAinvType;
986: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
988: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
989: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
990: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
991: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
992: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
993: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
994: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
995: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
996: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
997: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
998: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
999: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);
1001: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1002: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1003: PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1004: PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1005: PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
1006: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1007: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1008: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1009: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1010: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1011: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1013: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1014: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1015: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1016: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1017: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1018: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1019: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1020: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1021: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1022: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1023: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1024: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1025: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1026: PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *);
1027: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1028: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1029: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1030: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1031: PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1032: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1033: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1034: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
1036: /*E
1037: MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix
1039: Values:
1040: + `MAT_LMVM_MULT_RECURSIVE` - Use recursive formulas for products and solves
1041: . `MAT_LMVM_MULT_DENSE` - Use dense formulas for products and solves when possible
1042: - `MAT_LMVM_MULT_COMPACT_DENSE` - The same as `MATLMVM_MULT_DENSE`, but go further and ensure products and solves are computed in compact low-rank update form
1044: Level: advanced
1046: Options Database Keys:
1047: . -mat_lmvm_mult_algorithm - the algorithm to use for multiplication (recursive, dense, compact_dense)
1049: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()`
1050: E*/
1051: typedef enum {
1052: MAT_LMVM_MULT_RECURSIVE,
1053: MAT_LMVM_MULT_DENSE,
1054: MAT_LMVM_MULT_COMPACT_DENSE,
1055: } MatLMVMMultAlgorithm;
1057: PETSC_EXTERN const char *const MatLMVMMultAlgorithms[];
1059: PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm);
1060: PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *);
1062: /*E
1063: MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix.
1065: Values:
1066: + `MAT_LMVM_SYMBROYDEN_SCALE_NONE` - no rescaling
1067: . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR` - scalar rescaling
1068: . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling
1069: . `MAT_LMVM_SYMBROYDEN_SCALE_USER` - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE`
1070: - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE` - let PETSc decide rescaling
1072: Level: intermediate
1074: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1075: E*/
1076: typedef enum {
1077: MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0,
1078: MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1,
1079: MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1080: MAT_LMVM_SYMBROYDEN_SCALE_USER = 3,
1081: MAT_LMVM_SYMBROYDEN_SCALE_DECIDE = 4
1082: } MatLMVMSymBroydenScaleType;
1083: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
1085: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
1086: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *);
1087: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal);
1088: PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *);
1089: PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal);
1091: /*E
1092: MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.
1094: Values:
1095: + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1096: - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement
1098: Level: intermediate
1100: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1101: E*/
1102: typedef enum {
1103: MAT_LMVM_DENSE_REORDER,
1104: MAT_LMVM_DENSE_INPLACE
1105: } MatLMVMDenseType;
1106: PETSC_EXTERN const char *const MatLMVMDenseTypes[];
1108: PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);
1110: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1111: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1112: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1113: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1114: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
1116: /*S
1117: KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`
1119: Calling Sequence:
1120: + ksp - `ksp` context
1121: . b - output vector
1122: - ctx - [optional] user-defined function context
1124: Level: beginner
1126: .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
1127: S*/
1128: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeRHSFn)(KSP ksp, Vec b, void *ctx);
1130: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);
1132: /*S
1133: KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`
1135: Calling Sequence:
1136: + ksp - `KSP` context
1137: . A - the operator that defines the linear system
1138: . P - an operator from which to build the preconditioner (often the same as `A`)
1139: - ctx - [optional] user-defined function context
1141: Level: beginner
1143: .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
1144: S*/
1145: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeOperatorsFn)(KSP ksp, Mat A, Mat P, void *ctx);
1147: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);
1149: /*S
1150: KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`
1152: Calling Sequence:
1153: + ksp - `ksp` context
1154: . x - output vector
1155: - ctx - [optional] user-defined function context
1157: Level: beginner
1159: .seealso: [](ch_snes), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
1160: S*/
1161: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeInitialGuessFn)(KSP ksp, Vec x, void *ctx);
1163: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1164: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1165: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1166: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1167: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1168: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1169: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);
1171: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1172: 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);
1173: PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode);
1175: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1176: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
1178: PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1179: PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);
1181: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);