Home | History | Annotate | Download | only in Geometry
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 //
      6 // This Source Code Form is subject to the terms of the Mozilla
      7 // Public License v. 2.0. If a copy of the MPL was not distributed
      8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
      9 
     10 #ifndef EIGEN_HOMOGENEOUS_H
     11 #define EIGEN_HOMOGENEOUS_H
     12 
     13 namespace Eigen {
     14 
     15 /** \geometry_module \ingroup Geometry_Module
     16   *
     17   * \class Homogeneous
     18   *
     19   * \brief Expression of one (or a set of) homogeneous vector(s)
     20   *
     21   * \param MatrixType the type of the object in which we are making homogeneous
     22   *
     23   * This class represents an expression of one (or a set of) homogeneous vector(s).
     24   * It is the return type of MatrixBase::homogeneous() and most of the time
     25   * this is the only way it is used.
     26   *
     27   * \sa MatrixBase::homogeneous()
     28   */
     29 
     30 namespace internal {
     31 
     32 template<typename MatrixType,int Direction>
     33 struct traits<Homogeneous<MatrixType,Direction> >
     34  : traits<MatrixType>
     35 {
     36   typedef typename traits<MatrixType>::StorageKind StorageKind;
     37   typedef typename nested<MatrixType>::type MatrixTypeNested;
     38   typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
     39   enum {
     40     RowsPlusOne = (MatrixType::RowsAtCompileTime != Dynamic) ?
     41                   int(MatrixType::RowsAtCompileTime) + 1 : Dynamic,
     42     ColsPlusOne = (MatrixType::ColsAtCompileTime != Dynamic) ?
     43                   int(MatrixType::ColsAtCompileTime) + 1 : Dynamic,
     44     RowsAtCompileTime = Direction==Vertical  ?  RowsPlusOne : MatrixType::RowsAtCompileTime,
     45     ColsAtCompileTime = Direction==Horizontal ? ColsPlusOne : MatrixType::ColsAtCompileTime,
     46     MaxRowsAtCompileTime = RowsAtCompileTime,
     47     MaxColsAtCompileTime = ColsAtCompileTime,
     48     TmpFlags = _MatrixTypeNested::Flags & HereditaryBits,
     49     Flags = ColsAtCompileTime==1 ? (TmpFlags & ~RowMajorBit)
     50           : RowsAtCompileTime==1 ? (TmpFlags | RowMajorBit)
     51           : TmpFlags,
     52     CoeffReadCost = _MatrixTypeNested::CoeffReadCost
     53   };
     54 };
     55 
     56 template<typename MatrixType,typename Lhs> struct homogeneous_left_product_impl;
     57 template<typename MatrixType,typename Rhs> struct homogeneous_right_product_impl;
     58 
     59 } // end namespace internal
     60 
     61 template<typename MatrixType,int _Direction> class Homogeneous
     62   : public MatrixBase<Homogeneous<MatrixType,_Direction> >
     63 {
     64   public:
     65 
     66     enum { Direction = _Direction };
     67 
     68     typedef MatrixBase<Homogeneous> Base;
     69     EIGEN_DENSE_PUBLIC_INTERFACE(Homogeneous)
     70 
     71     inline Homogeneous(const MatrixType& matrix)
     72       : m_matrix(matrix)
     73     {}
     74 
     75     inline Index rows() const { return m_matrix.rows() + (int(Direction)==Vertical   ? 1 : 0); }
     76     inline Index cols() const { return m_matrix.cols() + (int(Direction)==Horizontal ? 1 : 0); }
     77 
     78     inline Scalar coeff(Index row, Index col) const
     79     {
     80       if(  (int(Direction)==Vertical   && row==m_matrix.rows())
     81         || (int(Direction)==Horizontal && col==m_matrix.cols()))
     82         return 1;
     83       return m_matrix.coeff(row, col);
     84     }
     85 
     86     template<typename Rhs>
     87     inline const internal::homogeneous_right_product_impl<Homogeneous,Rhs>
     88     operator* (const MatrixBase<Rhs>& rhs) const
     89     {
     90       eigen_assert(int(Direction)==Horizontal);
     91       return internal::homogeneous_right_product_impl<Homogeneous,Rhs>(m_matrix,rhs.derived());
     92     }
     93 
     94     template<typename Lhs> friend
     95     inline const internal::homogeneous_left_product_impl<Homogeneous,Lhs>
     96     operator* (const MatrixBase<Lhs>& lhs, const Homogeneous& rhs)
     97     {
     98       eigen_assert(int(Direction)==Vertical);
     99       return internal::homogeneous_left_product_impl<Homogeneous,Lhs>(lhs.derived(),rhs.m_matrix);
    100     }
    101 
    102     template<typename Scalar, int Dim, int Mode, int Options> friend
    103     inline const internal::homogeneous_left_product_impl<Homogeneous,Transform<Scalar,Dim,Mode,Options> >
    104     operator* (const Transform<Scalar,Dim,Mode,Options>& lhs, const Homogeneous& rhs)
    105     {
    106       eigen_assert(int(Direction)==Vertical);
    107       return internal::homogeneous_left_product_impl<Homogeneous,Transform<Scalar,Dim,Mode,Options> >(lhs,rhs.m_matrix);
    108     }
    109 
    110   protected:
    111     typename MatrixType::Nested m_matrix;
    112 };
    113 
    114 /** \geometry_module
    115   *
    116   * \return an expression of the equivalent homogeneous vector
    117   *
    118   * \only_for_vectors
    119   *
    120   * Example: \include MatrixBase_homogeneous.cpp
    121   * Output: \verbinclude MatrixBase_homogeneous.out
    122   *
    123   * \sa class Homogeneous
    124   */
    125 template<typename Derived>
    126 inline typename MatrixBase<Derived>::HomogeneousReturnType
    127 MatrixBase<Derived>::homogeneous() const
    128 {
    129   EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
    130   return derived();
    131 }
    132 
    133 /** \geometry_module
    134   *
    135   * \returns a matrix expression of homogeneous column (or row) vectors
    136   *
    137   * Example: \include VectorwiseOp_homogeneous.cpp
    138   * Output: \verbinclude VectorwiseOp_homogeneous.out
    139   *
    140   * \sa MatrixBase::homogeneous() */
    141 template<typename ExpressionType, int Direction>
    142 inline Homogeneous<ExpressionType,Direction>
    143 VectorwiseOp<ExpressionType,Direction>::homogeneous() const
    144 {
    145   return _expression();
    146 }
    147 
    148 /** \geometry_module
    149   *
    150   * \returns an expression of the homogeneous normalized vector of \c *this
    151   *
    152   * Example: \include MatrixBase_hnormalized.cpp
    153   * Output: \verbinclude MatrixBase_hnormalized.out
    154   *
    155   * \sa VectorwiseOp::hnormalized() */
    156 template<typename Derived>
    157 inline const typename MatrixBase<Derived>::HNormalizedReturnType
    158 MatrixBase<Derived>::hnormalized() const
    159 {
    160   EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
    161   return ConstStartMinusOne(derived(),0,0,
    162     ColsAtCompileTime==1?size()-1:1,
    163     ColsAtCompileTime==1?1:size()-1) / coeff(size()-1);
    164 }
    165 
    166 /** \geometry_module
    167   *
    168   * \returns an expression of the homogeneous normalized vector of \c *this
    169   *
    170   * Example: \include DirectionWise_hnormalized.cpp
    171   * Output: \verbinclude DirectionWise_hnormalized.out
    172   *
    173   * \sa MatrixBase::hnormalized() */
    174 template<typename ExpressionType, int Direction>
    175 inline const typename VectorwiseOp<ExpressionType,Direction>::HNormalizedReturnType
    176 VectorwiseOp<ExpressionType,Direction>::hnormalized() const
    177 {
    178   return HNormalized_Block(_expression(),0,0,
    179       Direction==Vertical   ? _expression().rows()-1 : _expression().rows(),
    180       Direction==Horizontal ? _expression().cols()-1 : _expression().cols()).cwiseQuotient(
    181       Replicate<HNormalized_Factors,
    182                 Direction==Vertical   ? HNormalized_SizeMinusOne : 1,
    183                 Direction==Horizontal ? HNormalized_SizeMinusOne : 1>
    184         (HNormalized_Factors(_expression(),
    185           Direction==Vertical    ? _expression().rows()-1:0,
    186           Direction==Horizontal  ? _expression().cols()-1:0,
    187           Direction==Vertical    ? 1 : _expression().rows(),
    188           Direction==Horizontal  ? 1 : _expression().cols()),
    189          Direction==Vertical   ? _expression().rows()-1 : 1,
    190          Direction==Horizontal ? _expression().cols()-1 : 1));
    191 }
    192 
    193 namespace internal {
    194 
    195 template<typename MatrixOrTransformType>
    196 struct take_matrix_for_product
    197 {
    198   typedef MatrixOrTransformType type;
    199   static const type& run(const type &x) { return x; }
    200 };
    201 
    202 template<typename Scalar, int Dim, int Mode,int Options>
    203 struct take_matrix_for_product<Transform<Scalar, Dim, Mode, Options> >
    204 {
    205   typedef Transform<Scalar, Dim, Mode, Options> TransformType;
    206   typedef typename internal::add_const<typename TransformType::ConstAffinePart>::type type;
    207   static type run (const TransformType& x) { return x.affine(); }
    208 };
    209 
    210 template<typename Scalar, int Dim, int Options>
    211 struct take_matrix_for_product<Transform<Scalar, Dim, Projective, Options> >
    212 {
    213   typedef Transform<Scalar, Dim, Projective, Options> TransformType;
    214   typedef typename TransformType::MatrixType type;
    215   static const type& run (const TransformType& x) { return x.matrix(); }
    216 };
    217 
    218 template<typename MatrixType,typename Lhs>
    219 struct traits<homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs> >
    220 {
    221   typedef typename take_matrix_for_product<Lhs>::type LhsMatrixType;
    222   typedef typename remove_all<MatrixType>::type MatrixTypeCleaned;
    223   typedef typename remove_all<LhsMatrixType>::type LhsMatrixTypeCleaned;
    224   typedef typename make_proper_matrix_type<
    225                  typename traits<MatrixTypeCleaned>::Scalar,
    226                  LhsMatrixTypeCleaned::RowsAtCompileTime,
    227                  MatrixTypeCleaned::ColsAtCompileTime,
    228                  MatrixTypeCleaned::PlainObject::Options,
    229                  LhsMatrixTypeCleaned::MaxRowsAtCompileTime,
    230                  MatrixTypeCleaned::MaxColsAtCompileTime>::type ReturnType;
    231 };
    232 
    233 template<typename MatrixType,typename Lhs>
    234 struct homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs>
    235   : public ReturnByValue<homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs> >
    236 {
    237   typedef typename traits<homogeneous_left_product_impl>::LhsMatrixType LhsMatrixType;
    238   typedef typename remove_all<LhsMatrixType>::type LhsMatrixTypeCleaned;
    239   typedef typename remove_all<typename LhsMatrixTypeCleaned::Nested>::type LhsMatrixTypeNested;
    240   typedef typename MatrixType::Index Index;
    241   homogeneous_left_product_impl(const Lhs& lhs, const MatrixType& rhs)
    242     : m_lhs(take_matrix_for_product<Lhs>::run(lhs)),
    243       m_rhs(rhs)
    244   {}
    245 
    246   inline Index rows() const { return m_lhs.rows(); }
    247   inline Index cols() const { return m_rhs.cols(); }
    248 
    249   template<typename Dest> void evalTo(Dest& dst) const
    250   {
    251     // FIXME investigate how to allow lazy evaluation of this product when possible
    252     dst = Block<const LhsMatrixTypeNested,
    253               LhsMatrixTypeNested::RowsAtCompileTime,
    254               LhsMatrixTypeNested::ColsAtCompileTime==Dynamic?Dynamic:LhsMatrixTypeNested::ColsAtCompileTime-1>
    255             (m_lhs,0,0,m_lhs.rows(),m_lhs.cols()-1) * m_rhs;
    256     dst += m_lhs.col(m_lhs.cols()-1).rowwise()
    257             .template replicate<MatrixType::ColsAtCompileTime>(m_rhs.cols());
    258   }
    259 
    260   typename LhsMatrixTypeCleaned::Nested m_lhs;
    261   typename MatrixType::Nested m_rhs;
    262 };
    263 
    264 template<typename MatrixType,typename Rhs>
    265 struct traits<homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs> >
    266 {
    267   typedef typename make_proper_matrix_type<typename traits<MatrixType>::Scalar,
    268                  MatrixType::RowsAtCompileTime,
    269                  Rhs::ColsAtCompileTime,
    270                  MatrixType::PlainObject::Options,
    271                  MatrixType::MaxRowsAtCompileTime,
    272                  Rhs::MaxColsAtCompileTime>::type ReturnType;
    273 };
    274 
    275 template<typename MatrixType,typename Rhs>
    276 struct homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs>
    277   : public ReturnByValue<homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs> >
    278 {
    279   typedef typename remove_all<typename Rhs::Nested>::type RhsNested;
    280   typedef typename MatrixType::Index Index;
    281   homogeneous_right_product_impl(const MatrixType& lhs, const Rhs& rhs)
    282     : m_lhs(lhs), m_rhs(rhs)
    283   {}
    284 
    285   inline Index rows() const { return m_lhs.rows(); }
    286   inline Index cols() const { return m_rhs.cols(); }
    287 
    288   template<typename Dest> void evalTo(Dest& dst) const
    289   {
    290     // FIXME investigate how to allow lazy evaluation of this product when possible
    291     dst = m_lhs * Block<const RhsNested,
    292                         RhsNested::RowsAtCompileTime==Dynamic?Dynamic:RhsNested::RowsAtCompileTime-1,
    293                         RhsNested::ColsAtCompileTime>
    294             (m_rhs,0,0,m_rhs.rows()-1,m_rhs.cols());
    295     dst += m_rhs.row(m_rhs.rows()-1).colwise()
    296             .template replicate<MatrixType::RowsAtCompileTime>(m_lhs.rows());
    297   }
    298 
    299   typename MatrixType::Nested m_lhs;
    300   typename Rhs::Nested m_rhs;
    301 };
    302 
    303 } // end namespace internal
    304 
    305 } // end namespace Eigen
    306 
    307 #endif // EIGEN_HOMOGENEOUS_H
    308