Clarify documentation of the tolerance and error returned in iterative solvers

This commit is contained in:
Gael Guennebaud 2015-06-25 13:51:13 +02:00
parent 84264ceebc
commit 973b0a90db
3 changed files with 15 additions and 3 deletions

View File

@ -136,6 +136,8 @@ struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
* and NumTraits<Scalar>::epsilon() for the tolerance.
*
* The tolerance is the relative residual error: |Ax-b|/|b|
*
* This class can be used as the direct solver classes. Here is a typical usage example:
* \include BiCGSTAB_simple.cpp
*

View File

@ -121,6 +121,8 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
* and NumTraits<Scalar>::epsilon() for the tolerance.
*
* The tolerance is the relative residual error: |Ax-b|/|b|
*
* This class can be used as the direct solver classes. Here is a typical usage example:
\code
int n = 10000;

View File

@ -126,10 +126,16 @@ public:
/** \internal */
Index cols() const { return mp_matrix.cols(); }
/** \returns the tolerance threshold used by the stopping criteria */
/** \returns the tolerance threshold used by the stopping criteria.
* \sa setTolerance()
*/
RealScalar tolerance() const { return m_tolerance; }
/** Sets the tolerance threshold used by the stopping criteria */
/** Sets the tolerance threshold used by the stopping criteria.
*
* This value is used as an upper bound to the relative residual error: |Ax-b|/|b|.
* The default value is the machine precision given by NumTraits<Scalar>::epsilon()
*/
Derived& setTolerance(const RealScalar& tolerance)
{
m_tolerance = tolerance;
@ -167,7 +173,9 @@ public:
return m_iterations;
}
/** \returns the tolerance error reached during the last solve */
/** \returns the tolerance error reached during the last solve.
* It is a close approximation of the true relative residual error |Ax-b|/|b|.
*/
RealScalar error() const
{
eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");