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