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