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