aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseCore/SparseFuzzy.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-09-22 23:33:28 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-09-22 23:33:28 +0200
commitff46ec0f240ef84e2293b33b265c703e9b765c2e (patch)
treed7be930881e5c9a4dd9d2e0e4109f02235b1e2e4 /Eigen/src/SparseCore/SparseFuzzy.h
parentae514ddfe561ef220bc9bbf8c0b7b5005b60d9c8 (diff)
bug #881: make SparseMatrixBase::isApprox(SparseMatrixBase) exploits sparse computations instead of converting the operands to dense matrices.
Diffstat (limited to 'Eigen/src/SparseCore/SparseFuzzy.h')
-rw-r--r--Eigen/src/SparseCore/SparseFuzzy.h30
1 files changed, 17 insertions, 13 deletions
diff --git a/Eigen/src/SparseCore/SparseFuzzy.h b/Eigen/src/SparseCore/SparseFuzzy.h
index 45f36e9eb..3e67cbf5f 100644
--- a/Eigen/src/SparseCore/SparseFuzzy.h
+++ b/Eigen/src/SparseCore/SparseFuzzy.h
@@ -1,7 +1,7 @@
// This file is part of Eigen, a lightweight C++ template library
// 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
// Public License v. 2.0. If a copy of the MPL was not distributed
@@ -10,17 +10,21 @@
#ifndef EIGEN_SPARSE_FUZZY_H
#define EIGEN_SPARSE_FUZZY_H
-// template<typename Derived>
-// template<typename OtherDerived>
-// bool SparseMatrixBase<Derived>::isApprox(
-// const OtherDerived& other,
-// typename NumTraits<Scalar>::Real prec
-// ) const
-// {
-// const typename internal::nested<Derived,2>::type nested(derived());
-// const typename internal::nested<OtherDerived,2>::type otherNested(other.derived());
-// return (nested - otherNested).cwise().abs2().sum()
-// <= prec * prec * (std::min)(nested.cwise().abs2().sum(), otherNested.cwise().abs2().sum());
-// }
+namespace Eigen {
+
+template<typename Derived>
+template<typename OtherDerived>
+bool SparseMatrixBase<Derived>::isApprox(const SparseMatrixBase<OtherDerived>& other, const RealScalar &prec) const
+{
+ using std::min;
+ const typename internal::nested_eval<Derived,2,PlainObject>::type actualA(derived());
+ typename internal::conditional<IsRowMajor==OtherDerived::IsRowMajor,
+ const typename internal::nested_eval<OtherDerived,2,PlainObject>::type,
+ 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