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