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 PETSCKSP_DLLEXPORT 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 KSPRICHARDSON "richardson" 32 #define KSPCHEBYCHEV "chebychev" 33 #define KSPCG "cg" 34 #define KSPCGNE "cgne" 35 #define KSPSTCG "stcg" 36 #define KSPGMRES "gmres" 37 #define KSPFGMRES "fgmres" 38 #define KSPLGMRES "lgmres" 39 #define KSPTCQMR "tcqmr" 40 #define KSPBCGS "bcgs" 41 #define KSPBCGSL "bcgsl" 42 #define KSPCGS "cgs" 43 #define KSPTFQMR "tfqmr" 44 #define KSPCR "cr" 45 #define KSPLSQR "lsqr" 46 #define KSPPREONLY "preonly" 47 #define KSPQCG "qcg" 48 #define KSPBICG "bicg" 49 #define KSPMINRES "minres" 50 #define KSPSYMMLQ "symmlq" 51 #define KSPLCD "lcd" 52 #define KSPType const char* 53 54 /* Logging support */ 55 extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE; 56 57 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *); 58 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType); 59 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP); 60 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP); 61 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec); 62 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec); 63 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP); 64 65 extern PetscFList KSPList; 66 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]); 67 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void); 68 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 69 70 /*MC 71 KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 72 73 Synopsis: 74 PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP)) 75 76 Not Collective 77 78 Input Parameters: 79 + name_solver - name of a new user-defined solver 80 . path - path (either absolute or relative) the library containing this solver 81 . name_create - name of routine to create method context 82 - routine_create - routine to create method context 83 84 Notes: 85 KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 86 87 If dynamic libraries are used, then the fourth input argument (routine_create) 88 is ignored. 89 90 Sample usage: 91 .vb 92 KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 93 "MySolverCreate",MySolverCreate); 94 .ve 95 96 Then, your solver can be chosen with the procedural interface via 97 $ KSPSetType(ksp,"my_solver") 98 or at runtime via the option 99 $ -ksp_type my_solver 100 101 Level: advanced 102 103 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 104 and others of the form ${any_environmental_variable} occuring in pathname will be 105 replaced with appropriate values. 106 If your function is not being put into a shared library then use KSPRegister() instead 107 108 .keywords: KSP, register 109 110 .seealso: KSPRegisterAll(), KSPRegisterDestroy() 111 112 M*/ 113 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 114 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 115 #else 116 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 117 #endif 118 119 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *); 120 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide); 121 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*); 122 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 123 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 124 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth); 125 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *); 126 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth); 127 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*); 128 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*); 129 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth); 130 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*); 131 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth); 132 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *); 133 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *); 134 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*); 135 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*); 136 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace); 137 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*); 138 139 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC); 140 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*); 141 142 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 143 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP); 144 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **); 145 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 146 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth); 147 148 /* not sure where to put this */ 149 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*); 150 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 151 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 152 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 153 154 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *); 155 156 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *); 157 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *); 158 159 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal); 160 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 161 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 162 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 163 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 164 165 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt); 166 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 167 168 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 169 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 170 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 171 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 172 173 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 174 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 175 176 /*E 177 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 178 179 Level: advanced 180 181 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 182 KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 183 184 E*/ 185 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 186 extern const char *KSPGMRESCGSRefinementTypes[]; 187 /*MC 188 KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 189 190 Level: advanced 191 192 Note: Possible unstable, but the fastest to compute 193 194 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 195 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 196 KSPGMRESModifiedGramSchmidtOrthogonalization() 197 M*/ 198 199 /*MC 200 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 201 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 202 poor orthogonality. 203 204 Level: advanced 205 206 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 207 estimate the orthogonality but is more stable. 208 209 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 210 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 211 KSPGMRESModifiedGramSchmidtOrthogonalization() 212 M*/ 213 214 /*MC 215 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 216 217 Level: advanced 218 219 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 220 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 221 222 You should only use this if you absolutely know that the iterative refinement is needed. 223 224 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 225 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 226 KSPGMRESModifiedGramSchmidtOrthogonalization() 227 M*/ 228 229 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 230 231 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 232 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 233 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 234 235 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 236 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 237 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 238 239 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal); 240 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth); 241 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int); 242 243 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 244 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 245 246 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *); 247 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *); 248 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *); 249 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *); 250 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *); 251 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *); 252 253 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 254 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 255 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 256 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 257 258 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 259 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 260 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 261 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 262 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 263 264 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth); 265 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*); 266 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth); 267 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*); 268 269 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 270 271 /*E 272 KSPNormType - Norm that is passed in the Krylov convergence 273 test routines. 274 275 Level: advanced 276 277 Notes: this must match finclude/petscksp.h 278 279 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 280 KSPSetConvergenceTest() 281 E*/ 282 typedef enum {KSP_NO_NORM = 0,KSP_PRECONDITIONED_NORM = 1,KSP_UNPRECONDITIONED_NORM = 2,KSP_NATURAL_NORM = 3} KSPNormType; 283 extern const char *KSPNormTypes[]; 284 /*MC 285 KSP_NO_NORM - Do not compute a norm during the Krylov process. This will 286 possibly save some computation but means the convergence test cannot 287 be based on a norm of a residual etc. 288 289 Level: advanced 290 291 Note: Some Krylov methods need to compute a residual norm and then this option is ignored 292 293 .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM 294 M*/ 295 296 /*MC 297 KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the 298 convergence test routine. 299 300 Level: advanced 301 302 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest() 303 M*/ 304 305 /*MC 306 KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the 307 convergence test routine. 308 309 Level: advanced 310 311 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest() 312 M*/ 313 314 /*MC 315 KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 316 convergence test routine. 317 318 Level: advanced 319 320 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest() 321 M*/ 322 323 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 324 325 /*E 326 KSPConvergedReason - reason a Krylov method was said to 327 have converged or diverged 328 329 Level: beginner 330 331 Notes: this must match finclude/petscksp.h 332 333 Developer note: The string versions of these are in 334 src/ksp/ksp/interface/itfunc.c called convergedreasons. 335 If these enums are changed you must change those. 336 337 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 338 E*/ 339 typedef enum {/* converged */ 340 KSP_CONVERGED_RTOL = 2, 341 KSP_CONVERGED_ATOL = 3, 342 KSP_CONVERGED_ITS = 4, 343 KSP_CONVERGED_STCG_NEG_CURVE = 5, 344 KSP_CONVERGED_STCG_CONSTRAINED = 6, 345 KSP_CONVERGED_STEP_LENGTH = 7, 346 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 347 /* diverged */ 348 KSP_DIVERGED_NULL = -2, 349 KSP_DIVERGED_ITS = -3, 350 KSP_DIVERGED_DTOL = -4, 351 KSP_DIVERGED_BREAKDOWN = -5, 352 KSP_DIVERGED_BREAKDOWN_BICG = -6, 353 KSP_DIVERGED_NONSYMMETRIC = -7, 354 KSP_DIVERGED_INDEFINITE_PC = -8, 355 KSP_DIVERGED_NAN = -9, 356 KSP_DIVERGED_INDEFINITE_MAT = -10, 357 358 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 359 extern const char **KSPConvergedReasons; 360 361 /*MC 362 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 363 364 Level: beginner 365 366 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 367 for left preconditioning it is the 2-norm of the preconditioned residual, and the 368 2-norm of the residual for right preconditioning 369 370 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 371 372 M*/ 373 374 /*MC 375 KSP_CONVERGED_ATOL - norm(r) <= atol 376 377 Level: beginner 378 379 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 380 for left preconditioning it is the 2-norm of the preconditioned residual, and the 381 2-norm of the residual for right preconditioning 382 383 Level: beginner 384 385 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 386 387 M*/ 388 389 /*MC 390 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 391 392 Level: beginner 393 394 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 395 for left preconditioning it is the 2-norm of the preconditioned residual, and the 396 2-norm of the residual for right preconditioning 397 398 Level: beginner 399 400 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 401 402 M*/ 403 404 /*MC 405 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 406 reached 407 408 Level: beginner 409 410 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 411 412 M*/ 413 414 /*MC 415 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the 416 preconditioner is applied. 417 418 419 Level: beginner 420 421 422 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 423 424 M*/ 425 426 /*MC 427 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 428 method could not continue to enlarge the Krylov space. 429 430 Level: beginner 431 432 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 433 434 M*/ 435 436 /*MC 437 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 438 method could not continue to enlarge the Krylov space. 439 440 441 Level: beginner 442 443 444 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 445 446 M*/ 447 448 /*MC 449 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 450 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 451 452 Level: beginner 453 454 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 455 456 M*/ 457 458 /*MC 459 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 460 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 461 be positive definite 462 463 Level: beginner 464 465 Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 466 the PCICC preconditioner to generate a positive definite preconditioner 467 468 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 469 470 M*/ 471 472 /*MC 473 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 474 while the KSPSolve() is still running. 475 476 Level: beginner 477 478 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 479 480 M*/ 481 482 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *); 483 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 484 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 485 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP); 486 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 487 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 488 489 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 490 491 /*E 492 KSPCGType - Determines what type of CG to use 493 494 Level: beginner 495 496 .seealso: KSPCGSetType() 497 E*/ 498 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 499 extern const char *KSPCGTypes[]; 500 501 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 502 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal); 503 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetQuadratic(KSP,PetscReal*); 504 505 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 506 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 507 508 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 509 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*); 510 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG); 511 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 512 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*); 513 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG); 514 515 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec)); 516 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec)); 517 518 PETSC_EXTERN_CXX_END 519 #endif 520