Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:56:17

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
0005 // Copyright (C) 2009 Ricard Marxer <email@ricardmarxer.com>
0006 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
0007 //
0008 // This Source Code Form is subject to the terms of the Mozilla
0009 // Public License v. 2.0. If a copy of the MPL was not distributed
0010 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0011 
0012 #ifndef EIGEN_REVERSE_H
0013 #define EIGEN_REVERSE_H
0014 
0015 namespace Eigen {
0016 
0017 namespace internal {
0018 
0019 template<typename MatrixType, int Direction>
0020 struct traits<Reverse<MatrixType, Direction> >
0021  : traits<MatrixType>
0022 {
0023   typedef typename MatrixType::Scalar Scalar;
0024   typedef typename traits<MatrixType>::StorageKind StorageKind;
0025   typedef typename traits<MatrixType>::XprKind XprKind;
0026   typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
0027   typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
0028   enum {
0029     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
0030     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
0031     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0032     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
0033     Flags = _MatrixTypeNested::Flags & (RowMajorBit | LvalueBit)
0034   };
0035 };
0036 
0037 template<typename PacketType, bool ReversePacket> struct reverse_packet_cond
0038 {
0039   static inline PacketType run(const PacketType& x) { return preverse(x); }
0040 };
0041 
0042 template<typename PacketType> struct reverse_packet_cond<PacketType,false>
0043 {
0044   static inline PacketType run(const PacketType& x) { return x; }
0045 };
0046 
0047 } // end namespace internal
0048 
0049 /** \class Reverse
0050   * \ingroup Core_Module
0051   *
0052   * \brief Expression of the reverse of a vector or matrix
0053   *
0054   * \tparam MatrixType the type of the object of which we are taking the reverse
0055   * \tparam Direction defines the direction of the reverse operation, can be Vertical, Horizontal, or BothDirections
0056   *
0057   * This class represents an expression of the reverse of a vector.
0058   * It is the return type of MatrixBase::reverse() and VectorwiseOp::reverse()
0059   * and most of the time this is the only way it is used.
0060   *
0061   * \sa MatrixBase::reverse(), VectorwiseOp::reverse()
0062   */
0063 template<typename MatrixType, int Direction> class Reverse
0064   : public internal::dense_xpr_base< Reverse<MatrixType, Direction> >::type
0065 {
0066   public:
0067 
0068     typedef typename internal::dense_xpr_base<Reverse>::type Base;
0069     EIGEN_DENSE_PUBLIC_INTERFACE(Reverse)
0070     typedef typename internal::remove_all<MatrixType>::type NestedExpression;
0071     using Base::IsRowMajor;
0072 
0073   protected:
0074     enum {
0075       PacketSize = internal::packet_traits<Scalar>::size,
0076       IsColMajor = !IsRowMajor,
0077       ReverseRow = (Direction == Vertical)   || (Direction == BothDirections),
0078       ReverseCol = (Direction == Horizontal) || (Direction == BothDirections),
0079       OffsetRow  = ReverseRow && IsColMajor ? PacketSize : 1,
0080       OffsetCol  = ReverseCol && IsRowMajor ? PacketSize : 1,
0081       ReversePacket = (Direction == BothDirections)
0082                     || ((Direction == Vertical)   && IsColMajor)
0083                     || ((Direction == Horizontal) && IsRowMajor)
0084     };
0085     typedef internal::reverse_packet_cond<PacketScalar,ReversePacket> reverse_packet;
0086   public:
0087 
0088     EIGEN_DEVICE_FUNC explicit inline Reverse(const MatrixType& matrix) : m_matrix(matrix) { }
0089 
0090     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Reverse)
0091 
0092     EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
0093     inline Index rows() const EIGEN_NOEXCEPT { return m_matrix.rows(); }
0094     EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
0095     inline Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
0096 
0097     EIGEN_DEVICE_FUNC inline Index innerStride() const
0098     {
0099       return -m_matrix.innerStride();
0100     }
0101 
0102     EIGEN_DEVICE_FUNC const typename internal::remove_all<typename MatrixType::Nested>::type&
0103     nestedExpression() const
0104     {
0105       return m_matrix;
0106     }
0107 
0108   protected:
0109     typename MatrixType::Nested m_matrix;
0110 };
0111 
0112 /** \returns an expression of the reverse of *this.
0113   *
0114   * Example: \include MatrixBase_reverse.cpp
0115   * Output: \verbinclude MatrixBase_reverse.out
0116   *
0117   */
0118 template<typename Derived>
0119 EIGEN_DEVICE_FUNC inline typename DenseBase<Derived>::ReverseReturnType
0120 DenseBase<Derived>::reverse()
0121 {
0122   return ReverseReturnType(derived());
0123 }
0124 
0125 
0126 //reverse const overload moved DenseBase.h due to a CUDA compiler bug
0127 
0128 /** This is the "in place" version of reverse: it reverses \c *this.
0129   *
0130   * In most cases it is probably better to simply use the reversed expression
0131   * of a matrix. However, when reversing the matrix data itself is really needed,
0132   * then this "in-place" version is probably the right choice because it provides
0133   * the following additional benefits:
0134   *  - less error prone: doing the same operation with .reverse() requires special care:
0135   *    \code m = m.reverse().eval(); \endcode
0136   *  - this API enables reverse operations without the need for a temporary
0137   *  - it allows future optimizations (cache friendliness, etc.)
0138   *
0139   * \sa VectorwiseOp::reverseInPlace(), reverse() */
0140 template<typename Derived>
0141 EIGEN_DEVICE_FUNC inline void DenseBase<Derived>::reverseInPlace()
0142 {
0143   if(cols()>rows())
0144   {
0145     Index half = cols()/2;
0146     leftCols(half).swap(rightCols(half).reverse());
0147     if((cols()%2)==1)
0148     {
0149       Index half2 = rows()/2;
0150       col(half).head(half2).swap(col(half).tail(half2).reverse());
0151     }
0152   }
0153   else
0154   {
0155     Index half = rows()/2;
0156     topRows(half).swap(bottomRows(half).reverse());
0157     if((rows()%2)==1)
0158     {
0159       Index half2 = cols()/2;
0160       row(half).head(half2).swap(row(half).tail(half2).reverse());
0161     }
0162   }
0163 }
0164 
0165 namespace internal {
0166 
0167 template<int Direction>
0168 struct vectorwise_reverse_inplace_impl;
0169 
0170 template<>
0171 struct vectorwise_reverse_inplace_impl<Vertical>
0172 {
0173   template<typename ExpressionType>
0174   static void run(ExpressionType &xpr)
0175   {
0176     const int HalfAtCompileTime = ExpressionType::RowsAtCompileTime==Dynamic?Dynamic:ExpressionType::RowsAtCompileTime/2;
0177     Index half = xpr.rows()/2;
0178     xpr.topRows(fix<HalfAtCompileTime>(half))
0179        .swap(xpr.bottomRows(fix<HalfAtCompileTime>(half)).colwise().reverse());
0180   }
0181 };
0182 
0183 template<>
0184 struct vectorwise_reverse_inplace_impl<Horizontal>
0185 {
0186   template<typename ExpressionType>
0187   static void run(ExpressionType &xpr)
0188   {
0189     const int HalfAtCompileTime = ExpressionType::ColsAtCompileTime==Dynamic?Dynamic:ExpressionType::ColsAtCompileTime/2;
0190     Index half = xpr.cols()/2;
0191     xpr.leftCols(fix<HalfAtCompileTime>(half))
0192        .swap(xpr.rightCols(fix<HalfAtCompileTime>(half)).rowwise().reverse());
0193   }
0194 };
0195 
0196 } // end namespace internal
0197 
0198 /** This is the "in place" version of VectorwiseOp::reverse: it reverses each column or row of \c *this.
0199   *
0200   * In most cases it is probably better to simply use the reversed expression
0201   * of a matrix. However, when reversing the matrix data itself is really needed,
0202   * then this "in-place" version is probably the right choice because it provides
0203   * the following additional benefits:
0204   *  - less error prone: doing the same operation with .reverse() requires special care:
0205   *    \code m = m.reverse().eval(); \endcode
0206   *  - this API enables reverse operations without the need for a temporary
0207   *
0208   * \sa DenseBase::reverseInPlace(), reverse() */
0209 template<typename ExpressionType, int Direction>
0210 EIGEN_DEVICE_FUNC void VectorwiseOp<ExpressionType,Direction>::reverseInPlace()
0211 {
0212   internal::vectorwise_reverse_inplace_impl<Direction>::run(m_matrix);
0213 }
0214 
0215 } // end namespace Eigen
0216 
0217 #endif // EIGEN_REVERSE_H