aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2012-02-03 23:18:26 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2012-02-03 23:18:26 +0100
commitfe85b7ebc69d98d65c0d70b189416c384408b6f7 (patch)
tree2f0e25b681150fc2e86ed3c802db47810bdd7a94 /Eigen/src/Core
parentbc7b251cd9f520097564ce611f314fb978f4f744 (diff)
fix several const qualifier issues: double ones, meaningless ones, some missing ones, etc.
(note that const qualifiers are set by internall::nested)
Diffstat (limited to 'Eigen/src/Core')
-rw-r--r--Eigen/src/Core/ArrayBase.h2
-rw-r--r--Eigen/src/Core/ArrayWrapper.h16
-rw-r--r--Eigen/src/Core/Block.h2
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h4
-rw-r--r--Eigen/src/Core/CwiseUnaryOp.h2
-rw-r--r--Eigen/src/Core/CwiseUnaryView.h2
-rw-r--r--Eigen/src/Core/DenseBase.h3
-rw-r--r--Eigen/src/Core/Diagonal.h3
-rw-r--r--Eigen/src/Core/DiagonalMatrix.h6
-rw-r--r--Eigen/src/Core/DiagonalProduct.h4
-rw-r--r--Eigen/src/Core/Fuzzy.h4
-rw-r--r--Eigen/src/Core/GeneralProduct.h12
-rw-r--r--Eigen/src/Core/IO.h2
-rw-r--r--Eigen/src/Core/MatrixBase.h2
-rw-r--r--Eigen/src/Core/PermutationMatrix.h4
-rw-r--r--Eigen/src/Core/ProductBase.h4
-rw-r--r--Eigen/src/Core/Replicate.h2
-rw-r--r--Eigen/src/Core/Reverse.h2
-rw-r--r--Eigen/src/Core/Select.h6
-rw-r--r--Eigen/src/Core/SelfAdjointView.h4
-rw-r--r--Eigen/src/Core/SolveTriangular.h4
-rw-r--r--Eigen/src/Core/Swap.h9
-rw-r--r--Eigen/src/Core/Transpose.h6
-rw-r--r--Eigen/src/Core/Transpositions.h2
-rw-r--r--Eigen/src/Core/TriangularMatrix.h16
-rw-r--r--Eigen/src/Core/VectorwiseOp.h2
-rw-r--r--Eigen/src/Core/products/CoeffBasedProduct.h4
-rw-r--r--Eigen/src/Core/products/GeneralMatrixMatrix.h4
-rw-r--r--Eigen/src/Core/products/GeneralMatrixMatrixTriangular.h4
-rw-r--r--Eigen/src/Core/products/SelfadjointMatrixMatrix.h4
-rw-r--r--Eigen/src/Core/products/SelfadjointMatrixVector.h4
-rw-r--r--Eigen/src/Core/products/SelfadjointProduct.h4
-rw-r--r--Eigen/src/Core/products/SelfadjointRank2Update.h4
-rw-r--r--Eigen/src/Core/products/TriangularMatrixMatrix.h4
-rw-r--r--Eigen/src/Core/products/TriangularMatrixVector.h4
-rw-r--r--Eigen/src/Core/util/BlasUtil.h12
-rw-r--r--Eigen/src/Core/util/XprHelper.h35
37 files changed, 109 insertions, 99 deletions
diff --git a/Eigen/src/Core/ArrayBase.h b/Eigen/src/Core/ArrayBase.h
index 2ddd9dea5..e16926ba4 100644
--- a/Eigen/src/Core/ArrayBase.h
+++ b/Eigen/src/Core/ArrayBase.h
@@ -159,7 +159,7 @@ template<typename Derived> class ArrayBase
/** \returns an \link MatrixBase Matrix \endlink expression of this array
* \sa MatrixBase::array() */
MatrixWrapper<Derived> matrix() { return derived(); }
- const MatrixWrapper<Derived> matrix() const { return derived(); }
+ const MatrixWrapper<const Derived> matrix() const { return derived(); }
// template<typename Dest>
// inline void evalTo(Dest& dst) const { dst = matrix(); }
diff --git a/Eigen/src/Core/ArrayWrapper.h b/Eigen/src/Core/ArrayWrapper.h
index db1518019..fd6e52a77 100644
--- a/Eigen/src/Core/ArrayWrapper.h
+++ b/Eigen/src/Core/ArrayWrapper.h
@@ -61,7 +61,7 @@ class ArrayWrapper : public ArrayBase<ArrayWrapper<ExpressionType> >
typedef typename internal::nested<ExpressionType>::type NestedExpressionType;
- inline ArrayWrapper(const ExpressionType& matrix) : m_expression(matrix) {}
+ inline ArrayWrapper(ExpressionType& matrix) : m_expression(matrix) {}
inline Index rows() const { return m_expression.rows(); }
inline Index cols() const { return m_expression.cols(); }
@@ -71,7 +71,7 @@ class ArrayWrapper : public ArrayBase<ArrayWrapper<ExpressionType> >
inline ScalarWithConstIfNotLvalue* data() { return m_expression.data(); }
inline const Scalar* data() const { return m_expression.data(); }
- inline const CoeffReturnType coeff(Index row, Index col) const
+ inline CoeffReturnType coeff(Index row, Index col) const
{
return m_expression.coeff(row, col);
}
@@ -86,7 +86,7 @@ class ArrayWrapper : public ArrayBase<ArrayWrapper<ExpressionType> >
return m_expression.const_cast_derived().coeffRef(row, col);
}
- inline const CoeffReturnType coeff(Index index) const
+ inline CoeffReturnType coeff(Index index) const
{
return m_expression.coeff(index);
}
@@ -135,7 +135,7 @@ class ArrayWrapper : public ArrayBase<ArrayWrapper<ExpressionType> >
}
protected:
- const NestedExpressionType m_expression;
+ NestedExpressionType m_expression;
};
/** \class MatrixWrapper
@@ -174,7 +174,7 @@ class MatrixWrapper : public MatrixBase<MatrixWrapper<ExpressionType> >
typedef typename internal::nested<ExpressionType>::type NestedExpressionType;
- inline MatrixWrapper(const ExpressionType& matrix) : m_expression(matrix) {}
+ inline MatrixWrapper(ExpressionType& matrix) : m_expression(matrix) {}
inline Index rows() const { return m_expression.rows(); }
inline Index cols() const { return m_expression.cols(); }
@@ -184,7 +184,7 @@ class MatrixWrapper : public MatrixBase<MatrixWrapper<ExpressionType> >
inline ScalarWithConstIfNotLvalue* data() { return m_expression.data(); }
inline const Scalar* data() const { return m_expression.data(); }
- inline const CoeffReturnType coeff(Index row, Index col) const
+ inline CoeffReturnType coeff(Index row, Index col) const
{
return m_expression.coeff(row, col);
}
@@ -199,7 +199,7 @@ class MatrixWrapper : public MatrixBase<MatrixWrapper<ExpressionType> >
return m_expression.derived().coeffRef(row, col);
}
- inline const CoeffReturnType coeff(Index index) const
+ inline CoeffReturnType coeff(Index index) const
{
return m_expression.coeff(index);
}
@@ -245,7 +245,7 @@ class MatrixWrapper : public MatrixBase<MatrixWrapper<ExpressionType> >
}
protected:
- const NestedExpressionType m_expression;
+ NestedExpressionType m_expression;
};
#endif // EIGEN_ARRAYWRAPPER_H
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index 427e6d5c9..dfbf21b73 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -361,7 +361,7 @@ class Block<XprType,BlockRows,BlockCols, InnerPanel,true>
: m_xpr.innerStride();
}
- const typename XprType::Nested m_xpr;
+ typename XprType::Nested m_xpr;
int m_outerStride;
};
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index 7386b2e18..6c891b5da 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -167,8 +167,8 @@ class CwiseBinaryOp : internal::no_assignment_operator,
const BinaryOp& functor() const { return m_functor; }
protected:
- const LhsNested m_lhs;
- const RhsNested m_rhs;
+ LhsNested m_lhs;
+ RhsNested m_rhs;
const BinaryOp m_functor;
};
diff --git a/Eigen/src/Core/CwiseUnaryOp.h b/Eigen/src/Core/CwiseUnaryOp.h
index 958571d64..30abc3096 100644
--- a/Eigen/src/Core/CwiseUnaryOp.h
+++ b/Eigen/src/Core/CwiseUnaryOp.h
@@ -95,7 +95,7 @@ class CwiseUnaryOp : internal::no_assignment_operator,
nestedExpression() { return m_xpr.const_cast_derived(); }
protected:
- const typename XprType::Nested m_xpr;
+ typename XprType::Nested m_xpr;
const UnaryOp m_functor;
};
diff --git a/Eigen/src/Core/CwiseUnaryView.h b/Eigen/src/Core/CwiseUnaryView.h
index d24ef0373..462e4e001 100644
--- a/Eigen/src/Core/CwiseUnaryView.h
+++ b/Eigen/src/Core/CwiseUnaryView.h
@@ -97,7 +97,7 @@ class CwiseUnaryView : internal::no_assignment_operator,
protected:
// FIXME changed from MatrixType::Nested because of a weird compilation error with sun CC
- const typename internal::nested<MatrixType>::type m_matrix;
+ typename internal::nested<MatrixType>::type m_matrix;
ViewOp m_functor;
};
diff --git a/Eigen/src/Core/DenseBase.h b/Eigen/src/Core/DenseBase.h
index 838fa4030..10131b2b7 100644
--- a/Eigen/src/Core/DenseBase.h
+++ b/Eigen/src/Core/DenseBase.h
@@ -376,12 +376,13 @@ template<typename Derived> class DenseBase
inline Derived& operator*=(const Scalar& other);
inline Derived& operator/=(const Scalar& other);
+ typedef typename internal::add_const_on_value_type<typename internal::eval<Derived>::type>::type EvalReturnType;
/** \returns the matrix or vector obtained by evaluating this expression.
*
* Notice that in the case of a plain matrix or vector (not an expression) this function just returns
* a const reference, in order to avoid a useless copy.
*/
- EIGEN_STRONG_INLINE const typename internal::eval<Derived>::type eval() const
+ EIGEN_STRONG_INLINE EvalReturnType eval() const
{
// Even though MSVC does not honor strong inlining when the return type
// is a dynamic matrix, we desperately need strong inlining for fixed
diff --git a/Eigen/src/Core/Diagonal.h b/Eigen/src/Core/Diagonal.h
index 8d9a423ef..58032dbcd 100644
--- a/Eigen/src/Core/Diagonal.h
+++ b/Eigen/src/Core/Diagonal.h
@@ -2,6 +2,7 @@
// for linear algebra.
//
// Copyright (C) 2007-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
+// Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
@@ -154,7 +155,7 @@ template<typename MatrixType, int DiagIndex> class Diagonal
}
protected:
- const typename MatrixType::Nested m_matrix;
+ typename MatrixType::Nested m_matrix;
const internal::variable_if_dynamic<Index, DiagIndex> m_index;
private:
diff --git a/Eigen/src/Core/DiagonalMatrix.h b/Eigen/src/Core/DiagonalMatrix.h
index f41a74bfa..1b7c90891 100644
--- a/Eigen/src/Core/DiagonalMatrix.h
+++ b/Eigen/src/Core/DiagonalMatrix.h
@@ -72,7 +72,7 @@ class DiagonalBase : public EigenBase<Derived>
const DiagonalProduct<MatrixDerived, Derived, OnTheLeft>
operator*(const MatrixBase<MatrixDerived> &matrix) const;
- inline const DiagonalWrapper<CwiseUnaryOp<internal::scalar_inverse_op<Scalar>, const DiagonalVectorType> >
+ inline const DiagonalWrapper<const CwiseUnaryOp<internal::scalar_inverse_op<Scalar>, const DiagonalVectorType> >
inverse() const
{
return diagonal().cwiseInverse();
@@ -251,13 +251,13 @@ class DiagonalWrapper
#endif
/** Constructor from expression of diagonal coefficients to wrap. */
- inline DiagonalWrapper(const DiagonalVectorType& diagonal) : m_diagonal(diagonal) {}
+ inline DiagonalWrapper(DiagonalVectorType& diagonal) : m_diagonal(diagonal) {}
/** \returns a const reference to the wrapped expression of diagonal coefficients. */
const DiagonalVectorType& diagonal() const { return m_diagonal; }
protected:
- const typename DiagonalVectorType::Nested m_diagonal;
+ typename DiagonalVectorType::Nested m_diagonal;
};
/** \returns a pseudo-expression of a diagonal matrix with *this as vector of diagonal coefficients
diff --git a/Eigen/src/Core/DiagonalProduct.h b/Eigen/src/Core/DiagonalProduct.h
index de0c6ed11..278e75a02 100644
--- a/Eigen/src/Core/DiagonalProduct.h
+++ b/Eigen/src/Core/DiagonalProduct.h
@@ -107,8 +107,8 @@ class DiagonalProduct : internal::no_assignment_operator,
m_diagonal.diagonal().template packet<DiagonalVectorPacketLoadMode>(id));
}
- const typename MatrixType::Nested m_matrix;
- const typename DiagonalType::Nested m_diagonal;
+ typename MatrixType::Nested m_matrix;
+ typename DiagonalType::Nested m_diagonal;
};
/** \returns the diagonal matrix product of \c *this by the diagonal matrix \a diagonal.
diff --git a/Eigen/src/Core/Fuzzy.h b/Eigen/src/Core/Fuzzy.h
index d266eed0a..6a802bb79 100644
--- a/Eigen/src/Core/Fuzzy.h
+++ b/Eigen/src/Core/Fuzzy.h
@@ -35,8 +35,8 @@ struct isApprox_selector
static bool run(const Derived& x, const OtherDerived& y, typename Derived::RealScalar prec)
{
using std::min;
- const typename internal::nested<Derived,2>::type nested(x);
- const typename internal::nested<OtherDerived,2>::type otherNested(y);
+ typename internal::nested<Derived,2>::type nested(x);
+ typename internal::nested<OtherDerived,2>::type otherNested(y);
return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * (min)(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
}
};
diff --git a/Eigen/src/Core/GeneralProduct.h b/Eigen/src/Core/GeneralProduct.h
index 61bc02a76..d6e3c7cb6 100644
--- a/Eigen/src/Core/GeneralProduct.h
+++ b/Eigen/src/Core/GeneralProduct.h
@@ -410,8 +410,8 @@ template<> struct gemv_selector<OnTheRight,ColMajor,true>
typedef typename ProductType::RhsBlasTraits RhsBlasTraits;
typedef Map<Matrix<ResScalar,Dynamic,1>, Aligned> MappedDest;
- const ActualLhsType actualLhs = LhsBlasTraits::extract(prod.lhs());
- const ActualRhsType actualRhs = RhsBlasTraits::extract(prod.rhs());
+ ActualLhsType actualLhs = LhsBlasTraits::extract(prod.lhs());
+ ActualRhsType actualRhs = RhsBlasTraits::extract(prod.rhs());
ResScalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs())
* RhsBlasTraits::extractScalarFactor(prod.rhs());
@@ -452,7 +452,7 @@ template<> struct gemv_selector<OnTheRight,ColMajor,true>
general_matrix_vector_product
<Index,LhsScalar,ColMajor,LhsBlasTraits::NeedToConjugate,RhsScalar,RhsBlasTraits::NeedToConjugate>::run(
actualLhs.rows(), actualLhs.cols(),
- &actualLhs.coeffRef(0,0), actualLhs.outerStride(),
+ actualLhs.data(), actualLhs.outerStride(),
actualRhs.data(), actualRhs.innerStride(),
actualDestPtr, 1,
compatibleAlpha);
@@ -511,9 +511,9 @@ template<> struct gemv_selector<OnTheRight,RowMajor,true>
general_matrix_vector_product
<Index,LhsScalar,RowMajor,LhsBlasTraits::NeedToConjugate,RhsScalar,RhsBlasTraits::NeedToConjugate>::run(
actualLhs.rows(), actualLhs.cols(),
- &actualLhs.coeffRef(0,0), actualLhs.outerStride(),
+ actualLhs.data(), actualLhs.outerStride(),
actualRhsPtr, 1,
- &dest.coeffRef(0,0), dest.innerStride(),
+ dest.data(), dest.innerStride(),
actualAlpha);
}
};
@@ -558,7 +558,7 @@ template<> struct gemv_selector<OnTheRight,RowMajor,false>
*/
template<typename Derived>
template<typename OtherDerived>
-inline const typename ProductReturnType<Derived,OtherDerived>::Type
+inline const typename ProductReturnType<Derived, OtherDerived>::Type
MatrixBase<Derived>::operator*(const MatrixBase<OtherDerived> &other) const
{
// A note regarding the function declaration: In MSVC, this function will sometimes
diff --git a/Eigen/src/Core/IO.h b/Eigen/src/Core/IO.h
index f3cfcdbf4..e2a692a80 100644
--- a/Eigen/src/Core/IO.h
+++ b/Eigen/src/Core/IO.h
@@ -171,7 +171,7 @@ std::ostream & print_matrix(std::ostream & s, const Derived& _m, const IOFormat&
return s;
}
- const typename Derived::Nested m = _m;
+ typename Derived::Nested m = _m;
typedef typename Derived::Scalar Scalar;
typedef typename Derived::Index Index;
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index 11d9970d5..9298f35eb 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -330,7 +330,7 @@ template<typename Derived> class MatrixBase
/** \returns an \link ArrayBase Array \endlink expression of this matrix
* \sa ArrayBase::matrix() */
ArrayWrapper<Derived> array() { return derived(); }
- const ArrayWrapper<Derived> array() const { return derived(); }
+ const ArrayWrapper<const Derived> array() const { return derived(); }
/////////// LU module ///////////
diff --git a/Eigen/src/Core/PermutationMatrix.h b/Eigen/src/Core/PermutationMatrix.h
index a064e053e..a56c4217a 100644
--- a/Eigen/src/Core/PermutationMatrix.h
+++ b/Eigen/src/Core/PermutationMatrix.h
@@ -511,7 +511,7 @@ class PermutationWrapper : public PermutationBase<PermutationWrapper<_IndicesTyp
protected:
- const typename IndicesType::Nested m_indices;
+ typename IndicesType::Nested m_indices;
};
/** \returns the matrix with the permutation applied to the columns.
@@ -608,7 +608,7 @@ struct permut_matrix_product_retval
protected:
const PermutationType& m_permutation;
- const typename MatrixType::Nested m_matrix;
+ typename MatrixType::Nested m_matrix;
};
/* Template partial specialization for transposed/inverse permutations */
diff --git a/Eigen/src/Core/ProductBase.h b/Eigen/src/Core/ProductBase.h
index e7e1c84c6..338f94111 100644
--- a/Eigen/src/Core/ProductBase.h
+++ b/Eigen/src/Core/ProductBase.h
@@ -179,8 +179,8 @@ class ProductBase : public MatrixBase<Derived>
protected:
- const LhsNested m_lhs;
- const RhsNested m_rhs;
+ LhsNested m_lhs;
+ RhsNested m_rhs;
mutable PlainObject m_result;
};
diff --git a/Eigen/src/Core/Replicate.h b/Eigen/src/Core/Replicate.h
index 2acb8ab77..0189c4072 100644
--- a/Eigen/src/Core/Replicate.h
+++ b/Eigen/src/Core/Replicate.h
@@ -128,7 +128,7 @@ template<typename MatrixType,int RowFactor,int ColFactor> class Replicate
}
protected:
- const typename MatrixType::Nested m_matrix;
+ typename MatrixType::Nested m_matrix;
const internal::variable_if_dynamic<Index, RowFactor> m_rowFactor;
const internal::variable_if_dynamic<Index, ColFactor> m_colFactor;
};
diff --git a/Eigen/src/Core/Reverse.h b/Eigen/src/Core/Reverse.h
index cef3114c7..a2fffe40b 100644
--- a/Eigen/src/Core/Reverse.h
+++ b/Eigen/src/Core/Reverse.h
@@ -190,7 +190,7 @@ template<typename MatrixType, int Direction> class Reverse
}
protected:
- const typename MatrixType::Nested m_matrix;
+ typename MatrixType::Nested m_matrix;
};
/** \returns an expression of the reverse of *this.
diff --git a/Eigen/src/Core/Select.h b/Eigen/src/Core/Select.h
index 87a071fc7..832ca2735 100644
--- a/Eigen/src/Core/Select.h
+++ b/Eigen/src/Core/Select.h
@@ -117,9 +117,9 @@ class Select : internal::no_assignment_operator,
}
protected:
- const typename ConditionMatrixType::Nested m_condition;
- const typename ThenMatrixType::Nested m_then;
- const typename ElseMatrixType::Nested m_else;
+ typename ConditionMatrixType::Nested m_condition;
+ typename ThenMatrixType::Nested m_then;
+ typename ElseMatrixType::Nested m_else;
};
diff --git a/Eigen/src/Core/SelfAdjointView.h b/Eigen/src/Core/SelfAdjointView.h
index 352f3536c..402c1f74f 100644
--- a/Eigen/src/Core/SelfAdjointView.h
+++ b/Eigen/src/Core/SelfAdjointView.h
@@ -82,7 +82,7 @@ template<typename MatrixType, unsigned int UpLo> class SelfAdjointView
};
typedef typename MatrixType::PlainObject PlainObject;
- inline SelfAdjointView(const MatrixType& matrix) : m_matrix(matrix)
+ inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix)
{}
inline Index rows() const { return m_matrix.rows(); }
@@ -199,7 +199,7 @@ template<typename MatrixType, unsigned int UpLo> class SelfAdjointView
#endif
protected:
- const MatrixTypeNested m_matrix;
+ MatrixTypeNested m_matrix;
};
diff --git a/Eigen/src/Core/SolveTriangular.h b/Eigen/src/Core/SolveTriangular.h
index 71e129c7f..bd3cab328 100644
--- a/Eigen/src/Core/SolveTriangular.h
+++ b/Eigen/src/Core/SolveTriangular.h
@@ -100,7 +100,7 @@ struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,Dynamic>
typedef typename LhsProductTraits::DirectLinearAccessType ActualLhsType;
static void run(const Lhs& lhs, Rhs& rhs)
{
- const ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
+ typename internal::add_const_on_value_type<ActualLhsType>::type actualLhs = LhsProductTraits::extract(lhs);
triangular_solve_matrix<Scalar,Index,Side,Mode,LhsProductTraits::NeedToConjugate,(int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor,
(Rhs::Flags&RowMajorBit) ? RowMajor : ColMajor>
::run(lhs.rows(), Side==OnTheLeft? rhs.cols() : rhs.rows(), &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &rhs.coeffRef(0,0), rhs.outerStride());
@@ -255,7 +255,7 @@ template<int Side, typename TriangularType, typename Rhs> struct triangular_solv
protected:
const TriangularType& m_triangularMatrix;
- const typename Rhs::Nested m_rhs;
+ typename Rhs::Nested m_rhs;
};
} // namespace internal
diff --git a/Eigen/src/Core/Swap.h b/Eigen/src/Core/Swap.h
index 5fdd36e3b..4bf23c8c1 100644
--- a/Eigen/src/Core/Swap.h
+++ b/Eigen/src/Core/Swap.h
@@ -52,6 +52,15 @@ template<typename ExpressionType> class SwapWrapper
inline Index cols() const { return m_expression.cols(); }
inline Index outerStride() const { return m_expression.outerStride(); }
inline Index innerStride() const { return m_expression.innerStride(); }
+
+ typedef typename internal::conditional<
+ internal::is_lvalue<ExpressionType>::value,
+ Scalar,
+ const Scalar
+ >::type ScalarWithConstIfNotLvalue;
+
+ inline ScalarWithConstIfNotLvalue* data() { return m_expression.data(); }
+ inline const Scalar* data() const { return m_expression.data(); }
inline Scalar& coeffRef(Index row, Index col)
{
diff --git a/Eigen/src/Core/Transpose.h b/Eigen/src/Core/Transpose.h
index 3f7c7df6e..6562e8409 100644
--- a/Eigen/src/Core/Transpose.h
+++ b/Eigen/src/Core/Transpose.h
@@ -91,7 +91,7 @@ template<typename MatrixType> class Transpose
nestedExpression() { return m_matrix.const_cast_derived(); }
protected:
- const typename MatrixType::Nested m_matrix;
+ typename MatrixType::Nested m_matrix;
};
namespace internal {
@@ -152,12 +152,12 @@ template<typename MatrixType> class TransposeImpl<MatrixType,Dense>
return derived().nestedExpression().coeffRef(index);
}
- inline const CoeffReturnType coeff(Index row, Index col) const
+ inline CoeffReturnType coeff(Index row, Index col) const
{
return derived().nestedExpression().coeff(col, row);
}
- inline const CoeffReturnType coeff(Index index) const
+ inline CoeffReturnType coeff(Index index) const
{
return derived().nestedExpression().coeff(index);
}
diff --git a/Eigen/src/Core/Transpositions.h b/Eigen/src/Core/Transpositions.h
index 88fdfb222..3200cce02 100644
--- a/Eigen/src/Core/Transpositions.h
+++ b/Eigen/src/Core/Transpositions.h
@@ -404,7 +404,7 @@ struct transposition_matrix_product_retval
protected:
const TranspositionType& m_transpositions;
- const typename MatrixType::Nested m_matrix;
+ typename MatrixType::Nested m_matrix;
};
} // end namespace internal
diff --git a/Eigen/src/Core/TriangularMatrix.h b/Eigen/src/Core/TriangularMatrix.h
index 94245101f..13c7f08af 100644
--- a/Eigen/src/Core/TriangularMatrix.h
+++ b/Eigen/src/Core/TriangularMatrix.h
@@ -273,11 +273,8 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
inline const TriangularView<MatrixConjugateReturnType,Mode> conjugate() const
{ return m_matrix.conjugate(); }
- /** \sa MatrixBase::adjoint() */
- inline TriangularView<typename MatrixType::AdjointReturnType,TransposeMode> adjoint()
- { return m_matrix.adjoint(); }
/** \sa MatrixBase::adjoint() const */
- inline const TriangularView<typename MatrixType::AdjointReturnType,TransposeMode> adjoint() const
+ inline const TriangularView<const typename MatrixType::AdjointReturnType,TransposeMode> adjoint() const
{ return m_matrix.adjoint(); }
/** \sa MatrixBase::transpose() */
@@ -288,11 +285,13 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
}
/** \sa MatrixBase::transpose() const */
inline const TriangularView<Transpose<MatrixType>,TransposeMode> transpose() const
- { return m_matrix.transpose(); }
+ {
+ return m_matrix.transpose();
+ }
/** Efficient triangular matrix times vector/matrix product */
template<typename OtherDerived>
- TriangularProduct<Mode,true,MatrixType,false,OtherDerived,OtherDerived::IsVectorAtCompileTime>
+ TriangularProduct<Mode,true,MatrixType,false,OtherDerived, OtherDerived::IsVectorAtCompileTime>
operator*(const MatrixBase<OtherDerived>& rhs) const
{
return TriangularProduct
@@ -375,7 +374,8 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
template<typename OtherDerived>
void swap(MatrixBase<OtherDerived> const & other)
{
- TriangularView<SwapWrapper<MatrixType>,Mode>(const_cast<MatrixType&>(m_matrix)).lazyAssign(other.derived());
+ SwapWrapper<MatrixType> swaper(const_cast<MatrixType&>(m_matrix));
+ TriangularView<SwapWrapper<MatrixType>,Mode>(swaper).lazyAssign(other.derived());
}
Scalar determinant() const
@@ -433,7 +433,7 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
template<typename ProductDerived, typename Lhs, typename Rhs>
EIGEN_STRONG_INLINE TriangularView& assignProduct(const ProductBase<ProductDerived, Lhs,Rhs>& prod, const Scalar& alpha);
- const MatrixTypeNested m_matrix;
+ MatrixTypeNested m_matrix;
};
/***************************************************************************
diff --git a/Eigen/src/Core/VectorwiseOp.h b/Eigen/src/Core/VectorwiseOp.h
index e8a905043..e31ce5f3a 100644
--- a/Eigen/src/Core/VectorwiseOp.h
+++ b/Eigen/src/Core/VectorwiseOp.h
@@ -110,7 +110,7 @@ class PartialReduxExpr : internal::no_assignment_operator,
}
protected:
- const MatrixTypeNested m_matrix;
+ MatrixTypeNested m_matrix;
const MemberOp m_functor;
};
diff --git a/Eigen/src/Core/products/CoeffBasedProduct.h b/Eigen/src/Core/products/CoeffBasedProduct.h
index e8c767f48..5ae99d7b9 100644
--- a/Eigen/src/Core/products/CoeffBasedProduct.h
+++ b/Eigen/src/Core/products/CoeffBasedProduct.h
@@ -224,8 +224,8 @@ class CoeffBasedProduct
{ return reinterpret_cast<const LazyCoeffBasedProductType&>(*this).diagonal(index); }
protected:
- const LhsNested m_lhs;
- const RhsNested m_rhs;
+ typename internal::add_const_on_value_type<LhsNested>::type m_lhs;
+ typename internal::add_const_on_value_type<RhsNested>::type m_rhs;
mutable PlainObject m_result;
};
diff --git a/Eigen/src/Core/products/GeneralMatrixMatrix.h b/Eigen/src/Core/products/GeneralMatrixMatrix.h
index ae94a2795..0e73c2499 100644
--- a/Eigen/src/Core/products/GeneralMatrixMatrix.h
+++ b/Eigen/src/Core/products/GeneralMatrixMatrix.h
@@ -412,8 +412,8 @@ class GeneralProduct<Lhs, Rhs, GemmProduct>
{
eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
- const ActualLhsType lhs = LhsBlasTraits::extract(m_lhs);
- const ActualRhsType rhs = RhsBlasTraits::extract(m_rhs);
+ typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(m_lhs);
+ typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(m_rhs);
Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(m_lhs)
* RhsBlasTraits::extractScalarFactor(m_rhs);
diff --git a/Eigen/src/Core/products/GeneralMatrixMatrixTriangular.h b/Eigen/src/Core/products/GeneralMatrixMatrixTriangular.h
index 0ed65f9b4..0499256b3 100644
--- a/Eigen/src/Core/products/GeneralMatrixMatrixTriangular.h
+++ b/Eigen/src/Core/products/GeneralMatrixMatrixTriangular.h
@@ -201,13 +201,13 @@ TriangularView<MatrixType,UpLo>& TriangularView<MatrixType,UpLo>::assignProduct(
typedef internal::blas_traits<Lhs> LhsBlasTraits;
typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhs;
typedef typename internal::remove_all<ActualLhs>::type _ActualLhs;
- const ActualLhs actualLhs = LhsBlasTraits::extract(prod.lhs());
+ typename internal::add_const_on_value_type<ActualLhs>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
typedef typename internal::remove_all<typename ProductDerived::RhsNested>::type Rhs;
typedef internal::blas_traits<Rhs> RhsBlasTraits;
typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhs;
typedef typename internal::remove_all<ActualRhs>::type _ActualRhs;
- const ActualRhs actualRhs = RhsBlasTraits::extract(prod.rhs());
+ typename internal::add_const_on_value_type<ActualRhs>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
typename ProductDerived::Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs().derived()) * RhsBlasTraits::extractScalarFactor(prod.rhs().derived());
diff --git a/Eigen/src/Core/products/SelfadjointMatrixMatrix.h b/Eigen/src/Core/products/SelfadjointMatrixMatrix.h
index ccd757cfa..d0a076e44 100644
--- a/Eigen/src/Core/products/SelfadjointMatrixMatrix.h
+++ b/Eigen/src/Core/products/SelfadjointMatrixMatrix.h
@@ -400,8 +400,8 @@ struct SelfadjointProductMatrix<Lhs,LhsMode,false,Rhs,RhsMode,false>
{
eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
- const ActualLhsType lhs = LhsBlasTraits::extract(m_lhs);
- const ActualRhsType rhs = RhsBlasTraits::extract(m_rhs);
+ typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(m_lhs);
+ typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(m_rhs);
Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(m_lhs)
* RhsBlasTraits::extractScalarFactor(m_rhs);
diff --git a/Eigen/src/Core/products/SelfadjointMatrixVector.h b/Eigen/src/Core/products/SelfadjointMatrixVector.h
index 8cceb566f..9daca4430 100644
--- a/Eigen/src/Core/products/SelfadjointMatrixVector.h
+++ b/Eigen/src/Core/products/SelfadjointMatrixVector.h
@@ -201,8 +201,8 @@ struct SelfadjointProductMatrix<Lhs,LhsMode,false,Rhs,0,true>
eigen_assert(dest.rows()==m_lhs.rows() && dest.cols()==m_rhs.cols());
- const ActualLhsType lhs = LhsBlasTraits::extract(m_lhs);
- const ActualRhsType rhs = RhsBlasTraits::extract(m_rhs);
+ typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(m_lhs);
+ typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(m_rhs);
Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(m_lhs)
* RhsBlasTraits::extractScalarFactor(m_rhs);
diff --git a/Eigen/src/Core/products/SelfadjointProduct.h b/Eigen/src/Core/products/SelfadjointProduct.h
index 1def96719..942bdefae 100644
--- a/Eigen/src/Core/products/SelfadjointProduct.h
+++ b/Eigen/src/Core/products/SelfadjointProduct.h
@@ -72,7 +72,7 @@ struct selfadjoint_product_selector<MatrixType,OtherType,UpLo,true>
typedef internal::blas_traits<OtherType> OtherBlasTraits;
typedef typename OtherBlasTraits::DirectLinearAccessType ActualOtherType;
typedef typename internal::remove_all<ActualOtherType>::type _ActualOtherType;
- const ActualOtherType actualOther = OtherBlasTraits::extract(other.derived());
+ typename internal::add_const_on_value_type<ActualOtherType>::type actualOther = OtherBlasTraits::extract(other.derived());
Scalar actualAlpha = alpha * OtherBlasTraits::extractScalarFactor(other.derived());
@@ -105,7 +105,7 @@ struct selfadjoint_product_selector<MatrixType,OtherType,UpLo,false>
typedef internal::blas_traits<OtherType> OtherBlasTraits;
typedef typename OtherBlasTraits::DirectLinearAccessType ActualOtherType;
typedef typename internal::remove_all<ActualOtherType>::type _ActualOtherType;
- const ActualOtherType actualOther = OtherBlasTraits::extract(other.derived());
+ typename internal::add_const_on_value_type<ActualOtherType>::type actualOther = OtherBlasTraits::extract(other.derived());
Scalar actualAlpha = alpha * OtherBlasTraits::extractScalarFactor(other.derived());
diff --git a/Eigen/src/Core/products/SelfadjointRank2Update.h b/Eigen/src/Core/products/SelfadjointRank2Update.h
index 9f8b8438a..4142b2a72 100644
--- a/Eigen/src/Core/products/SelfadjointRank2Update.h
+++ b/Eigen/src/Core/products/SelfadjointRank2Update.h
@@ -76,12 +76,12 @@ SelfAdjointView<MatrixType,UpLo>& SelfAdjointView<MatrixType,UpLo>
typedef internal::blas_traits<DerivedU> UBlasTraits;
typedef typename UBlasTraits::DirectLinearAccessType ActualUType;
typedef typename internal::remove_all<ActualUType>::type _ActualUType;
- const ActualUType actualU = UBlasTraits::extract(u.derived());
+ typename internal::add_const_on_value_type<ActualUType>::type actualU = UBlasTraits::extract(u.derived());
typedef internal::blas_traits<DerivedV> VBlasTraits;
typedef typename VBlasTraits::DirectLinearAccessType ActualVType;
typedef typename internal::remove_all<ActualVType>::type _ActualVType;
- const ActualVType actualV = VBlasTraits::extract(v.derived());
+ typename internal::add_const_on_value_type<ActualVType>::type actualV = VBlasTraits::extract(v.derived());
// If MatrixType is row major, then we use the routine for lower triangular in the upper triangular case and
// vice versa, and take the complex conjugate of all coefficients and vector entries.
diff --git a/Eigen/src/Core/products/TriangularMatrixMatrix.h b/Eigen/src/Core/products/TriangularMatrixMatrix.h
index 7632ba85a..8e5a409e2 100644
--- a/Eigen/src/Core/products/TriangularMatrixMatrix.h
+++ b/Eigen/src/Core/products/TriangularMatrixMatrix.h
@@ -378,8 +378,8 @@ struct TriangularProduct<Mode,LhsIsTriangular,Lhs,false,Rhs,false>
template<typename Dest> void scaleAndAddTo(Dest& dst, Scalar alpha) const
{
- const ActualLhsType lhs = LhsBlasTraits::extract(m_lhs);
- const ActualRhsType rhs = RhsBlasTraits::extract(m_rhs);
+ typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(m_lhs);
+ typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(m_rhs);
Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(m_lhs)
* RhsBlasTraits::extractScalarFactor(m_rhs);
diff --git a/Eigen/src/Core/products/TriangularMatrixVector.h b/Eigen/src/Core/products/TriangularMatrixVector.h
index ebda994f6..770641419 100644
--- a/Eigen/src/Core/products/TriangularMatrixVector.h
+++ b/Eigen/src/Core/products/TriangularMatrixVector.h
@@ -232,8 +232,8 @@ template<> struct trmv_selector<ColMajor>
typedef typename ProductType::RhsBlasTraits RhsBlasTraits;
typedef Map<Matrix<ResScalar,Dynamic,1>, Aligned> MappedDest;
- const ActualLhsType actualLhs = LhsBlasTraits::extract(prod.lhs());
- const ActualRhsType actualRhs = RhsBlasTraits::extract(prod.rhs());
+ typename internal::add_const_on_value_type<ActualLhsType>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
+ typename internal::add_const_on_value_type<ActualRhsType>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
ResScalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs())
* RhsBlasTraits::extractScalarFactor(prod.rhs());
diff --git a/Eigen/src/Core/util/BlasUtil.h b/Eigen/src/Core/util/BlasUtil.h
index 3b740b9c0..5dd2da049 100644
--- a/Eigen/src/Core/util/BlasUtil.h
+++ b/Eigen/src/Core/util/BlasUtil.h
@@ -175,7 +175,7 @@ template<typename XprType> struct blas_traits
ExtractType,
typename _ExtractType::PlainObject
>::type DirectLinearAccessType;
- static inline const ExtractType extract(const XprType& x) { return x; }
+ static inline ExtractType extract(const XprType& x) { return x; }
static inline const Scalar extractScalarFactor(const XprType&) { return Scalar(1); }
};
@@ -192,7 +192,7 @@ struct blas_traits<CwiseUnaryOp<scalar_conjugate_op<Scalar>, NestedXpr> >
IsComplex = NumTraits<Scalar>::IsComplex,
NeedToConjugate = Base::NeedToConjugate ? 0 : IsComplex
};
- static inline const ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
+ static inline ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
static inline Scalar extractScalarFactor(const XprType& x) { return conj(Base::extractScalarFactor(x.nestedExpression())); }
};
@@ -204,7 +204,7 @@ struct blas_traits<CwiseUnaryOp<scalar_multiple_op<Scalar>, NestedXpr> >
typedef blas_traits<NestedXpr> Base;
typedef CwiseUnaryOp<scalar_multiple_op<Scalar>, NestedXpr> XprType;
typedef typename Base::ExtractType ExtractType;
- static inline const ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
+ static inline ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
static inline Scalar extractScalarFactor(const XprType& x)
{ return x.functor().m_other * Base::extractScalarFactor(x.nestedExpression()); }
};
@@ -217,7 +217,7 @@ struct blas_traits<CwiseUnaryOp<scalar_opposite_op<Scalar>, NestedXpr> >
typedef blas_traits<NestedXpr> Base;
typedef CwiseUnaryOp<scalar_opposite_op<Scalar>, NestedXpr> XprType;
typedef typename Base::ExtractType ExtractType;
- static inline const ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
+ static inline ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
static inline Scalar extractScalarFactor(const XprType& x)
{ return - Base::extractScalarFactor(x.nestedExpression()); }
};
@@ -239,7 +239,7 @@ struct blas_traits<Transpose<NestedXpr> >
enum {
IsTransposed = Base::IsTransposed ? 0 : 1
};
- static inline const ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
+ static inline ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
static inline Scalar extractScalarFactor(const XprType& x) { return Base::extractScalarFactor(x.nestedExpression()); }
};
@@ -252,7 +252,7 @@ template<typename T, bool HasUsableDirectAccess=blas_traits<T>::HasUsableDirectA
struct extract_data_selector {
static const typename T::Scalar* run(const T& m)
{
- return const_cast<typename T::Scalar*>(&blas_traits<T>::extract(m).coeffRef(0,0)); // FIXME this should be .data()
+ return blas_traits<T>::extract(m).data();
}
};
diff --git a/Eigen/src/Core/util/XprHelper.h b/Eigen/src/Core/util/XprHelper.h
index c2078f137..ba1c4379d 100644
--- a/Eigen/src/Core/util/XprHelper.h
+++ b/Eigen/src/Core/util/XprHelper.h
@@ -260,30 +260,27 @@ template<typename T> struct plain_matrix_type_row_major
// we should be able to get rid of this one too
template<typename T> struct must_nest_by_value { enum { ret = false }; };
-template<class T>
-struct is_reference
-{
- enum { ret = false };
-};
-
-template<class T>
-struct is_reference<T&>
-{
- enum { ret = true };
-};
-
-/**
-* \internal The reference selector for template expressions. The idea is that we don't
-* need to use references for expressions since they are light weight proxy
-* objects which should generate no copying overhead.
-**/
+/** \internal The reference selector for template expressions. The idea is that we don't
+ * need to use references for expressions since they are light weight proxy
+ * objects which should generate no copying overhead. */
template <typename T>
struct ref_selector
{
typedef typename conditional<
bool(traits<T>::Flags & NestByRefBit),
T const&,
- T
+ const T
+ >::type type;
+};
+
+/** \internal Adds the const qualifier on the value-type of T2 if and only if T1 is a const type */
+template<typename T1, typename T2>
+struct transfer_constness
+{
+ typedef typename conditional<
+ bool(internal::is_const<T1>::value),
+ typename internal::add_const_on_value_type<T2>::type,
+ T2
>::type type;
};
@@ -297,6 +294,8 @@ struct ref_selector
* \param T the type of the expression being nested
* \param n the number of coefficient accesses in the nested expression for each coefficient access in the bigger expression.
*
+ * Note that if no evaluation occur, then the constness of T is preserved.
+ *
* Example. Suppose that a, b, and c are of type Matrix3d. The user forms the expression a*(b+c).
* b+c is an expression "sum of matrices", which we will denote by S. In order to determine how to nest it,
* the Product expression uses: nested<S, 3>::ret, which turns out to be Matrix3d because the internal logic of