aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2012-06-28 21:01:02 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2012-06-28 21:01:02 +0200
commit9629ba361a5d7b806053372c46f18c49ee971a10 (patch)
tree3d2f28735fefbb78df836ccdae2d5e894a7c1185 /Eigen/src
parent23184527faeae819480560592c7a2bdf6c272b82 (diff)
bug #482: pass scalar by const ref - pass on the sparse module
(also fix a compilation issue due to previous pass)
Diffstat (limited to 'Eigen/src')
-rw-r--r--Eigen/src/SparseCore/CompressedStorage.h8
-rw-r--r--Eigen/src/SparseCore/SparseDenseProduct.h14
-rw-r--r--Eigen/src/SparseCore/SparseMatrix.h9
-rw-r--r--Eigen/src/SparseCore/SparseMatrixBase.h4
-rw-r--r--Eigen/src/SparseCore/SparseProduct.h4
-rw-r--r--Eigen/src/SparseCore/SparseSelfAdjointView.h8
-rw-r--r--Eigen/src/SparseCore/SparseSparseProductWithPruning.h10
-rw-r--r--Eigen/src/SparseCore/SparseTranspose.h2
-rw-r--r--Eigen/src/SparseCore/SparseVector.h2
-rw-r--r--Eigen/src/SparseCore/SparseView.h2
10 files changed, 31 insertions, 32 deletions
diff --git a/Eigen/src/SparseCore/CompressedStorage.h b/Eigen/src/SparseCore/CompressedStorage.h
index fa2bfd763..84acbd48f 100644
--- a/Eigen/src/SparseCore/CompressedStorage.h
+++ b/Eigen/src/SparseCore/CompressedStorage.h
@@ -154,7 +154,7 @@ class CompressedStorage
/** \returns the stored value at index \a key
* If the value does not exist, then the value \a defaultValue is returned without any insertion. */
- inline Scalar at(Index key, Scalar defaultValue = Scalar(0)) const
+ inline Scalar at(Index key, const Scalar& defaultValue = Scalar(0)) const
{
if (m_size==0)
return defaultValue;
@@ -167,7 +167,7 @@ class CompressedStorage
}
/** Like at(), but the search is performed in the range [start,end) */
- inline Scalar atInRange(size_t start, size_t end, Index key, Scalar defaultValue = Scalar(0)) const
+ inline Scalar atInRange(size_t start, size_t end, Index key, const Scalar& defaultValue = Scalar(0)) const
{
if (start>=end)
return Scalar(0);
@@ -182,7 +182,7 @@ class CompressedStorage
/** \returns a reference to the value at index \a key
* If the value does not exist, then the value \a defaultValue is inserted
* such that the keys are sorted. */
- inline Scalar& atWithInsertion(Index key, Scalar defaultValue = Scalar(0))
+ inline Scalar& atWithInsertion(Index key, const Scalar& defaultValue = Scalar(0))
{
size_t id = searchLowerIndex(0,m_size,key);
if (id>=m_size || m_indices[id]!=key)
@@ -199,7 +199,7 @@ class CompressedStorage
return m_values[id];
}
- void prune(Scalar reference, RealScalar epsilon = NumTraits<RealScalar>::dummy_precision())
+ void prune(const Scalar& reference, const RealScalar& epsilon = NumTraits<RealScalar>::dummy_precision())
{
size_t k = 0;
size_t n = size();
diff --git a/Eigen/src/SparseCore/SparseDenseProduct.h b/Eigen/src/SparseCore/SparseDenseProduct.h
index 00ba606be..a2421025c 100644
--- a/Eigen/src/SparseCore/SparseDenseProduct.h
+++ b/Eigen/src/SparseCore/SparseDenseProduct.h
@@ -54,7 +54,7 @@ struct traits<SparseDenseOuterProduct<Lhs,Rhs,Tr> >
{
typedef Sparse StorageKind;
typedef typename scalar_product_traits<typename traits<Lhs>::Scalar,
- typename traits<Rhs>::Scalar>::ReturnType Scalar;
+ typename traits<Rhs>::Scalar>::ReturnType Scalar;
typedef typename Lhs::Index Index;
typedef typename Lhs::Nested LhsNested;
typedef typename Rhs::Nested RhsNested;
@@ -165,7 +165,7 @@ struct sparse_time_dense_product_impl<SparseLhsType,DenseRhsType,DenseResType, R
typedef typename internal::remove_all<DenseResType>::type Res;
typedef typename Lhs::Index Index;
typedef typename Lhs::InnerIterator LhsInnerIterator;
- static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, typename Res::Scalar alpha)
+ static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, const typename Res::Scalar& alpha)
{
for(Index c=0; c<rhs.cols(); ++c)
{
@@ -189,7 +189,7 @@ struct sparse_time_dense_product_impl<SparseLhsType,DenseRhsType,DenseResType, C
typedef typename internal::remove_all<DenseResType>::type Res;
typedef typename Lhs::InnerIterator LhsInnerIterator;
typedef typename Lhs::Index Index;
- static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, typename Res::Scalar alpha)
+ static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, const typename Res::Scalar& alpha)
{
for(Index c=0; c<rhs.cols(); ++c)
{
@@ -211,7 +211,7 @@ struct sparse_time_dense_product_impl<SparseLhsType,DenseRhsType,DenseResType, R
typedef typename internal::remove_all<DenseResType>::type Res;
typedef typename Lhs::InnerIterator LhsInnerIterator;
typedef typename Lhs::Index Index;
- static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, typename Res::Scalar alpha)
+ static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, const typename Res::Scalar& alpha)
{
for(Index j=0; j<lhs.outerSize(); ++j)
{
@@ -230,7 +230,7 @@ struct sparse_time_dense_product_impl<SparseLhsType,DenseRhsType,DenseResType, C
typedef typename internal::remove_all<DenseResType>::type Res;
typedef typename Lhs::InnerIterator LhsInnerIterator;
typedef typename Lhs::Index Index;
- static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, typename Res::Scalar alpha)
+ static void run(const SparseLhsType& lhs, const DenseRhsType& rhs, DenseResType& res, const typename Res::Scalar& alpha)
{
for(Index j=0; j<lhs.outerSize(); ++j)
{
@@ -259,7 +259,7 @@ class SparseTimeDenseProduct
SparseTimeDenseProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
{}
- template<typename Dest> void scaleAndAddTo(Dest& dest, Scalar alpha) const
+ template<typename Dest> void scaleAndAddTo(Dest& dest, const Scalar& alpha) const
{
internal::sparse_time_dense_product(m_lhs, m_rhs, dest, alpha);
}
@@ -289,7 +289,7 @@ class DenseTimeSparseProduct
DenseTimeSparseProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
{}
- template<typename Dest> void scaleAndAddTo(Dest& dest, Scalar alpha) const
+ template<typename Dest> void scaleAndAddTo(Dest& dest, const Scalar& alpha) const
{
Transpose<const _LhsNested> lhs_t(m_lhs);
Transpose<const _RhsNested> rhs_t(m_rhs);
diff --git a/Eigen/src/SparseCore/SparseMatrix.h b/Eigen/src/SparseCore/SparseMatrix.h
index 214f130f5..dbdfb1b5c 100644
--- a/Eigen/src/SparseCore/SparseMatrix.h
+++ b/Eigen/src/SparseCore/SparseMatrix.h
@@ -478,7 +478,7 @@ class SparseMatrix
}
/** Suppresses all nonzeros which are \b much \b smaller \b than \a reference under the tolerence \a epsilon */
- void prune(Scalar reference, RealScalar epsilon = NumTraits<RealScalar>::dummy_precision())
+ void prune(const Scalar& reference, const RealScalar& epsilon = NumTraits<RealScalar>::dummy_precision())
{
prune(default_prunning_func(reference,epsilon));
}
@@ -909,7 +909,7 @@ protected:
public:
/** \internal
* \sa insert(Index,Index) */
- inline Scalar& insertBackUncompressed(Index row, Index col)
+ EIGEN_STRONG_INLINE Scalar& insertBackUncompressed(Index row, Index col)
{
const Index outer = IsRowMajor ? row : col;
const Index inner = IsRowMajor ? col : row;
@@ -917,8 +917,7 @@ public:
eigen_assert(!isCompressed());
eigen_assert(m_innerNonZeros[outer]<=(m_outerIndex[outer+1] - m_outerIndex[outer]));
- Index p = m_outerIndex[outer] + m_innerNonZeros[outer];
- m_innerNonZeros[outer]++;
+ Index p = m_outerIndex[outer] + m_innerNonZeros[outer]++;
m_data.index(p) = inner;
return (m_data.value(p) = 0);
}
@@ -930,7 +929,7 @@ private:
}
struct default_prunning_func {
- default_prunning_func(Scalar ref, RealScalar eps) : reference(ref), epsilon(eps) {}
+ default_prunning_func(const Scalar& ref, const RealScalar& eps) : reference(ref), epsilon(eps) {}
inline bool operator() (const Index&, const Index&, const Scalar& value) const
{
return !internal::isMuchSmallerThan(value, reference, epsilon);
diff --git a/Eigen/src/SparseCore/SparseMatrixBase.h b/Eigen/src/SparseCore/SparseMatrixBase.h
index 9d7b83577..493988c0b 100644
--- a/Eigen/src/SparseCore/SparseMatrixBase.h
+++ b/Eigen/src/SparseCore/SparseMatrixBase.h
@@ -445,12 +445,12 @@ template<typename Derived> class SparseMatrixBase : public EigenBase<Derived>
template<typename OtherDerived>
bool isApprox(const SparseMatrixBase<OtherDerived>& other,
- RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
+ const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
{ return toDense().isApprox(other.toDense(),prec); }
template<typename OtherDerived>
bool isApprox(const MatrixBase<OtherDerived>& other,
- RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
+ const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
{ return toDense().isApprox(other,prec); }
/** \returns the matrix or vector obtained by evaluating this expression.
diff --git a/Eigen/src/SparseCore/SparseProduct.h b/Eigen/src/SparseCore/SparseProduct.h
index 813dbf624..ed974b9e2 100644
--- a/Eigen/src/SparseCore/SparseProduct.h
+++ b/Eigen/src/SparseCore/SparseProduct.h
@@ -114,13 +114,13 @@ class SparseSparseProduct : internal::no_assignment_operator,
}
template<typename Lhs, typename Rhs>
- EIGEN_STRONG_INLINE SparseSparseProduct(const Lhs& lhs, const Rhs& rhs, RealScalar tolerance)
+ EIGEN_STRONG_INLINE SparseSparseProduct(const Lhs& lhs, const Rhs& rhs, const RealScalar& tolerance)
: m_lhs(lhs), m_rhs(rhs), m_tolerance(tolerance), m_conservative(false)
{
init();
}
- SparseSparseProduct pruned(Scalar reference = 0, RealScalar epsilon = NumTraits<RealScalar>::dummy_precision()) const
+ SparseSparseProduct pruned(const Scalar& reference = 0, const RealScalar& epsilon = NumTraits<RealScalar>::dummy_precision()) const
{
return SparseSparseProduct(m_lhs,m_rhs,internal::abs(reference)*epsilon);
}
diff --git a/Eigen/src/SparseCore/SparseSelfAdjointView.h b/Eigen/src/SparseCore/SparseSelfAdjointView.h
index c925a894d..04537f338 100644
--- a/Eigen/src/SparseCore/SparseSelfAdjointView.h
+++ b/Eigen/src/SparseCore/SparseSelfAdjointView.h
@@ -109,7 +109,7 @@ template<typename MatrixType, unsigned int UpLo> class SparseSelfAdjointView
* call this function with u.adjoint().
*/
template<typename DerivedU>
- SparseSelfAdjointView& rankUpdate(const SparseMatrixBase<DerivedU>& u, Scalar alpha = Scalar(1));
+ SparseSelfAdjointView& rankUpdate(const SparseMatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1));
/** \internal triggered by sparse_matrix = SparseSelfadjointView; */
template<typename DestScalar,int StorageOrder> void evalTo(SparseMatrix<DestScalar,StorageOrder,Index>& _dest) const
@@ -188,7 +188,7 @@ SparseSelfAdjointView<Derived, UpLo> SparseMatrixBase<Derived>::selfadjointView(
template<typename MatrixType, unsigned int UpLo>
template<typename DerivedU>
SparseSelfAdjointView<MatrixType,UpLo>&
-SparseSelfAdjointView<MatrixType,UpLo>::rankUpdate(const SparseMatrixBase<DerivedU>& u, Scalar alpha)
+SparseSelfAdjointView<MatrixType,UpLo>::rankUpdate(const SparseMatrixBase<DerivedU>& u, const Scalar& alpha)
{
SparseMatrix<Scalar,MatrixType::Flags&RowMajorBit?RowMajor:ColMajor> tmp = u * u.adjoint();
if(alpha==Scalar(0))
@@ -222,7 +222,7 @@ class SparseSelfAdjointTimeDenseProduct
SparseSelfAdjointTimeDenseProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
{}
- template<typename Dest> void scaleAndAddTo(Dest& dest, Scalar alpha) const
+ template<typename Dest> void scaleAndAddTo(Dest& dest, const Scalar& alpha) const
{
// TODO use alpha
eigen_assert(alpha==Scalar(1) && "alpha != 1 is not implemented yet, sorry");
@@ -283,7 +283,7 @@ class DenseTimeSparseSelfAdjointProduct
DenseTimeSparseSelfAdjointProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
{}
- template<typename Dest> void scaleAndAddTo(Dest& /*dest*/, Scalar /*alpha*/) const
+ template<typename Dest> void scaleAndAddTo(Dest& /*dest*/, const Scalar& /*alpha*/) const
{
// TODO
}
diff --git a/Eigen/src/SparseCore/SparseSparseProductWithPruning.h b/Eigen/src/SparseCore/SparseSparseProductWithPruning.h
index abd4fda82..4d3b9ac51 100644
--- a/Eigen/src/SparseCore/SparseSparseProductWithPruning.h
+++ b/Eigen/src/SparseCore/SparseSparseProductWithPruning.h
@@ -32,7 +32,7 @@ namespace internal {
// perform a pseudo in-place sparse * sparse product assuming all matrices are col major
template<typename Lhs, typename Rhs, typename ResultType>
-static void sparse_sparse_product_with_pruning_impl(const Lhs& lhs, const Rhs& rhs, ResultType& res, typename ResultType::RealScalar tolerance)
+static void sparse_sparse_product_with_pruning_impl(const Lhs& lhs, const Rhs& rhs, ResultType& res, const typename ResultType::RealScalar& tolerance)
{
// return sparse_sparse_product_with_pruning_impl2(lhs,rhs,res);
@@ -100,7 +100,7 @@ struct sparse_sparse_product_with_pruning_selector<Lhs,Rhs,ResultType,ColMajor,C
typedef typename traits<typename remove_all<Lhs>::type>::Scalar Scalar;
typedef typename ResultType::RealScalar RealScalar;
- static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, RealScalar tolerance)
+ static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, const RealScalar& tolerance)
{
typename remove_all<ResultType>::type _res(res.rows(), res.cols());
internal::sparse_sparse_product_with_pruning_impl<Lhs,Rhs,ResultType>(lhs, rhs, _res, tolerance);
@@ -112,7 +112,7 @@ template<typename Lhs, typename Rhs, typename ResultType>
struct sparse_sparse_product_with_pruning_selector<Lhs,Rhs,ResultType,ColMajor,ColMajor,RowMajor>
{
typedef typename ResultType::RealScalar RealScalar;
- static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, RealScalar tolerance)
+ static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, const RealScalar& tolerance)
{
// we need a col-major matrix to hold the result
typedef SparseMatrix<typename ResultType::Scalar> SparseTemporaryType;
@@ -126,7 +126,7 @@ template<typename Lhs, typename Rhs, typename ResultType>
struct sparse_sparse_product_with_pruning_selector<Lhs,Rhs,ResultType,RowMajor,RowMajor,RowMajor>
{
typedef typename ResultType::RealScalar RealScalar;
- static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, RealScalar tolerance)
+ static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, const RealScalar& tolerance)
{
// let's transpose the product to get a column x column product
typename remove_all<ResultType>::type _res(res.rows(), res.cols());
@@ -139,7 +139,7 @@ template<typename Lhs, typename Rhs, typename ResultType>
struct sparse_sparse_product_with_pruning_selector<Lhs,Rhs,ResultType,RowMajor,RowMajor,ColMajor>
{
typedef typename ResultType::RealScalar RealScalar;
- static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, RealScalar tolerance)
+ static void run(const Lhs& lhs, const Rhs& rhs, ResultType& res, const RealScalar& tolerance)
{
typedef SparseMatrix<typename ResultType::Scalar,ColMajor> ColMajorMatrix;
ColMajorMatrix colLhs(lhs);
diff --git a/Eigen/src/SparseCore/SparseTranspose.h b/Eigen/src/SparseCore/SparseTranspose.h
index 07d9e0bbd..64763890c 100644
--- a/Eigen/src/SparseCore/SparseTranspose.h
+++ b/Eigen/src/SparseCore/SparseTranspose.h
@@ -33,7 +33,7 @@ template<typename MatrixType> class TransposeImpl<MatrixType,Sparse>
typedef typename internal::remove_all<typename MatrixType::Nested>::type _MatrixTypeNested;
public:
- EIGEN_SPARSE_PUBLIC_INTERFACE(Transpose<MatrixType>)
+ EIGEN_SPARSE_PUBLIC_INTERFACE(Transpose<MatrixType> )
class InnerIterator;
class ReverseInnerIterator;
diff --git a/Eigen/src/SparseCore/SparseVector.h b/Eigen/src/SparseCore/SparseVector.h
index e81347705..2a101aad6 100644
--- a/Eigen/src/SparseCore/SparseVector.h
+++ b/Eigen/src/SparseCore/SparseVector.h
@@ -184,7 +184,7 @@ class SparseVector
inline void finalize() {}
- void prune(Scalar reference, RealScalar epsilon = NumTraits<RealScalar>::dummy_precision())
+ void prune(const Scalar& reference, const RealScalar& epsilon = NumTraits<RealScalar>::dummy_precision())
{
m_data.prune(reference,epsilon);
}
diff --git a/Eigen/src/SparseCore/SparseView.h b/Eigen/src/SparseCore/SparseView.h
index 43a3adb24..cf3866d38 100644
--- a/Eigen/src/SparseCore/SparseView.h
+++ b/Eigen/src/SparseCore/SparseView.h
@@ -103,7 +103,7 @@ private:
template<typename Derived>
const SparseView<Derived> MatrixBase<Derived>::sparseView(const Scalar& m_reference,
- typename NumTraits<Scalar>::Real m_epsilon) const
+ const typename NumTraits<Scalar>::Real& m_epsilon) const
{
return SparseView<Derived>(derived(), m_reference, m_epsilon);
}