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