Home | History | Annotate | Download | only in Core
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2007-2009 Benoit Jacob <jacob.benoit.1 (at) gmail.com>
      5 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      6 //
      7 // This Source Code Form is subject to the terms of the Mozilla
      8 // Public License v. 2.0. If a copy of the MPL was not distributed
      9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
     10 
     11 #ifndef EIGEN_DIAGONAL_H
     12 #define EIGEN_DIAGONAL_H
     13 
     14 namespace Eigen {
     15 
     16 /** \class Diagonal
     17   * \ingroup Core_Module
     18   *
     19   * \brief Expression of a diagonal/subdiagonal/superdiagonal in a matrix
     20   *
     21   * \param MatrixType the type of the object in which we are taking a sub/main/super diagonal
     22   * \param DiagIndex the index of the sub/super diagonal. The default is 0 and it means the main diagonal.
     23   *              A positive value means a superdiagonal, a negative value means a subdiagonal.
     24   *              You can also use Dynamic so the index can be set at runtime.
     25   *
     26   * The matrix is not required to be square.
     27   *
     28   * This class represents an expression of the main diagonal, or any sub/super diagonal
     29   * of a square matrix. It is the return type of MatrixBase::diagonal() and MatrixBase::diagonal(Index) and most of the
     30   * time this is the only way it is used.
     31   *
     32   * \sa MatrixBase::diagonal(), MatrixBase::diagonal(Index)
     33   */
     34 
     35 namespace internal {
     36 template<typename MatrixType, int DiagIndex>
     37 struct traits<Diagonal<MatrixType,DiagIndex> >
     38  : traits<MatrixType>
     39 {
     40   typedef typename nested<MatrixType>::type MatrixTypeNested;
     41   typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
     42   typedef typename MatrixType::StorageKind StorageKind;
     43   enum {
     44     RowsAtCompileTime = (int(DiagIndex) == DynamicIndex || int(MatrixType::SizeAtCompileTime) == Dynamic) ? Dynamic
     45                       : (EIGEN_PLAIN_ENUM_MIN(MatrixType::RowsAtCompileTime - EIGEN_PLAIN_ENUM_MAX(-DiagIndex, 0),
     46                                               MatrixType::ColsAtCompileTime - EIGEN_PLAIN_ENUM_MAX( DiagIndex, 0))),
     47     ColsAtCompileTime = 1,
     48     MaxRowsAtCompileTime = int(MatrixType::MaxSizeAtCompileTime) == Dynamic ? Dynamic
     49                          : DiagIndex == DynamicIndex ? EIGEN_SIZE_MIN_PREFER_FIXED(MatrixType::MaxRowsAtCompileTime,
     50                                                                               MatrixType::MaxColsAtCompileTime)
     51                          : (EIGEN_PLAIN_ENUM_MIN(MatrixType::MaxRowsAtCompileTime - EIGEN_PLAIN_ENUM_MAX(-DiagIndex, 0),
     52                                                  MatrixType::MaxColsAtCompileTime - EIGEN_PLAIN_ENUM_MAX( DiagIndex, 0))),
     53     MaxColsAtCompileTime = 1,
     54     MaskLvalueBit = is_lvalue<MatrixType>::value ? LvalueBit : 0,
     55     Flags = (unsigned int)_MatrixTypeNested::Flags & (HereditaryBits | LinearAccessBit | MaskLvalueBit | DirectAccessBit) & ~RowMajorBit,
     56     CoeffReadCost = _MatrixTypeNested::CoeffReadCost,
     57     MatrixTypeOuterStride = outer_stride_at_compile_time<MatrixType>::ret,
     58     InnerStrideAtCompileTime = MatrixTypeOuterStride == Dynamic ? Dynamic : MatrixTypeOuterStride+1,
     59     OuterStrideAtCompileTime = 0
     60   };
     61 };
     62 }
     63 
     64 template<typename MatrixType, int _DiagIndex> class Diagonal
     65    : public internal::dense_xpr_base< Diagonal<MatrixType,_DiagIndex> >::type
     66 {
     67   public:
     68 
     69     enum { DiagIndex = _DiagIndex };
     70     typedef typename internal::dense_xpr_base<Diagonal>::type Base;
     71     EIGEN_DENSE_PUBLIC_INTERFACE(Diagonal)
     72 
     73     inline Diagonal(MatrixType& matrix, Index a_index = DiagIndex) : m_matrix(matrix), m_index(a_index) {}
     74 
     75     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Diagonal)
     76 
     77     inline Index rows() const
     78     { return m_index.value()<0 ? (std::min<Index>)(m_matrix.cols(),m_matrix.rows()+m_index.value()) : (std::min<Index>)(m_matrix.rows(),m_matrix.cols()-m_index.value()); }
     79 
     80     inline Index cols() const { return 1; }
     81 
     82     inline Index innerStride() const
     83     {
     84       return m_matrix.outerStride() + 1;
     85     }
     86 
     87     inline Index outerStride() const
     88     {
     89       return 0;
     90     }
     91 
     92     typedef typename internal::conditional<
     93                        internal::is_lvalue<MatrixType>::value,
     94                        Scalar,
     95                        const Scalar
     96                      >::type ScalarWithConstIfNotLvalue;
     97 
     98     inline ScalarWithConstIfNotLvalue* data() { return &(m_matrix.const_cast_derived().coeffRef(rowOffset(), colOffset())); }
     99     inline const Scalar* data() const { return &(m_matrix.const_cast_derived().coeffRef(rowOffset(), colOffset())); }
    100 
    101     inline Scalar& coeffRef(Index row, Index)
    102     {
    103       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
    104       return m_matrix.const_cast_derived().coeffRef(row+rowOffset(), row+colOffset());
    105     }
    106 
    107     inline const Scalar& coeffRef(Index row, Index) const
    108     {
    109       return m_matrix.const_cast_derived().coeffRef(row+rowOffset(), row+colOffset());
    110     }
    111 
    112     inline CoeffReturnType coeff(Index row, Index) const
    113     {
    114       return m_matrix.coeff(row+rowOffset(), row+colOffset());
    115     }
    116 
    117     inline Scalar& coeffRef(Index idx)
    118     {
    119       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
    120       return m_matrix.const_cast_derived().coeffRef(idx+rowOffset(), idx+colOffset());
    121     }
    122 
    123     inline const Scalar& coeffRef(Index idx) const
    124     {
    125       return m_matrix.const_cast_derived().coeffRef(idx+rowOffset(), idx+colOffset());
    126     }
    127 
    128     inline CoeffReturnType coeff(Index idx) const
    129     {
    130       return m_matrix.coeff(idx+rowOffset(), idx+colOffset());
    131     }
    132 
    133     const typename internal::remove_all<typename MatrixType::Nested>::type&
    134     nestedExpression() const
    135     {
    136       return m_matrix;
    137     }
    138 
    139     int index() const
    140     {
    141       return m_index.value();
    142     }
    143 
    144   protected:
    145     typename MatrixType::Nested m_matrix;
    146     const internal::variable_if_dynamicindex<Index, DiagIndex> m_index;
    147 
    148   private:
    149     // some compilers may fail to optimize std::max etc in case of compile-time constants...
    150     EIGEN_STRONG_INLINE Index absDiagIndex() const { return m_index.value()>0 ? m_index.value() : -m_index.value(); }
    151     EIGEN_STRONG_INLINE Index rowOffset() const { return m_index.value()>0 ? 0 : -m_index.value(); }
    152     EIGEN_STRONG_INLINE Index colOffset() const { return m_index.value()>0 ? m_index.value() : 0; }
    153     // triger a compile time error is someone try to call packet
    154     template<int LoadMode> typename MatrixType::PacketReturnType packet(Index) const;
    155     template<int LoadMode> typename MatrixType::PacketReturnType packet(Index,Index) const;
    156 };
    157 
    158 /** \returns an expression of the main diagonal of the matrix \c *this
    159   *
    160   * \c *this is not required to be square.
    161   *
    162   * Example: \include MatrixBase_diagonal.cpp
    163   * Output: \verbinclude MatrixBase_diagonal.out
    164   *
    165   * \sa class Diagonal */
    166 template<typename Derived>
    167 inline typename MatrixBase<Derived>::DiagonalReturnType
    168 MatrixBase<Derived>::diagonal()
    169 {
    170   return derived();
    171 }
    172 
    173 /** This is the const version of diagonal(). */
    174 template<typename Derived>
    175 inline typename MatrixBase<Derived>::ConstDiagonalReturnType
    176 MatrixBase<Derived>::diagonal() const
    177 {
    178   return ConstDiagonalReturnType(derived());
    179 }
    180 
    181 /** \returns an expression of the \a DiagIndex-th sub or super diagonal of the matrix \c *this
    182   *
    183   * \c *this is not required to be square.
    184   *
    185   * The template parameter \a DiagIndex represent a super diagonal if \a DiagIndex > 0
    186   * and a sub diagonal otherwise. \a DiagIndex == 0 is equivalent to the main diagonal.
    187   *
    188   * Example: \include MatrixBase_diagonal_int.cpp
    189   * Output: \verbinclude MatrixBase_diagonal_int.out
    190   *
    191   * \sa MatrixBase::diagonal(), class Diagonal */
    192 template<typename Derived>
    193 inline typename MatrixBase<Derived>::DiagonalDynamicIndexReturnType
    194 MatrixBase<Derived>::diagonal(Index index)
    195 {
    196   return DiagonalDynamicIndexReturnType(derived(), index);
    197 }
    198 
    199 /** This is the const version of diagonal(Index). */
    200 template<typename Derived>
    201 inline typename MatrixBase<Derived>::ConstDiagonalDynamicIndexReturnType
    202 MatrixBase<Derived>::diagonal(Index index) const
    203 {
    204   return ConstDiagonalDynamicIndexReturnType(derived(), index);
    205 }
    206 
    207 /** \returns an expression of the \a DiagIndex-th sub or super diagonal of the matrix \c *this
    208   *
    209   * \c *this is not required to be square.
    210   *
    211   * The template parameter \a DiagIndex represent a super diagonal if \a DiagIndex > 0
    212   * and a sub diagonal otherwise. \a DiagIndex == 0 is equivalent to the main diagonal.
    213   *
    214   * Example: \include MatrixBase_diagonal_template_int.cpp
    215   * Output: \verbinclude MatrixBase_diagonal_template_int.out
    216   *
    217   * \sa MatrixBase::diagonal(), class Diagonal */
    218 template<typename Derived>
    219 template<int Index>
    220 inline typename MatrixBase<Derived>::template DiagonalIndexReturnType<Index>::Type
    221 MatrixBase<Derived>::diagonal()
    222 {
    223   return derived();
    224 }
    225 
    226 /** This is the const version of diagonal<int>(). */
    227 template<typename Derived>
    228 template<int Index>
    229 inline typename MatrixBase<Derived>::template ConstDiagonalIndexReturnType<Index>::Type
    230 MatrixBase<Derived>::diagonal() const
    231 {
    232   return derived();
    233 }
    234 
    235 } // end namespace Eigen
    236 
    237 #endif // EIGEN_DIAGONAL_H
    238