aboutsummaryrefslogtreecommitdiff
path: root/Eigen/src/Core/PermutationMatrix.h
diff options
context:
space:
mode:
Diffstat (limited to 'Eigen/src/Core/PermutationMatrix.h')
-rw-r--r--Eigen/src/Core/PermutationMatrix.h46
1 files changed, 24 insertions, 22 deletions
diff --git a/Eigen/src/Core/PermutationMatrix.h b/Eigen/src/Core/PermutationMatrix.h
index bc29f8142..1297b8413 100644
--- a/Eigen/src/Core/PermutationMatrix.h
+++ b/Eigen/src/Core/PermutationMatrix.h
@@ -105,13 +105,13 @@ class PermutationBase : public EigenBase<Derived>
#endif
/** \returns the number of rows */
- inline Index rows() const { return indices().size(); }
+ inline Index rows() const { return Index(indices().size()); }
/** \returns the number of columns */
- inline Index cols() const { return indices().size(); }
+ inline Index cols() const { return Index(indices().size()); }
/** \returns the size of a side of the respective square matrix, i.e., the number of indices */
- inline Index size() const { return indices().size(); }
+ inline Index size() const { return Index(indices().size()); }
#ifndef EIGEN_PARSED_BY_DOXYGEN
template<typename DenseDerived>
@@ -139,9 +139,9 @@ class PermutationBase : public EigenBase<Derived>
/** Resizes to given size.
*/
- inline void resize(Index size)
+ inline void resize(Index newSize)
{
- indices().resize(size);
+ indices().resize(newSize);
}
/** Sets *this to be the identity permutation matrix */
@@ -153,9 +153,9 @@ class PermutationBase : public EigenBase<Derived>
/** Sets *this to be the identity permutation matrix of given size.
*/
- void setIdentity(Index size)
+ void setIdentity(Index newSize)
{
- resize(size);
+ resize(newSize);
setIdentity();
}
@@ -317,7 +317,7 @@ class PermutationMatrix : public PermutationBase<PermutationMatrix<SizeAtCompile
* array's size.
*/
template<typename Other>
- explicit inline PermutationMatrix(const MatrixBase<Other>& indices) : m_indices(indices)
+ explicit inline PermutationMatrix(const MatrixBase<Other>& a_indices) : m_indices(a_indices)
{}
/** Convert the Transpositions \a tr to a permutation matrix */
@@ -406,12 +406,12 @@ class Map<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime, IndexType>,
typedef typename IndicesType::Scalar Index;
#endif
- inline Map(const Index* indices)
- : m_indices(indices)
+ inline Map(const Index* indicesPtr)
+ : m_indices(indicesPtr)
{}
- inline Map(const Index* indices, Index size)
- : m_indices(indices,size)
+ inline Map(const Index* indicesPtr, Index size)
+ : m_indices(indicesPtr,size)
{}
/** Copies the other permutation into *this */
@@ -490,8 +490,8 @@ class PermutationWrapper : public PermutationBase<PermutationWrapper<_IndicesTyp
typedef typename Traits::IndicesType IndicesType;
#endif
- inline PermutationWrapper(const IndicesType& indices)
- : m_indices(indices)
+ inline PermutationWrapper(const IndicesType& a_indices)
+ : m_indices(a_indices)
{}
/** const version of indices(). */
@@ -541,24 +541,26 @@ struct permut_matrix_product_retval
: public ReturnByValue<permut_matrix_product_retval<PermutationType, MatrixType, Side, Transposed> >
{
typedef typename remove_all<typename MatrixType::Nested>::type MatrixTypeNestedCleaned;
+ typedef typename MatrixType::Index Index;
permut_matrix_product_retval(const PermutationType& perm, const MatrixType& matrix)
: m_permutation(perm), m_matrix(matrix)
{}
- inline int rows() const { return m_matrix.rows(); }
- inline int cols() const { return m_matrix.cols(); }
+ inline Index rows() const { return m_matrix.rows(); }
+ inline Index cols() const { return m_matrix.cols(); }
template<typename Dest> inline void evalTo(Dest& dst) const
{
- const int n = Side==OnTheLeft ? rows() : cols();
-
+ const Index n = Side==OnTheLeft ? rows() : cols();
+ // FIXME we need an is_same for expression that is not sensitive to constness. For instance
+ // is_same_xpr<Block<const Matrix>, Block<Matrix> >::value should be true.
if(is_same<MatrixTypeNestedCleaned,Dest>::value && extract_data(dst) == extract_data(m_matrix))
{
// apply the permutation inplace
Matrix<bool,PermutationType::RowsAtCompileTime,1,0,PermutationType::MaxRowsAtCompileTime> mask(m_permutation.size());
mask.fill(false);
- int r = 0;
+ Index r = 0;
while(r < m_permutation.size())
{
// search for the next seed
@@ -566,10 +568,10 @@ struct permut_matrix_product_retval
if(r>=m_permutation.size())
break;
// we got one, let's follow it until we are back to the seed
- int k0 = r++;
- int kPrev = k0;
+ Index k0 = r++;
+ Index kPrev = k0;
mask.coeffRef(k0) = true;
- for(int k=m_permutation.indices().coeff(k0); k!=k0; k=m_permutation.indices().coeff(k))
+ for(Index k=m_permutation.indices().coeff(k0); k!=k0; k=m_permutation.indices().coeff(k))
{
Block<Dest, Side==OnTheLeft ? 1 : Dest::RowsAtCompileTime, Side==OnTheRight ? 1 : Dest::ColsAtCompileTime>(dst, k)
.swap(Block<Dest, Side==OnTheLeft ? 1 : Dest::RowsAtCompileTime, Side==OnTheRight ? 1 : Dest::ColsAtCompileTime>