SelfAdjointView.h 13.9 KB
Newer Older
LM's avatar
LM committed
1 2 3 4 5
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
//
Don Gagne's avatar
Don Gagne committed
6 7 8
// 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/.
LM's avatar
LM committed
9 10 11 12

#ifndef EIGEN_SELFADJOINTMATRIX_H
#define EIGEN_SELFADJOINTMATRIX_H

Don Gagne's avatar
Don Gagne committed
13 14
namespace Eigen { 

LM's avatar
LM committed
15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34
/** \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<typename MatrixType, unsigned int UpLo>
struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
{
35
  typedef typename ref_selector<MatrixType>::non_const_type MatrixTypeNested;
LM's avatar
LM committed
36 37
  typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
  typedef MatrixType ExpressionType;
38
  typedef typename MatrixType::PlainObject FullMatrixType;
LM's avatar
LM committed
39 40
  enum {
    Mode = UpLo | SelfAdjoint,
41 42 43
    FlagsLvalueBit = is_lvalue<MatrixType>::value ? LvalueBit : 0,
    Flags =  MatrixTypeNestedCleaned::Flags & (HereditaryBits|FlagsLvalueBit)
           & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit)) // FIXME these flags should be preserved
LM's avatar
LM committed
44 45 46 47 48
  };
};
}


49 50
template<typename _MatrixType, unsigned int UpLo> class SelfAdjointView
  : public TriangularBase<SelfAdjointView<_MatrixType, UpLo> >
LM's avatar
LM committed
51 52 53
{
  public:

54
    typedef _MatrixType MatrixType;
LM's avatar
LM committed
55 56 57
    typedef TriangularBase<SelfAdjointView> Base;
    typedef typename internal::traits<SelfAdjointView>::MatrixTypeNested MatrixTypeNested;
    typedef typename internal::traits<SelfAdjointView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned;
58
    typedef MatrixTypeNestedCleaned NestedExpression;
LM's avatar
LM committed
59 60 61

    /** \brief The type of coefficients in this matrix */
    typedef typename internal::traits<SelfAdjointView>::Scalar Scalar; 
62 63
    typedef typename MatrixType::StorageIndex StorageIndex;
    typedef typename internal::remove_all<typename MatrixType::ConjugateReturnType>::type MatrixConjugateReturnType;
LM's avatar
LM committed
64 65

    enum {
66 67 68
      Mode = internal::traits<SelfAdjointView>::Mode,
      Flags = internal::traits<SelfAdjointView>::Flags,
      TransposeMode = ((Mode & Upper) ? Lower : 0) | ((Mode & Lower) ? Upper : 0)
LM's avatar
LM committed
69 70 71
    };
    typedef typename MatrixType::PlainObject PlainObject;

72 73 74 75 76
    EIGEN_DEVICE_FUNC
    explicit inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix)
    {
      EIGEN_STATIC_ASSERT(UpLo==Lower || UpLo==Upper,SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY);
    }
LM's avatar
LM committed
77

78
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
79
    inline Index rows() const { return m_matrix.rows(); }
80
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
81
    inline Index cols() const { return m_matrix.cols(); }
82
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
83
    inline Index outerStride() const { return m_matrix.outerStride(); }
84
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
85 86 87 88 89
    inline Index innerStride() const { return m_matrix.innerStride(); }

    /** \sa MatrixBase::coeff()
      * \warning the coordinates must fit into the referenced triangular part
      */
90
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
91 92 93 94 95 96 97 98 99
    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
      */
100
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
101 102
    inline Scalar& coeffRef(Index row, Index col)
    {
103
      EIGEN_STATIC_ASSERT_LVALUE(SelfAdjointView);
LM's avatar
LM committed
104
      Base::check_coordinates_internal(row, col);
105
      return m_matrix.coeffRef(row, col);
LM's avatar
LM committed
106 107 108
    }

    /** \internal */
109
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
110 111
    const MatrixTypeNestedCleaned& _expression() const { return m_matrix; }

112
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
113
    const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; }
114 115
    EIGEN_DEVICE_FUNC
    MatrixTypeNestedCleaned& nestedExpression() { return m_matrix; }
LM's avatar
LM committed
116

117
    /** Efficient triangular matrix times vector/matrix product */
LM's avatar
LM committed
118
    template<typename OtherDerived>
119 120
    EIGEN_DEVICE_FUNC
    const Product<SelfAdjointView,OtherDerived>
LM's avatar
LM committed
121 122
    operator*(const MatrixBase<OtherDerived>& rhs) const
    {
123
      return Product<SelfAdjointView,OtherDerived>(*this, rhs.derived());
LM's avatar
LM committed
124 125
    }

126
    /** Efficient vector/matrix times triangular matrix product */
LM's avatar
LM committed
127
    template<typename OtherDerived> friend
