mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-09 07:00:27 +08:00
Fix bug #611: diag * sparse * diag
This commit is contained in:
parent
9b9177f1ce
commit
4f14b3fa72
@ -78,7 +78,11 @@ class SparseDiagonalProduct
|
||||
EIGEN_SPARSE_PUBLIC_INTERFACE(SparseDiagonalProduct)
|
||||
|
||||
typedef internal::sparse_diagonal_product_inner_iterator_selector
|
||||
<_LhsNested,_RhsNested,SparseDiagonalProduct,LhsMode,RhsMode> InnerIterator;
|
||||
<_LhsNested,_RhsNested,SparseDiagonalProduct,LhsMode,RhsMode> InnerIterator;
|
||||
|
||||
// We do not want ReverseInnerIterator for diagonal-sparse products,
|
||||
// but this dummy declaration is needed to make diag * sparse * diag compile.
|
||||
class ReverseInnerIterator;
|
||||
|
||||
EIGEN_STRONG_INLINE SparseDiagonalProduct(const Lhs& lhs, const Rhs& rhs)
|
||||
: m_lhs(lhs), m_rhs(rhs)
|
||||
|
@ -161,6 +161,8 @@ template<typename SparseMatrixType> void sparse_product()
|
||||
VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2, refM3=v2.asDiagonal()*refM2);
|
||||
VERIFY_IS_APPROX(m3=v1.asDiagonal()*m2.transpose(), refM3=v1.asDiagonal()*refM2.transpose());
|
||||
|
||||
VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2*v1.asDiagonal(), refM3=v2.asDiagonal()*refM2*v1.asDiagonal());
|
||||
|
||||
// evaluate to a dense matrix to check the .row() and .col() iterator functions
|
||||
VERIFY_IS_APPROX(d3=m2*d1, refM3=refM2*d1);
|
||||
VERIFY_IS_APPROX(d3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
|
||||
|
Loading…
Reference in New Issue
Block a user