EPSGetErrorEstimate#
Returns the error estimate associated to the i-th computed eigenpair.
Synopsis#
#include "slepceps.h" 
PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)
Not Collective
Input Parameters#
- eps - the linear eigensolver context 
- i - index of eigenpair 
Output Parameter#
- errest - the error estimate 
Notes#
This is the error estimate used internally by the eigensolver. The actual error bound can be computed with EPSComputeError(). See also the users manual for details.
See Also#
Level#
advanced
Location#
Index of all EPS routines Table of Contents for all manual pages Index of all manual pages
