Fix sparse_extra_3, disable counting temporaries for testing DynamicSparseMatrix.

Multiplication of column-major `DynamicSparseMatrix`es involves three
temporaries:
- two for transposing twice to sort the coefficients
(`ConservativeSparseSparseProduct.h`, L160-161)
- one for a final copy assignment (`SparseAssign.h`, L108)
The latter is avoided in an optimization for `SparseMatrix`.

Since `DynamicSparseMatrix` is deprecated in favor of `SparseMatrix`, it's not
worth the effort to optimize further, so I simply disabled counting
temporaries via a macro.

Note that due to the inclusion of `sparse_product.cpp`, the `sparse_extra`
tests actually re-run all the original `sparse_product` tests as well.

We may want to simply drop the `DynamicSparseMatrix` tests altogether, which
would eliminate the test duplication.

Related to #2048
This commit is contained in:
Antonio Sanchez 2020-11-18 13:23:13 -08:00 committed by Rasmus Munk Larsen
parent 11e4056f6b
commit a8fdcae55d
3 changed files with 25 additions and 20 deletions

View File

@ -10,7 +10,7 @@
#ifndef EIGEN_CONSERVATIVESPARSESPARSEPRODUCT_H
#define EIGEN_CONSERVATIVESPARSESPARSEPRODUCT_H
namespace Eigen {
namespace Eigen {
namespace internal {
@ -25,16 +25,16 @@ static void conservative_sparse_sparse_product_impl(const Lhs& lhs, const Rhs& r
Index rows = lhs.innerSize();
Index cols = rhs.outerSize();
eigen_assert(lhs.outerSize() == rhs.innerSize());
ei_declare_aligned_stack_constructed_variable(bool, mask, rows, 0);
ei_declare_aligned_stack_constructed_variable(ResScalar, values, rows, 0);
ei_declare_aligned_stack_constructed_variable(Index, indices, rows, 0);
std::memset(mask,0,sizeof(bool)*rows);
evaluator<Lhs> lhsEval(lhs);
evaluator<Rhs> rhsEval(rhs);
// estimate the number of non zero entries
// given a rhs column containing Y non zeros, we assume that the respective Y columns
// of the lhs differs in average of one non zeros, thus the number of non zeros for
@ -141,7 +141,7 @@ struct conservative_sparse_sparse_product_selector<Lhs,Rhs,ResultType,ColMajor,C
typedef SparseMatrix<typename ResultType::Scalar,RowMajor,typename ResultType::StorageIndex> RowMajorMatrix;
typedef SparseMatrix<typename ResultType::Scalar,ColMajor,typename ResultType::StorageIndex> ColMajorMatrixAux;
typedef typename sparse_eval<ColMajorMatrixAux,ResultType::RowsAtCompileTime,ResultType::ColsAtCompileTime,ColMajorMatrixAux::Flags>::type ColMajorMatrix;
// If the result is tall and thin (in the extreme case a column vector)
// then it is faster to sort the coefficients inplace instead of transposing twice.
// FIXME, the following heuristic is probably not very good.
@ -155,7 +155,7 @@ struct conservative_sparse_sparse_product_selector<Lhs,Rhs,ResultType,ColMajor,C
else
{
ColMajorMatrixAux resCol(lhs.rows(),rhs.cols());
// ressort to transpose to sort the entries
// resort to transpose to sort the entries
internal::conservative_sparse_sparse_product_impl<Lhs,Rhs,ColMajorMatrixAux>(lhs, rhs, resCol, false);
RowMajorMatrix resRow(resCol);
res = resRow.markAsRValue();

View File

@ -100,6 +100,7 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_IS_APPROX(m4=(m2t.transpose()*m3t.transpose()).pruned(0), refMat4=refMat2t.transpose()*refMat3t.transpose());
VERIFY_IS_APPROX(m4=(m2*m3t.transpose()).pruned(0), refMat4=refMat2*refMat3t.transpose());
#ifndef EIGEN_SPARSE_PRODUCT_IGNORE_TEMPORARY_COUNT
// make sure the right product implementation is called:
if((!SparseMatrixType::IsRowMajor) && m2.rows()<=m3.cols())
{
@ -107,6 +108,7 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_EVALUATION_COUNT(m4 = (m2*m3).pruned(0), 1);
VERIFY_EVALUATION_COUNT(m4 = (m2*m3).eval().pruned(0), 4);
}
#endif
// and that pruning is effective:
{
@ -151,7 +153,7 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_IS_APPROX(dm4.noalias()-=m2*refMat3, refMat4-=refMat2*refMat3);
VERIFY_IS_APPROX(dm4=m2*(refMat3+refMat3), refMat4=refMat2*(refMat3+refMat3));
VERIFY_IS_APPROX(dm4=m2t.transpose()*(refMat3+refMat5)*0.5, refMat4=refMat2t.transpose()*(refMat3+refMat5)*0.5);
// sparse * dense vector
VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3.col(0), refMat4.col(0)=refMat2*refMat3.col(0));
VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3t.transpose().col(0), refMat4.col(0)=refMat2*refMat3t.transpose().col(0));
@ -182,7 +184,7 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_IS_APPROX( m4=m2.middleCols(c,1)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
VERIFY_IS_APPROX(dm4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.middleCols(c,1).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
@ -211,23 +213,23 @@ template<typename SparseMatrixType> void sparse_product()
}
VERIFY_IS_APPROX(m6=m6*m6, refMat6=refMat6*refMat6);
// sparse matrix * sparse vector
ColSpVector cv0(cols), cv1;
DenseVector dcv0(cols), dcv1;
initSparse(2*density,dcv0, cv0);
RowSpVector rv0(depth), rv1;
RowDenseVector drv0(depth), drv1(rv1);
initSparse(2*density,drv0, rv0);
VERIFY_IS_APPROX(cv1=m3*cv0, dcv1=refMat3*dcv0);
VERIFY_IS_APPROX(cv1=m3*cv0, dcv1=refMat3*dcv0);
VERIFY_IS_APPROX(rv1=rv0*m3, drv1=drv0*refMat3);
VERIFY_IS_APPROX(cv1=m3t.adjoint()*cv0, dcv1=refMat3t.adjoint()*dcv0);
VERIFY_IS_APPROX(cv1=rv0*m3, dcv1=drv0*refMat3);
VERIFY_IS_APPROX(rv1=m3*cv0, drv1=refMat3*dcv0);
}
// test matrix - diagonal product
{
DenseMatrix refM2 = DenseMatrix::Zero(rows, cols);
@ -243,7 +245,7 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_IS_APPROX(m3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
VERIFY_IS_APPROX(m3=d2*m2, refM3=d2*refM2);
VERIFY_IS_APPROX(m3=d1*m2.transpose(), refM3=d1*refM2.transpose());
// also check with a SparseWrapper:
DenseVector v1 = DenseVector::Random(cols);
DenseVector v2 = DenseVector::Random(rows);
@ -252,12 +254,12 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_IS_APPROX(m3=m2.transpose()*v2.asDiagonal(), refM3=refM2.transpose()*v2.asDiagonal());
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());
VERIFY_IS_APPROX(v2=m2*v1.asDiagonal()*v1, refM2*v1.asDiagonal()*v1);
VERIFY_IS_APPROX(v3=v2.asDiagonal()*m2*v1, v2.asDiagonal()*refM2*v1);
// 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);
@ -310,20 +312,20 @@ template<typename SparseMatrixType> void sparse_product()
VERIFY_IS_APPROX(x.noalias()+=mUp.template selfadjointView<Upper>()*b, refX+=refS*b);
VERIFY_IS_APPROX(x.noalias()-=mLo.template selfadjointView<Lower>()*b, refX-=refS*b);
VERIFY_IS_APPROX(x.noalias()+=mS.template selfadjointView<Upper|Lower>()*b, refX+=refS*b);
// sparse selfadjointView with sparse matrices
SparseMatrixType mSres(rows,rows);
VERIFY_IS_APPROX(mSres = mLo.template selfadjointView<Lower>()*mS,
refX = refLo.template selfadjointView<Lower>()*refS);
VERIFY_IS_APPROX(mSres = mS * mLo.template selfadjointView<Lower>(),
refX = refS * refLo.template selfadjointView<Lower>());
// sparse triangularView with dense matrices
VERIFY_IS_APPROX(x=mA.template triangularView<Upper>()*b, refX=refA.template triangularView<Upper>()*b);
VERIFY_IS_APPROX(x=mA.template triangularView<Lower>()*b, refX=refA.template triangularView<Lower>()*b);
VERIFY_IS_APPROX(x=b*mA.template triangularView<Upper>(), refX=b*refA.template triangularView<Upper>());
VERIFY_IS_APPROX(x=b*mA.template triangularView<Lower>(), refX=b*refA.template triangularView<Lower>());
// sparse triangularView with sparse matrices
VERIFY_IS_APPROX(mSres = mA.template triangularView<Lower>()*mS, refX = refA.template triangularView<Lower>()*refS);
VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Lower>(), refX = refS * refA.template triangularView<Lower>());
@ -368,9 +370,9 @@ void bug_942()
Vector d(1);
d[0] = 2;
double res = 2;
VERIFY_IS_APPROX( ( cmA*d.asDiagonal() ).eval().coeff(0,0), res );
VERIFY_IS_APPROX( ( d.asDiagonal()*rmA ).eval().coeff(0,0), res );
VERIFY_IS_APPROX( ( rmA*d.asDiagonal() ).eval().coeff(0,0), res );

View File

@ -22,6 +22,9 @@ static long g_dense_op_sparse_count = 0;
#endif
#define EIGEN_NO_DEPRECATED_WARNING
// Disable counting of temporaries, since sparse_product(DynamicSparseMatrix)
// has an extra copy-assignment.
#define EIGEN_SPARSE_PRODUCT_IGNORE_TEMPORARY_COUNT
#include "sparse_product.cpp"
#if 0 // sparse_basic(DynamicSparseMatrix) does not compile at all -> disabled