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