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) == Dynamic || 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 == Dynamic ? 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     typedef typename internal::dense_xpr_base<Diagonal>::type Base;
     70     EIGEN_DENSE_PUBLIC_INTERFACE(Diagonal)
     71 
     72     inline Diagonal(MatrixType& matrix, Index index = DiagIndex) : m_matrix(matrix), m_index(index) {}
     73 
     74     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Diagonal)
     75 
     76     inline Index rows() const
     77     { return m_index.value()<0 ? (std::min)(m_matrix.cols(),m_matrix.rows()+m_index.value()) : (std::min)(m_matrix.rows(),m_matrix.cols()-m_index.value()); }
     78 
     79     inline Index cols() const { return 1; }
     80 
     81     inline Index innerStride() const
     82     {
     83       return m_matrix.outerStride() + 1;
     84     }
     85 
     86     inline Index outerStride() const
     87     {
     88       return 0;
     89     }
     90 
     91     typedef typename internal::conditional<
     92                        internal::is_lvalue<MatrixType>::value,
     93                        Scalar,
     94                        const Scalar
     95                      >::type ScalarWithConstIfNotLvalue;
     96 
     97     inline ScalarWithConstIfNotLvalue* data() { return &(m_matrix.const_cast_derived().coeffRef(rowOffset(), colOffset())); }
     98     inline const Scalar* data() const { return &(m_matrix.const_cast_derived().coeffRef(rowOffset(), colOffset())); }
     99 
    100     inline Scalar& coeffRef(Index row, Index)
    101     {
    102       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
    103       return m_matrix.const_cast_derived().coeffRef(row+rowOffset(), row+colOffset());
    104     }
    105 
    106     inline const Scalar& coeffRef(Index row, Index) const
    107     {
    108       return m_matrix.const_cast_derived().coeffRef(row+rowOffset(), row+colOffset());
    109     }
    110 
    111     inline CoeffReturnType coeff(Index row, Index) const
    112     {
    113       return m_matrix.coeff(row+rowOffset(), row+colOffset());
    114     }
    115 
    116     inline Scalar& coeffRef(Index index)
    117     {
    118       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
    119       return m_matrix.const_cast_derived().coeffRef(index+rowOffset(), index+colOffset());
    120     }
    121 
    122     inline const Scalar& coeffRef(Index index) const
    123     {
    124       return m_matrix.const_cast_derived().coeffRef(index+rowOffset(), index+colOffset());
    125     }
    126 
    127     inline CoeffReturnType coeff(Index index) const
    128     {
    129       return m_matrix.coeff(index+rowOffset(), index+colOffset());
    130     }
    131 
    132     const typename internal::remove_all<typename MatrixType::Nested>::type&
    133     nestedExpression() const
    134     {
    135       return m_matrix;
    136     }
    137 
    138     int index() const
    139     {
    140       return m_index.value();
    141     }
    142 
    143   protected:
    144     typename MatrixType::Nested m_matrix;
    145     const internal::variable_if_dynamic<Index, DiagIndex> m_index;
    146 
    147   private:
    148     // some compilers may fail to optimize std::max etc in case of compile-time constants...
    149     EIGEN_STRONG_INLINE Index absDiagIndex() const { return m_index.value()>0 ? m_index.value() : -m_index.value(); }
    150     EIGEN_STRONG_INLINE Index rowOffset() const { return m_index.value()>0 ? 0 : -m_index.value(); }
    151     EIGEN_STRONG_INLINE Index colOffset() const { return m_index.value()>0 ? m_index.value() : 0; }
    152     // triger a compile time error is someone try to call packet
    153     template<int LoadMode> typename MatrixType::PacketReturnType packet(Index) const;
    154     template<int LoadMode> typename MatrixType::PacketReturnType packet(Index,Index) const;
    155 };
    156 
    157 /** \returns an expression of the main diagonal of the matrix \c *this
    158   *
    159   * \c *this is not required to be square.
    160   *
    161   * Example: \include MatrixBase_diagonal.cpp
    162   * Output: \verbinclude MatrixBase_diagonal.out
    163   *
    164   * \sa class Diagonal */
    165 template<typename Derived>
    166 inline typename MatrixBase<Derived>::DiagonalReturnType
    167 MatrixBase<Derived>::diagonal()
    168 {
    169   return derived();
    170 }
    171 
    172 /** This is the const version of diagonal(). */
    173 template<typename Derived>
    174 inline const typename MatrixBase<Derived>::ConstDiagonalReturnType
    175 MatrixBase<Derived>::diagonal() const
    176 {
    177   return ConstDiagonalReturnType(derived());
    178 }
    179 
    180 /** \returns an expression of the \a DiagIndex-th sub or super diagonal of the matrix \c *this
    181   *
    182   * \c *this is not required to be square.
    183   *
    184   * The template parameter \a DiagIndex represent a super diagonal if \a DiagIndex > 0
    185   * and a sub diagonal otherwise. \a DiagIndex == 0 is equivalent to the main diagonal.
    186   *
    187   * Example: \include MatrixBase_diagonal_int.cpp
    188   * Output: \verbinclude MatrixBase_diagonal_int.out
    189   *
    190   * \sa MatrixBase::diagonal(), class Diagonal */
    191 template<typename Derived>
    192 inline typename MatrixBase<Derived>::template DiagonalIndexReturnType<Dynamic>::Type
    193 MatrixBase<Derived>::diagonal(Index index)
    194 {
    195   return typename DiagonalIndexReturnType<Dynamic>::Type(derived(), index);
    196 }
    197 
    198 /** This is the const version of diagonal(Index). */
    199 template<typename Derived>
    200 inline typename MatrixBase<Derived>::template ConstDiagonalIndexReturnType<Dynamic>::Type
    201 MatrixBase<Derived>::diagonal(Index index) const
    202 {
    203   return typename ConstDiagonalIndexReturnType<Dynamic>::Type(derived(), index);
    204 }
    205 
    206 /** \returns an expression of the \a DiagIndex-th sub or super diagonal of the matrix \c *this
    207   *
    208   * \c *this is not required to be square.
    209   *
    210   * The template parameter \a DiagIndex represent a super diagonal if \a DiagIndex > 0
    211   * and a sub diagonal otherwise. \a DiagIndex == 0 is equivalent to the main diagonal.
    212   *
    213   * Example: \include MatrixBase_diagonal_template_int.cpp
    214   * Output: \verbinclude MatrixBase_diagonal_template_int.out
    215   *
    216   * \sa MatrixBase::diagonal(), class Diagonal */
    217 template<typename Derived>
    218 template<int Index>
    219 inline typename MatrixBase<Derived>::template DiagonalIndexReturnType<Index>::Type
    220 MatrixBase<Derived>::diagonal()
    221 {
    222   return derived();
    223 }
    224 
    225 /** This is the const version of diagonal<int>(). */
    226 template<typename Derived>
    227 template<int Index>
    228 inline typename MatrixBase<Derived>::template ConstDiagonalIndexReturnType<Index>::Type
    229 MatrixBase<Derived>::diagonal() const
    230 {
    231   return derived();
    232 }
    233 
    234 } // end namespace Eigen
    235 
    236 #endif // EIGEN_DIAGONAL_H
    237