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