xref: /petsc/include/petscksp.h (revision b90c6cbe8a28dbbfeb8633979a88a1769a41a59e)
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 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
170 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
171 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
172 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
173 
174 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
175 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
176 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
177 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
178 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
179 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
180 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
181 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
182 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
183 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
184 
185 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
186 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
187 
188 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
189 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
190 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
191 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
192 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
193 PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
194 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
195 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
196 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
197 
198 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
199 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
200 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
201 
202 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
203 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
204 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
205 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
206 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
207 
208 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
209 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
210 
211 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
212 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
213 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
214 
215 /*E
216     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
217 
218    Level: advanced
219 
220 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
221           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
222 
223 E*/
224 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
225 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
226 /*MC
227     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
228 
229    Level: advanced
230 
231    Note: Possible unstable, but the fastest to compute
232 
233 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
234           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
235           KSPGMRESModifiedGramSchmidtOrthogonalization()
236 M*/
237 
238 /*MC
239     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
240           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
241           poor orthogonality.
242 
243    Level: advanced
244 
245    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
246      estimate the orthogonality but is more stable.
247 
248 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
249           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
250           KSPGMRESModifiedGramSchmidtOrthogonalization()
251 M*/
252 
253 /*MC
254     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
255 
256    Level: advanced
257 
258    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
259      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
260 
261         You should only use this if you absolutely know that the iterative refinement is needed.
262 
263 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
264           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
265           KSPGMRESModifiedGramSchmidtOrthogonalization()
266 M*/
267 
268 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
269 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
270 
271 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
272 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
273 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
274 
275 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
276 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
277 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
278 
279 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
280 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
281 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
282 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
283 
284 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
285 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
286 
287 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
288 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
289 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
290 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
291 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
292 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
293 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
294 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
295 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
296 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
297 PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*);
298 PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**);
299 PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**);
300 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
301 
302 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
303 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
304 
305 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
306 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
307 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
308 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
309 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
310 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
311 PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
312 PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
313 
314 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
315 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
316 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
317 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
318 
319 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
320 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
321 
322 #define KSP_FILE_CLASSID 1211223
323 
324 PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
325 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
326 
327 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
328 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
329 
330 /*E
331     KSPNormType - Norm that is passed in the Krylov convergence
332        test routines.
333 
334    Level: advanced
335 
336    Each solver only supports a subset of these and some may support different ones
337    depending on left or right preconditioning, see KSPSetPCSide()
338 
339    Notes: this must match finclude/petscksp.h
340 
341 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
342           KSPSetConvergenceTest(), KSPSetPCSide()
343 E*/
344 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
345 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
346 PETSC_EXTERN const char *const*const KSPNormTypes;
347 
348 /*MC
349     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
350           possibly save some computation but means the convergence test cannot
351           be based on a norm of a residual etc.
352 
353    Level: advanced
354 
355     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
356 
357 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
358 M*/
359 
360 /*MC
361     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
362        convergence test routine.
363 
364    Level: advanced
365 
366 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
367 M*/
368 
369 /*MC
370     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
371        convergence test routine.
372 
373    Level: advanced
374 
375 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
376 M*/
377 
378 /*MC
379     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
380        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
381 
382    Level: advanced
383 
384 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
385 M*/
386 
387 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
388 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
389 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
390 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
391 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
392 
393 PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
394 
395 /*E
396     KSPConvergedReason - reason a Krylov method was said to
397          have converged or diverged
398 
399    Level: beginner
400 
401    Notes: See KSPGetConvergedReason() for explanation of each value
402 
403    Developer notes: this must match finclude/petscksp.h
404 
405       The string versions of these are KSPConvergedReasons; if you change
406       any of the values here also change them that array of names.
407 
408 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
409 E*/
410 typedef enum {/* converged */
411               KSP_CONVERGED_RTOL_NORMAL        =  1,
412               KSP_CONVERGED_ATOL_NORMAL        =  9,
413               KSP_CONVERGED_RTOL               =  2,
414               KSP_CONVERGED_ATOL               =  3,
415               KSP_CONVERGED_ITS                =  4,
416               KSP_CONVERGED_CG_NEG_CURVE       =  5,
417               KSP_CONVERGED_CG_CONSTRAINED     =  6,
418               KSP_CONVERGED_STEP_LENGTH        =  7,
419               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
420               /* diverged */
421               KSP_DIVERGED_NULL                = -2,
422               KSP_DIVERGED_ITS                 = -3,
423               KSP_DIVERGED_DTOL                = -4,
424               KSP_DIVERGED_BREAKDOWN           = -5,
425               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
426               KSP_DIVERGED_NONSYMMETRIC        = -7,
427               KSP_DIVERGED_INDEFINITE_PC       = -8,
428               KSP_DIVERGED_NANORINF            = -9,
429               KSP_DIVERGED_INDEFINITE_MAT      = -10,
430 
431               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
432 PETSC_EXTERN const char *const*KSPConvergedReasons;
433 
434 /*MC
435      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
436 
437    Level: beginner
438 
439    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
440        for left preconditioning it is the 2-norm of the preconditioned residual, and the
441        2-norm of the residual for right preconditioning
442 
443 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
444 
445 M*/
446 
447 /*MC
448      KSP_CONVERGED_ATOL - norm(r) <= atol
449 
450    Level: beginner
451 
452    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
453        for left preconditioning it is the 2-norm of the preconditioned residual, and the
454        2-norm of the residual for right preconditioning
455 
456    Level: beginner
457 
458 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
459 
460 M*/
461 
462 /*MC
463      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
464 
465    Level: beginner
466 
467    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
468        for left preconditioning it is the 2-norm of the preconditioned residual, and the
469        2-norm of the residual for right preconditioning
470 
471    Level: beginner
472 
473 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
474 
475 M*/
476 
477 /*MC
478      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
479       reached
480 
481    Level: beginner
482 
483 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
484 
485 M*/
486 
487 /*MC
488      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
489            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
490            test routine is set in KSP.
491 
492 
493    Level: beginner
494 
495 
496 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
497 
498 M*/
499 
500 /*MC
501      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
502           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
503           preconditioner.
504 
505    Level: beginner
506 
507 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
508 
509 M*/
510 
511 /*MC
512      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
513           method could not continue to enlarge the Krylov space.
514 
515 
516    Level: beginner
517 
518 
519 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
520 
521 M*/
522 
523 /*MC
524      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
525         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
526 
527    Level: beginner
528 
529 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
530 
531 M*/
532 
533 /*MC
534      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
535         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
536         be positive definite
537 
538    Level: beginner
539 
540      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
541   the PCICC preconditioner to generate a positive definite preconditioner
542 
543 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
544 
545 M*/
546 
547 /*MC
548      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
549         while the KSPSolve() is still running.
550 
551    Level: beginner
552 
553 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
554 
555 M*/
556 
557 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
558 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
559 PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
560 PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
561 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *);
562 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **);
563 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP);
564 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP);
565 PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
566 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
567 
568 PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
569 
570 /*E
571     KSPCGType - Determines what type of CG to use
572 
573    Level: beginner
574 
575 .seealso: KSPCGSetType()
576 E*/
577 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
578 PETSC_EXTERN const char *const KSPCGTypes[];
579 
580 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
581 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
582 
583 PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
584 PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
585 PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
586 
587 PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
588 PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
589 PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
590 
591 PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
592 PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
593 PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
594 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
595 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
596 
597 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
598 
599 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
600 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
601 
602 #include <petscdrawtypes.h>
603 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
604 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
605 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*);
606 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
607 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
608 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
609 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
610 
611 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
612 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
613 
614 /* see src/ksp/ksp/interface/iguess.c */
615 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
616 
617 PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
618 PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
619 PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
620 PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
621 PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
622 PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
623 
624 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
625 PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
626 PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
627 
628 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
629 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
630 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
631 PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat);
632 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
633 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
634 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
635 
636 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
637 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
638 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
639 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
640 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
641 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
642 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
643 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
644 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
645 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*);
646 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
647 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
648 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
649 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
650 
651 #endif
652