SparseDot.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_SPARSE_DOT_H
00011 #define EIGEN_SPARSE_DOT_H
00012 
00013 namespace Eigen { 
00014 
00015 template<typename Derived>
00016 template<typename OtherDerived>
00017 typename internal::traits<Derived>::Scalar
00018 SparseMatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
00019 {
00020   EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00021   EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
00022   EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
00023   EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
00024     YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
00025 
00026   eigen_assert(size() == other.size());
00027   eigen_assert(other.size()>0 && "you are using a non initialized vector");
00028 
00029   typename Derived::InnerIterator i(derived(),0);
00030   Scalar res(0);
00031   while (i)
00032   {
00033     res += internal::conj(i.value()) * other.coeff(i.index());
00034     ++i;
00035   }
00036   return res;
00037 }
00038 
00039 template<typename Derived>
00040 template<typename OtherDerived>
00041 typename internal::traits<Derived>::Scalar
00042 SparseMatrixBase<Derived>::dot(const SparseMatrixBase<OtherDerived>& other) const
00043 {
00044   EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00045   EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
00046   EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
00047   EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
00048     YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
00049 
00050   eigen_assert(size() == other.size());
00051 
00052   typedef typename Derived::Nested  Nested;
00053   typedef typename OtherDerived::Nested  OtherNested;
00054   typedef typename internal::remove_all<Nested>::type  NestedCleaned;
00055   typedef typename internal::remove_all<OtherNested>::type  OtherNestedCleaned;
00056 
00057   const Nested nthis(derived());
00058   const OtherNested nother(other.derived());
00059 
00060   typename NestedCleaned::InnerIterator i(nthis,0);
00061   typename OtherNestedCleaned::InnerIterator j(nother,0);
00062   Scalar res(0);
00063   while (i && j)
00064   {
00065     if (i.index()==j.index())
00066     {
00067       res += internal::conj(i.value()) * j.value();
00068       ++i; ++j;
00069     }
00070     else if (i.index()<j.index())
00071       ++i;
00072     else
00073       ++j;
00074   }
00075   return res;
00076 }
00077 
00078 template<typename Derived>
00079 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00080 SparseMatrixBase<Derived>::squaredNorm() const
00081 {
00082   return internal::real((*this).cwiseAbs2().sum());
00083 }
00084 
00085 template<typename Derived>
00086 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00087 SparseMatrixBase<Derived>::norm() const
00088 {
00089   return internal::sqrt(squaredNorm());
00090 }
00091 
00092 } // end namespace Eigen
00093 
00094 #endif // EIGEN_SPARSE_DOT_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:12:05