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