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