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