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