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