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