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