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