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