// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2009 Gael Guennebaud // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_SELFADJOINTMATRIX_H #define EIGEN_SELFADJOINTMATRIX_H namespace Eigen { /** \class SelfAdjointView * \ingroup Core_Module * * * \brief Expression of a selfadjoint matrix from a triangular part of a dense matrix * * \param MatrixType the type of the dense matrix storing the coefficients * \param TriangularPart can be either \c #Lower or \c #Upper * * This class is an expression of a sefladjoint matrix from a triangular part of a matrix * with given dense storage of the coefficients. It is the return type of MatrixBase::selfadjointView() * and most of the time this is the only way that it is used. * * \sa class TriangularBase, MatrixBase::selfadjointView() */ namespace internal { template struct traits > : traits { typedef typename nested::type MatrixTypeNested; typedef typename remove_all::type MatrixTypeNestedCleaned; typedef MatrixType ExpressionType; typedef typename MatrixType::PlainObject DenseMatrixType; enum { Mode = UpLo | SelfAdjoint, Flags = MatrixTypeNestedCleaned::Flags & (HereditaryBits) & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit)), // FIXME these flags should be preserved CoeffReadCost = MatrixTypeNestedCleaned::CoeffReadCost }; }; } template struct SelfadjointProductMatrix; // FIXME could also be called SelfAdjointWrapper to be consistent with DiagonalWrapper ?? template class SelfAdjointView : public TriangularBase > { public: typedef TriangularBase Base; typedef typename internal::traits::MatrixTypeNested MatrixTypeNested; typedef typename internal::traits::MatrixTypeNestedCleaned MatrixTypeNestedCleaned; /** \brief The type of coefficients in this matrix */ typedef typename internal::traits::Scalar Scalar; typedef typename MatrixType::Index Index; enum { Mode = internal::traits::Mode }; typedef typename MatrixType::PlainObject PlainObject; inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix) {} inline Index rows() const { return m_matrix.rows(); } inline Index cols() const { return m_matrix.cols(); } inline Index outerStride() const { return m_matrix.outerStride(); } inline Index innerStride() const { return m_matrix.innerStride(); } /** \sa MatrixBase::coeff() * \warning the coordinates must fit into the referenced triangular part */ inline Scalar coeff(Index row, Index col) const { Base::check_coordinates_internal(row, col); return m_matrix.coeff(row, col); } /** \sa MatrixBase::coeffRef() * \warning the coordinates must fit into the referenced triangular part */ inline Scalar& coeffRef(Index row, Index col) { Base::check_coordinates_internal(row, col); return m_matrix.const_cast_derived().coeffRef(row, col); } /** \internal */ const MatrixTypeNestedCleaned& _expression() const { return m_matrix; } const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; } MatrixTypeNestedCleaned& nestedExpression() { return *const_cast(&m_matrix); } /** Efficient self-adjoint matrix times vector/matrix product */ template SelfadjointProductMatrix operator*(const MatrixBase& rhs) const { return SelfadjointProductMatrix (m_matrix, rhs.derived()); } /** Efficient vector/matrix times self-adjoint matrix product */ template friend SelfadjointProductMatrix operator*(const MatrixBase& lhs, const SelfAdjointView& rhs) { return SelfadjointProductMatrix (lhs.derived(),rhs.m_matrix); } /** Perform a symmetric rank 2 update of the selfadjoint matrix \c *this: * \f$ this = this + \alpha u v^* + conj(\alpha) v u^* \f$ * \returns a reference to \c *this * * The vectors \a u and \c v \b must be column vectors, however they can be * a adjoint expression without any overhead. Only the meaningful triangular * part of the matrix is updated, the rest is left unchanged. * * \sa rankUpdate(const MatrixBase&, Scalar) */ template SelfAdjointView& rankUpdate(const MatrixBase& u, const MatrixBase& v, const Scalar& alpha = Scalar(1)); /** Perform a symmetric rank K update of the selfadjoint matrix \c *this: * \f$ this = this + \alpha ( u u^* ) \f$ where \a u is a vector or matrix. * * \returns a reference to \c *this * * Note that to perform \f$ this = this + \alpha ( u^* u ) \f$ you can simply * call this function with u.adjoint(). * * \sa rankUpdate(const MatrixBase&, const MatrixBase&, Scalar) */ template SelfAdjointView& rankUpdate(const MatrixBase& u, const Scalar& alpha = Scalar(1)); /////////// Cholesky module /////////// const LLT llt() const; const LDLT ldlt() const; /////////// Eigenvalue module /////////// /** Real part of #Scalar */ typedef typename NumTraits::Real RealScalar; /** Return type of eigenvalues() */ typedef Matrix::ColsAtCompileTime, 1> EigenvaluesReturnType; EigenvaluesReturnType eigenvalues() const; RealScalar operatorNorm() const; #ifdef EIGEN2_SUPPORT template SelfAdjointView& operator=(const MatrixBase& other) { enum { OtherPart = UpLo == Upper ? StrictlyLower : StrictlyUpper }; m_matrix.const_cast_derived().template triangularView() = other; m_matrix.const_cast_derived().template triangularView() = other.adjoint(); return *this; } template SelfAdjointView& operator=(const TriangularView& other) { enum { OtherPart = UpLo == Upper ? StrictlyLower : StrictlyUpper }; m_matrix.const_cast_derived().template triangularView() = other.toDenseMatrix(); m_matrix.const_cast_derived().template triangularView() = other.toDenseMatrix().adjoint(); return *this; } #endif protected: MatrixTypeNested m_matrix; }; // template // internal::selfadjoint_matrix_product_returntype > // operator*(const MatrixBase& lhs, const SelfAdjointView& rhs) // { // return internal::matrix_selfadjoint_product_returntype >(lhs.derived(),rhs); // } // selfadjoint to dense matrix namespace internal { template struct triangular_assignment_selector { enum { col = (UnrollCount-1) / Derived1::RowsAtCompileTime, row = (UnrollCount-1) % Derived1::RowsAtCompileTime }; static inline void run(Derived1 &dst, const Derived2 &src) { triangular_assignment_selector::run(dst, src); if(row == col) dst.coeffRef(row, col) = numext::real(src.coeff(row, col)); else if(row < col) dst.coeffRef(col, row) = numext::conj(dst.coeffRef(row, col) = src.coeff(row, col)); } }; template struct triangular_assignment_selector { static inline void run(Derived1 &, const Derived2 &) {} }; template struct triangular_assignment_selector { enum { col = (UnrollCount-1) / Derived1::RowsAtCompileTime, row = (UnrollCount-1) % Derived1::RowsAtCompileTime }; static inline void run(Derived1 &dst, const Derived2 &src) { triangular_assignment_selector::run(dst, src); if(row == col) dst.coeffRef(row, col) = numext::real(src.coeff(row, col)); else if(row > col) dst.coeffRef(col, row) = numext::conj(dst.coeffRef(row, col) = src.coeff(row, col)); } }; template struct triangular_assignment_selector { static inline void run(Derived1 &, const Derived2 &) {} }; template struct triangular_assignment_selector { typedef typename Derived1::Index Index; static inline void run(Derived1 &dst, const Derived2 &src) { for(Index j = 0; j < dst.cols(); ++j) { for(Index i = 0; i < j; ++i) { dst.copyCoeff(i, j, src); dst.coeffRef(j,i) = numext::conj(dst.coeff(i,j)); } dst.copyCoeff(j, j, src); } } }; template struct triangular_assignment_selector { static inline void run(Derived1 &dst, const Derived2 &src) { typedef typename Derived1::Index Index; for(Index i = 0; i < dst.rows(); ++i) { for(Index j = 0; j < i; ++j) { dst.copyCoeff(i, j, src); dst.coeffRef(j,i) = numext::conj(dst.coeff(i,j)); } dst.copyCoeff(i, i, src); } } }; } // end namespace internal /*************************************************************************** * Implementation of MatrixBase methods ***************************************************************************/ template template typename MatrixBase::template ConstSelfAdjointViewReturnType::Type MatrixBase::selfadjointView() const { return derived(); } template template typename MatrixBase::template SelfAdjointViewReturnType::Type MatrixBase::selfadjointView() { return derived(); } } // end namespace Eigen #endif // EIGEN_SELFADJOINTMATRIX_H