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