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