aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-02-09 10:23:45 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-02-09 10:23:45 +0100
commit3af29caae87b00ed8f002605573d227ad1b629a4 (patch)
tree9ac73034204b7fd7446a88a21cc7d1d59a209c8c
parentf2ff8c091e02b4aab8c7568807c09e43f51d1156 (diff)
Cleaning and add more unit tests for Ref<SparseMatrix> and Map<SparseMatrix>
-rw-r--r--Eigen/src/SparseCore/MappedSparseMatrix.h172
-rw-r--r--Eigen/src/SparseCore/SparseMap.h34
-rw-r--r--Eigen/src/SparseCore/SparseRef.h12
-rw-r--r--test/sparse_basic.cpp20
-rw-r--r--test/sparse_ref.cpp5
5 files changed, 72 insertions, 171 deletions
diff --git a/Eigen/src/SparseCore/MappedSparseMatrix.h b/Eigen/src/SparseCore/MappedSparseMatrix.h
index 2852c669a..533479fd0 100644
--- a/Eigen/src/SparseCore/MappedSparseMatrix.h
+++ b/Eigen/src/SparseCore/MappedSparseMatrix.h
@@ -1,7 +1,7 @@
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
-// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
+// Copyright (C) 2008-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
@@ -10,9 +10,10 @@
#ifndef EIGEN_MAPPED_SPARSEMATRIX_H
#define EIGEN_MAPPED_SPARSEMATRIX_H
-namespace Eigen {
+namespace Eigen {
-/** \class MappedSparseMatrix
+/** \deprecated Use Map<SparseMatrix<> >
+ * \class MappedSparseMatrix
*
* \brief Sparse matrix
*
@@ -25,179 +26,38 @@ namespace internal {
template<typename _Scalar, int _Flags, typename _Index>
struct traits<MappedSparseMatrix<_Scalar, _Flags, _Index> > : traits<SparseMatrix<_Scalar, _Flags, _Index> >
{};
-}
+} // end namespace internal
template<typename _Scalar, int _Flags, typename _Index>
class MappedSparseMatrix
- : public SparseMatrixBase<MappedSparseMatrix<_Scalar, _Flags, _Index> >
+ : public Map<SparseMatrix<_Scalar, _Flags, _Index> >
{
- public:
- EIGEN_SPARSE_PUBLIC_INTERFACE(MappedSparseMatrix)
- enum { IsRowMajor = Base::IsRowMajor };
-
- protected:
-
- Index m_outerSize;
- Index m_innerSize;
- Index m_nnz;
- Index* m_outerIndex;
- Index* m_innerIndices;
- Scalar* m_values;
+ typedef Map<SparseMatrix<_Scalar, _Flags, _Index> > Base;
public:
-
- inline Index rows() const { return IsRowMajor ? m_outerSize : m_innerSize; }
- inline Index cols() const { return IsRowMajor ? m_innerSize : m_outerSize; }
- inline Index innerSize() const { return m_innerSize; }
- inline Index outerSize() const { return m_outerSize; }
- bool isCompressed() const { return true; }
-
- //----------------------------------------
- // direct access interface
- inline const Scalar* valuePtr() const { return m_values; }
- inline Scalar* valuePtr() { return m_values; }
-
- inline const Index* innerIndexPtr() const { return m_innerIndices; }
- inline Index* innerIndexPtr() { return m_innerIndices; }
-
- inline const Index* outerIndexPtr() const { return m_outerIndex; }
- inline Index* outerIndexPtr() { return m_outerIndex; }
- //----------------------------------------
-
- inline Scalar coeff(Index row, Index col) const
- {
- const Index outer = IsRowMajor ? row : col;
- const Index inner = IsRowMajor ? col : row;
-
- Index start = m_outerIndex[outer];
- Index end = m_outerIndex[outer+1];
- if (start==end)
- return Scalar(0);
- else if (end>0 && inner==m_innerIndices[end-1])
- return m_values[end-1];
- // ^^ optimization: let's first check if it is the last coefficient
- // (very common in high level algorithms)
+ typedef typename Base::Index Index;
+ typedef typename Base::Scalar Scalar;
- const Index* r = std::lower_bound(&m_innerIndices[start],&m_innerIndices[end-1],inner);
- const Index id = r-&m_innerIndices[0];
- return ((*r==inner) && (id<end)) ? m_values[id] : Scalar(0);
- }
-
- inline Scalar& coeffRef(Index row, Index col)
- {
- const Index outer = IsRowMajor ? row : col;
- const Index inner = IsRowMajor ? col : row;
-
- Index start = m_outerIndex[outer];
- Index end = m_outerIndex[outer+1];
- eigen_assert(end>=start && "you probably called coeffRef on a non finalized matrix");
- eigen_assert(end>start && "coeffRef cannot be called on a zero coefficient");
- Index* r = std::lower_bound(&m_innerIndices[start],&m_innerIndices[end],inner);
- const Index id = r-&m_innerIndices[0];
- eigen_assert((*r==inner) && (id<end) && "coeffRef cannot be called on a zero coefficient");
- return m_values[id];
- }
-
- class InnerIterator;
- class ReverseInnerIterator;
-
- /** \returns the number of non zero coefficients */
- inline Index nonZeros() const { return m_nnz; }
-
- inline MappedSparseMatrix(Index rows, Index cols, Index nnz, Index* outerIndexPtr, Index* innerIndexPtr, Scalar* valuePtr)
- : m_outerSize(IsRowMajor?rows:cols), m_innerSize(IsRowMajor?cols:rows), m_nnz(nnz), m_outerIndex(outerIndexPtr),
- m_innerIndices(innerIndexPtr), m_values(valuePtr)
+ inline MappedSparseMatrix(Index rows, Index cols, Index nnz, Index* outerIndexPtr, Index* innerIndexPtr, Scalar* valuePtr, Index* innerNonZeroPtr = 0)
+ : Base(rows, cols, nnz, outerIndexPtr, innerIndexPtr, valuePtr, innerNonZeroPtr)
{}
/** Empty destructor */
inline ~MappedSparseMatrix() {}
};
-template<typename Scalar, int _Flags, typename _Index>
-class MappedSparseMatrix<Scalar,_Flags,_Index>::InnerIterator
-{
- public:
- InnerIterator(const MappedSparseMatrix& mat, Index outer)
- : m_matrix(mat),
- m_outer(outer),
- m_id(mat.outerIndexPtr()[outer]),
- m_start(m_id),
- m_end(mat.outerIndexPtr()[outer+1])
- {}
-
- inline InnerIterator& operator++() { m_id++; return *this; }
-
- inline Scalar value() const { return m_matrix.valuePtr()[m_id]; }
- inline Scalar& valueRef() { return const_cast<Scalar&>(m_matrix.valuePtr()[m_id]); }
-
- inline Index index() const { return m_matrix.innerIndexPtr()[m_id]; }
- inline Index row() const { return IsRowMajor ? m_outer : index(); }
- inline Index col() const { return IsRowMajor ? index() : m_outer; }
-
- inline operator bool() const { return (m_id < m_end) && (m_id>=m_start); }
-
- protected:
- const MappedSparseMatrix& m_matrix;
- const Index m_outer;
- Index m_id;
- const Index m_start;
- const Index m_end;
-};
-
-template<typename Scalar, int _Flags, typename _Index>
-class MappedSparseMatrix<Scalar,_Flags,_Index>::ReverseInnerIterator
-{
- public:
- ReverseInnerIterator(const MappedSparseMatrix& mat, Index outer)
- : m_matrix(mat),
- m_outer(outer),
- m_id(mat.outerIndexPtr()[outer+1]),
- m_start(mat.outerIndexPtr()[outer]),
- m_end(m_id)
- {}
-
- inline ReverseInnerIterator& operator--() { m_id--; return *this; }
-
- inline Scalar value() const { return m_matrix.valuePtr()[m_id-1]; }
- inline Scalar& valueRef() { return const_cast<Scalar&>(m_matrix.valuePtr()[m_id-1]); }
-
- inline Index index() const { return m_matrix.innerIndexPtr()[m_id-1]; }
- inline Index row() const { return IsRowMajor ? m_outer : index(); }
- inline Index col() const { return IsRowMajor ? index() : m_outer; }
-
- inline operator bool() const { return (m_id <= m_end) && (m_id>m_start); }
-
- protected:
- const MappedSparseMatrix& m_matrix;
- const Index m_outer;
- Index m_id;
- const Index m_start;
- const Index m_end;
-};
-
namespace internal {
template<typename _Scalar, int _Options, typename _Index>
struct evaluator<MappedSparseMatrix<_Scalar,_Options,_Index> >
- : evaluator_base<MappedSparseMatrix<_Scalar,_Options,_Index> >
+ : evaluator<SparseCompressedBase<MappedSparseMatrix<_Scalar,_Options,_Index> > >
{
- typedef MappedSparseMatrix<_Scalar,_Options,_Index> MappedSparseMatrixType;
- typedef typename MappedSparseMatrixType::InnerIterator InnerIterator;
- typedef typename MappedSparseMatrixType::ReverseInnerIterator ReverseInnerIterator;
-
- enum {
- CoeffReadCost = NumTraits<_Scalar>::ReadCost,
- Flags = MappedSparseMatrixType::Flags
- };
-
- evaluator() : m_matrix(0) {}
- explicit evaluator(const MappedSparseMatrixType &mat) : m_matrix(&mat) {}
-
- operator MappedSparseMatrixType&() { return m_matrix->const_cast_derived(); }
- operator const MappedSparseMatrixType&() const { return *m_matrix; }
+ typedef MappedSparseMatrix<_Scalar,_Options,_Index> XprType;
+ typedef evaluator<SparseCompressedBase<XprType> > Base;
- const MappedSparseMatrixType *m_matrix;
+ evaluator() : Base() {}
+ explicit evaluator(const XprType &mat) : Base(mat) {}
};
}
diff --git a/Eigen/src/SparseCore/SparseMap.h b/Eigen/src/SparseCore/SparseMap.h
index 91e8f7480..72dedb1ec 100644
--- a/Eigen/src/SparseCore/SparseMap.h
+++ b/Eigen/src/SparseCore/SparseMap.h
@@ -12,6 +12,32 @@
namespace Eigen {
+namespace internal {
+
+template<typename MatScalar, int MatOptions, typename MatIndex, int Options, typename StrideType>
+struct traits<Map<SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType> >
+ : public traits<SparseMatrix<MatScalar,MatOptions,MatIndex> >
+{
+ typedef SparseMatrix<MatScalar,MatOptions,MatIndex> PlainObjectType;
+ typedef traits<PlainObjectType> TraitsBase;
+ enum {
+ Flags = TraitsBase::Flags & (~NestByRefBit)
+ };
+};
+
+template<typename MatScalar, int MatOptions, typename MatIndex, int Options, typename StrideType>
+struct traits<Map<const SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType> >
+ : public traits<SparseMatrix<MatScalar,MatOptions,MatIndex> >
+{
+ typedef SparseMatrix<MatScalar,MatOptions,MatIndex> PlainObjectType;
+ typedef traits<PlainObjectType> TraitsBase;
+ enum {
+ Flags = TraitsBase::Flags & (~ (NestByRefBit | LvalueBit))
+ };
+};
+
+} // end namespace internal
+
template<typename Derived,
int Level = internal::accessors_level<Derived>::has_write_access ? WriteAccessors : ReadOnlyAccessors
> class SparseMapBase;
@@ -25,7 +51,7 @@ class SparseMapBase<Derived,ReadOnlyAccessors>
typedef typename Base::Scalar Scalar;
typedef typename Base::Index Index;
enum { IsRowMajor = Base::IsRowMajor };
-
+ using Base::operator=;
protected:
typedef typename internal::conditional<
@@ -103,6 +129,8 @@ class SparseMapBase<Derived,WriteAccessors>
typedef typename Base::Scalar Scalar;
typedef typename Base::Index Index;
enum { IsRowMajor = Base::IsRowMajor };
+
+ using Base::operator=;
public:
@@ -147,7 +175,7 @@ class Map<SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType>
: public SparseMapBase<Map<SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType> >
{
public:
- typedef SparseMapBase<Map<MatScalar, MatOptions, MatIndex> > Base;
+ typedef SparseMapBase<Map> Base;
_EIGEN_SPARSE_PUBLIC_INTERFACE(Map)
enum { IsRowMajor = Base::IsRowMajor };
@@ -167,7 +195,7 @@ class Map<const SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType
: public SparseMapBase<Map<const SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType> >
{
public:
- typedef SparseMapBase<Map<MatScalar, MatOptions, MatIndex> > Base;
+ typedef SparseMapBase<Map> Base;
_EIGEN_SPARSE_PUBLIC_INTERFACE(Map)
enum { IsRowMajor = Base::IsRowMajor };
diff --git a/Eigen/src/SparseCore/SparseRef.h b/Eigen/src/SparseCore/SparseRef.h
index cfea6ce8a..2ca039323 100644
--- a/Eigen/src/SparseCore/SparseRef.h
+++ b/Eigen/src/SparseCore/SparseRef.h
@@ -23,7 +23,7 @@ struct traits<Ref<SparseMatrix<MatScalar,MatOptions,MatIndex>, _Options, _Stride
typedef SparseMatrix<MatScalar,MatOptions,MatIndex> PlainObjectType;
enum {
Options = _Options,
- Flags = traits<MappedSparseMatrix<MatScalar,MatOptions,MatIndex> >::Flags | CompressedAccessBit | NestByRefBit
+ Flags = traits<SparseMatrix<MatScalar,MatOptions,MatIndex> >::Flags | CompressedAccessBit | NestByRefBit
};
template<typename Derived> struct match {
@@ -41,7 +41,7 @@ struct traits<Ref<const SparseMatrix<MatScalar,MatOptions,MatIndex>, _Options, _
: public traits<Ref<SparseMatrix<MatScalar,MatOptions,MatIndex>, _Options, _StrideType> >
{
enum {
- Flags = (traits<MappedSparseMatrix<MatScalar,MatOptions,MatIndex> >::Flags | CompressedAccessBit | NestByRefBit) & ~LvalueBit
+ Flags = (traits<SparseMatrix<MatScalar,MatOptions,MatIndex> >::Flags | CompressedAccessBit | NestByRefBit) & ~LvalueBit
};
};
@@ -49,11 +49,8 @@ template<typename Derived>
struct traits<SparseRefBase<Derived> > : public traits<Derived> {};
template<typename Derived> class SparseRefBase
-// : public MappedSparseMatrix<MatScalar,MatOptions,MatIndex>
: public SparseMapBase<Derived>
{
-// typedef typename internal::traits<Derived>::PlainObjectType PlainObjectType;
-
public:
typedef SparseMapBase<Derived> Base;
@@ -63,8 +60,6 @@ public:
: Base(RowsAtCompileTime==Dynamic?0:RowsAtCompileTime,ColsAtCompileTime==Dynamic?0:ColsAtCompileTime, 0, 0, 0, 0, 0)
{}
- EIGEN_INHERIT_ASSIGNMENT_OPERATORS(SparseRefBase)
-
protected:
@@ -119,9 +114,6 @@ class Ref<SparseMatrix<MatScalar,MatOptions,MatIndex>, Options, StrideType >
EIGEN_STATIC_ASSERT(bool(Traits::template match<Derived>::MatchAtCompileTime), STORAGE_LAYOUT_DOES_NOT_MATCH);
Base::construct(expr.const_cast_derived());
}
-
- EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Ref)
-
};
// this is the const ref version
diff --git a/test/sparse_basic.cpp b/test/sparse_basic.cpp
index 097959c84..a19de296a 100644
--- a/test/sparse_basic.cpp
+++ b/test/sparse_basic.cpp
@@ -408,6 +408,26 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
m.setFromTriplets(triplets.begin(), triplets.end());
VERIFY_IS_APPROX(m, refMat);
}
+
+ // test Map
+ {
+ DenseMatrix refMat2(rows, cols), refMat3(rows, cols);
+ SparseMatrixType m2(rows, cols), m3(rows, cols);
+ initSparse<Scalar>(density, refMat2, m2);
+ initSparse<Scalar>(density, refMat3, m3);
+ {
+ Map<SparseMatrixType> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
+ Map<SparseMatrixType> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr());
+ VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
+ VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
+ }
+ {
+ MappedSparseMatrix<Scalar,SparseMatrixType::Options,Index> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
+ MappedSparseMatrix<Scalar,SparseMatrixType::Options,Index> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr());
+ VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
+ VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
+ }
+ }
// test triangularView
{
diff --git a/test/sparse_ref.cpp b/test/sparse_ref.cpp
index b261cecf3..27700f827 100644
--- a/test/sparse_ref.cpp
+++ b/test/sparse_ref.cpp
@@ -69,8 +69,9 @@ void call_ref()
VERIFY_EVALUATION_COUNT( call_ref_2(A*A, A*A), 1);
Ref<SparseMatrix<float> > Ar(A);
- VERIFY_EVALUATION_COUNT( call_ref_1(Ar, Ar), 0);
- VERIFY_EVALUATION_COUNT( call_ref_2(Ar, Ar), 0);
+ VERIFY_IS_APPROX(Ar+Ar, A+A);
+ VERIFY_EVALUATION_COUNT( call_ref_1(Ar, A), 0);
+ VERIFY_EVALUATION_COUNT( call_ref_2(Ar, A), 0);
Ref<SparseMatrix<float,RowMajor> > Br(B);
VERIFY_EVALUATION_COUNT( call_ref_1(Br.transpose(), Br.transpose()), 0);