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) 2006-2008 Benoit Jacob <jacob.benoit.1 (at) gmail.com>
      5 // Copyright (C) 2009 Ricard Marxer <email (at) ricardmarxer.com>
      6 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      7 //
      8 // This Source Code Form is subject to the terms of the Mozilla
      9 // Public License v. 2.0. If a copy of the MPL was not distributed
     10 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
     11 
     12 #ifndef EIGEN_REVERSE_H
     13 #define EIGEN_REVERSE_H
     14 
     15 namespace Eigen {
     16 
     17 /** \class Reverse
     18   * \ingroup Core_Module
     19   *
     20   * \brief Expression of the reverse of a vector or matrix
     21   *
     22   * \param MatrixType the type of the object of which we are taking the reverse
     23   *
     24   * This class represents an expression of the reverse of a vector.
     25   * It is the return type of MatrixBase::reverse() and VectorwiseOp::reverse()
     26   * and most of the time this is the only way it is used.
     27   *
     28   * \sa MatrixBase::reverse(), VectorwiseOp::reverse()
     29   */
     30 
     31 namespace internal {
     32 
     33 template<typename MatrixType, int Direction>
     34 struct traits<Reverse<MatrixType, Direction> >
     35  : traits<MatrixType>
     36 {
     37   typedef typename MatrixType::Scalar Scalar;
     38   typedef typename traits<MatrixType>::StorageKind StorageKind;
     39   typedef typename traits<MatrixType>::XprKind XprKind;
     40   typedef typename nested<MatrixType>::type MatrixTypeNested;
     41   typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
     42   enum {
     43     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
     44     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
     45     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
     46     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
     47 
     48     // let's enable LinearAccess only with vectorization because of the product overhead
     49     LinearAccess = ( (Direction==BothDirections) && (int(_MatrixTypeNested::Flags)&PacketAccessBit) )
     50                  ? LinearAccessBit : 0,
     51 
     52     Flags = int(_MatrixTypeNested::Flags) & (HereditaryBits | LvalueBit | PacketAccessBit | LinearAccess),
     53 
     54     CoeffReadCost = _MatrixTypeNested::CoeffReadCost
     55   };
     56 };
     57 
     58 template<typename PacketScalar, bool ReversePacket> struct reverse_packet_cond
     59 {
     60   static inline PacketScalar run(const PacketScalar& x) { return preverse(x); }
     61 };
     62 
     63 template<typename PacketScalar> struct reverse_packet_cond<PacketScalar,false>
     64 {
     65   static inline PacketScalar run(const PacketScalar& x) { return x; }
     66 };
     67 
     68 } // end namespace internal
     69 
     70 template<typename MatrixType, int Direction> class Reverse
     71   : public internal::dense_xpr_base< Reverse<MatrixType, Direction> >::type
     72 {
     73   public:
     74 
     75     typedef typename internal::dense_xpr_base<Reverse>::type Base;
     76     EIGEN_DENSE_PUBLIC_INTERFACE(Reverse)
     77     using Base::IsRowMajor;
     78 
     79     // next line is necessary because otherwise const version of operator()
     80     // is hidden by non-const version defined in this file
     81     using Base::operator();
     82 
     83   protected:
     84     enum {
     85       PacketSize = internal::packet_traits<Scalar>::size,
     86       IsColMajor = !IsRowMajor,
     87       ReverseRow = (Direction == Vertical)   || (Direction == BothDirections),
     88       ReverseCol = (Direction == Horizontal) || (Direction == BothDirections),
     89       OffsetRow  = ReverseRow && IsColMajor ? PacketSize : 1,
     90       OffsetCol  = ReverseCol && IsRowMajor ? PacketSize : 1,
     91       ReversePacket = (Direction == BothDirections)
     92                     || ((Direction == Vertical)   && IsColMajor)
     93                     || ((Direction == Horizontal) && IsRowMajor)
     94     };
     95     typedef internal::reverse_packet_cond<PacketScalar,ReversePacket> reverse_packet;
     96   public:
     97 
     98     inline Reverse(const MatrixType& matrix) : m_matrix(matrix) { }
     99 
    100     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Reverse)
    101 
    102     inline Index rows() const { return m_matrix.rows(); }
    103     inline Index cols() const { return m_matrix.cols(); }
    104 
    105     inline Index innerStride() const
    106     {
    107       return -m_matrix.innerStride();
    108     }
    109 
    110     inline Scalar& operator()(Index row, Index col)
    111     {
    112       eigen_assert(row >= 0 && row < rows() && col >= 0 && col < cols());
    113       return coeffRef(row, col);
    114     }
    115 
    116     inline Scalar& coeffRef(Index row, Index col)
    117     {
    118       return m_matrix.const_cast_derived().coeffRef(ReverseRow ? m_matrix.rows() - row - 1 : row,
    119                                                     ReverseCol ? m_matrix.cols() - col - 1 : col);
    120     }
    121 
    122     inline CoeffReturnType coeff(Index row, Index col) const
    123     {
    124       return m_matrix.coeff(ReverseRow ? m_matrix.rows() - row - 1 : row,
    125                             ReverseCol ? m_matrix.cols() - col - 1 : col);
    126     }
    127 
    128     inline CoeffReturnType coeff(Index index) const
    129     {
    130       return m_matrix.coeff(m_matrix.size() - index - 1);
    131     }
    132 
    133     inline Scalar& coeffRef(Index index)
    134     {
    135       return m_matrix.const_cast_derived().coeffRef(m_matrix.size() - index - 1);
    136     }
    137 
    138     inline Scalar& operator()(Index index)
    139     {
    140       eigen_assert(index >= 0 && index < m_matrix.size());
    141       return coeffRef(index);
    142     }
    143 
    144     template<int LoadMode>
    145     inline const PacketScalar packet(Index row, Index col) const
    146     {
    147       return reverse_packet::run(m_matrix.template packet<LoadMode>(
    148                                     ReverseRow ? m_matrix.rows() - row - OffsetRow : row,
    149                                     ReverseCol ? m_matrix.cols() - col - OffsetCol : col));
    150     }
    151 
    152     template<int LoadMode>
    153     inline void writePacket(Index row, Index col, const PacketScalar& x)
    154     {
    155       m_matrix.const_cast_derived().template writePacket<LoadMode>(
    156                                       ReverseRow ? m_matrix.rows() - row - OffsetRow : row,
    157                                       ReverseCol ? m_matrix.cols() - col - OffsetCol : col,
    158                                       reverse_packet::run(x));
    159     }
    160 
    161     template<int LoadMode>
    162     inline const PacketScalar packet(Index index) const
    163     {
    164       return internal::preverse(m_matrix.template packet<LoadMode>( m_matrix.size() - index - PacketSize ));
    165     }
    166 
    167     template<int LoadMode>
    168     inline void writePacket(Index index, const PacketScalar& x)
    169     {
    170       m_matrix.const_cast_derived().template writePacket<LoadMode>(m_matrix.size() - index - PacketSize, internal::preverse(x));
    171     }
    172 
    173     const typename internal::remove_all<typename MatrixType::Nested>::type&
    174     nestedExpression() const
    175     {
    176       return m_matrix;
    177     }
    178 
    179   protected:
    180     typename MatrixType::Nested m_matrix;
    181 };
    182 
    183 /** \returns an expression of the reverse of *this.
    184   *
    185   * Example: \include MatrixBase_reverse.cpp
    186   * Output: \verbinclude MatrixBase_reverse.out
    187   *
    188   */
    189 template<typename Derived>
    190 inline typename DenseBase<Derived>::ReverseReturnType
    191 DenseBase<Derived>::reverse()
    192 {
    193   return derived();
    194 }
    195 
    196 /** This is the const version of reverse(). */
    197 template<typename Derived>
    198 inline const typename DenseBase<Derived>::ConstReverseReturnType
    199 DenseBase<Derived>::reverse() const
    200 {
    201   return derived();
    202 }
    203 
    204 /** This is the "in place" version of reverse: it reverses \c *this.
    205   *
    206   * In most cases it is probably better to simply use the reversed expression
    207   * of a matrix. However, when reversing the matrix data itself is really needed,
    208   * then this "in-place" version is probably the right choice because it provides
    209   * the following additional features:
    210   *  - less error prone: doing the same operation with .reverse() requires special care:
    211   *    \code m = m.reverse().eval(); \endcode
    212   *  - this API allows to avoid creating a temporary (the current implementation creates a temporary, but that could be avoided using swap)
    213   *  - it allows future optimizations (cache friendliness, etc.)
    214   *
    215   * \sa reverse() */
    216 template<typename Derived>
    217 inline void DenseBase<Derived>::reverseInPlace()
    218 {
    219   derived() = derived().reverse().eval();
    220 }
    221 
    222 } // end namespace Eigen
    223 
    224 #endif // EIGEN_REVERSE_H
    225