aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/IterativeLinearSolvers
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-02-16 13:19:05 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-02-16 13:19:05 +0100
commitaa6c516ec17fb44dff85b1abf3a1b05d58d3bc01 (patch)
treecd0b6ce4d023ff6734beecc167749d00f301c017 /Eigen/src/IterativeLinearSolvers
parentfc202bab398ed9b78ef8452f8e4ef35e233965f6 (diff)
Fix many long to int conversion warnings:
- fix usage of Index (API) versus StorageIndex (when multiple indexes are stored) - use StorageIndex(val) when the input has already been check - use internal::convert_index<StorageIndex>(val) when val is potentially unsafe (directly comes from user input)
Diffstat (limited to 'Eigen/src/IterativeLinearSolvers')
-rw-r--r--Eigen/src/IterativeLinearSolvers/BiCGSTAB.h12
-rw-r--r--Eigen/src/IterativeLinearSolvers/ConjugateGradient.h10
-rw-r--r--Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h16
3 files changed, 19 insertions, 19 deletions
diff --git a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
index a715c7285..e67f09184 100644
--- a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
+++ b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
@@ -27,7 +27,7 @@ namespace internal {
*/
template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
bool bicgstab(const MatrixType& mat, const Rhs& rhs, Dest& x,
- const Preconditioner& precond, int& iters,
+ const Preconditioner& precond, Index& iters,
typename Dest::RealScalar& tol_error)
{
using std::sqrt;
@@ -36,9 +36,9 @@ bool bicgstab(const MatrixType& mat, const Rhs& rhs, Dest& x,
typedef typename Dest::Scalar Scalar;
typedef Matrix<Scalar,Dynamic,1> VectorType;
RealScalar tol = tol_error;
- int maxIters = iters;
+ Index maxIters = iters;
- int n = mat.cols();
+ Index n = mat.cols();
VectorType r = rhs - mat * x;
VectorType r0 = r;
@@ -61,8 +61,8 @@ bool bicgstab(const MatrixType& mat, const Rhs& rhs, Dest& x,
RealScalar tol2 = tol*tol;
RealScalar eps2 = NumTraits<Scalar>::epsilon()*NumTraits<Scalar>::epsilon();
- int i = 0;
- int restarts = 0;
+ Index i = 0;
+ Index restarts = 0;
while ( r.squaredNorm()/rhs_sqnorm > tol2 && i<maxIters )
{
@@ -182,7 +182,7 @@ public:
void _solve_with_guess_impl(const Rhs& b, Dest& x) const
{
bool failed = false;
- for(int j=0; j<b.cols(); ++j)
+ for(Index j=0; j<b.cols(); ++j)
{
m_iterations = Base::maxIterations();
m_error = Base::m_tolerance;
diff --git a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
index d8bc13f5f..a799c3ef5 100644
--- a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
@@ -26,7 +26,7 @@ namespace internal {
template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
EIGEN_DONT_INLINE
void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
- const Preconditioner& precond, int& iters,
+ const Preconditioner& precond, Index& iters,
typename Dest::RealScalar& tol_error)
{
using std::sqrt;
@@ -36,9 +36,9 @@ void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
typedef Matrix<Scalar,Dynamic,1> VectorType;
RealScalar tol = tol_error;
- int maxIters = iters;
+ Index maxIters = iters;
- int n = mat.cols();
+ Index n = mat.cols();
VectorType residual = rhs - mat * x; //initial residual
@@ -64,7 +64,7 @@ void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
VectorType z(n), tmp(n);
RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
- int i = 0;
+ Index i = 0;
while(i < maxIters)
{
tmp.noalias() = mat * p; // the bottleneck of the algorithm
@@ -190,7 +190,7 @@ public:
m_iterations = Base::maxIterations();
m_error = Base::m_tolerance;
- for(int j=0; j<b.cols(); ++j)
+ for(Index j=0; j<b.cols(); ++j)
{
m_iterations = Base::maxIterations();
m_error = Base::m_tolerance;
diff --git a/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h b/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h
index 2afd80f4d..38334028a 100644
--- a/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h
+++ b/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h
@@ -145,7 +145,7 @@ public:
* It is either the value setted by setMaxIterations or, by default,
* twice the number of columns of the matrix.
*/
- int maxIterations() const
+ Index maxIterations() const
{
return (m_maxIterations<0) ? 2*mp_matrix.cols() : m_maxIterations;
}
@@ -153,14 +153,14 @@ public:
/** Sets the max number of iterations.
* Default is twice the number of columns of the matrix.
*/
- Derived& setMaxIterations(int maxIters)
+ Derived& setMaxIterations(Index maxIters)
{
m_maxIterations = maxIters;
return derived();
}
/** \returns the number of iterations performed during the last solve */
- int iterations() const
+ Index iterations() const
{
eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
return m_iterations;
@@ -200,11 +200,11 @@ public:
{
eigen_assert(rows()==b.rows());
- int rhsCols = b.cols();
- int size = b.rows();
+ Index rhsCols = b.cols();
+ Index size = b.rows();
Eigen::Matrix<DestScalar,Dynamic,1> tb(size);
Eigen::Matrix<DestScalar,Dynamic,1> tx(size);
- for(int k=0; k<rhsCols; ++k)
+ for(Index k=0; k<rhsCols; ++k)
{
tb = b.col(k);
tx = derived().solve(tb);
@@ -233,11 +233,11 @@ protected:
Ref<const MatrixType> mp_matrix;
Preconditioner m_preconditioner;
- int m_maxIterations;
+ Index m_maxIterations;
RealScalar m_tolerance;
mutable RealScalar m_error;
- mutable int m_iterations;
+ mutable Index m_iterations;
mutable ComputationInfo m_info;
mutable bool m_analysisIsOk, m_factorizationIsOk;
};