xref: /petsc/include/petscksp.h (revision ff218e97a57ed641f3ebc93f697e38ef0f3aa217)
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_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
318 extern const char *KSPNormTypes[];
319 
320 /*MC
321     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
322           possibly save some computation but means the convergence test cannot
323           be based on a norm of a residual etc.
324 
325    Level: advanced
326 
327     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
328 
329 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
330 M*/
331 
332 /*MC
333     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
334        convergence test routine.
335 
336    Level: advanced
337 
338 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
339 M*/
340 
341 /*MC
342     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
343        convergence test routine.
344 
345    Level: advanced
346 
347 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
348 M*/
349 
350 /*MC
351     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
352        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
353 
354    Level: advanced
355 
356 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
357 M*/
358 
359 extern PetscErrorCode  KSPSetNormType(KSP,KSPNormType);
360 extern PetscErrorCode  KSPGetNormType(KSP,KSPNormType*);
361 extern PetscErrorCode  KSPSetCheckNormIteration(KSP,PetscInt);
362 extern PetscErrorCode  KSPSetLagNorm(KSP,PetscBool );
363 
364 /*E
365     KSPConvergedReason - reason a Krylov method was said to
366          have converged or diverged
367 
368    Level: beginner
369 
370    Notes: See KSPGetConvergedReason() for explanation of each value
371 
372    Developer notes: this must match finclude/petscksp.h
373 
374       The string versions of these are KSPConvergedReasons; if you change
375       any of the values here also change them that array of names.
376 
377 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
378 E*/
379 typedef enum {/* converged */
380               KSP_CONVERGED_RTOL_NORMAL        =  1,
381               KSP_CONVERGED_ATOL_NORMAL        =  9,
382               KSP_CONVERGED_RTOL               =  2,
383               KSP_CONVERGED_ATOL               =  3,
384               KSP_CONVERGED_ITS                =  4,
385               KSP_CONVERGED_CG_NEG_CURVE       =  5,
386               KSP_CONVERGED_CG_CONSTRAINED     =  6,
387               KSP_CONVERGED_STEP_LENGTH        =  7,
388               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
389               /* diverged */
390               KSP_DIVERGED_NULL                = -2,
391               KSP_DIVERGED_ITS                 = -3,
392               KSP_DIVERGED_DTOL                = -4,
393               KSP_DIVERGED_BREAKDOWN           = -5,
394               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
395               KSP_DIVERGED_NONSYMMETRIC        = -7,
396               KSP_DIVERGED_INDEFINITE_PC       = -8,
397               KSP_DIVERGED_NAN                 = -9,
398               KSP_DIVERGED_INDEFINITE_MAT      = -10,
399 
400               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
401 extern const char *const*KSPConvergedReasons;
402 
403 /*MC
404      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
405 
406    Level: beginner
407 
408    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
409        for left preconditioning it is the 2-norm of the preconditioned residual, and the
410        2-norm of the residual for right preconditioning
411 
412 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
413 
414 M*/
415 
416 /*MC
417      KSP_CONVERGED_ATOL - norm(r) <= atol
418 
419    Level: beginner
420 
421    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
422        for left preconditioning it is the 2-norm of the preconditioned residual, and the
423        2-norm of the residual for right preconditioning
424 
425    Level: beginner
426 
427 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
428 
429 M*/
430 
431 /*MC
432      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
433 
434    Level: beginner
435 
436    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
437        for left preconditioning it is the 2-norm of the preconditioned residual, and the
438        2-norm of the residual for right preconditioning
439 
440    Level: beginner
441 
442 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
443 
444 M*/
445 
446 /*MC
447      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
448       reached
449 
450    Level: beginner
451 
452 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
453 
454 M*/
455 
456 /*MC
457      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
458            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
459            test routine is set in KSP.
460 
461 
462    Level: beginner
463 
464 
465 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
466 
467 M*/
468 
469 /*MC
470      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
471           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
472           preconditioner.
473 
474    Level: beginner
475 
476 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
477 
478 M*/
479 
480 /*MC
481      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
482           method could not continue to enlarge the Krylov space.
483 
484 
485    Level: beginner
486 
487 
488 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
489 
490 M*/
491 
492 /*MC
493      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
494         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
495 
496    Level: beginner
497 
498 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
499 
500 M*/
501 
502 /*MC
503      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
504         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
505         be positive definite
506 
507    Level: beginner
508 
509      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
510   the PCICC preconditioner to generate a positive definite preconditioner
511 
512 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
513 
514 M*/
515 
516 /*MC
517      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
518         while the KSPSolve() is still running.
519 
520    Level: beginner
521 
522 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
523 
524 M*/
525 
526 extern PetscErrorCode  KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
527 extern PetscErrorCode  KSPGetConvergenceContext(KSP,void **);
528 extern PetscErrorCode  KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
529 extern PetscErrorCode  KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
530 extern PetscErrorCode  KSPDefaultConvergedDestroy(void *);
531 extern PetscErrorCode  KSPDefaultConvergedCreate(void **);
532 extern PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP);
533 extern PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP);
534 extern PetscErrorCode  KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
535 extern PetscErrorCode  KSPGetConvergedReason(KSP,KSPConvergedReason *);
536 
537 extern PetscErrorCode  KSPComputeExplicitOperator(KSP,Mat *);
538 
539 /*E
540     KSPCGType - Determines what type of CG to use
541 
542    Level: beginner
543 
544 .seealso: KSPCGSetType()
545 E*/
546 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
547 extern const char *KSPCGTypes[];
548 
549 extern PetscErrorCode  KSPCGSetType(KSP,KSPCGType);
550 extern PetscErrorCode  KSPCGUseSingleReduction(KSP,PetscBool );
551 
552 extern PetscErrorCode  KSPNASHSetRadius(KSP,PetscReal);
553 extern PetscErrorCode  KSPNASHGetNormD(KSP,PetscReal *);
554 extern PetscErrorCode  KSPNASHGetObjFcn(KSP,PetscReal *);
555 
556 extern PetscErrorCode  KSPSTCGSetRadius(KSP,PetscReal);
557 extern PetscErrorCode  KSPSTCGGetNormD(KSP,PetscReal *);
558 extern PetscErrorCode  KSPSTCGGetObjFcn(KSP,PetscReal *);
559 
560 extern PetscErrorCode  KSPGLTRSetRadius(KSP,PetscReal);
561 extern PetscErrorCode  KSPGLTRGetNormD(KSP,PetscReal *);
562 extern PetscErrorCode  KSPGLTRGetObjFcn(KSP,PetscReal *);
563 extern PetscErrorCode  KSPGLTRGetMinEig(KSP,PetscReal *);
564 extern PetscErrorCode  KSPGLTRGetLambda(KSP,PetscReal *);
565 
566 extern PetscErrorCode  KSPPythonSetType(KSP,const char[]);
567 
568 extern PetscErrorCode  PCPreSolve(PC,KSP);
569 extern PetscErrorCode  PCPostSolve(PC,KSP);
570 
571 extern PetscErrorCode  KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
572 extern PetscErrorCode  KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
573 extern PetscErrorCode  KSPMonitorLGDestroy(PetscDrawLG*);
574 extern PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
575 extern PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
576 extern PetscErrorCode  KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
577 extern PetscErrorCode  KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
578 extern PetscErrorCode  KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
579 extern PetscErrorCode  KSPMonitorLGRangeDestroy(PetscDrawLG*);
580 
581 extern PetscErrorCode  PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
582 extern PetscErrorCode  PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
583 
584 /* see src/ksp/ksp/interface/iguess.c */
585 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
586 
587 extern PetscErrorCode  KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
588 extern PetscErrorCode  KSPFischerGuessDestroy(KSPFischerGuess*);
589 extern PetscErrorCode  KSPFischerGuessReset(KSPFischerGuess);
590 extern PetscErrorCode  KSPFischerGuessUpdate(KSPFischerGuess,Vec);
591 extern PetscErrorCode  KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
592 extern PetscErrorCode  KSPFischerGuessSetFromOptions(KSPFischerGuess);
593 
594 extern PetscErrorCode  KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
595 extern PetscErrorCode  KSPSetFischerGuess(KSP,KSPFischerGuess);
596 extern PetscErrorCode  KSPGetFischerGuess(KSP,KSPFischerGuess*);
597 
598 extern PetscErrorCode  MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
599 extern PetscErrorCode  MatSchurComplementGetKSP(Mat,KSP*);
600 extern PetscErrorCode  MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
601 extern PetscErrorCode  MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
602 extern PetscErrorCode  MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
603 
604 extern PetscErrorCode  MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat);
605 
606 extern PetscErrorCode  KSPSetDM(KSP,DM);
607 extern PetscErrorCode  KSPSetDMActive(KSP,PetscBool );
608 extern PetscErrorCode  KSPGetDM(KSP,DM*);
609 extern PetscErrorCode  KSPSetApplicationContext(KSP,void*);
610 extern PetscErrorCode  KSPGetApplicationContext(KSP,void*);
611 
612 PETSC_EXTERN_CXX_END
613 #endif
614