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