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