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