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