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 KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 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 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 174 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 175 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 176 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 177 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 178 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 179 180 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 181 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 182 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 183 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 184 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 185 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 186 187 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 188 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 189 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 190 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 191 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 192 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 193 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 194 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 195 196 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 197 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 198 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 199 200 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 201 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 202 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 203 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 204 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 205 206 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 207 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 208 209 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 210 211 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 212 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 213 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 214 215 /*E 216 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 217 218 Level: advanced 219 220 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 221 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 222 223 E*/ 224 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 225 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 226 /*MC 227 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 228 229 Level: advanced 230 231 Note: Possible unstable, but the fastest to compute 232 233 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 234 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 235 KSPGMRESModifiedGramSchmidtOrthogonalization() 236 M*/ 237 238 /*MC 239 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 240 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 241 poor orthogonality. 242 243 Level: advanced 244 245 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 246 estimate the orthogonality but is more stable. 247 248 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 249 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 250 KSPGMRESModifiedGramSchmidtOrthogonalization() 251 M*/ 252 253 /*MC 254 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 255 256 Level: advanced 257 258 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 259 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 260 261 You should only use this if you absolutely know that the iterative refinement is needed. 262 263 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 264 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 265 KSPGMRESModifiedGramSchmidtOrthogonalization() 266 M*/ 267 268 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 269 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 270 271 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 272 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 273 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 274 275 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 276 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 277 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 278 279 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 280 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 281 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 282 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 283 284 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 285 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 286 287 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 288 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 289 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 290 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 291 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 292 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 293 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 294 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 295 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 296 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 297 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 298 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 299 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 300 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 301 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*)); 302 303 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 304 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 305 306 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 307 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 308 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 309 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 310 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 311 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 312 PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 313 PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 314 315 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 316 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 317 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 318 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 319 320 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 321 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 322 PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 323 PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 324 PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 325 326 #define KSP_FILE_CLASSID 1211223 327 328 PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 329 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 330 331 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 332 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 333 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 334 335 /*E 336 KSPNormType - Norm that is passed in the Krylov convergence 337 test routines. 338 339 Level: advanced 340 341 Each solver only supports a subset of these and some may support different ones 342 depending on left or right preconditioning, see KSPSetPCSide() 343 344 Notes: this must match petsc/finclude/petscksp.h 345 346 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 347 KSPSetConvergenceTest(), KSPSetPCSide() 348 E*/ 349 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 350 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 351 PETSC_EXTERN const char *const*const KSPNormTypes; 352 353 /*MC 354 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 355 possibly save some computation but means the convergence test cannot 356 be based on a norm of a residual etc. 357 358 Level: advanced 359 360 Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 361 362 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 363 M*/ 364 365 /*MC 366 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 367 convergence test routine. 368 369 Level: advanced 370 371 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 372 M*/ 373 374 /*MC 375 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 376 convergence test routine. 377 378 Level: advanced 379 380 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 381 M*/ 382 383 /*MC 384 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 385 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 386 387 Level: advanced 388 389 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 390 M*/ 391 392 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 393 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 394 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 395 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 396 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 397 398 /*E 399 KSPConvergedReason - reason a Krylov method was said to have converged or diverged 400 401 Level: beginner 402 403 Notes: See KSPGetConvergedReason() for explanation of each value 404 405 Developer notes: this must match petsc/finclude/petscksp.h 406 407 The string versions of these are KSPConvergedReasons; if you change 408 any of the values here also change them that array of names. 409 410 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 411 E*/ 412 typedef enum {/* converged */ 413 KSP_CONVERGED_RTOL_NORMAL = 1, 414 KSP_CONVERGED_ATOL_NORMAL = 9, 415 KSP_CONVERGED_RTOL = 2, 416 KSP_CONVERGED_ATOL = 3, 417 KSP_CONVERGED_ITS = 4, 418 KSP_CONVERGED_CG_NEG_CURVE = 5, 419 KSP_CONVERGED_CG_CONSTRAINED = 6, 420 KSP_CONVERGED_STEP_LENGTH = 7, 421 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 422 /* diverged */ 423 KSP_DIVERGED_NULL = -2, 424 KSP_DIVERGED_ITS = -3, 425 KSP_DIVERGED_DTOL = -4, 426 KSP_DIVERGED_BREAKDOWN = -5, 427 KSP_DIVERGED_BREAKDOWN_BICG = -6, 428 KSP_DIVERGED_NONSYMMETRIC = -7, 429 KSP_DIVERGED_INDEFINITE_PC = -8, 430 KSP_DIVERGED_NANORINF = -9, 431 KSP_DIVERGED_INDEFINITE_MAT = -10, 432 KSP_DIVERGED_PCSETUP_FAILED = -11, 433 434 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 435 PETSC_EXTERN const char *const*KSPConvergedReasons; 436 437 /*MC 438 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 439 440 Level: beginner 441 442 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 443 for left preconditioning it is the 2-norm of the preconditioned residual, and the 444 2-norm of the residual for right preconditioning 445 446 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 447 448 M*/ 449 450 /*MC 451 KSP_CONVERGED_ATOL - norm(r) <= atol 452 453 Level: beginner 454 455 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 456 for left preconditioning it is the 2-norm of the preconditioned residual, and the 457 2-norm of the residual for right preconditioning 458 459 Level: beginner 460 461 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 462 463 M*/ 464 465 /*MC 466 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 467 468 Level: beginner 469 470 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 471 for left preconditioning it is the 2-norm of the preconditioned residual, and the 472 2-norm of the residual for right preconditioning 473 474 Level: beginner 475 476 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 477 478 M*/ 479 480 /*MC 481 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 482 reached 483 484 Level: beginner 485 486 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 487 488 M*/ 489 490 /*MC 491 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 492 the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 493 test routine is set in KSP. 494 495 Level: beginner 496 497 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 498 499 M*/ 500 501 /*MC 502 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 503 method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 504 preconditioner. 505 506 Level: beginner 507 508 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 509 510 M*/ 511 512 /*MC 513 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 514 method could not continue to enlarge the Krylov space. 515 516 Level: beginner 517 518 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 519 520 M*/ 521 522 /*MC 523 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 524 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 525 526 Level: beginner 527 528 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 529 530 M*/ 531 532 /*MC 533 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 534 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 535 be positive definite 536 537 Level: beginner 538 539 Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 540 the PCICC preconditioner to generate a positive definite preconditioner 541 542 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 543 544 M*/ 545 546 /*MC 547 KSP_DIVERGED_PCSETUP_FAILED - It was not possible to build the requested preconditioner. This is usually due to a 548 zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 549 such as PCFIELDSPLIT. 550 551 Level: beginner 552 553 Notes: Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 554 555 556 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 557 558 M*/ 559 560 /*MC 561 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 562 while the KSPSolve() is still running. 563 564 Level: beginner 565 566 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 567 568 M*/ 569 570 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 571 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 572 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 573 PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 574 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 575 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 576 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 577 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 578 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 579 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 580 581 PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 582 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 583 PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 584 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 585 PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 586 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 587 PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 588 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 589 PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 590 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 591 PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 592 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 593 594 PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 595 596 /*E 597 KSPCGType - Determines what type of CG to use 598 599 Level: beginner 600 601 .seealso: KSPCGSetType() 602 E*/ 603 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 604 PETSC_EXTERN const char *const KSPCGTypes[]; 605 606 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 607 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 608 609 PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 610 PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 611 PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 612 613 PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 614 PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 615 PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 616 617 PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 618 PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 619 PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 620 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 621 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 622 623 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 624 625 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 626 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 627 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*); 634 635 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 636 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 637 638 /* see src/ksp/ksp/interface/iguess.c */ 639 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 640 641 PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 642 PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 643 PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 644 PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 645 PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 646 PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 647 648 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 649 PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 650 PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 651 652 /*E 653 MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 654 655 Level: intermediate 656 657 .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 658 E*/ 659 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 660 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 661 662 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 663 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 664 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 665 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 666 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 667 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 668 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 669 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 670 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 671 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 672 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 673 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 674 675 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 676 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 677 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 678 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 679 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 680 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 681 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 682 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 683 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 684 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 685 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 686 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 687 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 688 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 689 690 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 691 PETSC_EXTERN PetscErrorCode DMProjectField(DM, Vec, 692 void (**)(PetscInt, PetscInt, PetscInt, 693 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 694 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 695 PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec); 696 #endif 697