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