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