mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-09-17 20:03:17 +08:00
set the default number of iteration to the size of the problem
This commit is contained in:
parent
15ea999f84
commit
bdee0c9baa
@ -125,8 +125,8 @@ struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
|
||||
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
|
||||
*
|
||||
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
|
||||
* and setTolerance() methods. The default are 1000 max iterations and NumTraits<Scalar>::epsilon()
|
||||
* for the tolerance.
|
||||
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
|
||||
* and NumTraits<Scalar>::epsilon() for the tolerance.
|
||||
*
|
||||
* This class can be used as the direct solver classes. Here is a typical usage example:
|
||||
* \code
|
||||
@ -218,7 +218,7 @@ public:
|
||||
{
|
||||
for(int j=0; j<b.cols(); ++j)
|
||||
{
|
||||
m_iterations = Base::m_maxIterations;
|
||||
m_iterations = Base::maxIterations();
|
||||
m_error = Base::m_tolerance;
|
||||
|
||||
typename Dest::ColXpr xj(x,j);
|
||||
|
@ -111,8 +111,8 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
|
||||
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
|
||||
*
|
||||
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
|
||||
* and setTolerance() methods. The default are 1000 max iterations and NumTraits<Scalar>::epsilon()
|
||||
* for the tolerance.
|
||||
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
|
||||
* and NumTraits<Scalar>::epsilon() for the tolerance.
|
||||
*
|
||||
* This class can be used as the direct solver classes. Here is a typical usage example:
|
||||
* \code
|
||||
@ -206,12 +206,12 @@ public:
|
||||
template<typename Rhs,typename Dest>
|
||||
void _solveWithGuess(const Rhs& b, Dest& x) const
|
||||
{
|
||||
m_iterations = Base::m_maxIterations;
|
||||
m_iterations = Base::maxIterations();
|
||||
m_error = Base::m_tolerance;
|
||||
|
||||
for(int j=0; j<b.cols(); ++j)
|
||||
{
|
||||
m_iterations = Base::m_maxIterations;
|
||||
m_iterations = Base::maxIterations();
|
||||
m_error = Base::m_tolerance;
|
||||
|
||||
typename Dest::ColXpr xj(x,j);
|
||||
|
@ -91,9 +91,9 @@ public:
|
||||
}
|
||||
|
||||
/** \internal */
|
||||
Index rows() const { return mp_matrix->rows(); }
|
||||
Index rows() const { return mp_matrix ? mp_matrix->rows() : 0; }
|
||||
/** \internal */
|
||||
Index cols() const { return mp_matrix->cols(); }
|
||||
Index cols() const { return mp_matrix ? mp_matrix->cols() : 0; }
|
||||
|
||||
/** \returns the tolerance threshold used by the stopping criteria */
|
||||
RealScalar tolerance() const { return m_tolerance; }
|
||||
@ -112,7 +112,10 @@ public:
|
||||
const Preconditioner& preconditioner() const { return m_preconditioner; }
|
||||
|
||||
/** \returns the max number of iterations */
|
||||
int maxIterations() const { return m_maxIterations; }
|
||||
int maxIterations() const
|
||||
{
|
||||
return (mp_matrix && m_maxIterations<0) ? mp_matrix->cols() : m_maxIterations;
|
||||
}
|
||||
|
||||
/** Sets the max number of iterations */
|
||||
Derived& setMaxIterations(int maxIters)
|
||||
@ -191,7 +194,7 @@ protected:
|
||||
void init()
|
||||
{
|
||||
m_isInitialized = false;
|
||||
m_maxIterations = 1000;
|
||||
m_maxIterations = -1;
|
||||
m_tolerance = NumTraits<Scalar>::epsilon();
|
||||
}
|
||||
const MatrixType* mp_matrix;
|
||||
|
Loading…
x
Reference in New Issue
Block a user