xref: /petsc/include/petscksp.h (revision b698fc57f0bea7237255b29c1b77df0acc362ffd)
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         Notes:
17     When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
18        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).
19 
20 .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES
21 S*/
22 typedef struct _p_KSP*     KSP;
23 
24 /*J
25     KSPType - String with the name of a PETSc Krylov method.
26 
27    Level: beginner
28 
29 .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
30 J*/
31 typedef const char* KSPType;
32 #define KSPRICHARDSON "richardson"
33 #define KSPCHEBYSHEV  "chebyshev"
34 #define KSPCG         "cg"
35 #define KSPGROPPCG    "groppcg"
36 #define KSPPIPECG     "pipecg"
37 #define KSPPIPECGRR   "pipecgrr"
38 #define KSPPIPELCG     "pipelcg"
39 #define KSPPIPEPRCG    "pipeprcg"
40 #define KSPPIPECG2     "pipecg2"
41 #define   KSPCGNE       "cgne"
42 #define   KSPNASH       "nash"
43 #define   KSPSTCG       "stcg"
44 #define   KSPGLTR       "gltr"
45 #define     KSPCGNASH  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"")  "nash"
46 #define     KSPCGSTCG  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"")  "stcg"
47 #define     KSPCGGLTR  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
48 #define KSPFCG        "fcg"
49 #define KSPPIPEFCG    "pipefcg"
50 #define KSPGMRES      "gmres"
51 #define KSPPIPEFGMRES "pipefgmres"
52 #define   KSPFGMRES     "fgmres"
53 #define   KSPLGMRES     "lgmres"
54 #define   KSPDGMRES     "dgmres"
55 #define   KSPPGMRES     "pgmres"
56 #define KSPTCQMR      "tcqmr"
57 #define KSPBCGS       "bcgs"
58 #define   KSPIBCGS      "ibcgs"
59 #define   KSPFBCGS      "fbcgs"
60 #define   KSPFBCGSR     "fbcgsr"
61 #define   KSPBCGSL      "bcgsl"
62 #define   KSPPIPEBCGS   "pipebcgs"
63 #define KSPCGS        "cgs"
64 #define KSPTFQMR      "tfqmr"
65 #define KSPCR         "cr"
66 #define KSPPIPECR     "pipecr"
67 #define KSPLSQR       "lsqr"
68 #define KSPPREONLY    "preonly"
69 #define KSPQCG        "qcg"
70 #define KSPBICG       "bicg"
71 #define KSPMINRES     "minres"
72 #define KSPSYMMLQ     "symmlq"
73 #define KSPLCD        "lcd"
74 #define KSPPYTHON     "python"
75 #define KSPGCR        "gcr"
76 #define KSPPIPEGCR    "pipegcr"
77 #define KSPTSIRM      "tsirm"
78 #define KSPCGLS       "cgls"
79 #define KSPFETIDP     "fetidp"
80 #define KSPHPDDM      "hpddm"
81 
82 /* Logging support */
83 PETSC_EXTERN PetscClassId KSP_CLASSID;
84 PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
85 PETSC_EXTERN PetscClassId DMKSP_CLASSID;
86 
87 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
88 PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
89 PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
90 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
91 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
92 PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
93 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
94 PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat);
95 PETSC_EXTERN PetscErrorCode KSPSetMatSolveBlockSize(KSP,PetscInt);
96 PETSC_EXTERN PetscErrorCode KSPGetMatSolveBlockSize(KSP,PetscInt*);
97 PETSC_EXTERN PetscErrorCode KSPReset(KSP);
98 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
99 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
100 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
101 PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
102 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
103 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
104 
105 PETSC_EXTERN PetscFunctionList KSPList;
106 PETSC_EXTERN PetscFunctionList KSPGuessList;
107 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
108 
109 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
110 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
111 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
112 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
113 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
114 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
115 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
116 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
117 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
118 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
119 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
120 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
121 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
122 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
123 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
124 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
125 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
126 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
127 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
128 PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
129 
130 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
131 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
132 
133 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
134 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
135 
136 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
137 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
138 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
139 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
140 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*);
141 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);
142 PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*);
143 PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool);
144 
145 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
146 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
147 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
148 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
149 
150 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
151 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
152 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
153 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
154 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
155 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
156 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
157 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
158 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
159 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
160 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
161 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
162 PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
163 /*
164   PCMGCoarseList contains the list of coarse space constructor currently registered
165   These are added with PCMGRegisterCoarseSpaceConstructor()
166 */
167 PETSC_EXTERN PetscFunctionList PCMGCoarseList;
168 PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
169 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
170 
171 
172 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
173 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
174 
175 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
176 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
177 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
178 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
179 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
180 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
181 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
182 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
183 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
184 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
185 
186 /*E
187 
188   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
189 
190   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
191   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
192 
193    Level: intermediate
194 .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
195 
196 E*/
197 typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
198 PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
199 
200 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
201 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
202 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
203 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
204 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
205 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
206 
207 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
208 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
209 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
210 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
211 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
212 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
213 
214 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
215 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
216 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
217 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
218 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
219 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
220 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
221 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
222 
223 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
224 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
225 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
226 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal);
227 
228 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
229 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
230 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
231 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
232 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
233 
234 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
235 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
236 
237 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
238 
239 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
240 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
241 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
242 
243 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
244 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
245 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
246 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
247 
248 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
249 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
250 PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
251 /*E
252     KSPHPDDMType - Type of Krylov method used by KSPHPDDM
253 
254     Level: intermediate
255 
256     Values:
257 +   KSP_HPDDM_TYPE_GMRES (default)
258 .   KSP_HPDDM_TYPE_BGMRES
259 .   KSP_HPDDM_TYPE_CG
260 .   KSP_HPDDM_TYPE_BCG
261 .   KSP_HPDDM_TYPE_GCRODR
262 .   KSP_HPDDM_TYPE_BGCRODR
263 .   KSP_HPDDM_TYPE_BFBCG
264 -   KSP_HPDDM_TYPE_PREONLY
265 
266 .seealso: KSPHPDDM, KSPHPDDMSetType()
267 E*/
268 typedef enum {
269   KSP_HPDDM_TYPE_GMRES = 0,
270   KSP_HPDDM_TYPE_BGMRES = 1,
271   KSP_HPDDM_TYPE_CG = 2,
272   KSP_HPDDM_TYPE_BCG = 3,
273   KSP_HPDDM_TYPE_GCRODR = 4,
274   KSP_HPDDM_TYPE_BGCRODR = 5,
275   KSP_HPDDM_TYPE_BFBCG = 6,
276   KSP_HPDDM_TYPE_PREONLY = 7
277 } KSPHPDDMType;
278 PETSC_EXTERN const char *const KSPHPDDMTypes[];
279 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
280 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
281 
282 /*E
283     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
284 
285    Level: advanced
286 
287 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
288           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
289 
290 E*/
291 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
292 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
293 /*MC
294     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
295 
296    Level: advanced
297 
298    Note: Possible unstable, but the fastest to compute
299 
300 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
301           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
302           KSPGMRESModifiedGramSchmidtOrthogonalization()
303 M*/
304 
305 /*MC
306     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
307           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
308           poor orthogonality.
309 
310    Level: advanced
311 
312    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
313      estimate the orthogonality but is more stable.
314 
315 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
316           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
317           KSPGMRESModifiedGramSchmidtOrthogonalization()
318 M*/
319 
320 /*MC
321     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
322 
323    Level: advanced
324 
325    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
326      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
327 
328         You should only use this if you absolutely know that the iterative refinement is needed.
329 
330 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
331           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
332           KSPGMRESModifiedGramSchmidtOrthogonalization()
333 M*/
334 
335 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
336 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
337 
338 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
339 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
340 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
341 
342 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
343 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
344 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
345 
346 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
347 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
348 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
349 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
350 
351 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
352 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
353 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
354 
355 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
356 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
357 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
358 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
359 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
360 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
361 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
362 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
363 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
364 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
365 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
366 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
367 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
368 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
369 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
370 
371 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
372 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
373 
374 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
375 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
376 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
377 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
378 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
379 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
380 
381 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
382 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
383 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
384 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
385 
386 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
387 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
388 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
389 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer);
390 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
391 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer);
392 
393 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);}
394 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}
395 
396 #define KSP_FILE_CLASSID 1211223
397 
398 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
399 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
400 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
401 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
402 
403 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
404 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
405 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
406 
407 /*E
408     KSPNormType - Norm that is passed in the Krylov convergence
409        test routines.
410 
411    Level: advanced
412 
413    Each solver only supports a subset of these and some may support different ones
414    depending on left or right preconditioning, see KSPSetPCSide()
415 
416    Notes:
417     this must match petsc/finclude/petscksp.h
418 
419 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
420           KSPSetConvergenceTest(), KSPSetPCSide()
421 E*/
422 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
423 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
424 PETSC_EXTERN const char *const*const KSPNormTypes;
425 
426 /*MC
427     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
428           possibly save some computation but means the convergence test cannot
429           be based on a norm of a residual etc.
430 
431    Level: advanced
432 
433     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
434 
435 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
436 M*/
437 
438 /*MC
439     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
440        convergence test routine.
441 
442    Level: advanced
443 
444 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
445 M*/
446 
447 /*MC
448     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
449        convergence test routine.
450 
451    Level: advanced
452 
453 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
454 M*/
455 
456 /*MC
457     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
458        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
459 
460    Level: advanced
461 
462 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
463 M*/
464 
465 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
466 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
467 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
468 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
469 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
470 
471 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
472 /*E
473     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
474 
475    Level: beginner
476 
477    Notes:
478     See KSPGetConvergedReason() for explanation of each value
479 
480    Developer Notes:
481     this must match petsc/finclude/petscksp.h
482 
483       The string versions of these are KSPConvergedReasons; if you change
484       any of the values here also change them that array of names.
485 
486 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView()
487 E*/
488 typedef enum {/* converged */
489               KSP_CONVERGED_RTOL_NORMAL        =  1,
490               KSP_CONVERGED_ATOL_NORMAL        =  9,
491               KSP_CONVERGED_RTOL               =  2,
492               KSP_CONVERGED_ATOL               =  3,
493               KSP_CONVERGED_ITS                =  4,
494               KSP_CONVERGED_CG_NEG_CURVE       =  5,
495               KSP_CONVERGED_CG_CONSTRAINED     =  6,
496               KSP_CONVERGED_STEP_LENGTH        =  7,
497               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
498               /* diverged */
499               KSP_DIVERGED_NULL                = -2,
500               KSP_DIVERGED_ITS                 = -3,
501               KSP_DIVERGED_DTOL                = -4,
502               KSP_DIVERGED_BREAKDOWN           = -5,
503               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
504               KSP_DIVERGED_NONSYMMETRIC        = -7,
505               KSP_DIVERGED_INDEFINITE_PC       = -8,
506               KSP_DIVERGED_NANORINF            = -9,
507               KSP_DIVERGED_INDEFINITE_MAT      = -10,
508               KSP_DIVERGED_PC_FAILED           = -11,
509               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
510 
511               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
512 PETSC_EXTERN const char *const*KSPConvergedReasons;
513 
514 /*MC
515      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
516 
517    Level: beginner
518 
519    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
520        for left preconditioning it is the 2-norm of the preconditioned residual, and the
521        2-norm of the residual for right preconditioning
522 
523    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
524 
525 .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
526 
527 M*/
528 
529 /*MC
530      KSP_CONVERGED_ATOL - norm(r) <= atol
531 
532    Level: beginner
533 
534    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
535        for left preconditioning it is the 2-norm of the preconditioned residual, and the
536        2-norm of the residual for right preconditioning
537 
538    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
539 
540 .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
541 
542 M*/
543 
544 /*MC
545      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
546 
547    Level: beginner
548 
549    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
550        for left preconditioning it is the 2-norm of the preconditioned residual, and the
551        2-norm of the residual for right preconditioning
552 
553    Level: beginner
554 
555 .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
556 
557 M*/
558 
559 /*MC
560      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
561       reached
562 
563    Level: beginner
564 
565 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
566 
567 M*/
568 
569 /*MC
570      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
571            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
572            test routine is set in KSP.
573 
574    Level: beginner
575 
576 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
577 
578 M*/
579 
580 /*MC
581      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
582           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
583           preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
584           are colinear.
585 
586    Level: beginner
587 
588 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
589 
590 M*/
591 
592 /*MC
593      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
594           method could not continue to enlarge the Krylov space.
595 
596    Level: beginner
597 
598 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
599 
600 M*/
601 
602 /*MC
603      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
604         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
605 
606    Level: beginner
607 
608 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
609 
610 M*/
611 
612 /*MC
613      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
614         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
615         be positive definite
616 
617    Level: beginner
618 
619      Notes:
620     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
621   the PCICC preconditioner to generate a positive definite preconditioner
622 
623 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
624 
625 M*/
626 
627 /*MC
628      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
629      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
630      such as PCFIELDSPLIT.
631 
632    Level: beginner
633 
634     Notes:
635     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
636 
637 
638 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
639 
640 M*/
641 
642 /*MC
643      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
644         while the KSPSolve() is still running.
645 
646    Level: beginner
647 
648 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
649 
650 M*/
651 
652 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
653 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
654 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
655 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
656 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
657 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
658 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
659 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
660 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
661 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
662 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
663 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
664 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
665 PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
666 
667 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
668 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
669 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
670 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
671 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
672 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
673 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
674 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
675 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
676 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
677 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
678 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
679 
680 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
681 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
682 
683 /*E
684     KSPCGType - Determines what type of CG to use
685 
686    Level: beginner
687 
688 .seealso: KSPCGSetType()
689 E*/
690 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
691 PETSC_EXTERN const char *const KSPCGTypes[];
692 
693 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
694 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);
695 
696 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
697 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
698 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
699 
700 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
701 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
702 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
703 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
704 
705 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
706 
707 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
708 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
709 
710 #include <petscdrawtypes.h>
711 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
712 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
713 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
714 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
715 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
716 
717 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
718 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
719 
720 /*S
721      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
722 
723    Level: beginner
724 
725 .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
726 S*/
727 typedef struct _p_KSPGuess* KSPGuess;
728 /*J
729     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
730 
731    Level: beginner
732 
733 .seealso: KSPGuess
734 J*/
735 typedef const char* KSPGuessType;
736 #define KSPGUESSFISCHER "fischer"
737 #define KSPGUESSPOD     "pod"
738 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
739 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
740 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
741 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
742 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
743 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
744 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
745 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
746 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
747 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
748 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
749 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
750 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
751 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
752 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
753 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
754 
755 /*E
756     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
757 
758     Level: intermediate
759 
760 .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
761 E*/
762 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
763 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
764 
765 typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
766               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
767               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
768               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
769 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
770 
771 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
772 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
773 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
774 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
775 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
776 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
777 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
778 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
779 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
780 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
781 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
782 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
783 
784 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
785 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
786 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
787 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
788 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
789 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
790 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
791 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
792 
793 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
794 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
795 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
796 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
797 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
798 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
799 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
800 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
801 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
802 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
803 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
804 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
805 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
806 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
807 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
808 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
809 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
810 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
811 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
812 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
813 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
814 
815 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
816 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
817 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
818 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
819 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
820 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
821 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
822 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
823 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
824 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
825 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
826 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
827 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
828 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
829 
830 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
831 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
832                                            void (**)(PetscInt, PetscInt, PetscInt,
833                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
834                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
835                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
836 
837 
838 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *);
839 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal);
840 #endif
841