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