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