SelfAdjointView.h 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 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. #ifndef EIGEN_SELFADJOINTMATRIX_H
  10. #define EIGEN_SELFADJOINTMATRIX_H
  11. namespace Eigen {
  12. /** \class SelfAdjointView
  13. * \ingroup Core_Module
  14. *
  15. *
  16. * \brief Expression of a selfadjoint matrix from a triangular part of a dense matrix
  17. *
  18. * \param MatrixType the type of the dense matrix storing the coefficients
  19. * \param TriangularPart can be either \c #Lower or \c #Upper
  20. *
  21. * This class is an expression of a sefladjoint matrix from a triangular part of a matrix
  22. * with given dense storage of the coefficients. It is the return type of MatrixBase::selfadjointView()
  23. * and most of the time this is the only way that it is used.
  24. *
  25. * \sa class TriangularBase, MatrixBase::selfadjointView()
  26. */
  27. namespace internal {
  28. template<typename MatrixType, unsigned int UpLo>
  29. struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
  30. {
  31. typedef typename ref_selector<MatrixType>::non_const_type MatrixTypeNested;
  32. typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
  33. typedef MatrixType ExpressionType;
  34. typedef typename MatrixType::PlainObject FullMatrixType;
  35. enum {
  36. Mode = UpLo | SelfAdjoint,
  37. FlagsLvalueBit = is_lvalue<MatrixType>::value ? LvalueBit : 0,
  38. Flags = MatrixTypeNestedCleaned::Flags & (HereditaryBits|FlagsLvalueBit)
  39. & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit)) // FIXME these flags should be preserved
  40. };
  41. };
  42. }
  43. template<typename _MatrixType, unsigned int UpLo> class SelfAdjointView
  44. : public TriangularBase<SelfAdjointView<_MatrixType, UpLo> >
  45. {
  46. public:
  47. typedef _MatrixType MatrixType;
  48. typedef TriangularBase<SelfAdjointView> Base;
  49. typedef typename internal::traits<SelfAdjointView>::MatrixTypeNested MatrixTypeNested;
  50. typedef typename internal::traits<SelfAdjointView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned;
  51. typedef MatrixTypeNestedCleaned NestedExpression;
  52. /** \brief The type of coefficients in this matrix */
  53. typedef typename internal::traits<SelfAdjointView>::Scalar Scalar;
  54. typedef typename MatrixType::StorageIndex StorageIndex;
  55. typedef typename internal::remove_all<typename MatrixType::ConjugateReturnType>::type MatrixConjugateReturnType;
  56. enum {
  57. Mode = internal::traits<SelfAdjointView>::Mode,
  58. Flags = internal::traits<SelfAdjointView>::Flags,
  59. TransposeMode = ((Mode & Upper) ? Lower : 0) | ((Mode & Lower) ? Upper : 0)
  60. };
  61. typedef typename MatrixType::PlainObject PlainObject;
  62. EIGEN_DEVICE_FUNC
  63. explicit inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix)
  64. {
  65. EIGEN_STATIC_ASSERT(UpLo==Lower || UpLo==Upper,SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY);
  66. }
  67. EIGEN_DEVICE_FUNC
  68. inline Index rows() const { return m_matrix.rows(); }
  69. EIGEN_DEVICE_FUNC
  70. inline Index cols() const { return m_matrix.cols(); }
  71. EIGEN_DEVICE_FUNC
  72. inline Index outerStride() const { return m_matrix.outerStride(); }
  73. EIGEN_DEVICE_FUNC
  74. inline Index innerStride() const { return m_matrix.innerStride(); }
  75. /** \sa MatrixBase::coeff()
  76. * \warning the coordinates must fit into the referenced triangular part
  77. */
  78. EIGEN_DEVICE_FUNC
  79. inline Scalar coeff(Index row, Index col) const
  80. {
  81. Base::check_coordinates_internal(row, col);
  82. return m_matrix.coeff(row, col);
  83. }
  84. /** \sa MatrixBase::coeffRef()
  85. * \warning the coordinates must fit into the referenced triangular part
  86. */
  87. EIGEN_DEVICE_FUNC
  88. inline Scalar& coeffRef(Index row, Index col)
  89. {
  90. EIGEN_STATIC_ASSERT_LVALUE(SelfAdjointView);
  91. Base::check_coordinates_internal(row, col);
  92. return m_matrix.coeffRef(row, col);
  93. }
  94. /** \internal */
  95. EIGEN_DEVICE_FUNC
  96. const MatrixTypeNestedCleaned& _expression() const { return m_matrix; }
  97. EIGEN_DEVICE_FUNC
  98. const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; }
  99. EIGEN_DEVICE_FUNC
  100. MatrixTypeNestedCleaned& nestedExpression() { return m_matrix; }
  101. /** Efficient triangular matrix times vector/matrix product */
  102. template<typename OtherDerived>
  103. EIGEN_DEVICE_FUNC
  104. const Product<SelfAdjointView,OtherDerived>
  105. operator*(const MatrixBase<OtherDerived>& rhs) const
  106. {
  107. return Product<SelfAdjointView,OtherDerived>(*this, rhs.derived());
  108. }
  109. /** Efficient vector/matrix times triangular matrix product */
  110. template<typename OtherDerived> friend
  111. EIGEN_DEVICE_FUNC
  112. const Product<OtherDerived,SelfAdjointView>
  113. operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView& rhs)
  114. {
  115. return Product<OtherDerived,SelfAdjointView>(lhs.derived(),rhs);
  116. }
  117. friend EIGEN_DEVICE_FUNC
  118. const SelfAdjointView<const EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(Scalar,MatrixType,product),UpLo>
  119. operator*(const Scalar& s, const SelfAdjointView& mat)
  120. {
  121. return (s*mat.nestedExpression()).template selfadjointView<UpLo>();
  122. }
  123. /** Perform a symmetric rank 2 update of the selfadjoint matrix \c *this:
  124. * \f$ this = this + \alpha u v^* + conj(\alpha) v u^* \f$
  125. * \returns a reference to \c *this
  126. *
  127. * The vectors \a u and \c v \b must be column vectors, however they can be
  128. * a adjoint expression without any overhead. Only the meaningful triangular
  129. * part of the matrix is updated, the rest is left unchanged.
  130. *
  131. * \sa rankUpdate(const MatrixBase<DerivedU>&, Scalar)
  132. */
  133. template<typename DerivedU, typename DerivedV>
  134. EIGEN_DEVICE_FUNC
  135. SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, const Scalar& alpha = Scalar(1));
  136. /** Perform a symmetric rank K update of the selfadjoint matrix \c *this:
  137. * \f$ this = this + \alpha ( u u^* ) \f$ where \a u is a vector or matrix.
  138. *
  139. * \returns a reference to \c *this
  140. *
  141. * Note that to perform \f$ this = this + \alpha ( u^* u ) \f$ you can simply
  142. * call this function with u.adjoint().
  143. *
  144. * \sa rankUpdate(const MatrixBase<DerivedU>&, const MatrixBase<DerivedV>&, Scalar)
  145. */
  146. template<typename DerivedU>
  147. EIGEN_DEVICE_FUNC
  148. SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1));
  149. /** \returns an expression of a triangular view extracted from the current selfadjoint view of a given triangular part
  150. *
  151. * The parameter \a TriMode can have the following values: \c #Upper, \c #StrictlyUpper, \c #UnitUpper,
  152. * \c #Lower, \c #StrictlyLower, \c #UnitLower.
  153. *
  154. * If \c TriMode references the same triangular part than \c *this, then this method simply return a \c TriangularView of the nested expression,
  155. * otherwise, the nested expression is first transposed, thus returning a \c TriangularView<Transpose<MatrixType>> object.
  156. *
  157. * \sa MatrixBase::triangularView(), class TriangularView
  158. */
  159. template<unsigned int TriMode>
  160. EIGEN_DEVICE_FUNC
  161. typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
  162. TriangularView<MatrixType,TriMode>,
  163. TriangularView<typename MatrixType::AdjointReturnType,TriMode> >::type
  164. triangularView() const
  165. {
  166. typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::ConstTransposeReturnType>::type tmp1(m_matrix);
  167. typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::AdjointReturnType>::type tmp2(tmp1);
  168. return typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
  169. TriangularView<MatrixType,TriMode>,
  170. TriangularView<typename MatrixType::AdjointReturnType,TriMode> >::type(tmp2);
  171. }
  172. typedef SelfAdjointView<const MatrixConjugateReturnType,UpLo> ConjugateReturnType;
  173. /** \sa MatrixBase::conjugate() const */
  174. EIGEN_DEVICE_FUNC
  175. inline const ConjugateReturnType conjugate() const
  176. { return ConjugateReturnType(m_matrix.conjugate()); }
  177. typedef SelfAdjointView<const typename MatrixType::AdjointReturnType,TransposeMode> AdjointReturnType;
  178. /** \sa MatrixBase::adjoint() const */
  179. EIGEN_DEVICE_FUNC
  180. inline const AdjointReturnType adjoint() const
  181. { return AdjointReturnType(m_matrix.adjoint()); }
  182. typedef SelfAdjointView<typename MatrixType::TransposeReturnType,TransposeMode> TransposeReturnType;
  183. /** \sa MatrixBase::transpose() */
  184. EIGEN_DEVICE_FUNC
  185. inline TransposeReturnType transpose()
  186. {
  187. EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
  188. typename MatrixType::TransposeReturnType tmp(m_matrix);
  189. return TransposeReturnType(tmp);
  190. }
  191. typedef SelfAdjointView<const typename MatrixType::ConstTransposeReturnType,TransposeMode> ConstTransposeReturnType;
  192. /** \sa MatrixBase::transpose() const */
  193. EIGEN_DEVICE_FUNC
  194. inline const ConstTransposeReturnType transpose() const
  195. {
  196. return ConstTransposeReturnType(m_matrix.transpose());
  197. }
  198. /** \returns a const expression of the main diagonal of the matrix \c *this
  199. *
  200. * This method simply returns the diagonal of the nested expression, thus by-passing the SelfAdjointView decorator.
  201. *
  202. * \sa MatrixBase::diagonal(), class Diagonal */
  203. EIGEN_DEVICE_FUNC
  204. typename MatrixType::ConstDiagonalReturnType diagonal() const
  205. {
  206. return typename MatrixType::ConstDiagonalReturnType(m_matrix);
  207. }
  208. /////////// Cholesky module ///////////
  209. const LLT<PlainObject, UpLo> llt() const;
  210. const LDLT<PlainObject, UpLo> ldlt() const;
  211. /////////// Eigenvalue module ///////////
  212. /** Real part of #Scalar */
  213. typedef typename NumTraits<Scalar>::Real RealScalar;
  214. /** Return type of eigenvalues() */
  215. typedef Matrix<RealScalar, internal::traits<MatrixType>::ColsAtCompileTime, 1> EigenvaluesReturnType;
  216. EIGEN_DEVICE_FUNC
  217. EigenvaluesReturnType eigenvalues() const;
  218. EIGEN_DEVICE_FUNC
  219. RealScalar operatorNorm() const;
  220. protected:
  221. MatrixTypeNested m_matrix;
  222. };
  223. // template<typename OtherDerived, typename MatrixType, unsigned int UpLo>
  224. // internal::selfadjoint_matrix_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >
  225. // operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView<MatrixType,UpLo>& rhs)
  226. // {
  227. // return internal::matrix_selfadjoint_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >(lhs.derived(),rhs);
  228. // }
  229. // selfadjoint to dense matrix
  230. namespace internal {
  231. // TODO currently a selfadjoint expression has the form SelfAdjointView<.,.>
  232. // in the future selfadjoint-ness should be defined by the expression traits
  233. // such that Transpose<SelfAdjointView<.,.> > is valid. (currently TriangularBase::transpose() is overloaded to make it work)
  234. template<typename MatrixType, unsigned int Mode>
  235. struct evaluator_traits<SelfAdjointView<MatrixType,Mode> >
  236. {
  237. typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
  238. typedef SelfAdjointShape Shape;
  239. };
  240. template<int UpLo, int SetOpposite, typename DstEvaluatorTypeT, typename SrcEvaluatorTypeT, typename Functor, int Version>
  241. class triangular_dense_assignment_kernel<UpLo,SelfAdjoint,SetOpposite,DstEvaluatorTypeT,SrcEvaluatorTypeT,Functor,Version>
  242. : public generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version>
  243. {
  244. protected:
  245. typedef generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version> Base;
  246. typedef typename Base::DstXprType DstXprType;
  247. typedef typename Base::SrcXprType SrcXprType;
  248. using Base::m_dst;
  249. using Base::m_src;
  250. using Base::m_functor;
  251. public:
  252. typedef typename Base::DstEvaluatorType DstEvaluatorType;
  253. typedef typename Base::SrcEvaluatorType SrcEvaluatorType;
  254. typedef typename Base::Scalar Scalar;
  255. typedef typename Base::AssignmentTraits AssignmentTraits;
  256. EIGEN_DEVICE_FUNC triangular_dense_assignment_kernel(DstEvaluatorType &dst, const SrcEvaluatorType &src, const Functor &func, DstXprType& dstExpr)
  257. : Base(dst, src, func, dstExpr)
  258. {}
  259. EIGEN_DEVICE_FUNC void assignCoeff(Index row, Index col)
  260. {
  261. eigen_internal_assert(row!=col);
  262. Scalar tmp = m_src.coeff(row,col);
  263. m_functor.assignCoeff(m_dst.coeffRef(row,col), tmp);
  264. m_functor.assignCoeff(m_dst.coeffRef(col,row), numext::conj(tmp));
  265. }
  266. EIGEN_DEVICE_FUNC void assignDiagonalCoeff(Index id)
  267. {
  268. Base::assignCoeff(id,id);
  269. }
  270. EIGEN_DEVICE_FUNC void assignOppositeCoeff(Index, Index)
  271. { eigen_internal_assert(false && "should never be called"); }
  272. };
  273. } // end namespace internal
  274. /***************************************************************************
  275. * Implementation of MatrixBase methods
  276. ***************************************************************************/
  277. /** This is the const version of MatrixBase::selfadjointView() */
  278. template<typename Derived>
  279. template<unsigned int UpLo>
  280. typename MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
  281. MatrixBase<Derived>::selfadjointView() const
  282. {
  283. return typename ConstSelfAdjointViewReturnType<UpLo>::Type(derived());
  284. }
  285. /** \returns an expression of a symmetric/self-adjoint view extracted from the upper or lower triangular part of the current matrix
  286. *
  287. * The parameter \a UpLo can be either \c #Upper or \c #Lower
  288. *
  289. * Example: \include MatrixBase_selfadjointView.cpp
  290. * Output: \verbinclude MatrixBase_selfadjointView.out
  291. *
  292. * \sa class SelfAdjointView
  293. */
  294. template<typename Derived>
  295. template<unsigned int UpLo>
  296. typename MatrixBase<Derived>::template SelfAdjointViewReturnType<UpLo>::Type
  297. MatrixBase<Derived>::selfadjointView()
  298. {
  299. return typename SelfAdjointViewReturnType<UpLo>::Type(derived());
  300. }
  301. } // end namespace Eigen
  302. #endif // EIGEN_SELFADJOINTMATRIX_H