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