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