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