xref: /petsc/include/petscksp.h (revision ef0bb6c736604ce380bf8bea4ebd4a7bda431d97)
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 
231 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
232 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
233 
234 /*E
235     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
236 
237    Level: advanced
238 
239 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
240           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
241 
242 E*/
243 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
244 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
245 /*MC
246     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
247 
248    Level: advanced
249 
250    Note: Possible unstable, but the fastest to compute
251 
252 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
253           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
254           KSPGMRESModifiedGramSchmidtOrthogonalization()
255 M*/
256 
257 /*MC
258     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
259           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
260           poor orthogonality.
261 
262    Level: advanced
263 
264    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
265      estimate the orthogonality but is more stable.
266 
267 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
268           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
269           KSPGMRESModifiedGramSchmidtOrthogonalization()
270 M*/
271 
272 /*MC
273     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
274 
275    Level: advanced
276 
277    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
278      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
279 
280         You should only use this if you absolutely know that the iterative refinement is needed.
281 
282 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
283           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
284           KSPGMRESModifiedGramSchmidtOrthogonalization()
285 M*/
286 
287 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
288 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
289 
290 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
291 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
292 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
293 
294 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
295 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
296 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
297 
298 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
299 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
300 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
301 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
302 
303 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
304 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
305 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
306 
307 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
308 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
309 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
310 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
311 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
312 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
313 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
314 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
315 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
316 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
317 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
318 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
319 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
320 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
321 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
322 
323 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
324 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
325 
326 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
327 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
328 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
329 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
330 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
331 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
332 
333 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
334 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
335 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
336 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
337 
338 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
339 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
340 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
341 PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
342 PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
343 
344 #define KSP_FILE_CLASSID 1211223
345 
346 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
347 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
348 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
349 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
350 
351 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
352 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
353 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
354 
355 /*E
356     KSPNormType - Norm that is passed in the Krylov convergence
357        test routines.
358 
359    Level: advanced
360 
361    Each solver only supports a subset of these and some may support different ones
362    depending on left or right preconditioning, see KSPSetPCSide()
363 
364    Notes:
365     this must match petsc/finclude/petscksp.h
366 
367 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
368           KSPSetConvergenceTest(), KSPSetPCSide()
369 E*/
370 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
371 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
372 PETSC_EXTERN const char *const*const KSPNormTypes;
373 
374 /*MC
375     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
376           possibly save some computation but means the convergence test cannot
377           be based on a norm of a residual etc.
378 
379    Level: advanced
380 
381     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
382 
383 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
384 M*/
385 
386 /*MC
387     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
388        convergence test routine.
389 
390    Level: advanced
391 
392 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
393 M*/
394 
395 /*MC
396     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
397        convergence test routine.
398 
399    Level: advanced
400 
401 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
402 M*/
403 
404 /*MC
405     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
406        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
407 
408    Level: advanced
409 
410 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
411 M*/
412 
413 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
414 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
415 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
416 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
417 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
418 
419 /*E
420     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
421 
422    Level: beginner
423 
424    Notes:
425     See KSPGetConvergedReason() for explanation of each value
426 
427    Developer Notes:
428     this must match petsc/finclude/petscksp.h
429 
430       The string versions of these are KSPConvergedReasons; if you change
431       any of the values here also change them that array of names.
432 
433 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
434 E*/
435 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
436 typedef enum {/* converged */
437               KSP_CONVERGED_RTOL_NORMAL        =  1,
438               KSP_CONVERGED_ATOL_NORMAL        =  9,
439               KSP_CONVERGED_RTOL               =  2,
440               KSP_CONVERGED_ATOL               =  3,
441               KSP_CONVERGED_ITS                =  4,
442               KSP_CONVERGED_CG_NEG_CURVE       =  5,
443               KSP_CONVERGED_CG_CONSTRAINED     =  6,
444               KSP_CONVERGED_STEP_LENGTH        =  7,
445               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
446               /* diverged */
447               KSP_DIVERGED_NULL                = -2,
448               KSP_DIVERGED_ITS                 = -3,
449               KSP_DIVERGED_DTOL                = -4,
450               KSP_DIVERGED_BREAKDOWN           = -5,
451               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
452               KSP_DIVERGED_NONSYMMETRIC        = -7,
453               KSP_DIVERGED_INDEFINITE_PC       = -8,
454               KSP_DIVERGED_NANORINF            = -9,
455               KSP_DIVERGED_INDEFINITE_MAT      = -10,
456               KSP_DIVERGED_PC_FAILED           = -11,
457               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
458 
459               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
460 PETSC_EXTERN const char *const*KSPConvergedReasons;
461 
462 /*MC
463      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
464 
465    Level: beginner
466 
467    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
468        for left preconditioning it is the 2-norm of the preconditioned residual, and the
469        2-norm of the residual for right preconditioning
470 
471    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
472 
473 .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
474 
475 M*/
476 
477 /*MC
478      KSP_CONVERGED_ATOL - norm(r) <= atol
479 
480    Level: beginner
481 
482    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
483        for left preconditioning it is the 2-norm of the preconditioned residual, and the
484        2-norm of the residual for right preconditioning
485 
486    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
487 
488 .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
489 
490 M*/
491 
492 /*MC
493      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
494 
495    Level: beginner
496 
497    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
498        for left preconditioning it is the 2-norm of the preconditioned residual, and the
499        2-norm of the residual for right preconditioning
500 
501    Level: beginner
502 
503 .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
504 
505 M*/
506 
507 /*MC
508      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
509       reached
510 
511    Level: beginner
512 
513 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
514 
515 M*/
516 
517 /*MC
518      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
519            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
520            test routine is set in KSP.
521 
522    Level: beginner
523 
524 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
525 
526 M*/
527 
528 /*MC
529      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
530           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
531           preconditioner.
532 
533    Level: beginner
534 
535 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
536 
537 M*/
538 
539 /*MC
540      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
541           method could not continue to enlarge the Krylov space.
542 
543    Level: beginner
544 
545 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
546 
547 M*/
548 
549 /*MC
550      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
551         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
552 
553    Level: beginner
554 
555 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
556 
557 M*/
558 
559 /*MC
560      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
561         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
562         be positive definite
563 
564    Level: beginner
565 
566      Notes:
567     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
568   the PCICC preconditioner to generate a positive definite preconditioner
569 
570 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
571 
572 M*/
573 
574 /*MC
575      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
576      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
577      such as PCFIELDSPLIT.
578 
579    Level: beginner
580 
581     Notes:
582     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
583 
584 
585 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
586 
587 M*/
588 
589 /*MC
590      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
591         while the KSPSolve() is still running.
592 
593    Level: beginner
594 
595 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
596 
597 M*/
598 
599 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
600 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
601 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
602 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
603 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
604 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
605 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
606 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
607 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
608 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
609 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
610 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
611 
612 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
613 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
614 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
615 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
616 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
617 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
618 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
619 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
620 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
621 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
622 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
623 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
624 
625 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
626 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
627 
628 /*E
629     KSPCGType - Determines what type of CG to use
630 
631    Level: beginner
632 
633 .seealso: KSPCGSetType()
634 E*/
635 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
636 PETSC_EXTERN const char *const KSPCGTypes[];
637 
638 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
639 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
640 
641 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
642 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
643 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
644 
645 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
646 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
647 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
648 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
649 
650 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
651 
652 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
653 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
654 
655 #include <petscdrawtypes.h>
656 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
657 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
658 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
659 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
660 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
661 
662 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
663 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
664 
665 /*S
666      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
667 
668    Level: beginner
669 
670 .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
671 S*/
672 typedef struct _p_KSPGuess* KSPGuess;
673 /*J
674     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
675 
676    Level: beginner
677 
678 .seealso: KSPGuess
679 J*/
680 typedef const char* KSPGuessType;
681 #define KSPGUESSFISCHER "fischer"
682 #define KSPGUESSPOD     "pod"
683 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
684 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
685 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
686 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
687 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
688 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
689 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
690 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
691 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
692 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
693 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
694 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
695 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
696 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
697 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
698 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
699 
700 /*E
701     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
702 
703     Level: intermediate
704 
705 .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
706 E*/
707 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
708 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
709 
710 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
711 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
712 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
713 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
714 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
715 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
716 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
717 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
718 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
719 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
720 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
721 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
722 
723 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
724 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
725 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
726 PETSC_EXTERN PetscErrorCode MatCreateLMVMBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
727 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
728 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
729 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
730 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
731 
732 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
733 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
734 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
735 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
736 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
737 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
738 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
739 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
740 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
741 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
742 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
743 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
744 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
745 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
746 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
747 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
748 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
749 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
750 PETSC_EXTERN PetscErrorCode MatSymBrdnSetDelta(Mat, PetscScalar);
751 
752 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
753 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
754 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
755 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
756 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
757 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
758 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
759 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
760 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
761 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
762 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
763 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
764 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
765 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
766 
767 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
768 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
769                                            void (**)(PetscInt, PetscInt, PetscInt,
770                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
771                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
772                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
773 #endif
774