xref: /petsc/include/petscksp.h (revision 6ce1633cb736e3bd2a11b0bc146401a5bd4cb96c)
1 /*
2    Defines the interface functions for the Krylov subspace accelerators.
3 */
4 #ifndef __PETSCKSP_H
5 #define __PETSCKSP_H
6 #include "petscpc.h"
7 PETSC_EXTERN_CXX_BEGIN
8 
9 extern PetscErrorCode  KSPInitializePackage(const char[]);
10 
11 /*S
12      KSP - Abstract PETSc object that manages all Krylov methods
13 
14    Level: beginner
15 
16   Concepts: Krylov methods
17 
18 .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
19 S*/
20 typedef struct _p_KSP*     KSP;
21 
22 /*E
23     KSPType - String with the name of a PETSc Krylov method or the creation function
24        with an optional dynamic library name, for example
25        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
26 
27    Level: beginner
28 
29 .seealso: KSPSetType(), KSP
30 E*/
31 #define KSPType char*
32 #define KSPRICHARDSON "richardson"
33 #define KSPCHEBYCHEV  "chebychev"
34 #define KSPCG         "cg"
35 #define   KSPCGNE       "cgne"
36 #define   KSPNASH       "nash"
37 #define   KSPSTCG       "stcg"
38 #define   KSPGLTR       "gltr"
39 #define KSPGMRES      "gmres"
40 #define   KSPFGMRES     "fgmres"
41 #define   KSPLGMRES     "lgmres"
42 #define KSPTCQMR      "tcqmr"
43 #define KSPBCGS       "bcgs"
44 #define KSPIBCGS        "ibcgs"
45 #define KSPBCGSL        "bcgsl"
46 #define KSPCGS        "cgs"
47 #define KSPTFQMR      "tfqmr"
48 #define KSPCR         "cr"
49 #define KSPLSQR       "lsqr"
50 #define KSPPREONLY    "preonly"
51 #define KSPQCG        "qcg"
52 #define KSPBICG       "bicg"
53 #define KSPMINRES     "minres"
54 #define KSPSYMMLQ     "symmlq"
55 #define KSPLCD        "lcd"
56 #define KSPPYTHON     "python"
57 #define KSPBROYDEN    "broyden"
58 #define KSPGCR        "gcr"
59 #define KSPNGMRES     "ngmres"
60 #define KSPSPECEST    "specest"
61 
62 /* Logging support */
63 extern PetscClassId  KSP_CLASSID;
64 
65 extern PetscErrorCode  KSPCreate(MPI_Comm,KSP *);
66 extern PetscErrorCode  KSPSetType(KSP,const KSPType);
67 extern PetscErrorCode  KSPSetUp(KSP);
68 extern PetscErrorCode  KSPSetUpOnBlocks(KSP);
69 extern PetscErrorCode  KSPSolve(KSP,Vec,Vec);
70 extern PetscErrorCode  KSPSolveTranspose(KSP,Vec,Vec);
71 extern PetscErrorCode  KSPDestroy(KSP);
72 
73 extern PetscFList KSPList;
74 extern PetscBool  KSPRegisterAllCalled;
75 extern PetscErrorCode  KSPRegisterAll(const char[]);
76 extern PetscErrorCode  KSPRegisterDestroy(void);
77 extern PetscErrorCode  KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
78 
79 /*MC
80    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
81 
82    Synopsis:
83    PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP))
84 
85    Not Collective
86 
87    Input Parameters:
88 +  name_solver - name of a new user-defined solver
89 .  path - path (either absolute or relative) the library containing this solver
90 .  name_create - name of routine to create method context
91 -  routine_create - routine to create method context
92 
93    Notes:
94    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
95 
96    If dynamic libraries are used, then the fourth input argument (routine_create)
97    is ignored.
98 
99    Sample usage:
100 .vb
101    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
102                "MySolverCreate",MySolverCreate);
103 .ve
104 
105    Then, your solver can be chosen with the procedural interface via
106 $     KSPSetType(ksp,"my_solver")
107    or at runtime via the option
108 $     -ksp_type my_solver
109 
110    Level: advanced
111 
112    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
113           and others of the form ${any_environmental_variable} occuring in pathname will be
114           replaced with appropriate values.
115          If your function is not being put into a shared library then use KSPRegister() instead
116 
117 .keywords: KSP, register
118 
119 .seealso: KSPRegisterAll(), KSPRegisterDestroy()
120 
121 M*/
122 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
123 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
124 #else
125 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
126 #endif
127 
128 extern PetscErrorCode  KSPGetType(KSP,const KSPType *);
129 extern PetscErrorCode  KSPSetPCSide(KSP,PCSide);
130 extern PetscErrorCode  KSPGetPCSide(KSP,PCSide*);
131 extern PetscErrorCode  KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
132 extern PetscErrorCode  KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
133 extern PetscErrorCode  KSPSetInitialGuessNonzero(KSP,PetscBool );
134 extern PetscErrorCode  KSPGetInitialGuessNonzero(KSP,PetscBool  *);
135 extern PetscErrorCode  KSPSetInitialGuessKnoll(KSP,PetscBool );
136 extern PetscErrorCode  KSPGetInitialGuessKnoll(KSP,PetscBool *);
137 extern PetscErrorCode  KSPSetErrorIfNotConverged(KSP,PetscBool );
138 extern PetscErrorCode  KSPGetErrorIfNotConverged(KSP,PetscBool  *);
139 extern PetscErrorCode  KSPGetComputeEigenvalues(KSP,PetscBool *);
140 extern PetscErrorCode  KSPSetComputeEigenvalues(KSP,PetscBool );
141 extern PetscErrorCode  KSPGetComputeSingularValues(KSP,PetscBool *);
142 extern PetscErrorCode  KSPSetComputeSingularValues(KSP,PetscBool );
143 extern PetscErrorCode  KSPGetRhs(KSP,Vec *);
144 extern PetscErrorCode  KSPGetSolution(KSP,Vec *);
145 extern PetscErrorCode  KSPGetResidualNorm(KSP,PetscReal*);
146 extern PetscErrorCode  KSPGetIterationNumber(KSP,PetscInt*);
147 extern PetscErrorCode  KSPSetNullSpace(KSP,MatNullSpace);
148 extern PetscErrorCode  KSPGetNullSpace(KSP,MatNullSpace*);
149 extern PetscErrorCode  KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
150 
151 extern PetscErrorCode  KSPSetPC(KSP,PC);
152 extern PetscErrorCode  KSPGetPC(KSP,PC*);
153 
154 extern PetscErrorCode  KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
155 extern PetscErrorCode  KSPMonitorCancel(KSP);
156 extern PetscErrorCode  KSPGetMonitorContext(KSP,void **);
157 extern PetscErrorCode  KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
158 extern PetscErrorCode  KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
159 
160 /* not sure where to put this */
161 extern PetscErrorCode  PCKSPGetKSP(PC,KSP*);
162 extern PetscErrorCode  PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
163 extern PetscErrorCode  PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
164 extern PetscErrorCode  PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
165 extern PetscErrorCode  PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
166 
167 extern PetscErrorCode  PCGalerkinGetKSP(PC,KSP *);
168 
169 extern PetscErrorCode  KSPBuildSolution(KSP,Vec,Vec *);
170 extern PetscErrorCode  KSPBuildResidual(KSP,Vec,Vec,Vec *);
171 
172 extern PetscErrorCode  KSPRichardsonSetScale(KSP,PetscReal);
173 extern PetscErrorCode  KSPRichardsonSetSelfScale(KSP,PetscBool );
174 extern PetscErrorCode  KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
175 extern PetscErrorCode  KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
176 extern PetscErrorCode  KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
177 extern PetscErrorCode  KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
178 
179 extern PetscErrorCode  KSPGMRESSetRestart(KSP, PetscInt);
180 extern PetscErrorCode  KSPGMRESGetRestart(KSP, PetscInt*);
181 extern PetscErrorCode  KSPGMRESSetHapTol(KSP,PetscReal);
182 
183 extern PetscErrorCode  KSPGMRESSetPreAllocateVectors(KSP);
184 extern PetscErrorCode  KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
185 extern PetscErrorCode  KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
186 extern PetscErrorCode  KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
187 extern PetscErrorCode  KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
188 
189 extern PetscErrorCode  KSPLGMRESSetAugDim(KSP,PetscInt);
190 extern PetscErrorCode  KSPLGMRESSetConstant(KSP);
191 
192 extern PetscErrorCode  KSPGCRSetRestart(KSP,PetscInt);
193 extern PetscErrorCode  KSPGCRGetRestart(KSP,PetscInt*);
194 extern PetscErrorCode  KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
195 
196 /*E
197     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
198 
199    Level: advanced
200 
201 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
202           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
203 
204 E*/
205 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
206 extern const char *KSPGMRESCGSRefinementTypes[];
207 /*MC
208     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
209 
210    Level: advanced
211 
212    Note: Possible unstable, but the fastest to compute
213 
214 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
215           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
216           KSPGMRESModifiedGramSchmidtOrthogonalization()
217 M*/
218 
219 /*MC
220     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
221           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
222           poor orthogonality.
223 
224    Level: advanced
225 
226    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
227      estimate the orthogonality but is more stable.
228 
229 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
230           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
231           KSPGMRESModifiedGramSchmidtOrthogonalization()
232 M*/
233 
234 /*MC
235     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
236 
237    Level: advanced
238 
239    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
240      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
241 
242         You should only use this if you absolutely know that the iterative refinement is needed.
243 
244 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
245           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
246           KSPGMRESModifiedGramSchmidtOrthogonalization()
247 M*/
248 
249 extern PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
250 extern PetscErrorCode  KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
251 
252 extern PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
253 extern PetscErrorCode  KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
254 extern PetscErrorCode  KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
255 
256 extern PetscErrorCode  KSPQCGSetTrustRegionRadius(KSP,PetscReal);
257 extern PetscErrorCode  KSPQCGGetQuadratic(KSP,PetscReal*);
258 extern PetscErrorCode  KSPQCGGetTrialStepNorm(KSP,PetscReal*);
259 
260 extern PetscErrorCode  KSPBCGSLSetXRes(KSP,PetscReal);
261 extern PetscErrorCode  KSPBCGSLSetPol(KSP,PetscBool );
262 extern PetscErrorCode  KSPBCGSLSetEll(KSP,PetscInt);
263 
264 extern PetscErrorCode  KSPSetFromOptions(KSP);
265 extern PetscErrorCode  KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
266 
267 extern PetscErrorCode  KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
268 extern PetscErrorCode  KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
269 extern PetscErrorCode  KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *);
270 extern PetscErrorCode  KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
271 extern PetscErrorCode  KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
272 extern PetscErrorCode  KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
273 extern PetscErrorCode  KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
274 extern PetscErrorCode  KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
275 
276 extern PetscErrorCode  KSPUnwindPreconditioner(KSP,Vec,Vec);
277 extern PetscErrorCode  KSPDefaultBuildSolution(KSP,Vec,Vec*);
278 extern PetscErrorCode  KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
279 extern PetscErrorCode  KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
280 
281 extern PetscErrorCode  KSPSetOperators(KSP,Mat,Mat,MatStructure);
282 extern PetscErrorCode  KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
283 extern PetscErrorCode  KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
284 extern PetscErrorCode  KSPSetOptionsPrefix(KSP,const char[]);
285 extern PetscErrorCode  KSPAppendOptionsPrefix(KSP,const char[]);
286 extern PetscErrorCode  KSPGetOptionsPrefix(KSP,const char*[]);
287 
288 extern PetscErrorCode  KSPSetDiagonalScale(KSP,PetscBool );
289 extern PetscErrorCode  KSPGetDiagonalScale(KSP,PetscBool *);
290 extern PetscErrorCode  KSPSetDiagonalScaleFix(KSP,PetscBool );
291 extern PetscErrorCode  KSPGetDiagonalScaleFix(KSP,PetscBool *);
292 
293 extern PetscErrorCode  KSPView(KSP,PetscViewer);
294 
295 extern PetscErrorCode  KSPLSQRSetStandardErrorVec(KSP,Vec);
296 extern PetscErrorCode  KSPLSQRGetStandardErrorVec(KSP,Vec*);
297 
298 extern PetscErrorCode  PCRedundantGetKSP(PC,KSP*);
299 extern PetscErrorCode  PCRedistributeGetKSP(PC,KSP*);
300 
301 /*E
302     KSPNormType - Norm that is passed in the Krylov convergence
303        test routines.
304 
305    Level: advanced
306 
307    Each solver only supports a subset of these and some may support different ones
308    depending on left or right preconditioning, see KSPSetPCSide()
309 
310    Notes: this must match finclude/petscksp.h
311 
312 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
313           KSPSetConvergenceTest(), KSPSetPCSide()
314 E*/
315 typedef enum {KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
316 extern const char *KSPNormTypes[];
317 /*MC
318     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
319           possibly save some computation but means the convergence test cannot
320           be based on a norm of a residual etc.
321 
322    Level: advanced
323 
324     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
325 
326 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
327 M*/
328 
329 /*MC
330     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
331        convergence test routine.
332 
333    Level: advanced
334 
335 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
336 M*/
337 
338 /*MC
339     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
340        convergence test routine.
341 
342    Level: advanced
343 
344 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
345 M*/
346 
347 /*MC
348     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
349        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
350 
351    Level: advanced
352 
353 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
354 M*/
355 
356 extern PetscErrorCode  KSPSetNormType(KSP,KSPNormType);
357 extern PetscErrorCode  KSPGetNormType(KSP,KSPNormType*);
358 extern PetscErrorCode  KSPSetCheckNormIteration(KSP,PetscInt);
359 extern PetscErrorCode  KSPSetLagNorm(KSP,PetscBool );
360 
361 /*E
362     KSPConvergedReason - reason a Krylov method was said to
363          have converged or diverged
364 
365    Level: beginner
366 
367    Notes: See KSPGetConvergedReason() for explanation of each value
368 
369    Developer notes: this must match finclude/petscksp.h
370 
371       The string versions of these are KSPConvergedReasons; if you change
372       any of the values here also change them that array of names.
373 
374 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
375 E*/
376 typedef enum {/* converged */
377               KSP_CONVERGED_RTOL_NORMAL        =  1,
378               KSP_CONVERGED_ATOL_NORMAL        =  9,
379               KSP_CONVERGED_RTOL               =  2,
380               KSP_CONVERGED_ATOL               =  3,
381               KSP_CONVERGED_ITS                =  4,
382               KSP_CONVERGED_CG_NEG_CURVE       =  5,
383               KSP_CONVERGED_CG_CONSTRAINED     =  6,
384               KSP_CONVERGED_STEP_LENGTH        =  7,
385               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
386               /* diverged */
387               KSP_DIVERGED_NULL                = -2,
388               KSP_DIVERGED_ITS                 = -3,
389               KSP_DIVERGED_DTOL                = -4,
390               KSP_DIVERGED_BREAKDOWN           = -5,
391               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
392               KSP_DIVERGED_NONSYMMETRIC        = -7,
393               KSP_DIVERGED_INDEFINITE_PC       = -8,
394               KSP_DIVERGED_NAN                 = -9,
395               KSP_DIVERGED_INDEFINITE_MAT      = -10,
396 
397               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
398 extern const char **KSPConvergedReasons;
399 
400 /*MC
401      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
402 
403    Level: beginner
404 
405    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
406        for left preconditioning it is the 2-norm of the preconditioned residual, and the
407        2-norm of the residual for right preconditioning
408 
409 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
410 
411 M*/
412 
413 /*MC
414      KSP_CONVERGED_ATOL - norm(r) <= atol
415 
416    Level: beginner
417 
418    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
419        for left preconditioning it is the 2-norm of the preconditioned residual, and the
420        2-norm of the residual for right preconditioning
421 
422    Level: beginner
423 
424 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
425 
426 M*/
427 
428 /*MC
429      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
430 
431    Level: beginner
432 
433    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
434        for left preconditioning it is the 2-norm of the preconditioned residual, and the
435        2-norm of the residual for right preconditioning
436 
437    Level: beginner
438 
439 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
440 
441 M*/
442 
443 /*MC
444      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
445       reached
446 
447    Level: beginner
448 
449 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
450 
451 M*/
452 
453 /*MC
454      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
455            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
456            test routine is set in KSP.
457 
458 
459    Level: beginner
460 
461 
462 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
463 
464 M*/
465 
466 /*MC
467      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
468           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
469           preconditioner.
470 
471    Level: beginner
472 
473 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
474 
475 M*/
476 
477 /*MC
478      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
479           method could not continue to enlarge the Krylov space.
480 
481 
482    Level: beginner
483 
484 
485 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
486 
487 M*/
488 
489 /*MC
490      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
491         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
492 
493    Level: beginner
494 
495 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
496 
497 M*/
498 
499 /*MC
500      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
501         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
502         be positive definite
503 
504    Level: beginner
505 
506      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
507   the PCICC preconditioner to generate a positive definite preconditioner
508 
509 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
510 
511 M*/
512 
513 /*MC
514      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
515         while the KSPSolve() is still running.
516 
517    Level: beginner
518 
519 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
520 
521 M*/
522 
523 extern PetscErrorCode  KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
524 extern PetscErrorCode  KSPGetConvergenceContext(KSP,void **);
525 extern PetscErrorCode  KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
526 extern PetscErrorCode  KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
527 extern PetscErrorCode  KSPDefaultConvergedDestroy(void *);
528 extern PetscErrorCode  KSPDefaultConvergedCreate(void **);
529 extern PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP);
530 extern PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP);
531 extern PetscErrorCode  KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
532 extern PetscErrorCode  KSPGetConvergedReason(KSP,KSPConvergedReason *);
533 
534 extern PetscErrorCode  KSPComputeExplicitOperator(KSP,Mat *);
535 
536 /*E
537     KSPCGType - Determines what type of CG to use
538 
539    Level: beginner
540 
541 .seealso: KSPCGSetType()
542 E*/
543 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
544 extern const char *KSPCGTypes[];
545 
546 extern PetscErrorCode  KSPCGSetType(KSP,KSPCGType);
547 extern PetscErrorCode  KSPCGUseSingleReduction(KSP,PetscBool );
548 
549 extern PetscErrorCode  KSPNASHSetRadius(KSP,PetscReal);
550 extern PetscErrorCode  KSPNASHGetNormD(KSP,PetscReal *);
551 extern PetscErrorCode  KSPNASHGetObjFcn(KSP,PetscReal *);
552 
553 extern PetscErrorCode  KSPSTCGSetRadius(KSP,PetscReal);
554 extern PetscErrorCode  KSPSTCGGetNormD(KSP,PetscReal *);
555 extern PetscErrorCode  KSPSTCGGetObjFcn(KSP,PetscReal *);
556 
557 extern PetscErrorCode  KSPGLTRSetRadius(KSP,PetscReal);
558 extern PetscErrorCode  KSPGLTRGetNormD(KSP,PetscReal *);
559 extern PetscErrorCode  KSPGLTRGetObjFcn(KSP,PetscReal *);
560 extern PetscErrorCode  KSPGLTRGetMinEig(KSP,PetscReal *);
561 extern PetscErrorCode  KSPGLTRGetLambda(KSP,PetscReal *);
562 
563 extern PetscErrorCode  KSPPythonSetType(KSP,const char[]);
564 
565 extern PetscErrorCode  PCPreSolve(PC,KSP);
566 extern PetscErrorCode  PCPostSolve(PC,KSP);
567 
568 extern PetscErrorCode  KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
569 extern PetscErrorCode  KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
570 extern PetscErrorCode  KSPMonitorLGDestroy(PetscDrawLG);
571 extern PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
572 extern PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
573 extern PetscErrorCode  KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
574 extern PetscErrorCode  KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
575 extern PetscErrorCode  KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
576 extern PetscErrorCode  KSPMonitorLGRangeDestroy(PetscDrawLG);
577 
578 extern PetscErrorCode  PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
579 extern PetscErrorCode  PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
580 
581 /* see src/ksp/ksp/interface/iguess.c */
582 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
583 
584 extern PetscErrorCode  KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
585 extern PetscErrorCode  KSPFischerGuessDestroy(KSPFischerGuess);
586 extern PetscErrorCode  KSPFischerGuessReset(KSPFischerGuess);
587 extern PetscErrorCode  KSPFischerGuessUpdate(KSPFischerGuess,Vec);
588 extern PetscErrorCode  KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
589 extern PetscErrorCode  KSPFischerGuessSetFromOptions(KSPFischerGuess);
590 
591 extern PetscErrorCode  KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
592 extern PetscErrorCode  KSPSetFischerGuess(KSP,KSPFischerGuess);
593 extern PetscErrorCode  KSPGetFischerGuess(KSP,KSPFischerGuess*);
594 
595 extern PetscErrorCode  MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
596 extern PetscErrorCode  MatSchurComplementGetKSP(Mat,KSP*);
597 extern PetscErrorCode  MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
598 extern PetscErrorCode  MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
599 extern PetscErrorCode  MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
600 
601 extern PetscErrorCode  KSPSetDM(KSP,DM);
602 extern PetscErrorCode  KSPSetDMActive(KSP,PetscBool );
603 extern PetscErrorCode  KSPGetDM(KSP,DM*);
604 
605 PETSC_EXTERN_CXX_END
606 #endif
607