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 KSPTCQMR "tcqmr" 43 #define KSPBCGS "bcgs" 44 #define KSPIBCGS "ibcgs" 45 #define KSPBCGSL "bcgsl" 46 #define KSPCGS "cgs" 47 #define KSPTFQMR "tfqmr" 48 #define KSPCR "cr" 49 #define KSPLSQR "lsqr" 50 #define KSPPREONLY "preonly" 51 #define KSPQCG "qcg" 52 #define KSPBICG "bicg" 53 #define KSPMINRES "minres" 54 #define KSPSYMMLQ "symmlq" 55 #define KSPLCD "lcd" 56 #define KSPPYTHON "python" 57 #define KSPBROYDEN "broyden" 58 #define KSPGCR "gcr" 59 #define KSPNGMRES "ngmres" 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 KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 156 extern PetscErrorCode KSPMonitorCancel(KSP); 157 extern PetscErrorCode KSPGetMonitorContext(KSP,void **); 158 extern PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 159 extern PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 160 161 /* not sure where to put this */ 162 extern PetscErrorCode PCKSPGetKSP(PC,KSP*); 163 extern PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 164 extern PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 165 extern PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 166 extern PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 167 168 extern PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 169 170 extern PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 171 extern PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 172 173 extern PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 174 extern PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 175 extern PetscErrorCode KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 176 extern PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 177 extern PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 178 extern PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 179 180 extern PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 181 extern PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 182 extern PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 183 184 extern PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 185 extern PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 186 extern PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 187 extern PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 188 extern PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 189 190 extern PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 191 extern PetscErrorCode KSPLGMRESSetConstant(KSP); 192 193 extern PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 194 extern PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 195 extern PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 196 197 /*E 198 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 199 200 Level: advanced 201 202 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 203 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 204 205 E*/ 206 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 207 extern const char *KSPGMRESCGSRefinementTypes[]; 208 /*MC 209 KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 210 211 Level: advanced 212 213 Note: Possible unstable, but the fastest to compute 214 215 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 216 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 217 KSPGMRESModifiedGramSchmidtOrthogonalization() 218 M*/ 219 220 /*MC 221 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 222 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 223 poor orthogonality. 224 225 Level: advanced 226 227 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 228 estimate the orthogonality but is more stable. 229 230 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 231 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 232 KSPGMRESModifiedGramSchmidtOrthogonalization() 233 M*/ 234 235 /*MC 236 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 237 238 Level: advanced 239 240 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 241 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 242 243 You should only use this if you absolutely know that the iterative refinement is needed. 244 245 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 246 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 247 KSPGMRESModifiedGramSchmidtOrthogonalization() 248 M*/ 249 250 extern PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 251 extern PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 252 253 extern PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 254 extern PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 255 extern PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 256 257 extern PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 258 extern PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 259 extern PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 260 261 extern PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 262 extern PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 263 extern PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 264 265 extern PetscErrorCode KSPSetFromOptions(KSP); 266 extern PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 267 268 extern PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 269 extern PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 270 extern PetscErrorCode KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *); 271 extern PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 272 extern PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 273 extern PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 274 extern PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 275 extern PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 276 277 extern PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 278 extern PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 279 extern PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 280 extern PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 281 282 extern PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 283 extern PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 284 extern PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 285 extern PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 286 extern PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 287 extern PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 288 289 extern PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 290 extern PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 291 extern PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 292 extern PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 293 294 extern PetscErrorCode KSPView(KSP,PetscViewer); 295 296 extern PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 297 extern PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 298 299 extern PetscErrorCode PCRedundantGetKSP(PC,KSP*); 300 extern PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 301 302 /*E 303 KSPNormType - Norm that is passed in the Krylov convergence 304 test routines. 305 306 Level: advanced 307 308 Each solver only supports a subset of these and some may support different ones 309 depending on left or right preconditioning, see KSPSetPCSide() 310 311 Notes: this must match finclude/petscksp.h 312 313 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 314 KSPSetConvergenceTest(), KSPSetPCSide() 315 E*/ 316 typedef enum {KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 317 extern const char *KSPNormTypes[]; 318 /*MC 319 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 320 possibly save some computation but means the convergence test cannot 321 be based on a norm of a residual etc. 322 323 Level: advanced 324 325 Note: Some Krylov methods need to compute a residual norm and then this option is ignored 326 327 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 328 M*/ 329 330 /*MC 331 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 332 convergence test routine. 333 334 Level: advanced 335 336 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 337 M*/ 338 339 /*MC 340 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 341 convergence test routine. 342 343 Level: advanced 344 345 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 346 M*/ 347 348 /*MC 349 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 350 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 351 352 Level: advanced 353 354 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 355 M*/ 356 357 extern PetscErrorCode KSPSetNormType(KSP,KSPNormType); 358 extern PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 359 extern PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 360 extern PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 361 362 /*E 363 KSPConvergedReason - reason a Krylov method was said to 364 have converged or diverged 365 366 Level: beginner 367 368 Notes: See KSPGetConvergedReason() for explanation of each value 369 370 Developer notes: this must match finclude/petscksp.h 371 372 The string versions of these are KSPConvergedReasons; if you change 373 any of the values here also change them that array of names. 374 375 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 376 E*/ 377 typedef enum {/* converged */ 378 KSP_CONVERGED_RTOL_NORMAL = 1, 379 KSP_CONVERGED_ATOL_NORMAL = 9, 380 KSP_CONVERGED_RTOL = 2, 381 KSP_CONVERGED_ATOL = 3, 382 KSP_CONVERGED_ITS = 4, 383 KSP_CONVERGED_CG_NEG_CURVE = 5, 384 KSP_CONVERGED_CG_CONSTRAINED = 6, 385 KSP_CONVERGED_STEP_LENGTH = 7, 386 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 387 /* diverged */ 388 KSP_DIVERGED_NULL = -2, 389 KSP_DIVERGED_ITS = -3, 390 KSP_DIVERGED_DTOL = -4, 391 KSP_DIVERGED_BREAKDOWN = -5, 392 KSP_DIVERGED_BREAKDOWN_BICG = -6, 393 KSP_DIVERGED_NONSYMMETRIC = -7, 394 KSP_DIVERGED_INDEFINITE_PC = -8, 395 KSP_DIVERGED_NAN = -9, 396 KSP_DIVERGED_INDEFINITE_MAT = -10, 397 398 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 399 extern const char *const*KSPConvergedReasons; 400 401 /*MC 402 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 403 404 Level: beginner 405 406 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 407 for left preconditioning it is the 2-norm of the preconditioned residual, and the 408 2-norm of the residual for right preconditioning 409 410 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 411 412 M*/ 413 414 /*MC 415 KSP_CONVERGED_ATOL - norm(r) <= atol 416 417 Level: beginner 418 419 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 420 for left preconditioning it is the 2-norm of the preconditioned residual, and the 421 2-norm of the residual for right preconditioning 422 423 Level: beginner 424 425 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 426 427 M*/ 428 429 /*MC 430 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 431 432 Level: beginner 433 434 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 435 for left preconditioning it is the 2-norm of the preconditioned residual, and the 436 2-norm of the residual for right preconditioning 437 438 Level: beginner 439 440 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 441 442 M*/ 443 444 /*MC 445 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 446 reached 447 448 Level: beginner 449 450 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 451 452 M*/ 453 454 /*MC 455 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 456 the preconditioner is applied. Also used when the KSPSkipConverged() convergence 457 test routine is set in KSP. 458 459 460 Level: beginner 461 462 463 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 464 465 M*/ 466 467 /*MC 468 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 469 method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 470 preconditioner. 471 472 Level: beginner 473 474 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 475 476 M*/ 477 478 /*MC 479 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 480 method could not continue to enlarge the Krylov space. 481 482 483 Level: beginner 484 485 486 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 487 488 M*/ 489 490 /*MC 491 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 492 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 493 494 Level: beginner 495 496 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 497 498 M*/ 499 500 /*MC 501 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 502 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 503 be positive definite 504 505 Level: beginner 506 507 Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 508 the PCICC preconditioner to generate a positive definite preconditioner 509 510 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 511 512 M*/ 513 514 /*MC 515 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 516 while the KSPSolve() is still running. 517 518 Level: beginner 519 520 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 521 522 M*/ 523 524 extern PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 525 extern PetscErrorCode KSPGetConvergenceContext(KSP,void **); 526 extern PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 527 extern PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 528 extern PetscErrorCode KSPDefaultConvergedDestroy(void *); 529 extern PetscErrorCode KSPDefaultConvergedCreate(void **); 530 extern PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 531 extern PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 532 extern PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 533 extern PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 534 535 extern PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 536 537 /*E 538 KSPCGType - Determines what type of CG to use 539 540 Level: beginner 541 542 .seealso: KSPCGSetType() 543 E*/ 544 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 545 extern const char *KSPCGTypes[]; 546 547 extern PetscErrorCode KSPCGSetType(KSP,KSPCGType); 548 extern PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 549 550 extern PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 551 extern PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 552 extern PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 553 554 extern PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 555 extern PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 556 extern PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 557 558 extern PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 559 extern PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 560 extern PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 561 extern PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 562 extern PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 563 564 extern PetscErrorCode KSPPythonSetType(KSP,const char[]); 565 566 extern PetscErrorCode PCPreSolve(PC,KSP); 567 extern PetscErrorCode PCPostSolve(PC,KSP); 568 569 extern PetscErrorCode KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 570 extern PetscErrorCode KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 571 extern PetscErrorCode KSPMonitorLGDestroy(PetscDrawLG); 572 extern PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 573 extern PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 574 extern PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 575 extern PetscErrorCode KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 576 extern PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 577 extern PetscErrorCode KSPMonitorLGRangeDestroy(PetscDrawLG); 578 579 extern PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 580 extern PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 581 582 /* see src/ksp/ksp/interface/iguess.c */ 583 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 584 585 extern PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 586 extern PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess); 587 extern PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 588 extern PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 589 extern PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 590 extern PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 591 592 extern PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 593 extern PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 594 extern PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 595 596 extern PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 597 extern PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 598 extern PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 599 extern PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 600 extern PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 601 602 extern PetscErrorCode KSPSetDM(KSP,DM); 603 extern PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 604 extern PetscErrorCode KSPGetDM(KSP,DM*); 605 606 PETSC_EXTERN_CXX_END 607 #endif 608