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