#include "petscksp.h" PetscErrorCode KSPConvergedDefault(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)Collective on ksp
ksp | - iterative context | |
n | - iteration number | |
rnorm | - residual norm (may be estimated, depending on the method may be the preconditioned residual norm) | |
ctx | - convergence context which must be created by KSPConvergedDefaultCreate() |
positive | - if the iteration has converged | |
negative | - if the iteration has diverged | |
KSP_CONVERGED_ITERATING | - otherwise. |
rtol = relative tolerance, | ||
abstol = absolute tolerance. | ||
dtol = divergence tolerance, | ||
rnorm_0 is the two norm of the right hand side (or the preconditioned norm, depending on what was set with | - KSPSetNormType(). When initial guess is non-zero you can call KSPConvergedDefaultSetUIRNorm() to use the norm of (b - A*(initial guess)) as the starting point for relative norm convergence testing, that is as rnorm_0 |
Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.
Use KSPSetNormType() (or -ksp_norm_type <none,preconditioned,unpreconditioned,natural>) to change the norm used for computing rnorm
The precise values of reason are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.
This routine is used by KSP by default so the user generally never needs call it directly.
Use KSPSetConvergenceTest() to provide your own test instead of using this one.