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