slepc-3.10.0 2018-09-18
Report Typos and Errors

EPSSetConvergenceTestFunction

Sets a function to compute the error estimate used in the convergence test.

Synopsis

#include "slepceps.h" 
PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
Logically Collective on EPS

Input Parameters

eps  - eigensolver context obtained from EPSCreate()
func  - a pointer to the convergence test 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 func

  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 EPSSetConvergenceTestFunction()

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

EPSSetConvergenceTest(), EPSSetTolerances()

Location: src/eps/interface/epsopts.c
Index of all EPS routines
Table of Contents for all manual pages
Index of all manual pages