#include "slepceps.h" PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))Logically Collective on EPS
eps | - eigensolver context obtained from EPSCreate() | |
func | - a pointer to the convergence test function | |
ctx | - [optional] context for private data for the convergence routine | |
destroy | - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran) |
func(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
eps | - eigensolver context obtained from EPSCreate() | |
eigr | - real part of the eigenvalue | |
eigi | - imaginary part of the eigenvalue | |
res | - residual norm associated to the eigenpair | |
errest | - (output) computed error estimate | |
ctx | - optional context, as set by EPSSetConvergenceTest() |
Location: src/eps/interface/epsopts.c
Index of all EPS routines
Table of Contents for all manual pages
Index of all manual pages