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