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