aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--CMakeLists.txt12
-rw-r--r--CTestConfig.cmake6
-rw-r--r--Eigen/src/LU/PartialLU.h4
-rw-r--r--Eigen/src/Sparse/DynamicSparseMatrix.h30
-rw-r--r--Eigen/src/Sparse/SparseBlock.h27
-rw-r--r--Eigen/src/Sparse/SparseCwiseBinaryOp.h4
-rw-r--r--Eigen/src/Sparse/SparseMatrix.h32
-rw-r--r--Eigen/src/Sparse/SparseMatrixBase.h35
-rw-r--r--Eigen/src/Sparse/SparseRedux.h16
-rw-r--r--Eigen/src/Sparse/SparseTranspose.h10
-rw-r--r--Eigen/src/Sparse/SparseVector.h27
-rw-r--r--eigen2.pc.in7
-rw-r--r--test/testsuite.cmake6
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h21
14 files changed, 160 insertions, 77 deletions
diff --git a/CMakeLists.txt b/CMakeLists.txt
index 1d84e9107..dc26e71d1 100644
--- a/CMakeLists.txt
+++ b/CMakeLists.txt
@@ -37,6 +37,9 @@ if(NOT WIN32)
option(EIGEN_BUILD_LIB "Build the binary shared library" OFF)
endif(NOT WIN32)
option(EIGEN_BUILD_BTL "Build benchmark suite" OFF)
+if(NOT WIN32)
+ option(EIGEN_BUILD_PKGCONFIG "Build pkg-config .pc file for Eigen" ON)
+endif(NOT WIN32)
if(EIGEN_BUILD_LIB)
option(EIGEN_TEST_LIB "Build the unit tests using the library (disable -pedantic)" OFF)
@@ -108,6 +111,13 @@ set(INCLUDE_INSTALL_DIR
"The directory where we install the header files"
FORCE)
+if(EIGEN_BUILD_PKGCONFIG)
+ configure_file(eigen2.pc.in eigen2.pc)
+ install(FILES eigen2.pc
+ DESTINATION lib/pkgconfig
+ )
+endif(EIGEN_BUILD_PKGCONFIG)
+
add_subdirectory(Eigen)
add_subdirectory(doc)
@@ -129,4 +139,4 @@ endif(EIGEN_BUILD_BTL)
if(EIGEN_BUILD_TESTS)
ei_testing_print_summary()
-endif(EIGEN_BUILD_TESTS) \ No newline at end of file
+endif(EIGEN_BUILD_TESTS)
diff --git a/CTestConfig.cmake b/CTestConfig.cmake
index be43daca8..263cda480 100644
--- a/CTestConfig.cmake
+++ b/CTestConfig.cmake
@@ -3,11 +3,11 @@
## project to incorporate the testing dashboard.
## # The following are required to uses Dart and the Cdash dashboard
## ENABLE_TESTING()
-## INCLUDE(Dart)
+## INCLUDE(CTest)
set(CTEST_PROJECT_NAME "Eigen")
set(CTEST_NIGHTLY_START_TIME "05:00:00 UTC")
set(CTEST_DROP_METHOD "http")
-set(CTEST_DROP_SITE "www.cdash.org")
-set(CTEST_DROP_LOCATION "/CDashPublic/submit.php?project=Eigen")
+set(CTEST_DROP_SITE "my.cdash.org")
+set(CTEST_DROP_LOCATION "/submit.php?project=Eigen")
set(CTEST_DROP_SITE_CDASH TRUE)
diff --git a/Eigen/src/LU/PartialLU.h b/Eigen/src/LU/PartialLU.h
index 7fdbeac38..3a6e9f286 100644
--- a/Eigen/src/LU/PartialLU.h
+++ b/Eigen/src/LU/PartialLU.h
@@ -225,8 +225,8 @@ void PartialLU<MatrixType>::solve(
/* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
* So we proceed as follows:
* Step 1: compute c = Pb.
- * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
- * Step 3: replace c by the solution x to Ux = c. Check if a solution really exists.
+ * Step 2: replace c by the solution x to Lx = c.
+ * Step 3: replace c by the solution x to Ux = c.
*/
const int size = m_lu.rows();
diff --git a/Eigen/src/Sparse/DynamicSparseMatrix.h b/Eigen/src/Sparse/DynamicSparseMatrix.h
index 2927cd583..ec730dd3f 100644
--- a/Eigen/src/Sparse/DynamicSparseMatrix.h
+++ b/Eigen/src/Sparse/DynamicSparseMatrix.h
@@ -35,7 +35,7 @@
* random read/write accesses in log(rho*outer_size) where \c rho is the probability that a coefficient is
* nonzero and outer_size is the number of columns if the matrix is column-major and the number of rows
* otherwise.
- *
+ *
* Internally, the data are stored as a std::vector of compressed vector. The performances of random writes might
* decrease as the number of nonzeros per inner-vector increase. In practice, we observed very good performance
* till about 100 nonzeros/vector, and the performance remains relatively good till 500 nonzeros/vectors.
@@ -67,23 +67,23 @@ class DynamicSparseMatrix
// EIGEN_SPARSE_INHERIT_ASSIGNMENT_OPERATOR(DynamicSparseMatrix, +=)
// EIGEN_SPARSE_INHERIT_ASSIGNMENT_OPERATOR(DynamicSparseMatrix, -=)
typedef MappedSparseMatrix<Scalar,Flags> Map;
+ using Base::IsRowMajor;
protected:
- enum { IsRowMajor = Base::IsRowMajor };
typedef DynamicSparseMatrix<Scalar,(Flags&~RowMajorBit)|(IsRowMajor?RowMajorBit:0)> TransposedSparseMatrix;
int m_innerSize;
std::vector<CompressedStorage<Scalar> > m_data;
public:
-
+
inline int rows() const { return IsRowMajor ? outerSize() : m_innerSize; }
inline int cols() const { return IsRowMajor ? m_innerSize : outerSize(); }
inline int innerSize() const { return m_innerSize; }
inline int outerSize() const { return m_data.size(); }
inline int innerNonZeros(int j) const { return m_data[j].size(); }
-
+
std::vector<CompressedStorage<Scalar> >& _data() { return m_data; }
const std::vector<CompressedStorage<Scalar> >& _data() const { return m_data; }
@@ -132,7 +132,7 @@ class DynamicSparseMatrix
setZero();
reserve(reserveSize);
}
-
+
void reserve(int reserveSize = 1000)
{
if (outerSize()>0)
@@ -144,13 +144,13 @@ class DynamicSparseMatrix
}
}
}
-
+
inline void startVec(int /*outer*/) {}
-
+
inline Scalar& insertBack(int outer, int inner)
{
ei_assert(outer<int(m_data.size()) && inner<m_innerSize && "out of range");
- ei_assert(((m_data[outer].size()==0) || (m_data[outer].index(m_data[outer].size()-1)<inner))
+ ei_assert(((m_data[outer].size()==0) || (m_data[outer].index(m_data[outer].size()-1)<inner))
&& "wrong sorted insertion");
m_data[outer].append(0, inner);
return m_data[outer].value(m_data[outer].size()-1);
@@ -162,7 +162,7 @@ class DynamicSparseMatrix
* 2 - this the coefficient with greater inner coordinate for the given outer coordinate.
* In other words, assuming \c *this is column-major, then there must not exists any nonzero coefficient of coordinates
* \c i \c x \a col such that \c i >= \a row. Otherwise the matrix is invalid.
- *
+ *
* \see fillrand(), coeffRef()
*/
EIGEN_DEPRECATED Scalar& fill(int row, int col)
@@ -181,12 +181,12 @@ class DynamicSparseMatrix
{
return insert(row,col);
}
-
+
inline Scalar& insert(int row, int col)
{
const int outer = IsRowMajor ? row : col;
const int inner = IsRowMajor ? col : row;
-
+
int startId = 0;
int id = m_data[outer].size() - 1;
m_data[outer].resize(id+2,1);
@@ -205,9 +205,9 @@ class DynamicSparseMatrix
/** \deprecated use finalize()
* Does nothing. Provided for compatibility with SparseMatrix. */
EIGEN_DEPRECATED void endFill() {}
-
+
inline void finalize() {}
-
+
void prune(Scalar reference, RealScalar epsilon = precision<RealScalar>())
{
for (int j=0; j<outerSize(); ++j)
@@ -226,7 +226,7 @@ class DynamicSparseMatrix
m_data.resize(outerSize);
}
}
-
+
void resizeAndKeepData(int rows, int cols)
{
const int outerSize = IsRowMajor ? rows : cols;
@@ -309,7 +309,7 @@ class DynamicSparseMatrix<Scalar,_Flags>::InnerIterator : public SparseVector<Sc
InnerIterator(const DynamicSparseMatrix& mat, int outer)
: Base(mat.m_data[outer]), m_outer(outer)
{}
-
+
inline int row() const { return IsRowMajor ? m_outer : Base::index(); }
inline int col() const { return IsRowMajor ? Base::index() : m_outer; }
diff --git a/Eigen/src/Sparse/SparseBlock.h b/Eigen/src/Sparse/SparseBlock.h
index 72589d0c0..655ba6d93 100644
--- a/Eigen/src/Sparse/SparseBlock.h
+++ b/Eigen/src/Sparse/SparseBlock.h
@@ -43,16 +43,21 @@ template<typename MatrixType, int Size>
class SparseInnerVectorSet : ei_no_assignment_operator,
public SparseMatrixBase<SparseInnerVectorSet<MatrixType, Size> >
{
- enum { IsRowMajor = ei_traits<SparseInnerVectorSet>::IsRowMajor };
public:
+ enum { IsRowMajor = ei_traits<SparseInnerVectorSet>::IsRowMajor };
+
EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseInnerVectorSet)
class InnerIterator: public MatrixType::InnerIterator
{
public:
inline InnerIterator(const SparseInnerVectorSet& xpr, int outer)
- : MatrixType::InnerIterator(xpr.m_matrix, xpr.m_outerStart + outer)
+ : MatrixType::InnerIterator(xpr.m_matrix, xpr.m_outerStart + outer), m_outer(outer)
{}
+ inline int row() const { return IsRowMajor ? m_outer : this->index(); }
+ inline int col() const { return IsRowMajor ? this->index() : m_outer; }
+ protected:
+ int m_outer;
};
inline SparseInnerVectorSet(const MatrixType& matrix, int outerStart, int outerSize)
@@ -100,16 +105,21 @@ class SparseInnerVectorSet<DynamicSparseMatrix<_Scalar, _Options>, Size>
: public SparseMatrixBase<SparseInnerVectorSet<DynamicSparseMatrix<_Scalar, _Options>, Size> >
{
typedef DynamicSparseMatrix<_Scalar, _Options> MatrixType;
- enum { IsRowMajor = ei_traits<SparseInnerVectorSet>::IsRowMajor };
public:
+ enum { IsRowMajor = ei_traits<SparseInnerVectorSet>::IsRowMajor };
+
EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseInnerVectorSet)
class InnerIterator: public MatrixType::InnerIterator
{
public:
inline InnerIterator(const SparseInnerVectorSet& xpr, int outer)
- : MatrixType::InnerIterator(xpr.m_matrix, xpr.m_outerStart + outer)
+ : MatrixType::InnerIterator(xpr.m_matrix, xpr.m_outerStart + outer), m_outer(outer)
{}
+ inline int row() const { return IsRowMajor ? m_outer : this->index(); }
+ inline int col() const { return IsRowMajor ? this->index() : m_outer; }
+ protected:
+ int m_outer;
};
inline SparseInnerVectorSet(const MatrixType& matrix, int outerStart, int outerSize)
@@ -193,16 +203,21 @@ class SparseInnerVectorSet<SparseMatrix<_Scalar, _Options>, Size>
: public SparseMatrixBase<SparseInnerVectorSet<SparseMatrix<_Scalar, _Options>, Size> >
{
typedef SparseMatrix<_Scalar, _Options> MatrixType;
- enum { IsRowMajor = ei_traits<SparseInnerVectorSet>::IsRowMajor };
public:
+ enum { IsRowMajor = ei_traits<SparseInnerVectorSet>::IsRowMajor };
+
EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseInnerVectorSet)
class InnerIterator: public MatrixType::InnerIterator
{
public:
inline InnerIterator(const SparseInnerVectorSet& xpr, int outer)
- : MatrixType::InnerIterator(xpr.m_matrix, xpr.m_outerStart + outer)
+ : MatrixType::InnerIterator(xpr.m_matrix, xpr.m_outerStart + outer), m_outer(outer)
{}
+ inline int row() const { return IsRowMajor ? m_outer : this->index(); }
+ inline int col() const { return IsRowMajor ? this->index() : m_outer; }
+ protected:
+ int m_outer;
};
inline SparseInnerVectorSet(const MatrixType& matrix, int outerStart, int outerSize)
diff --git a/Eigen/src/Sparse/SparseCwiseBinaryOp.h b/Eigen/src/Sparse/SparseCwiseBinaryOp.h
index d19970efc..96fd1d36a 100644
--- a/Eigen/src/Sparse/SparseCwiseBinaryOp.h
+++ b/Eigen/src/Sparse/SparseCwiseBinaryOp.h
@@ -186,8 +186,8 @@ class ei_sparse_cwise_binary_op_inner_iterator_selector<BinaryOp, Lhs, Rhs, Deri
EIGEN_STRONG_INLINE Scalar value() const { return m_value; }
EIGEN_STRONG_INLINE int index() const { return m_id; }
- EIGEN_STRONG_INLINE int row() const { return m_lhsIter.row(); }
- EIGEN_STRONG_INLINE int col() const { return m_lhsIter.col(); }
+ EIGEN_STRONG_INLINE int row() const { return Lhs::IsRowMajor ? m_lhsIter.row() : index(); }
+ EIGEN_STRONG_INLINE int col() const { return Lhs::IsRowMajor ? index() : m_lhsIter.col(); }
EIGEN_STRONG_INLINE operator bool() const { return m_id>=0; }
diff --git a/Eigen/src/Sparse/SparseMatrix.h b/Eigen/src/Sparse/SparseMatrix.h
index 8bcf3c8b9..65cee69cf 100644
--- a/Eigen/src/Sparse/SparseMatrix.h
+++ b/Eigen/src/Sparse/SparseMatrix.h
@@ -25,17 +25,24 @@
#ifndef EIGEN_SPARSEMATRIX_H
#define EIGEN_SPARSEMATRIX_H
-/** \class SparseMatrix
+/** \ingroup Sparse_Module
*
- * \brief Sparse matrix
+ * \class SparseMatrix
+ *
+ * \brief The main sparse matrix class
+ *
+ * This class implements a sparse matrix using the very common compressed row/column storage
+ * scheme.
*
* \param _Scalar the scalar type, i.e. the type of the coefficients
+ * \param _Options Union of bit flags controlling the storage scheme. Currently the only possibility
+ * is RowMajor. The default is 0 which means column-major.
*
* See http://www.netlib.org/linalg/html_templates/node91.html for details on the storage scheme.
*
*/
-template<typename _Scalar, int _Flags>
-struct ei_traits<SparseMatrix<_Scalar, _Flags> >
+template<typename _Scalar, int _Options>
+struct ei_traits<SparseMatrix<_Scalar, _Options> >
{
typedef _Scalar Scalar;
enum {
@@ -43,17 +50,15 @@ struct ei_traits<SparseMatrix<_Scalar, _Flags> >
ColsAtCompileTime = Dynamic,
MaxRowsAtCompileTime = Dynamic,
MaxColsAtCompileTime = Dynamic,
- Flags = SparseBit | _Flags,
+ Flags = SparseBit | _Options,
CoeffReadCost = NumTraits<Scalar>::ReadCost,
SupportedAccessPatterns = InnerRandomAccessPattern
};
};
-
-
-template<typename _Scalar, int _Flags>
+template<typename _Scalar, int _Options>
class SparseMatrix
- : public SparseMatrixBase<SparseMatrix<_Scalar, _Flags> >
+ : public SparseMatrixBase<SparseMatrix<_Scalar, _Options> >
{
public:
EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseMatrix)
@@ -64,10 +69,10 @@ class SparseMatrix
// EIGEN_SPARSE_INHERIT_SCALAR_ASSIGNMENT_OPERATOR(SparseMatrix, /=)
typedef MappedSparseMatrix<Scalar,Flags> Map;
+ using Base::IsRowMajor;
protected:
- enum { IsRowMajor = Base::IsRowMajor };
typedef SparseMatrix<Scalar,(Flags&~RowMajorBit)|(IsRowMajor?RowMajorBit:0)> TransposedSparseMatrix;
int m_outerSize;
@@ -508,10 +513,13 @@ class SparseMatrix
{
delete[] m_outerIndex;
}
+
+ /** Overloaded for performance */
+ Scalar sum() const;
};
-template<typename Scalar, int _Flags>
-class SparseMatrix<Scalar,_Flags>::InnerIterator
+template<typename Scalar, int _Options>
+class SparseMatrix<Scalar,_Options>::InnerIterator
{
public:
InnerIterator(const SparseMatrix& mat, int outer)
diff --git a/Eigen/src/Sparse/SparseMatrixBase.h b/Eigen/src/Sparse/SparseMatrixBase.h
index 948424f90..34779d9e6 100644
--- a/Eigen/src/Sparse/SparseMatrixBase.h
+++ b/Eigen/src/Sparse/SparseMatrixBase.h
@@ -25,6 +25,17 @@
#ifndef EIGEN_SPARSEMATRIXBASE_H
#define EIGEN_SPARSEMATRIXBASE_H
+/** \ingroup Sparse_Module
+ *
+ * \class SparseMatrixBase
+ *
+ * \brief Base class of any sparse matrices or sparse expressions
+ *
+ * \param Derived
+ *
+ *
+ *
+ */
template<typename Derived> class SparseMatrixBase
{
public:
@@ -69,7 +80,7 @@ template<typename Derived> class SparseMatrixBase
/**< This stores expression \ref flags flags which may or may not be inherited by new expressions
* constructed from this one. See the \ref flags "list of flags".
*/
-
+
CoeffReadCost = ei_traits<Derived>::CoeffReadCost,
/**< This is a rough measure of how expensive it is to read one coefficient from
* this expression.
@@ -156,9 +167,9 @@ template<typename Derived> class SparseMatrixBase
ei_assert(( ((ei_traits<Derived>::SupportedAccessPatterns&OuterRandomAccessPattern)==OuterRandomAccessPattern) ||
(!((Flags & RowMajorBit) != (OtherDerived::Flags & RowMajorBit)))) &&
"the transpose operation is supposed to be handled in SparseMatrix::operator=");
-
+
enum { Flip = (Flags & RowMajorBit) != (OtherDerived::Flags & RowMajorBit) };
-
+
const int outerSize = other.outerSize();
//typedef typename ei_meta_if<transpose, LinkedVectorMatrix<Scalar,Flags&RowMajorBit>, Derived>::ret TempType;
// thanks to shallow copies, we always eval to a tempary
@@ -293,13 +304,13 @@ template<typename Derived> class SparseMatrixBase
template<typename OtherDerived>
const typename SparseProductReturnType<Derived,OtherDerived>::Type
operator*(const SparseMatrixBase<OtherDerived> &other) const;
-
+
// dense * sparse (return a dense object)
- template<typename OtherDerived> friend
+ template<typename OtherDerived> friend
const typename SparseProductReturnType<OtherDerived,Derived>::Type
operator*(const MatrixBase<OtherDerived>& lhs, const Derived& rhs)
{ return typename SparseProductReturnType<OtherDerived,Derived>::Type(lhs.derived(),rhs); }
-
+
template<typename OtherDerived>
const typename SparseProductReturnType<Derived,OtherDerived>::Type
operator*(const MatrixBase<OtherDerived> &other) const;
@@ -340,7 +351,7 @@ template<typename Derived> class SparseMatrixBase
const SparseInnerVectorSet<Derived,1> col(int j) const;
SparseInnerVectorSet<Derived,1> innerVector(int outer);
const SparseInnerVectorSet<Derived,1> innerVector(int outer) const;
-
+
// set of sub-vectors
SparseInnerVectorSet<Derived,Dynamic> subrows(int start, int size);
const SparseInnerVectorSet<Derived,Dynamic> subrows(int start, int size) const;
@@ -432,10 +443,16 @@ template<typename Derived> class SparseMatrixBase
for (int j=0; j<outerSize(); ++j)
{
for (typename Derived::InnerIterator i(derived(),j); i; ++i)
+ {
if(IsRowMajor)
- res.coeffRef(j,i.index()) = i.value();
+ std::cerr << i.row() << "," << i.col() << " == " << j << "," << i.index() << "\n";
else
- res.coeffRef(i.index(),j) = i.value();
+ std::cerr << i.row() << "," << i.col() << " == " << i.index() << "," << j << "\n";
+// if(IsRowMajor)
+ res.coeffRef(i.row(),i.col()) = i.value();
+// else
+// res.coeffRef(i.index(),j) = i.value();
+ }
}
return res;
}
diff --git a/Eigen/src/Sparse/SparseRedux.h b/Eigen/src/Sparse/SparseRedux.h
index f0d370548..e223c37ec 100644
--- a/Eigen/src/Sparse/SparseRedux.h
+++ b/Eigen/src/Sparse/SparseRedux.h
@@ -37,4 +37,20 @@ SparseMatrixBase<Derived>::sum() const
return res;
}
+template<typename _Scalar, int _Options>
+typename ei_traits<SparseMatrix<_Scalar,_Options> >::Scalar
+SparseMatrix<_Scalar,_Options>::sum() const
+{
+ ei_assert(rows()>0 && cols()>0 && "you are using a non initialized matrix");
+ return Matrix<Scalar,1,Dynamic>::Map(m_data.value(0), m_data.size()).sum();
+}
+
+template<typename _Scalar, int _Options>
+typename ei_traits<SparseVector<_Scalar,_Options> >::Scalar
+SparseVector<_Scalar,_Options>::sum() const
+{
+ ei_assert(rows()>0 && cols()>0 && "you are using a non initialized matrix");
+ return Matrix<Scalar,1,Dynamic>::Map(m_data.value(0), m_data.size()).sum();
+}
+
#endif // EIGEN_SPARSEREDUX_H
diff --git a/Eigen/src/Sparse/SparseTranspose.h b/Eigen/src/Sparse/SparseTranspose.h
index 89a14d707..879c377a9 100644
--- a/Eigen/src/Sparse/SparseTranspose.h
+++ b/Eigen/src/Sparse/SparseTranspose.h
@@ -66,20 +66,26 @@ template<typename MatrixType> class SparseTranspose
template<typename MatrixType> class SparseTranspose<MatrixType>::InnerIterator : public MatrixType::InnerIterator
{
+ typedef typename MatrixType::InnerIterator Base;
public:
EIGEN_STRONG_INLINE InnerIterator(const SparseTranspose& trans, int outer)
- : MatrixType::InnerIterator(trans.m_matrix, outer)
+ : Base(trans.m_matrix, outer)
{}
+ inline int row() const { return Base::col(); }
+ inline int col() const { return Base::row(); }
};
template<typename MatrixType> class SparseTranspose<MatrixType>::ReverseInnerIterator : public MatrixType::ReverseInnerIterator
{
+ typedef typename MatrixType::ReverseInnerIterator Base;
public:
EIGEN_STRONG_INLINE ReverseInnerIterator(const SparseTranspose& xpr, int outer)
- : MatrixType::ReverseInnerIterator(xpr.m_matrix, outer)
+ : Base(xpr.m_matrix, outer)
{}
+ inline int row() const { return Base::col(); }
+ inline int col() const { return Base::row(); }
};
#endif // EIGEN_SPARSETRANSPOSE_H
diff --git a/Eigen/src/Sparse/SparseVector.h b/Eigen/src/Sparse/SparseVector.h
index ac815e698..05356942b 100644
--- a/Eigen/src/Sparse/SparseVector.h
+++ b/Eigen/src/Sparse/SparseVector.h
@@ -34,26 +34,26 @@
* See http://www.netlib.org/linalg/html_templates/node91.html for details on the storage scheme.
*
*/
-template<typename _Scalar, int _Flags>
-struct ei_traits<SparseVector<_Scalar, _Flags> >
+template<typename _Scalar, int _Options>
+struct ei_traits<SparseVector<_Scalar, _Options> >
{
typedef _Scalar Scalar;
enum {
- IsColVector = _Flags & RowMajorBit ? 0 : 1,
+ IsColVector = _Options & RowMajorBit ? 0 : 1,
RowsAtCompileTime = IsColVector ? Dynamic : 1,
ColsAtCompileTime = IsColVector ? 1 : Dynamic,
MaxRowsAtCompileTime = RowsAtCompileTime,
MaxColsAtCompileTime = ColsAtCompileTime,
- Flags = SparseBit | _Flags,
+ Flags = SparseBit | _Options,
CoeffReadCost = NumTraits<Scalar>::ReadCost,
SupportedAccessPatterns = InnerRandomAccessPattern
};
};
-template<typename _Scalar, int _Flags>
+template<typename _Scalar, int _Options>
class SparseVector
- : public SparseMatrixBase<SparseVector<_Scalar, _Flags> >
+ : public SparseMatrixBase<SparseVector<_Scalar, _Options> >
{
public:
EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseVector)
@@ -124,7 +124,7 @@ class SparseVector
{
ei_assert(outer==0);
}
-
+
inline Scalar& insertBack(int outer, int inner)
{
ei_assert(outer==0);
@@ -183,7 +183,7 @@ class SparseVector
m_data.append(0, i);
return m_data.value(m_data.size()-1);
}
-
+
/** \deprecated use insert(int,int) */
EIGEN_DEPRECATED Scalar& fillrand(int r, int c)
{
@@ -196,11 +196,11 @@ class SparseVector
{
return insert(i);
}
-
+
/** \deprecated use finalize() */
EIGEN_DEPRECATED void endFill() {}
inline void finalize() {}
-
+
void prune(Scalar reference, RealScalar epsilon = precision<RealScalar>())
{
m_data.prune(reference,epsilon);
@@ -357,10 +357,13 @@ class SparseVector
/** Destructor */
inline ~SparseVector() {}
+
+ /** Overloaded for performance */
+ Scalar sum() const;
};
-template<typename Scalar, int _Flags>
-class SparseVector<Scalar,_Flags>::InnerIterator
+template<typename Scalar, int _Options>
+class SparseVector<Scalar,_Options>::InnerIterator
{
public:
InnerIterator(const SparseVector& vec, int outer=0)
diff --git a/eigen2.pc.in b/eigen2.pc.in
new file mode 100644
index 000000000..b508a58d3
--- /dev/null
+++ b/eigen2.pc.in
@@ -0,0 +1,7 @@
+
+Name: Eigen2
+Description: A C++ template library for linear algebra: vectors, matrices, and related algorithms
+Requires:
+Version: ${EIGEN_VERSION_NUMBER}
+Libs:
+Cflags: -I${INCLUDE_INSTALL_DIR}
diff --git a/test/testsuite.cmake b/test/testsuite.cmake
index 9bbe01dcb..f3b0e4172 100644
--- a/test/testsuite.cmake
+++ b/test/testsuite.cmake
@@ -39,7 +39,7 @@
# VERSION=opensuse-11.1
# WORK_DIR=/home/gael/Coding/eigen2/cdash
# # get the last version of the script
-# svn cat svn://anonsvn.kde.org/home/kde/trunk/kdesupport/eigen2/test/testsuite.cmake > $WORK_DIR/testsuite.cmake
+# wget http://bitbucket.org/eigen/eigen2/raw/tip/test/testsuite.cmake -o $WORK_DIR/testsuite.cmake
# COMMON="ctest -S $WORK_DIR/testsuite.cmake,EIGEN_WORK_DIR=$WORK_DIR,EIGEN_SITE=$SITE,EIGEN_MODE=$1,EIGEN_BUILD_STRING=$OS_VERSION-$ARCH"
# $COMMON-gcc-3.4.6,EIGEN_CXX=g++-3.4
# $COMMON-gcc-4.0.1,EIGEN_CXX=g++-4.0.1
@@ -132,8 +132,8 @@ endif(NOT EIGEN_MODE)
## mandatory variables (the default should be ok in most cases):
-SET (CTEST_CVS_COMMAND "svn")
-SET (CTEST_CVS_CHECKOUT "${CTEST_CVS_COMMAND} co svn://anonsvn.kde.org/home/kde/trunk/kdesupport/eigen2 \"${CTEST_SOURCE_DIRECTORY}\"")
+SET (CTEST_CVS_COMMAND "hg")
+SET (CTEST_CVS_CHECKOUT "${CTEST_CVS_COMMAND} clone http://bitbucket.org/eigen/eigen2 \"${CTEST_SOURCE_DIRECTORY}\"")
# which ctest command to use for running the dashboard
SET (CTEST_COMMAND "${EIGEN_CMAKE_DIR}ctest -D ${EIGEN_MODE}")
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h b/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
index bc1bcceb2..7950d22b4 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
@@ -72,7 +72,8 @@ void ei_matrix_exponential(const MatrixBase<Derived> &M, typename ei_plain_matri
PlainMatrixType num, den, U, V;
PlainMatrixType Id = PlainMatrixType::Identity(M.rows(), M.cols());
- RealScalar l1norm = M.cwise().abs().colwise().sum().maxCoeff();
+ typename ei_eval<Derived>::type Meval = M.eval();
+ RealScalar l1norm = Meval.cwise().abs().colwise().sum().maxCoeff();
int squarings = 0;
// Choose degree of Pade approximant, depending on norm of M
@@ -81,9 +82,9 @@ void ei_matrix_exponential(const MatrixBase<Derived> &M, typename ei_plain_matri
// Use (3,3)-Pade
const Scalar b[] = {120., 60., 12., 1.};
PlainMatrixType M2;
- M2 = (M * M).lazy();
+ M2 = (Meval * Meval).lazy();
num = b[3]*M2 + b[1]*Id;
- U = (M * num).lazy();
+ U = (Meval * num).lazy();
V = b[2]*M2 + b[0]*Id;
} else if (l1norm < 2.539398330063230e-001) {
@@ -91,10 +92,10 @@ void ei_matrix_exponential(const MatrixBase<Derived> &M, typename ei_plain_matri
// Use (5,5)-Pade
const Scalar b[] = {30240., 15120., 3360., 420., 30., 1.};
PlainMatrixType M2, M4;
- M2 = (M * M).lazy();
+ M2 = (Meval * Meval).lazy();
M4 = (M2 * M2).lazy();
num = b[5]*M4 + b[3]*M2 + b[1]*Id;
- U = (M * num).lazy();
+ U = (Meval * num).lazy();
V = b[4]*M4 + b[2]*M2 + b[0]*Id;
} else if (l1norm < 9.504178996162932e-001) {
@@ -102,11 +103,11 @@ void ei_matrix_exponential(const MatrixBase<Derived> &M, typename ei_plain_matri
// Use (7,7)-Pade
const Scalar b[] = {17297280., 8648640., 1995840., 277200., 25200., 1512., 56., 1.};
PlainMatrixType M2, M4, M6;
- M2 = (M * M).lazy();
+ M2 = (Meval * Meval).lazy();
M4 = (M2 * M2).lazy();
M6 = (M4 * M2).lazy();
num = b[7]*M6 + b[5]*M4 + b[3]*M2 + b[1]*Id;
- U = (M * num).lazy();
+ U = (Meval * num).lazy();
V = b[6]*M6 + b[4]*M4 + b[2]*M2 + b[0]*Id;
} else if (l1norm < 2.097847961257068e+000) {
@@ -115,12 +116,12 @@ void ei_matrix_exponential(const MatrixBase<Derived> &M, typename ei_plain_matri
const Scalar b[] = {17643225600., 8821612800., 2075673600., 302702400., 30270240.,
2162160., 110880., 3960., 90., 1.};
PlainMatrixType M2, M4, M6, M8;
- M2 = (M * M).lazy();
+ M2 = (Meval * Meval).lazy();
M4 = (M2 * M2).lazy();
M6 = (M4 * M2).lazy();
M8 = (M6 * M2).lazy();
num = b[9]*M8 + b[7]*M6 + b[5]*M4 + b[3]*M2 + b[1]*Id;
- U = (M * num).lazy();
+ U = (Meval * num).lazy();
V = b[8]*M8 + b[6]*M6 + b[4]*M4 + b[2]*M2 + b[0]*Id;
} else {
@@ -135,7 +136,7 @@ void ei_matrix_exponential(const MatrixBase<Derived> &M, typename ei_plain_matri
squarings = std::max(0, (int)ceil(log2(l1norm / maxnorm)));
PlainMatrixType A, A2, A4, A6;
- A = M / pow(Scalar(2), squarings);
+ A = Meval / pow(Scalar(2), squarings);
A2 = (A * A).lazy();
A4 = (A2 * A2).lazy();
A6 = (A4 * A2).lazy();