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