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 /*E 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 E*/ 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 KSPNGMRES "ngmres" 61 #define KSPSPECEST "specest" 62 63 /* Logging support */ 64 extern PetscClassId KSP_CLASSID; 65 66 extern PetscErrorCode KSPCreate(MPI_Comm,KSP *); 67 extern PetscErrorCode KSPSetType(KSP,const KSPType); 68 extern PetscErrorCode KSPSetUp(KSP); 69 extern PetscErrorCode KSPSetUpOnBlocks(KSP); 70 extern PetscErrorCode KSPSolve(KSP,Vec,Vec); 71 extern PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 72 extern PetscErrorCode KSPReset(KSP); 73 extern PetscErrorCode KSPDestroy_(KSP); 74 #define KSPDestroy(a) (KSPDestroy_(a) || (((a) = 0),0)) 75 76 extern PetscFList KSPList; 77 extern PetscBool KSPRegisterAllCalled; 78 extern PetscErrorCode KSPRegisterAll(const char[]); 79 extern PetscErrorCode KSPRegisterDestroy(void); 80 extern PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 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 extern PetscErrorCode KSPGetType(KSP,const KSPType *); 132 extern PetscErrorCode KSPSetPCSide(KSP,PCSide); 133 extern PetscErrorCode KSPGetPCSide(KSP,PCSide*); 134 extern PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 135 extern PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 136 extern PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 137 extern PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 138 extern PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 139 extern PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 140 extern PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 141 extern PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 142 extern PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 143 extern PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 144 extern PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 145 extern PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 146 extern PetscErrorCode KSPGetRhs(KSP,Vec *); 147 extern PetscErrorCode KSPGetSolution(KSP,Vec *); 148 extern PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 149 extern PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 150 extern PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 151 extern PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 152 extern PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 153 154 extern PetscErrorCode KSPSetPC(KSP,PC); 155 extern PetscErrorCode KSPGetPC(KSP,PC*); 156 157 extern PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 158 extern PetscErrorCode KSPMonitorCancel(KSP); 159 extern PetscErrorCode KSPGetMonitorContext(KSP,void **); 160 extern PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 161 extern PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 162 163 /* not sure where to put this */ 164 extern PetscErrorCode PCKSPGetKSP(PC,KSP*); 165 extern PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 166 extern PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 167 extern PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 168 extern PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 169 170 extern PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 171 172 extern PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 173 extern PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 174 175 extern PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 176 extern PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 177 extern PetscErrorCode KSPChebychevSetEigenvalues(KSP,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 KSPMonitorDefaultLSQR(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 KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 278 279 extern PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 280 extern PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 281 extern PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 282 extern PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 283 284 extern PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 285 extern PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 286 extern PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 287 extern PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 288 extern PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 289 extern PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 290 291 extern PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 292 extern PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 293 extern PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 294 extern PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 295 296 extern PetscErrorCode KSPView(KSP,PetscViewer); 297 298 extern PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 299 extern PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 300 301 extern PetscErrorCode PCRedundantGetKSP(PC,KSP*); 302 extern PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 303 304 /*E 305 KSPNormType - Norm that is passed in the Krylov convergence 306 test routines. 307 308 Level: advanced 309 310 Each solver only supports a subset of these and some may support different ones 311 depending on left or right preconditioning, see KSPSetPCSide() 312 313 Notes: this must match finclude/petscksp.h 314 315 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 316 KSPSetConvergenceTest(), KSPSetPCSide() 317 E*/ 318 typedef enum {KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 319 extern const char *KSPNormTypes[]; 320 /*MC 321 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 322 possibly save some computation but means the convergence test cannot 323 be based on a norm of a residual etc. 324 325 Level: advanced 326 327 Note: Some Krylov methods need to compute a residual norm and then this option is ignored 328 329 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 330 M*/ 331 332 /*MC 333 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 334 convergence test routine. 335 336 Level: advanced 337 338 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 339 M*/ 340 341 /*MC 342 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 343 convergence test routine. 344 345 Level: advanced 346 347 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 348 M*/ 349 350 /*MC 351 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 352 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 353 354 Level: advanced 355 356 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 357 M*/ 358 359 extern PetscErrorCode KSPSetNormType(KSP,KSPNormType); 360 extern PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 361 extern PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 362 extern PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 363 364 /*E 365 KSPConvergedReason - reason a Krylov method was said to 366 have converged or diverged 367 368 Level: beginner 369 370 Notes: See KSPGetConvergedReason() for explanation of each value 371 372 Developer notes: this must match finclude/petscksp.h 373 374 The string versions of these are KSPConvergedReasons; if you change 375 any of the values here also change them that array of names. 376 377 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 378 E*/ 379 typedef enum {/* converged */ 380 KSP_CONVERGED_RTOL_NORMAL = 1, 381 KSP_CONVERGED_ATOL_NORMAL = 9, 382 KSP_CONVERGED_RTOL = 2, 383 KSP_CONVERGED_ATOL = 3, 384 KSP_CONVERGED_ITS = 4, 385 KSP_CONVERGED_CG_NEG_CURVE = 5, 386 KSP_CONVERGED_CG_CONSTRAINED = 6, 387 KSP_CONVERGED_STEP_LENGTH = 7, 388 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 389 /* diverged */ 390 KSP_DIVERGED_NULL = -2, 391 KSP_DIVERGED_ITS = -3, 392 KSP_DIVERGED_DTOL = -4, 393 KSP_DIVERGED_BREAKDOWN = -5, 394 KSP_DIVERGED_BREAKDOWN_BICG = -6, 395 KSP_DIVERGED_NONSYMMETRIC = -7, 396 KSP_DIVERGED_INDEFINITE_PC = -8, 397 KSP_DIVERGED_NAN = -9, 398 KSP_DIVERGED_INDEFINITE_MAT = -10, 399 400 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 401 extern const char *const*KSPConvergedReasons; 402 403 /*MC 404 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 405 406 Level: beginner 407 408 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 409 for left preconditioning it is the 2-norm of the preconditioned residual, and the 410 2-norm of the residual for right preconditioning 411 412 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 413 414 M*/ 415 416 /*MC 417 KSP_CONVERGED_ATOL - norm(r) <= atol 418 419 Level: beginner 420 421 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 422 for left preconditioning it is the 2-norm of the preconditioned residual, and the 423 2-norm of the residual for right preconditioning 424 425 Level: beginner 426 427 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 428 429 M*/ 430 431 /*MC 432 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 433 434 Level: beginner 435 436 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 437 for left preconditioning it is the 2-norm of the preconditioned residual, and the 438 2-norm of the residual for right preconditioning 439 440 Level: beginner 441 442 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 443 444 M*/ 445 446 /*MC 447 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 448 reached 449 450 Level: beginner 451 452 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 453 454 M*/ 455 456 /*MC 457 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 458 the preconditioner is applied. Also used when the KSPSkipConverged() convergence 459 test routine is set in KSP. 460 461 462 Level: beginner 463 464 465 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 466 467 M*/ 468 469 /*MC 470 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 471 method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 472 preconditioner. 473 474 Level: beginner 475 476 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 477 478 M*/ 479 480 /*MC 481 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 482 method could not continue to enlarge the Krylov space. 483 484 485 Level: beginner 486 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 extern PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 527 extern PetscErrorCode KSPGetConvergenceContext(KSP,void **); 528 extern PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 529 extern PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 530 extern PetscErrorCode KSPDefaultConvergedDestroy(void *); 531 extern PetscErrorCode KSPDefaultConvergedCreate(void **); 532 extern PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 533 extern PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 534 extern PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 535 extern PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 536 537 extern PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 538 539 /*E 540 KSPCGType - Determines what type of CG to use 541 542 Level: beginner 543 544 .seealso: KSPCGSetType() 545 E*/ 546 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 547 extern const char *KSPCGTypes[]; 548 549 extern PetscErrorCode KSPCGSetType(KSP,KSPCGType); 550 extern PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 551 552 extern PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 553 extern PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 554 extern PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 555 556 extern PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 557 extern PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 558 extern PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 559 560 extern PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 561 extern PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 562 extern PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 563 extern PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 564 extern PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 565 566 extern PetscErrorCode KSPPythonSetType(KSP,const char[]); 567 568 extern PetscErrorCode PCPreSolve(PC,KSP); 569 extern PetscErrorCode PCPostSolve(PC,KSP); 570 571 extern PetscErrorCode KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 572 extern PetscErrorCode KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 573 extern PetscErrorCode KSPMonitorLGDestroy(PetscDrawLG); 574 extern PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 575 extern PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 576 extern PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 577 extern PetscErrorCode KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 578 extern PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 579 extern PetscErrorCode KSPMonitorLGRangeDestroy(PetscDrawLG); 580 581 extern PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 582 extern PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 583 584 /* see src/ksp/ksp/interface/iguess.c */ 585 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 586 587 extern PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 588 extern PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 589 extern PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 590 extern PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 591 extern PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 592 extern PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 593 594 extern PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 595 extern PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 596 extern PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 597 598 extern PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 599 extern PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 600 extern PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 601 extern PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 602 extern PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 603 604 extern PetscErrorCode KSPSetDM(KSP,DM); 605 extern PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 606 extern PetscErrorCode KSPGetDM(KSP,DM*); 607 608 PETSC_EXTERN_CXX_END 609 #endif 610