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