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-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_TRANSPOSE_H
     12 #define EIGEN_TRANSPOSE_H
     13 
     14 namespace Eigen {
     15 
     16 /** \class Transpose
     17   * \ingroup Core_Module
     18   *
     19   * \brief Expression of the transpose of a matrix
     20   *
     21   * \param MatrixType the type of the object of which we are taking the transpose
     22   *
     23   * This class represents an expression of the transpose of a matrix.
     24   * It is the return type of MatrixBase::transpose() and MatrixBase::adjoint()
     25   * and most of the time this is the only way it is used.
     26   *
     27   * \sa MatrixBase::transpose(), MatrixBase::adjoint()
     28   */
     29 
     30 namespace internal {
     31 template<typename MatrixType>
     32 struct traits<Transpose<MatrixType> > : traits<MatrixType>
     33 {
     34   typedef typename MatrixType::Scalar Scalar;
     35   typedef typename nested<MatrixType>::type MatrixTypeNested;
     36   typedef typename remove_reference<MatrixTypeNested>::type MatrixTypeNestedPlain;
     37   typedef typename traits<MatrixType>::StorageKind StorageKind;
     38   typedef typename traits<MatrixType>::XprKind XprKind;
     39   enum {
     40     RowsAtCompileTime = MatrixType::ColsAtCompileTime,
     41     ColsAtCompileTime = MatrixType::RowsAtCompileTime,
     42     MaxRowsAtCompileTime = MatrixType::MaxColsAtCompileTime,
     43     MaxColsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
     44     FlagsLvalueBit = is_lvalue<MatrixType>::value ? LvalueBit : 0,
     45     Flags0 = MatrixTypeNestedPlain::Flags & ~(LvalueBit | NestByRefBit),
     46     Flags1 = Flags0 | FlagsLvalueBit,
     47     Flags = Flags1 ^ RowMajorBit,
     48     CoeffReadCost = MatrixTypeNestedPlain::CoeffReadCost,
     49     InnerStrideAtCompileTime = inner_stride_at_compile_time<MatrixType>::ret,
     50     OuterStrideAtCompileTime = outer_stride_at_compile_time<MatrixType>::ret
     51   };
     52 };
     53 }
     54 
     55 template<typename MatrixType, typename StorageKind> class TransposeImpl;
     56 
     57 template<typename MatrixType> class Transpose
     58   : public TransposeImpl<MatrixType,typename internal::traits<MatrixType>::StorageKind>
     59 {
     60   public:
     61 
     62     typedef typename TransposeImpl<MatrixType,typename internal::traits<MatrixType>::StorageKind>::Base Base;
     63     EIGEN_GENERIC_PUBLIC_INTERFACE(Transpose)
     64 
     65     inline Transpose(MatrixType& matrix) : m_matrix(matrix) {}
     66 
     67     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Transpose)
     68 
     69     inline Index rows() const { return m_matrix.cols(); }
     70     inline Index cols() const { return m_matrix.rows(); }
     71 
     72     /** \returns the nested expression */
     73     const typename internal::remove_all<typename MatrixType::Nested>::type&
     74     nestedExpression() const { return m_matrix; }
     75 
     76     /** \returns the nested expression */
     77     typename internal::remove_all<typename MatrixType::Nested>::type&
     78     nestedExpression() { return m_matrix.const_cast_derived(); }
     79 
     80   protected:
     81     typename MatrixType::Nested m_matrix;
     82 };
     83 
     84 namespace internal {
     85 
     86 template<typename MatrixType, bool HasDirectAccess = has_direct_access<MatrixType>::ret>
     87 struct TransposeImpl_base
     88 {
     89   typedef typename dense_xpr_base<Transpose<MatrixType> >::type type;
     90 };
     91 
     92 template<typename MatrixType>
     93 struct TransposeImpl_base<MatrixType, false>
     94 {
     95   typedef typename dense_xpr_base<Transpose<MatrixType> >::type type;
     96 };
     97 
     98 } // end namespace internal
     99 
    100 template<typename MatrixType> class TransposeImpl<MatrixType,Dense>
    101   : public internal::TransposeImpl_base<MatrixType>::type
    102 {
    103   public:
    104 
    105     typedef typename internal::TransposeImpl_base<MatrixType>::type Base;
    106     EIGEN_DENSE_PUBLIC_INTERFACE(Transpose<MatrixType>)
    107 
    108     inline Index innerStride() const { return derived().nestedExpression().innerStride(); }
    109     inline Index outerStride() const { return derived().nestedExpression().outerStride(); }
    110 
    111     typedef typename internal::conditional<
    112                        internal::is_lvalue<MatrixType>::value,
    113                        Scalar,
    114                        const Scalar
    115                      >::type ScalarWithConstIfNotLvalue;
    116 
    117     inline ScalarWithConstIfNotLvalue* data() { return derived().nestedExpression().data(); }
    118     inline const Scalar* data() const { return derived().nestedExpression().data(); }
    119 
    120     inline ScalarWithConstIfNotLvalue& coeffRef(Index row, Index col)
    121     {
    122       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
    123       return derived().nestedExpression().const_cast_derived().coeffRef(col, row);
    124     }
    125 
    126     inline ScalarWithConstIfNotLvalue& coeffRef(Index index)
    127     {
    128       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
    129       return derived().nestedExpression().const_cast_derived().coeffRef(index);
    130     }
    131 
    132     inline const Scalar& coeffRef(Index row, Index col) const
    133     {
    134       return derived().nestedExpression().coeffRef(col, row);
    135     }
    136 
    137     inline const Scalar& coeffRef(Index index) const
    138     {
    139       return derived().nestedExpression().coeffRef(index);
    140     }
    141 
    142     inline CoeffReturnType coeff(Index row, Index col) const
    143     {
    144       return derived().nestedExpression().coeff(col, row);
    145     }
    146 
    147     inline CoeffReturnType coeff(Index index) const
    148     {
    149       return derived().nestedExpression().coeff(index);
    150     }
    151 
    152     template<int LoadMode>
    153     inline const PacketScalar packet(Index row, Index col) const
    154     {
    155       return derived().nestedExpression().template packet<LoadMode>(col, row);
    156     }
    157 
    158     template<int LoadMode>
    159     inline void writePacket(Index row, Index col, const PacketScalar& x)
    160     {
    161       derived().nestedExpression().const_cast_derived().template writePacket<LoadMode>(col, row, x);
    162     }
    163 
    164     template<int LoadMode>
    165     inline const PacketScalar packet(Index index) const
    166     {
    167       return derived().nestedExpression().template packet<LoadMode>(index);
    168     }
    169 
    170     template<int LoadMode>
    171     inline void writePacket(Index index, const PacketScalar& x)
    172     {
    173       derived().nestedExpression().const_cast_derived().template writePacket<LoadMode>(index, x);
    174     }
    175 };
    176 
    177 /** \returns an expression of the transpose of *this.
    178   *
    179   * Example: \include MatrixBase_transpose.cpp
    180   * Output: \verbinclude MatrixBase_transpose.out
    181   *
    182   * \warning If you want to replace a matrix by its own transpose, do \b NOT do this:
    183   * \code
    184   * m = m.transpose(); // bug!!! caused by aliasing effect
    185   * \endcode
    186   * Instead, use the transposeInPlace() method:
    187   * \code
    188   * m.transposeInPlace();
    189   * \endcode
    190   * which gives Eigen good opportunities for optimization, or alternatively you can also do:
    191   * \code
    192   * m = m.transpose().eval();
    193   * \endcode
    194   *
    195   * \sa transposeInPlace(), adjoint() */
    196 template<typename Derived>
    197 inline Transpose<Derived>
    198 DenseBase<Derived>::transpose()
    199 {
    200   return derived();
    201 }
    202 
    203 /** This is the const version of transpose().
    204   *
    205   * Make sure you read the warning for transpose() !
    206   *
    207   * \sa transposeInPlace(), adjoint() */
    208 template<typename Derived>
    209 inline const typename DenseBase<Derived>::ConstTransposeReturnType
    210 DenseBase<Derived>::transpose() const
    211 {
    212   return ConstTransposeReturnType(derived());
    213 }
    214 
    215 /** \returns an expression of the adjoint (i.e. conjugate transpose) of *this.
    216   *
    217   * Example: \include MatrixBase_adjoint.cpp
    218   * Output: \verbinclude MatrixBase_adjoint.out
    219   *
    220   * \warning If you want to replace a matrix by its own adjoint, do \b NOT do this:
    221   * \code
    222   * m = m.adjoint(); // bug!!! caused by aliasing effect
    223   * \endcode
    224   * Instead, use the adjointInPlace() method:
    225   * \code
    226   * m.adjointInPlace();
    227   * \endcode
    228   * which gives Eigen good opportunities for optimization, or alternatively you can also do:
    229   * \code
    230   * m = m.adjoint().eval();
    231   * \endcode
    232   *
    233   * \sa adjointInPlace(), transpose(), conjugate(), class Transpose, class internal::scalar_conjugate_op */
    234 template<typename Derived>
    235 inline const typename MatrixBase<Derived>::AdjointReturnType
    236 MatrixBase<Derived>::adjoint() const
    237 {
    238   return this->transpose(); // in the complex case, the .conjugate() is be implicit here
    239                             // due to implicit conversion to return type
    240 }
    241 
    242 /***************************************************************************
    243 * "in place" transpose implementation
    244 ***************************************************************************/
    245 
    246 namespace internal {
    247 
    248 template<typename MatrixType,
    249   bool IsSquare = (MatrixType::RowsAtCompileTime == MatrixType::ColsAtCompileTime) && MatrixType::RowsAtCompileTime!=Dynamic>
    250 struct inplace_transpose_selector;
    251 
    252 template<typename MatrixType>
    253 struct inplace_transpose_selector<MatrixType,true> { // square matrix
    254   static void run(MatrixType& m) {
    255     m.template triangularView<StrictlyUpper>().swap(m.transpose());
    256   }
    257 };
    258 
    259 template<typename MatrixType>
    260 struct inplace_transpose_selector<MatrixType,false> { // non square matrix
    261   static void run(MatrixType& m) {
    262     if (m.rows()==m.cols())
    263       m.template triangularView<StrictlyUpper>().swap(m.transpose());
    264     else
    265       m = m.transpose().eval();
    266   }
    267 };
    268 
    269 } // end namespace internal
    270 
    271 /** This is the "in place" version of transpose(): it replaces \c *this by its own transpose.
    272   * Thus, doing
    273   * \code
    274   * m.transposeInPlace();
    275   * \endcode
    276   * has the same effect on m as doing
    277   * \code
    278   * m = m.transpose().eval();
    279   * \endcode
    280   * and is faster and also safer because in the latter line of code, forgetting the eval() results
    281   * in a bug caused by aliasing.
    282   *
    283   * Notice however that this method is only useful if you want to replace a matrix by its own transpose.
    284   * If you just need the transpose of a matrix, use transpose().
    285   *
    286   * \note if the matrix is not square, then \c *this must be a resizable matrix.
    287   *
    288   * \sa transpose(), adjoint(), adjointInPlace() */
    289 template<typename Derived>
    290 inline void DenseBase<Derived>::transposeInPlace()
    291 {
    292   internal::inplace_transpose_selector<Derived>::run(derived());
    293 }
    294 
    295 /***************************************************************************
    296 * "in place" adjoint implementation
    297 ***************************************************************************/
    298 
    299 /** This is the "in place" version of adjoint(): it replaces \c *this by its own transpose.
    300   * Thus, doing
    301   * \code
    302   * m.adjointInPlace();
    303   * \endcode
    304   * has the same effect on m as doing
    305   * \code
    306   * m = m.adjoint().eval();
    307   * \endcode
    308   * and is faster and also safer because in the latter line of code, forgetting the eval() results
    309   * in a bug caused by aliasing.
    310   *
    311   * Notice however that this method is only useful if you want to replace a matrix by its own adjoint.
    312   * If you just need the adjoint of a matrix, use adjoint().
    313   *
    314   * \note if the matrix is not square, then \c *this must be a resizable matrix.
    315   *
    316   * \sa transpose(), adjoint(), transposeInPlace() */
    317 template<typename Derived>
    318 inline void MatrixBase<Derived>::adjointInPlace()
    319 {
    320   derived() = adjoint().eval();
    321 }
    322 
    323 #ifndef EIGEN_NO_DEBUG
    324 
    325 // The following is to detect aliasing problems in most common cases.
    326 
    327 namespace internal {
    328 
    329 template<typename BinOp,typename NestedXpr,typename Rhs>
    330 struct blas_traits<SelfCwiseBinaryOp<BinOp,NestedXpr,Rhs> >
    331  : blas_traits<NestedXpr>
    332 {
    333   typedef SelfCwiseBinaryOp<BinOp,NestedXpr,Rhs> XprType;
    334   static inline const XprType extract(const XprType& x) { return x; }
    335 };
    336 
    337 template<bool DestIsTransposed, typename OtherDerived>
    338 struct check_transpose_aliasing_compile_time_selector
    339 {
    340   enum { ret = bool(blas_traits<OtherDerived>::IsTransposed) != DestIsTransposed };
    341 };
    342 
    343 template<bool DestIsTransposed, typename BinOp, typename DerivedA, typename DerivedB>
    344 struct check_transpose_aliasing_compile_time_selector<DestIsTransposed,CwiseBinaryOp<BinOp,DerivedA,DerivedB> >
    345 {
    346   enum { ret =    bool(blas_traits<DerivedA>::IsTransposed) != DestIsTransposed
    347                || bool(blas_traits<DerivedB>::IsTransposed) != DestIsTransposed
    348   };
    349 };
    350 
    351 template<typename Scalar, bool DestIsTransposed, typename OtherDerived>
    352 struct check_transpose_aliasing_run_time_selector
    353 {
    354   static bool run(const Scalar* dest, const OtherDerived& src)
    355   {
    356     return (bool(blas_traits<OtherDerived>::IsTransposed) != DestIsTransposed) && (dest!=0 && dest==(Scalar*)extract_data(src));
    357   }
    358 };
    359 
    360 template<typename Scalar, bool DestIsTransposed, typename BinOp, typename DerivedA, typename DerivedB>
    361 struct check_transpose_aliasing_run_time_selector<Scalar,DestIsTransposed,CwiseBinaryOp<BinOp,DerivedA,DerivedB> >
    362 {
    363   static bool run(const Scalar* dest, const CwiseBinaryOp<BinOp,DerivedA,DerivedB>& src)
    364   {
    365     return ((blas_traits<DerivedA>::IsTransposed != DestIsTransposed) && (dest!=0 && dest==(Scalar*)extract_data(src.lhs())))
    366         || ((blas_traits<DerivedB>::IsTransposed != DestIsTransposed) && (dest!=0 && dest==(Scalar*)extract_data(src.rhs())));
    367   }
    368 };
    369 
    370 // the following selector, checkTransposeAliasing_impl, based on MightHaveTransposeAliasing,
    371 // is because when the condition controlling the assert is known at compile time, ICC emits a warning.
    372 // This is actually a good warning: in expressions that don't have any transposing, the condition is
    373 // known at compile time to be false, and using that, we can avoid generating the code of the assert again
    374 // and again for all these expressions that don't need it.
    375 
    376 template<typename Derived, typename OtherDerived,
    377          bool MightHaveTransposeAliasing
    378                  = check_transpose_aliasing_compile_time_selector
    379                      <blas_traits<Derived>::IsTransposed,OtherDerived>::ret
    380         >
    381 struct checkTransposeAliasing_impl
    382 {
    383     static void run(const Derived& dst, const OtherDerived& other)
    384     {
    385         eigen_assert((!check_transpose_aliasing_run_time_selector
    386                       <typename Derived::Scalar,blas_traits<Derived>::IsTransposed,OtherDerived>
    387                       ::run(extract_data(dst), other))
    388           && "aliasing detected during tranposition, use transposeInPlace() "
    389              "or evaluate the rhs into a temporary using .eval()");
    390 
    391     }
    392 };
    393 
    394 template<typename Derived, typename OtherDerived>
    395 struct checkTransposeAliasing_impl<Derived, OtherDerived, false>
    396 {
    397     static void run(const Derived&, const OtherDerived&)
    398     {
    399     }
    400 };
    401 
    402 } // end namespace internal
    403 
    404 template<typename Derived>
    405 template<typename OtherDerived>
    406 void DenseBase<Derived>::checkTransposeAliasing(const OtherDerived& other) const
    407 {
    408     internal::checkTransposeAliasing_impl<Derived, OtherDerived>::run(derived(), other);
    409 }
    410 #endif
    411 
    412 } // end namespace Eigen
    413 
    414 #endif // EIGEN_TRANSPOSE_H
    415