// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008 Gael Guennebaud // // 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 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_SPARSE_DOT_H #define EIGEN_SPARSE_DOT_H namespace Eigen { template template typename internal::traits::Scalar SparseMatrixBase::dot(const MatrixBase& other) const { EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived) EIGEN_STATIC_ASSERT((internal::is_same::value), YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) eigen_assert(size() == other.size()); eigen_assert(other.size()>0 && "you are using a non initialized vector"); #ifndef EIGEN_TEST_EVALUATORS typename Derived::InnerIterator i(derived(),0); #else typename internal::evaluator::type thisEval(derived()); typename internal::evaluator::InnerIterator i(thisEval, 0); #endif Scalar res(0); while (i) { res += numext::conj(i.value()) * other.coeff(i.index()); ++i; } return res; } template template typename internal::traits::Scalar SparseMatrixBase::dot(const SparseMatrixBase& other) const { EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived) EIGEN_STATIC_ASSERT((internal::is_same::value), YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) eigen_assert(size() == other.size()); #ifndef EIGEN_TEST_EVALUATORS typedef typename Derived::Nested Nested; typedef typename OtherDerived::Nested OtherNested; typedef typename internal::remove_all::type NestedCleaned; typedef typename internal::remove_all::type OtherNestedCleaned; Nested nthis(derived()); OtherNested nother(other.derived()); typename NestedCleaned::InnerIterator i(nthis,0); typename OtherNestedCleaned::InnerIterator j(nother,0); #else typename internal::evaluator::type thisEval(derived()); typename internal::evaluator::InnerIterator i(thisEval, 0); typename internal::evaluator::type otherEval(other.derived()); typename internal::evaluator::InnerIterator j(otherEval, 0); #endif Scalar res(0); while (i && j) { if (i.index()==j.index()) { res += numext::conj(i.value()) * j.value(); ++i; ++j; } else if (i.index() inline typename NumTraits::Scalar>::Real SparseMatrixBase::squaredNorm() const { return numext::real((*this).cwiseAbs2().sum()); } template inline typename NumTraits::Scalar>::Real SparseMatrixBase::norm() const { using std::sqrt; return sqrt(squaredNorm()); } template inline typename NumTraits::Scalar>::Real SparseMatrixBase::blueNorm() const { return internal::blueNorm_impl(*this); } } // end namespace Eigen #endif // EIGEN_SPARSE_DOT_H