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