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