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