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