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