PEPSetConvergenceTestFunction#
Sets a function to compute the error estimate used in the convergence test.
Synopsis#
#include "slepcpep.h"
PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PEPConvergenceTestFn *conv,void *ctx,PetscCtxDestroyFn *destroy)
Logically Collective
Input Parameters#
pep - eigensolver context obtained from PEPCreate()
conv - convergence test function, see PEPConvergenceTestFn for the calling sequence
ctx - context for private data for the convergence routine (may be NULL)
destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
Note#
If the error estimate returned by the convergence test function is less than the tolerance, then the eigenvalue is accepted as converged.
See Also#
PEPSetConvergenceTest(), PEPSetTolerances()
Level#
advanced
Location#
Index of all PEP routines Table of Contents for all manual pages Index of all manual pages