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