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