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