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