// This file is part of Eigen, a lightweight C++ template library // for linear algebra. Eigen itself is part of the KDE project. // // Copyright (C) 2008-2009 Gael Guennebaud // // Eigen is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public // License as published by the Free Software Foundation; either // version 3 of the License, or (at your option) any later version. // // Alternatively, you can redistribute it and/or // modify it under the terms of the GNU General Public License as // published by the Free Software Foundation; either version 2 of // the License, or (at your option) any later version. // // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the // GNU General Public License for more details. // // You should have received a copy of the GNU Lesser General Public // License and a copy of the GNU General Public License along with // Eigen. If not, see . #ifndef EIGEN_SPARSEMATRIX_H #define EIGEN_SPARSEMATRIX_H /** \class SparseMatrix * * \brief Sparse matrix * * \param _Scalar the scalar type, i.e. the type of the coefficients * * See http://www.netlib.org/linalg/html_templates/node91.html for details on the storage scheme. * */ template struct ei_traits > { typedef _Scalar Scalar; enum { RowsAtCompileTime = Dynamic, ColsAtCompileTime = Dynamic, MaxRowsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic, Flags = SparseBit | _Flags, CoeffReadCost = NumTraits::ReadCost, SupportedAccessPatterns = InnerRandomAccessPattern }; }; template class SparseMatrix : public SparseMatrixBase > { public: EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseMatrix) EIGEN_SPARSE_INHERIT_ASSIGNMENT_OPERATOR(SparseMatrix, +=) EIGEN_SPARSE_INHERIT_ASSIGNMENT_OPERATOR(SparseMatrix, -=) // FIXME: why are these operator already alvailable ??? // EIGEN_SPARSE_INHERIT_SCALAR_ASSIGNMENT_OPERATOR(SparseMatrix, *=) // EIGEN_SPARSE_INHERIT_SCALAR_ASSIGNMENT_OPERATOR(SparseMatrix, /=) typedef MappedSparseMatrix Map; protected: enum { IsRowMajor = Base::IsRowMajor }; typedef SparseMatrix TransposedSparseMatrix; int m_outerSize; int m_innerSize; int* m_outerIndex; CompressedStorage m_data; public: inline int rows() const { return IsRowMajor ? m_outerSize : m_innerSize; } inline int cols() const { return IsRowMajor ? m_innerSize : m_outerSize; } inline int innerSize() const { return m_innerSize; } inline int outerSize() const { return m_outerSize; } inline int innerNonZeros(int j) const { return m_outerIndex[j+1]-m_outerIndex[j]; } inline const Scalar* _valuePtr() const { return &m_data.value(0); } inline Scalar* _valuePtr() { return &m_data.value(0); } inline const int* _innerIndexPtr() const { return &m_data.index(0); } inline int* _innerIndexPtr() { return &m_data.index(0); } inline const int* _outerIndexPtr() const { return m_outerIndex; } inline int* _outerIndexPtr() { return m_outerIndex; } inline Scalar coeff(int row, int col) const { const int outer = IsRowMajor ? row : col; const int inner = IsRowMajor ? col : row; return m_data.atInRange(m_outerIndex[outer], m_outerIndex[outer+1], inner); } inline Scalar& coeffRef(int row, int col) { const int outer = IsRowMajor ? row : col; const int inner = IsRowMajor ? col : row; int start = m_outerIndex[outer]; int end = m_outerIndex[outer+1]; ei_assert(end>=start && "you probably called coeffRef on a non finalized matrix"); ei_assert(end>start && "coeffRef cannot be called on a zero coefficient"); const int id = m_data.searchLowerIndex(start,end-1,inner); ei_assert((id=0 && m_outerIndex[i]==0) { m_outerIndex[i] = m_data.size(); --i; } m_outerIndex[outer+1] = m_outerIndex[outer]; } assert(size_t(m_outerIndex[outer+1]) == m_data.size()); int id = m_outerIndex[outer+1]; ++m_outerIndex[outer+1]; m_data.append(0, inner); return m_data.value(id); } /** Like fill() but with random inner coordinates. */ inline Scalar& fillrand(int row, int col) { const int outer = IsRowMajor ? row : col; const int inner = IsRowMajor ? col : row; if (m_outerIndex[outer+1]==0) { // we start a new inner vector // nothing special to do here int i = outer; while (i>=0 && m_outerIndex[i]==0) { m_outerIndex[i] = m_data.size(); --i; } m_outerIndex[outer+1] = m_outerIndex[outer]; } assert(size_t(m_outerIndex[outer+1]) == m_data.size() && "invalid outer index"); size_t startId = m_outerIndex[outer]; // FIXME let's make sure sizeof(long int) == sizeof(size_t) size_t id = m_outerIndex[outer+1]; ++m_outerIndex[outer+1]; float reallocRatio = 1; if (m_data.allocatedSize() startId) && (m_data.index(id-1) > inner) ) { m_data.index(id) = m_data.index(id-1); m_data.value(id) = m_data.value(id-1); --id; } m_data.index(id) = inner; return (m_data.value(id) = 0); } inline void endFill() { int size = m_data.size(); int i = m_outerSize; // find the last filled column while (i>=0 && m_outerIndex[i]==0) --i; ++i; while (i<=m_outerSize) { m_outerIndex[i] = size; ++i; } } void prune(Scalar reference, RealScalar epsilon = precision()) { int k = 0; for (int j=0; j inline SparseMatrix(const SparseMatrixBase& other) : m_outerSize(0), m_innerSize(0), m_outerIndex(0) { *this = other.derived(); } inline SparseMatrix(const SparseMatrix& other) : Base(), m_outerSize(0), m_innerSize(0), m_outerIndex(0) { *this = other.derived(); } inline void swap(SparseMatrix& other) { //EIGEN_DBG_SPARSE(std::cout << "SparseMatrix:: swap\n"); std::swap(m_outerIndex, other.m_outerIndex); std::swap(m_innerSize, other.m_innerSize); std::swap(m_outerSize, other.m_outerSize); m_data.swap(other.m_data); } inline SparseMatrix& operator=(const SparseMatrix& other) { // std::cout << "SparseMatrix& operator=(const SparseMatrix& other)\n"; if (other.isRValue()) { swap(other.const_cast_derived()); } else { resize(other.rows(), other.cols()); memcpy(m_outerIndex, other.m_outerIndex, (m_outerSize+1)*sizeof(int)); m_data = other.m_data; } return *this; } template inline SparseMatrix& operator=(const SparseMatrixBase& other) { const bool needToTranspose = (Flags & RowMajorBit) != (OtherDerived::Flags & RowMajorBit); if (needToTranspose) { // two passes algorithm: // 1 - compute the number of coeffs per dest inner vector // 2 - do the actual copy/eval // Since each coeff of the rhs has to be evaluated twice, let's evauluate it if needed //typedef typename ei_nested::type OtherCopy; typedef typename ei_eval::type OtherCopy; typedef typename ei_cleantype::type _OtherCopy; OtherCopy otherCopy(other.derived()); resize(other.rows(), other.cols()); Eigen::Map(m_outerIndex,outerSize()).setZero(); // pass 1 // FIXME the above copy could be merged with that pass for (int j=0; j::operator=(other.derived()); } } friend std::ostream & operator << (std::ostream & s, const SparseMatrix& m) { EIGEN_DBG_SPARSE( s << "Nonzero entries:\n"; for (int i=0; i&>(m); return s; } /** Destructor */ inline ~SparseMatrix() { delete[] m_outerIndex; } }; template class SparseMatrix::InnerIterator { public: InnerIterator(const SparseMatrix& mat, int outer) : m_matrix(mat), m_outer(outer), m_id(mat.m_outerIndex[outer]), m_start(m_id), m_end(mat.m_outerIndex[outer+1]) {} template InnerIterator(const Flagged& mat, int outer) : m_matrix(mat._expression()), m_outer(outer), m_id(m_matrix.m_outerIndex[outer]), m_start(m_id), m_end(m_matrix.m_outerIndex[outer+1]) {} inline InnerIterator& operator++() { m_id++; return *this; } inline Scalar value() const { return m_matrix.m_data.value(m_id); } inline Scalar& valueRef() { return const_cast(m_matrix.m_data.value(m_id)); } inline int index() const { return m_matrix.m_data.index(m_id); } inline int row() const { return IsRowMajor ? m_outer : index(); } inline int col() const { return IsRowMajor ? index() : m_outer; } inline operator bool() const { return (m_id < m_end) && (m_id>=m_start); } protected: const SparseMatrix& m_matrix; const int m_outer; int m_id; const int m_start; const int m_end; }; #endif // EIGEN_SPARSEMATRIX_H