KSPSetConvergenceTest#
Sets the function to be used to determine convergence.
Synopsis#
#include "petscksp.h"
#include "petscmat.h"
PetscErrorCode KSPSetConvergenceTest(KSP ksp, PetscErrorCode (*converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void *ctx, PetscErrorCode (*destroy)(void *ctx))
Logically Collective
Input Parameters#
ksp - iterative context obtained from
KSPCreate()
converge - pointer to the function
ctx - context for private data for the convergence routine (may be
NULL
)destroy - a routine for destroying the context (may be
NULL
)
Calling sequence of converge
#
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
ctx - optional convergence context, as set by
KSPSetConvergenceTest()
Calling sequence of destroy
#
ctx - the context
Notes#
Must be called after the KSP
type has been set so put this after
a call to KSPSetType()
, or KSPSetFromOptions()
.
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.
See Also#
KSP: Linear System Solvers, KSP
, KSPConvergedDefault()
, KSPGetConvergenceContext()
, KSPSetTolerances()
, KSPGetConvergenceTest()
, KSPGetAndClearConvergenceTest()
Level#
advanced
Location#
Examples#
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages