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