xref: /petsc/include/petscksp.h (revision dfbbaf821b4c49d07b4ce746493b0d955783fdf9)
1 /*
2    Defines the interface functions for the Krylov subspace accelerators.
3 */
4 #ifndef PETSCKSP_H
5 #define PETSCKSP_H
6 
7 #include <petscpc.h>
8 
9 /* SUBMANSEC = KSP */
10 
11 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
12 
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).
16 
17    Level: beginner
18 
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.
22 
23 .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
24 S*/
25 typedef struct _p_KSP*     KSP;
26 
27 /*J
28     KSPType - String with the name of a PETSc Krylov method.
29 
30    Level: beginner
31 
32 .seealso: `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"
86 
87 /* Logging support */
88 PETSC_EXTERN PetscClassId KSP_CLASSID;
89 PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
90 PETSC_EXTERN PetscClassId DMKSP_CLASSID;
91 
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) {return KSPSetMatSolveBatchSize(ksp,n);}
103 PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP,PetscInt*);
104 PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp,PetscInt *n) {return KSPGetMatSolveBatchSize(ksp,n);}
105 PETSC_EXTERN PetscErrorCode KSPReset(KSP);
106 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
107 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
108 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
109 PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
110 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
111 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
112 
113 PETSC_EXTERN PetscFunctionList KSPList;
114 PETSC_EXTERN PetscFunctionList KSPGuessList;
115 PETSC_EXTERN PetscFunctionList KSPMonitorList;
116 PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
117 PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
118 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
119 PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
120 
121 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
122 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
123 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
124 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
125 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
126 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
127 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
128 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
129 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
130 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
131 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
132 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
133 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
134 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
135 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
136 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
137 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
138 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
139 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
140 PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
141 
142 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
143 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
144 
145 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
146 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
147 
148 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
149 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
150 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
151 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void*);
152 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*);
153 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);
154 PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*);
155 PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool);
156 
157 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
158 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
159 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
160 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
161 
162 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
163 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
164 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
165 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
166 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
167 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
168 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
169 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
170 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
171 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
172 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
173 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
174 PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
175 /*
176   PCMGCoarseList contains the list of coarse space constructor currently registered
177   These are added with PCMGRegisterCoarseSpaceConstructor()
178 */
179 PETSC_EXTERN PetscFunctionList PCMGCoarseList;
180 PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*));
181 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*));
182 
183 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
184 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
185 
186 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
187 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
188 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
189 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
190 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
191 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
192 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
193 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
194 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
195 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
196 
197 /*E
198 
199   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
200 
201   Values:
202 $ `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
203 $ `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
204 
205    Level: intermediate
206 
207 .seealso `:` `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
208 E*/
209 typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
210 PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
211 
212 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
213 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
214 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
215 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
216 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
217 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
218 
219 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
220 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
221 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
222 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
223 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
224 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
225 
226 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
227 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
228 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
229 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
230 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
231 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
232 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
233 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
234 
235 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
236 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
237 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
238 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal);
239 
240 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
241 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
242 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
243 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
244 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
245 
246 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
247 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
248 
249 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
250 
251 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
252 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
253 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
254 
255 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
256 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
257 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
258 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
259 
260 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP,Mat);
261 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP,Mat*);
262 #if PetscDefined(HAVE_HPDDM)
263 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) { return KSPHPDDMSetDeflationMat(ksp, U); }
264 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) { return KSPHPDDMGetDeflationMat(ksp, U); }
265 #endif
266 PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
267 /*E
268     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`
269 
270     Level: intermediate
271 
272     Values:
273 $   `KSP_HPDDM_TYPE_GMRES` (default)
274 $   `KSP_HPDDM_TYPE_BGMRES`
275 $   `KSP_HPDDM_TYPE_CG`
276 $   `KSP_HPDDM_TYPE_BCG`
277 $   `KSP_HPDDM_TYPE_GCRODR`
278 $   `KSP_HPDDM_TYPE_BGCRODR`
279 $   `KSP_HPDDM_TYPE_BFBCG`
280 $   `KSP_HPDDM_TYPE_PREONLY`
281 
282 .seealso: `KSPHPDDM`, `KSPHPDDMSetType()`
283 E*/
284 typedef enum {
285   KSP_HPDDM_TYPE_GMRES = 0,
286   KSP_HPDDM_TYPE_BGMRES = 1,
287   KSP_HPDDM_TYPE_CG = 2,
288   KSP_HPDDM_TYPE_BCG = 3,
289   KSP_HPDDM_TYPE_GCRODR = 4,
290   KSP_HPDDM_TYPE_BGCRODR = 5,
291   KSP_HPDDM_TYPE_BFBCG = 6,
292   KSP_HPDDM_TYPE_PREONLY = 7
293 } KSPHPDDMType;
294 PETSC_EXTERN const char *const KSPHPDDMTypes[];
295 /*E
296     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`
297 
298     Level: intermediate
299 
300     Values:
301 $   `KSP_HPDDM_PRECISION_HALF` (currently unsupported)
302 $   `KSP_HPDDM_PRECISION_SINGLE` (default when PETSc is configured --with-precision=single)
303 $   `KSP_HPDDM_PRECISION_DOUBLE` (default when PETSc is configured --with-precision=double)
304 $   `KSP_HPDDM_PRECISION_QUADRUPLE` (default when PETSc is configured --with-precision=__float128)
305 
306 .seealso: `KSPHPDDM`
307 E*/
308 typedef enum {
309   KSP_HPDDM_PRECISION_HALF = 0,
310   KSP_HPDDM_PRECISION_SINGLE = 1,
311   KSP_HPDDM_PRECISION_DOUBLE = 2,
312   KSP_HPDDM_PRECISION_QUADRUPLE = 3
313 } KSPHPDDMPrecision;
314 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
315 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
316 
317 /*E
318     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
319 
320    Level: advanced
321 
322 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
323           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
324 
325 E*/
326 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
327 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
328 /*MC
329     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
330 
331    Level: advanced
332 
333    Note:
334    Possibly unstable, but the fastest to compute
335 
336 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
337           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
338           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
339 M*/
340 
341 /*MC
342     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
343           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
344           poor orthogonality.
345 
346    Level: advanced
347 
348    Note: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
349      estimate the orthogonality but is more stable.
350 
351 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
352           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
353           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
354 M*/
355 
356 /*MC
357     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
358 
359    Level: advanced
360 
361    Notes:
362    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
363    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.
364 
365    You should only use this if you absolutely know that the iterative refinement is needed.
366 
367 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
368           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
369           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
370 M*/
371 
372 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
373 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
374 
375 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
376 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
377 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
378 
379 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
380 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
381 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
382 
383 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
384 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
385 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
386 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
387 
388 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
389 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
390 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
391 
392 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*);
393 PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *);
394 PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
395 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
396 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
397 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
398 PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
399 PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
400 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
401 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
402 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
403 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
404 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
405 PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
406 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
407 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
408 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
409 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
410 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
411 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
412 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
413 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
414 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
415 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);}
416 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);}
417 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);}
418 
419 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
420 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
421 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
422 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
423 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
424 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
425 
426 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
427 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
428 
429 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
430 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
431 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
432 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
433 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
434 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
435 
436 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
437 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
438 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
439 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
440 
441 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
442 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
443 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
444 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer);
445 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP,PetscErrorCode (*)(KSP,void*),void *vctx,PetscErrorCode (*)(void**));
446 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
447 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
448 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer);
449 
450 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);}
451 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}
452 
453 #define KSP_FILE_CLASSID 1211223
454 
455 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
456 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
457 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
458 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
459 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
460 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
461 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
462 
463 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
464 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
465 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
466 
467 /*E
468     KSPNormType - Norm that is passed in the Krylov convergence
469        test routines.
470 
471    Level: advanced
472 
473    Each solver only supports a subset of these and some may support different ones
474    depending on left or right preconditioning, see `KSPSetPCSide()`
475 
476    Developer Note:
477     this must match petsc/finclude/petscksp.h
478 
479 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
480           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
481 E*/
482 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
483 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
484 PETSC_EXTERN const char *const*const KSPNormTypes;
485 
486 /*MC
487     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
488           possibly save some computation but means the convergence test cannot
489           be based on a norm of a residual etc.
490 
491    Level: advanced
492 
493     Note: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored
494 
495 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
496 M*/
497 
498 /*MC
499     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
500        convergence test routine.
501 
502    Level: advanced
503 
504 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
505 M*/
506 
507 /*MC
508     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
509        convergence test routine.
510 
511    Level: advanced
512 
513 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
514 M*/
515 
516 /*MC
517     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
518        convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`
519 
520    Level: advanced
521 
522 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
523 M*/
524 
525 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
526 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
527 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
528 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
529 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
530 
531 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
532 /*E
533     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
534 
535    Level: beginner
536 
537    Notes:
538     See `KSPGetConvergedReason()` for explanation of each value
539 
540    Developer Notes:
541    This must match petsc/finclude/petscksp.h
542 
543    The string versions of these are `KSPConvergedReasons`; if you change
544    any of the values here also change them that array of names.
545 
546 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
547 E*/
548 typedef enum {/* converged */
549               KSP_CONVERGED_RTOL_NORMAL        =  1,
550               KSP_CONVERGED_ATOL_NORMAL        =  9,
551               KSP_CONVERGED_RTOL               =  2,
552               KSP_CONVERGED_ATOL               =  3,
553               KSP_CONVERGED_ITS                =  4,
554               KSP_CONVERGED_CG_NEG_CURVE       =  5,
555               KSP_CONVERGED_CG_CONSTRAINED     =  6,
556               KSP_CONVERGED_STEP_LENGTH        =  7,
557               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
558               /* diverged */
559               KSP_DIVERGED_NULL                = -2,
560               KSP_DIVERGED_ITS                 = -3,
561               KSP_DIVERGED_DTOL                = -4,
562               KSP_DIVERGED_BREAKDOWN           = -5,
563               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
564               KSP_DIVERGED_NONSYMMETRIC        = -7,
565               KSP_DIVERGED_INDEFINITE_PC       = -8,
566               KSP_DIVERGED_NANORINF            = -9,
567               KSP_DIVERGED_INDEFINITE_MAT      = -10,
568               KSP_DIVERGED_PC_FAILED           = -11,
569               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
570 
571               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
572 PETSC_EXTERN const char *const*KSPConvergedReasons;
573 
574 /*MC
575      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called
576 
577    Level: beginner
578 
579    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
580    for left preconditioning it is the 2-norm of the preconditioned residual, and the
581    2-norm of the residual for right preconditioning
582 
583    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.
584 
585 .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
586 
587 M*/
588 
589 /*MC
590      KSP_CONVERGED_ATOL - norm(r) <= atol
591 
592    Level: beginner
593 
594    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
595        for left preconditioning it is the 2-norm of the preconditioned residual, and the
596        2-norm of the residual for right preconditioning
597 
598    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.
599 
600 .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
601 
602 M*/
603 
604 /*MC
605      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
606 
607    Level: beginner
608 
609    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
610        for left preconditioning it is the 2-norm of the preconditioned residual, and the
611        2-norm of the residual for right preconditioning
612 
613    Level: beginner
614 
615 .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
616 
617 M*/
618 
619 /*MC
620      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
621       reached
622 
623    Level: beginner
624 
625 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
626 
627 M*/
628 
629 /*MC
630      KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
631            the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
632            test routine is set in KSP.
633 
634    Level: beginner
635 
636 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
637 
638 M*/
639 
640 /*MC
641      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
642           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
643           preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
644           are colinear.
645 
646    Level: beginner
647 
648 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
649 
650 M*/
651 
652 /*MC
653      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
654           method could not continue to enlarge the Krylov space.
655 
656    Level: beginner
657 
658 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
659 
660 M*/
661 
662 /*MC
663      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
664         symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry
665 
666    Level: beginner
667 
668 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
669 
670 M*/
671 
672 /*MC
673      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
674         positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
675         be positive definite
676 
677    Level: beginner
678 
679      Note:
680     This can happen with the `PCICC` preconditioner, use -pc_factor_shift_positive_definite to force
681   the `PCICC` preconditioner to generate a positive definite preconditioner
682 
683 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
684 
685 M*/
686 
687 /*MC
688      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
689      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
690      such as `PCFIELDSPLIT`.
691 
692    Level: beginner
693 
694     Notes:
695     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
696 
697 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
698 
699 M*/
700 
701 /*MC
702      KSP_CONVERGED_ITERATING - This flag is returned if you call `KSPGetConvergedReason()`
703         while the `KSPSolve()` is still running.
704 
705    Level: beginner
706 
707 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
708 
709 M*/
710 
711 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
712 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
713 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
714 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void*);
715 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
716 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
717 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
718 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
719 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
720 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
721 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
722 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
723 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
724 PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP,const char**);
725 PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
726 
727 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) { /* never called */ }
728 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
729 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) { /* never called */ }
730 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
731 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) { /* never called */ }
732 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
733 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
734 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
735 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
736 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
737 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) { /* never called */ }
738 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
739 
740 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
741 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
742 
743 /*E
744     KSPCGType - Determines what type of CG to use
745 
746    Level: beginner
747 
748 .seealso: `KSPCGSetType()`
749 E*/
750 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
751 PETSC_EXTERN const char *const KSPCGTypes[];
752 
753 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
754 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);
755 
756 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
757 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
758 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
759 
760 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
761 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
762 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
763 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
764 
765 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
766 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP,const char*[]);
767 
768 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC,PetscErrorCode (*)(PC,KSP));
769 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
770 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
771 
772 #include <petscdrawtypes.h>
773 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
774 
775 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
776 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
777 
778 /*S
779      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
780 
781    Level: beginner
782 
783 .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
784 S*/
785 typedef struct _p_KSPGuess* KSPGuess;
786 /*J
787     KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.
788 
789    Level: beginner
790 
791 .seealso: `KSPGuess`
792 J*/
793 typedef const char* KSPGuessType;
794 #define KSPGUESSFISCHER "fischer"
795 #define KSPGUESSPOD     "pod"
796 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
797 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
798 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
799 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
800 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
801 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
802 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
803 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
804 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess,PetscReal);
805 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
806 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
807 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
808 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
809 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
810 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
811 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
812 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
813 
814 /*E
815     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
816 
817     Level: intermediate
818 
819 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()`
820 E*/
821 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, MAT_SCHUR_COMPLEMENT_AINV_FULL} MatSchurComplementAinvType;
822 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
823 
824 typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
825               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
826               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
827               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
828 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
829 
830 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
831 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
832 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
833 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
834 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
835 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
836 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
837 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
838 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
839 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
840 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
841 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
842 
843 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
844 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
845 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
846 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
847 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
848 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
849 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
850 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
851 
852 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
853 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
854 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
855 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
856 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
857 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
858 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
859 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
860 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
861 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
862 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
863 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
864 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
865 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
866 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
867 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
868 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
869 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
870 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
871 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
872 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
873 
874 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
875 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
876 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
877 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
878 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
879 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
880 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
881 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
882 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
883 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
884 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
885 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
886 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
887 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
888 
889 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
890 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
891                                            void (**)(PetscInt, PetscInt, PetscInt,
892                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
893                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
894                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
895 
896 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
897 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
898 #endif
899