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