1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10 #ifndef EIGEN_SPARSE_DOT_H
11 #define EIGEN_SPARSE_DOT_H
12
13 namespace Eigen {
14
15 template<typename Derived>
16 template<typename OtherDerived>
17 typename internal::traits<Derived>::Scalar
dot(const MatrixBase<OtherDerived> & other)18 SparseMatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
19 {
20 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
21 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
22 EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
23 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
24 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
25
26 eigen_assert(size() == other.size());
27 eigen_assert(other.size()>0 && "you are using a non initialized vector");
28
29 internal::evaluator<Derived> thisEval(derived());
30 typename internal::evaluator<Derived>::InnerIterator i(thisEval, 0);
31 Scalar res(0);
32 while (i)
33 {
34 res += numext::conj(i.value()) * other.coeff(i.index());
35 ++i;
36 }
37 return res;
38 }
39
40 template<typename Derived>
41 template<typename OtherDerived>
42 typename internal::traits<Derived>::Scalar
dot(const SparseMatrixBase<OtherDerived> & other)43 SparseMatrixBase<Derived>::dot(const SparseMatrixBase<OtherDerived>& other) const
44 {
45 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
46 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
47 EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
48 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
49 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
50
51 eigen_assert(size() == other.size());
52
53 internal::evaluator<Derived> thisEval(derived());
54 typename internal::evaluator<Derived>::InnerIterator i(thisEval, 0);
55
56 internal::evaluator<OtherDerived> otherEval(other.derived());
57 typename internal::evaluator<OtherDerived>::InnerIterator j(otherEval, 0);
58
59 Scalar res(0);
60 while (i && j)
61 {
62 if (i.index()==j.index())
63 {
64 res += numext::conj(i.value()) * j.value();
65 ++i; ++j;
66 }
67 else if (i.index()<j.index())
68 ++i;
69 else
70 ++j;
71 }
72 return res;
73 }
74
75 template<typename Derived>
76 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
squaredNorm()77 SparseMatrixBase<Derived>::squaredNorm() const
78 {
79 return numext::real((*this).cwiseAbs2().sum());
80 }
81
82 template<typename Derived>
83 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
norm()84 SparseMatrixBase<Derived>::norm() const
85 {
86 using std::sqrt;
87 return sqrt(squaredNorm());
88 }
89
90 template<typename Derived>
91 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
blueNorm()92 SparseMatrixBase<Derived>::blueNorm() const
93 {
94 return internal::blueNorm_impl(*this);
95 }
96 } // end namespace Eigen
97
98 #endif // EIGEN_SPARSE_DOT_H
99