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