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