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