#include "petscksp.h" PetscErrorCode KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))Logically Collective on ksp
ksp | - iterative context obtained from KSPCreate() | |
converge | - pointer to the function | |
cctx | - context for private data for the convergence routine (may be null) | |
destroy | - a routine for destroying the context (may be null) |
converge (KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)
ksp | - iterative context obtained from KSPCreate() | |
it | - iteration number | |
rnorm | - (estimated) 2-norm of (preconditioned) residual | |
reason | - the reason why it has converged or diverged | |
cctx | - optional convergence context, as set by KSPSetConvergenceTest() |
The default convergence test, KSPConvergedDefault(), aborts if the residual grows to more than 10000 times the initial residual.
The default is a combination of relative and absolute tolerances. The residual value that is tested may be an approximation; routines that need exact values should compute them.
In the default PETSc convergence test, the precise values of reason are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.