128 129
    EIGEN_DEVICE_FUNC
    const Product<OtherDerived,SelfAdjointView>
LM's avatar
LM committed
130 131
    operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView& rhs)
    {
132 133 134 135 136 137 138 139
      return Product<OtherDerived,SelfAdjointView>(lhs.derived(),rhs);
    }
    
    friend EIGEN_DEVICE_FUNC
    const SelfAdjointView<const EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(Scalar,MatrixType,product),UpLo>
    operator*(const Scalar& s, const SelfAdjointView& mat)
    {
      return (s*mat.nestedExpression()).template selfadjointView<UpLo>();
LM's avatar
LM committed
140 141 142 143 144 145 146 147 148 149 150 151 152
    }

    /** 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<DerivedU>&, Scalar)
      */
    template<typename DerivedU, typename DerivedV>
153
    EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
154
    SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, const Scalar& alpha = Scalar(1));
LM's avatar
LM committed
155 156 157 158 159 160 161 162 163 164 165 166

    /** 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<DerivedU>&, const MatrixBase<DerivedV>&, Scalar)
      */
    template<typename DerivedU>
167
    EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
168
    SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1));
LM's avatar
LM committed
169

170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234
    /** \returns an expression of a triangular view extracted from the current selfadjoint view of a given triangular part
      *
      * The parameter \a TriMode can have the following values: \c #Upper, \c #StrictlyUpper, \c #UnitUpper,
      * \c #Lower, \c #StrictlyLower, \c #UnitLower.
      *
      * If \c TriMode references the same triangular part than \c *this, then this method simply return a \c TriangularView of the nested expression,
      * otherwise, the nested expression is first transposed, thus returning a \c TriangularView<Transpose<MatrixType>> object.
      *
      * \sa MatrixBase::triangularView(), class TriangularView
      */
    template<unsigned int TriMode>
    EIGEN_DEVICE_FUNC
    typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
                                   TriangularView<MatrixType,TriMode>,
                                   TriangularView<typename MatrixType::AdjointReturnType,TriMode> >::type
    triangularView() const
    {
      typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::ConstTransposeReturnType>::type tmp1(m_matrix);
      typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::AdjointReturnType>::type tmp2(tmp1);
      return typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
                                   TriangularView<MatrixType,TriMode>,
                                   TriangularView<typename MatrixType::AdjointReturnType,TriMode> >::type(tmp2);
    }

    typedef SelfAdjointView<const MatrixConjugateReturnType,UpLo> ConjugateReturnType;
    /** \sa MatrixBase::conjugate() const */
    EIGEN_DEVICE_FUNC
    inline const ConjugateReturnType conjugate() const
    { return ConjugateReturnType(m_matrix.conjugate()); }

    typedef SelfAdjointView<const typename MatrixType::AdjointReturnType,TransposeMode> AdjointReturnType;
    /** \sa MatrixBase::adjoint() const */
    EIGEN_DEVICE_FUNC
    inline const AdjointReturnType adjoint() const
    { return AdjointReturnType(m_matrix.adjoint()); }

    typedef SelfAdjointView<typename MatrixType::TransposeReturnType,TransposeMode> TransposeReturnType;
     /** \sa MatrixBase::transpose() */
    EIGEN_DEVICE_FUNC
    inline TransposeReturnType transpose()
    {
      EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
      typename MatrixType::TransposeReturnType tmp(m_matrix);
      return TransposeReturnType(tmp);
    }

    typedef SelfAdjointView<const typename MatrixType::ConstTransposeReturnType,TransposeMode> ConstTransposeReturnType;
    /** \sa MatrixBase::transpose() const */
    EIGEN_DEVICE_FUNC
    inline const ConstTransposeReturnType transpose() const
    {
      return ConstTransposeReturnType(m_matrix.transpose());
    }

    /** \returns a const expression of the main diagonal of the matrix \c *this
      *
      * This method simply returns the diagonal of the nested expression, thus by-passing the SelfAdjointView decorator.
      *
      * \sa MatrixBase::diagonal(), class Diagonal */
    EIGEN_DEVICE_FUNC
    typename MatrixType::ConstDiagonalReturnType diagonal() const
    {
      return typename MatrixType::ConstDiagonalReturnType(m_matrix);
    }

LM's avatar
LM committed
235 236 237 238 239 240 241 242 243 244 245 246
/////////// Cholesky module ///////////

    const LLT<PlainObject, UpLo> llt() const;
    const LDLT<PlainObject, UpLo> ldlt() const;

/////////// Eigenvalue module ///////////

    /** Real part of #Scalar */
    typedef typename NumTraits<Scalar>::Real RealScalar;
    /** Return type of eigenvalues() */
    typedef Matrix<RealScalar, internal::traits<MatrixType>::ColsAtCompileTime, 1> EigenvaluesReturnType;

247
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
248
    EigenvaluesReturnType eigenvalues() const;
