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