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