249
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
250 251 252
    RealScalar operatorNorm() const;

  protected:
Don Gagne's avatar
Don Gagne committed
253
    MatrixTypeNested m_matrix;
LM's avatar
LM committed
254 255 256 257 258 259 260 261 262 263 264 265 266 267
};


// template<typename OtherDerived, typename MatrixType, unsigned int UpLo>
// internal::selfadjoint_matrix_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >
// operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView<MatrixType,UpLo>& rhs)
// {
//   return internal::matrix_selfadjoint_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >(lhs.derived(),rhs);
// }

// selfadjoint to dense matrix

namespace internal {

268 269 270 271 272
// TODO currently a selfadjoint expression has the form SelfAdjointView<.,.>
//      in the future selfadjoint-ness should be defined by the expression traits
//      such that Transpose<SelfAdjointView<.,.> > is valid. (currently TriangularBase::transpose() is overloaded to make it work)
template<typename MatrixType, unsigned int Mode>
struct evaluator_traits<SelfAdjointView<MatrixType,Mode> >
LM's avatar
LM committed
273
{
274 275
  typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
  typedef SelfAdjointShape Shape;
LM's avatar
LM committed
276 277
};

278 279 280
template<int UpLo, int SetOpposite, typename DstEvaluatorTypeT, typename SrcEvaluatorTypeT, typename Functor, int Version>
class triangular_dense_assignment_kernel<UpLo,SelfAdjoint,SetOpposite,DstEvaluatorTypeT,SrcEvaluatorTypeT,Functor,Version>
  : public generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version>
LM's avatar
LM committed
281
{
282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301
protected:
  typedef generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version> Base;
  typedef typename Base::DstXprType DstXprType;
  typedef typename Base::SrcXprType SrcXprType;
  using Base::m_dst;
  using Base::m_src;
  using Base::m_functor;
public:
  
  typedef typename Base::DstEvaluatorType DstEvaluatorType;
  typedef typename Base::SrcEvaluatorType SrcEvaluatorType;
  typedef typename Base::Scalar Scalar;
  typedef typename Base::AssignmentTraits AssignmentTraits;
  
  
  EIGEN_DEVICE_FUNC triangular_dense_assignment_kernel(DstEvaluatorType &dst, const SrcEvaluatorType &src, const Functor &func, DstXprType& dstExpr)
    : Base(dst, src, func, dstExpr)
  {}
  
  EIGEN_DEVICE_FUNC void assignCoeff(Index row, Index col)
LM's avatar
LM committed
302
  {
303 304 305 306
    eigen_internal_assert(row!=col);
    Scalar tmp = m_src.coeff(row,col);
    m_functor.assignCoeff(m_dst.coeffRef(row,col), tmp);
    m_functor.assignCoeff(m_dst.coeffRef(col,row), numext::conj(tmp));
LM's avatar
LM committed
307
  }
308 309
  
  EIGEN_DEVICE_FUNC void assignDiagonalCoeff(Index id)
LM's avatar
LM committed
310
  {
311
    Base::assignCoeff(id,id);
LM's avatar
LM committed
312
  }
313 314 315
  
  EIGEN_DEVICE_FUNC void assignOppositeCoeff(Index, Index)
  { eigen_internal_assert(false && "should never be called"); }
LM's avatar
LM committed
316 317 318 319 320 321 322 323
};

} // end namespace internal

/***************************************************************************
* Implementation of MatrixBase methods
***************************************************************************/

324
/** This is the const version of MatrixBase::selfadjointView() */
LM's avatar
LM committed
325 326 327 328 329
template<typename Derived>
template<unsigned int UpLo>
typename MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
MatrixBase<Derived>::selfadjointView() const
{
330
  return typename ConstSelfAdjointViewReturnType<UpLo>::Type(derived());
LM's avatar
LM committed
331 332
}

333 334 335 336 337 338 339 340 341
/** \returns an expression of a symmetric/self-adjoint view extracted from the upper or lower triangular part of the current matrix
  *
  * The parameter \a UpLo can be either \c #Upper or \c #Lower
  *
  * Example: \include MatrixBase_selfadjointView.cpp
  * Output: \verbinclude MatrixBase_selfadjointView.out
  *
  * \sa class SelfAdjointView
  */
LM's avatar
LM committed
342 343 344 345 346
template<typename Derived>
template<unsigned int UpLo>
typename MatrixBase<Derived>::template SelfAdjointViewReturnType<UpLo>::Type
MatrixBase<Derived>::selfadjointView()
{
347
  return typename SelfAdjointViewReturnType<UpLo>::Type(derived());
LM's avatar
LM committed
348 349
}

Don Gagne's avatar
Don Gagne committed
350 351
} // end namespace Eigen

LM's avatar
LM committed
352
#endif // EIGEN_SELFADJOINTMATRIX_H