bug #881: make SparseMatrixBase::isApprox(SparseMatrixBase) exploits sparse computations instead of converting the operands to dense matrices.

This commit is contained in:
Gael Guennebaud 2014-09-22 23:33:28 +02:00
parent ae514ddfe5
commit ff46ec0f24
3 changed files with 22 additions and 15 deletions

View File

@ -1,7 +1,7 @@
// This file is part of Eigen, a lightweight C++ template library // This file is part of Eigen, a lightweight C++ template library
// for linear algebra. // 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 // 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 // Public License v. 2.0. If a copy of the MPL was not distributed
@ -10,17 +10,21 @@
#ifndef EIGEN_SPARSE_FUZZY_H #ifndef EIGEN_SPARSE_FUZZY_H
#define EIGEN_SPARSE_FUZZY_H #define EIGEN_SPARSE_FUZZY_H
// template<typename Derived> namespace Eigen {
// template<typename OtherDerived>
// bool SparseMatrixBase<Derived>::isApprox( template<typename Derived>
// const OtherDerived& other, template<typename OtherDerived>
// typename NumTraits<Scalar>::Real prec bool SparseMatrixBase<Derived>::isApprox(const SparseMatrixBase<OtherDerived>& other, const RealScalar &prec) const
// ) const {
// { using std::min;
// const typename internal::nested<Derived,2>::type nested(derived()); const typename internal::nested_eval<Derived,2,PlainObject>::type actualA(derived());
// const typename internal::nested<OtherDerived,2>::type otherNested(other.derived()); typename internal::conditional<IsRowMajor==OtherDerived::IsRowMajor,
// return (nested - otherNested).cwise().abs2().sum() const typename internal::nested_eval<OtherDerived,2,PlainObject>::type,
// <= prec * prec * (std::min)(nested.cwise().abs2().sum(), otherNested.cwise().abs2().sum()); const PlainObject>::type actualB(other.derived());
// }
return (actualA - actualB).squaredNorm() <= prec * prec * (min)(actualA.squaredNorm(), actualB.squaredNorm());
}
} // end namespace Eigen
#endif // EIGEN_SPARSE_FUZZY_H #endif // EIGEN_SPARSE_FUZZY_H

View File

@ -324,8 +324,7 @@ template<typename Derived> class SparseMatrixBase : public EigenBase<Derived>
template<typename OtherDerived> template<typename OtherDerived>
bool isApprox(const SparseMatrixBase<OtherDerived>& other, bool isApprox(const SparseMatrixBase<OtherDerived>& other,
const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
{ return toDense().isApprox(other.toDense(),prec); }
template<typename OtherDerived> template<typename OtherDerived>
bool isApprox(const MatrixBase<OtherDerived>& other, bool isApprox(const MatrixBase<OtherDerived>& other,

View File

@ -304,6 +304,10 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose()); VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint()); VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
// check isApprox handles opposite storage order
typename Transpose<SparseMatrixType>::PlainObject m3(m2);
VERIFY(m2.isApprox(m3));
} }