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 Parameter

eps - 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 EPSComputeRelativeError(). See also the user's manual for details.

See Also

EPSComputeRelativeError()

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