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