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 ref_selector<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   };
     53 };
     54 
     55 template<typename MatrixType,typename Lhs> struct homogeneous_left_product_impl;
     56 template<typename MatrixType,typename Rhs> struct homogeneous_right_product_impl;
     57 
     58 } // end namespace internal
     59 
     60 template<typename MatrixType,int _Direction> class Homogeneous
     61   : public MatrixBase<Homogeneous<MatrixType,_Direction> >, internal::no_assignment_operator
     62 {
     63   public:
     64 
     65     typedef MatrixType NestedExpression;
     66     enum { Direction = _Direction };
     67 
     68     typedef MatrixBase<Homogeneous> Base;
     69     EIGEN_DENSE_PUBLIC_INTERFACE(Homogeneous)
     70 
     71     EIGEN_DEVICE_FUNC explicit inline Homogeneous(const MatrixType& matrix)
     72       : m_matrix(matrix)
     73     {}
     74 
     75     EIGEN_DEVICE_FUNC inline Index rows() const { return m_matrix.rows() + (int(Direction)==Vertical   ? 1 : 0); }
     76     EIGEN_DEVICE_FUNC inline Index cols() const { return m_matrix.cols() + (int(Direction)==Horizontal ? 1 : 0); }
     77 
     78     EIGEN_DEVICE_FUNC const NestedExpression& nestedExpression() const { return m_matrix; }
     79 
     80     template<typename Rhs>
     81     EIGEN_DEVICE_FUNC inline const Product<Homogeneous,Rhs>
     82     operator* (const MatrixBase<Rhs>& rhs) const
     83     {
     84       eigen_assert(int(Direction)==Horizontal);
     85       return Product<Homogeneous,Rhs>(*this,rhs.derived());
     86     }
     87 
     88     template<typename Lhs> friend
     89     EIGEN_DEVICE_FUNC inline const Product<Lhs,Homogeneous>
     90     operator* (const MatrixBase<Lhs>& lhs, const Homogeneous& rhs)
     91     {
     92       eigen_assert(int(Direction)==Vertical);
     93       return Product<Lhs,Homogeneous>(lhs.derived(),rhs);
     94     }
     95 
     96     template<typename Scalar, int Dim, int Mode, int Options> friend
     97     EIGEN_DEVICE_FUNC inline const Product<Transform<Scalar,Dim,Mode,Options>, Homogeneous >
     98     operator* (const Transform<Scalar,Dim,Mode,Options>& lhs, const Homogeneous& rhs)
     99     {
    100       eigen_assert(int(Direction)==Vertical);
    101       return Product<Transform<Scalar,Dim,Mode,Options>, Homogeneous>(lhs,rhs);
    102     }
    103 
    104     template<typename Func>
    105     EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename internal::result_of<Func(Scalar,Scalar)>::type
    106     redux(const Func& func) const
    107     {
    108       return func(m_matrix.redux(func), Scalar(1));
    109     }
    110 
    111   protected:
    112     typename MatrixType::Nested m_matrix;
    113 };
    114 
    115 /** \geometry_module \ingroup Geometry_Module
    116   *
    117   * \returns a vector expression that is one longer than the vector argument, with the value 1 symbolically appended as the last coefficient.
    118   *
    119   * This can be used to convert affine coordinates to homogeneous coordinates.
    120   *
    121   * \only_for_vectors
    122   *
    123   * Example: \include MatrixBase_homogeneous.cpp
    124   * Output: \verbinclude MatrixBase_homogeneous.out
    125   *
    126   * \sa VectorwiseOp::homogeneous(), class Homogeneous
    127   */
    128 template<typename Derived>
    129 EIGEN_DEVICE_FUNC inline typename MatrixBase<Derived>::HomogeneousReturnType
    130 MatrixBase<Derived>::homogeneous() const
    131 {
    132   EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
    133   return HomogeneousReturnType(derived());
    134 }
    135 
    136 /** \geometry_module \ingroup Geometry_Module
    137   *
    138   * \returns an expression where the value 1 is symbolically appended as the final coefficient to each column (or row) of the matrix.
    139   *
    140   * This can be used to convert affine coordinates to homogeneous coordinates.
    141   *
    142   * Example: \include VectorwiseOp_homogeneous.cpp
    143   * Output: \verbinclude VectorwiseOp_homogeneous.out
    144   *
    145   * \sa MatrixBase::homogeneous(), class Homogeneous */
    146 template<typename ExpressionType, int Direction>
    147 EIGEN_DEVICE_FUNC inline Homogeneous<ExpressionType,Direction>
    148 VectorwiseOp<ExpressionType,Direction>::homogeneous() const
    149 {
    150   return HomogeneousReturnType(_expression());
    151 }
    152 
    153 /** \geometry_module \ingroup Geometry_Module
    154   *
    155   * \brief homogeneous normalization
    156   *
    157   * \returns a vector expression of the N-1 first coefficients of \c *this divided by that last coefficient.
    158   *
    159   * This can be used to convert homogeneous coordinates to affine coordinates.
    160   *
    161   * It is essentially a shortcut for:
    162   * \code
    163     this->head(this->size()-1)/this->coeff(this->size()-1);
    164     \endcode
    165   *
    166   * Example: \include MatrixBase_hnormalized.cpp
    167   * Output: \verbinclude MatrixBase_hnormalized.out
    168   *
    169   * \sa VectorwiseOp::hnormalized() */
    170 template<typename Derived>
    171 EIGEN_DEVICE_FUNC inline const typename MatrixBase<Derived>::HNormalizedReturnType
    172 MatrixBase<Derived>::hnormalized() const
    173 {
    174   EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
    175   return ConstStartMinusOne(derived(),0,0,
    176     ColsAtCompileTime==1?size()-1:1,
    177     ColsAtCompileTime==1?1:size()-1) / coeff(size()-1);
    178 }
    179 
    180 /** \geometry_module \ingroup Geometry_Module
    181   *
    182   * \brief column or row-wise homogeneous normalization
    183   *
    184   * \returns an expression of the first N-1 coefficients of each column (or row) of \c *this divided by the last coefficient of each column (or row).
    185   *
    186   * This can be used to convert homogeneous coordinates to affine coordinates.
    187   *
    188   * It is conceptually equivalent to calling MatrixBase::hnormalized() to each column (or row) of \c *this.
    189   *
    190   * Example: \include DirectionWise_hnormalized.cpp
    191   * Output: \verbinclude DirectionWise_hnormalized.out
    192   *
    193   * \sa MatrixBase::hnormalized() */
    194 template<typename ExpressionType, int Direction>
    195 EIGEN_DEVICE_FUNC inline const typename VectorwiseOp<ExpressionType,Direction>::HNormalizedReturnType
    196 VectorwiseOp<ExpressionType,Direction>::hnormalized() const
    197 {
    198   return HNormalized_Block(_expression(),0,0,
    199       Direction==Vertical   ? _expression().rows()-1 : _expression().rows(),
    200       Direction==Horizontal ? _expression().cols()-1 : _expression().cols()).cwiseQuotient(
    201       Replicate<HNormalized_Factors,
    202                 Direction==Vertical   ? HNormalized_SizeMinusOne : 1,
    203                 Direction==Horizontal ? HNormalized_SizeMinusOne : 1>
    204         (HNormalized_Factors(_expression(),
    205           Direction==Vertical    ? _expression().rows()-1:0,
    206           Direction==Horizontal  ? _expression().cols()-1:0,
    207           Direction==Vertical    ? 1 : _expression().rows(),
    208           Direction==Horizontal  ? 1 : _expression().cols()),
    209          Direction==Vertical   ? _expression().rows()-1 : 1,
    210          Direction==Horizontal ? _expression().cols()-1 : 1));
    211 }
    212 
    213 namespace internal {
    214 
    215 template<typename MatrixOrTransformType>
    216 struct take_matrix_for_product
    217 {
    218   typedef MatrixOrTransformType type;
    219   EIGEN_DEVICE_FUNC static const type& run(const type &x) { return x; }
    220 };
    221 
    222 template<typename Scalar, int Dim, int Mode,int Options>
    223 struct take_matrix_for_product<Transform<Scalar, Dim, Mode, Options> >
    224 {
    225   typedef Transform<Scalar, Dim, Mode, Options> TransformType;
    226   typedef typename internal::add_const<typename TransformType::ConstAffinePart>::type type;
    227   EIGEN_DEVICE_FUNC static type run (const TransformType& x) { return x.affine(); }
    228 };
    229 
    230 template<typename Scalar, int Dim, int Options>
    231 struct take_matrix_for_product<Transform<Scalar, Dim, Projective, Options> >
    232 {
    233   typedef Transform<Scalar, Dim, Projective, Options> TransformType;
    234   typedef typename TransformType::MatrixType type;
    235   EIGEN_DEVICE_FUNC static const type& run (const TransformType& x) { return x.matrix(); }
    236 };
    237 
    238 template<typename MatrixType,typename Lhs>
    239 struct traits<homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs> >
    240 {
    241   typedef typename take_matrix_for_product<Lhs>::type LhsMatrixType;
    242   typedef typename remove_all<MatrixType>::type MatrixTypeCleaned;
    243   typedef typename remove_all<LhsMatrixType>::type LhsMatrixTypeCleaned;
    244   typedef typename make_proper_matrix_type<
    245                  typename traits<MatrixTypeCleaned>::Scalar,
    246                  LhsMatrixTypeCleaned::RowsAtCompileTime,
    247                  MatrixTypeCleaned::ColsAtCompileTime,
    248                  MatrixTypeCleaned::PlainObject::Options,
    249                  LhsMatrixTypeCleaned::MaxRowsAtCompileTime,
    250                  MatrixTypeCleaned::MaxColsAtCompileTime>::type ReturnType;
    251 };
    252 
    253 template<typename MatrixType,typename Lhs>
    254 struct homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs>
    255   : public ReturnByValue<homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs> >
    256 {
    257   typedef typename traits<homogeneous_left_product_impl>::LhsMatrixType LhsMatrixType;
    258   typedef typename remove_all<LhsMatrixType>::type LhsMatrixTypeCleaned;
    259   typedef typename remove_all<typename LhsMatrixTypeCleaned::Nested>::type LhsMatrixTypeNested;
    260   EIGEN_DEVICE_FUNC homogeneous_left_product_impl(const Lhs& lhs, const MatrixType& rhs)
    261     : m_lhs(take_matrix_for_product<Lhs>::run(lhs)),
    262       m_rhs(rhs)
    263   {}
    264 
    265   EIGEN_DEVICE_FUNC inline Index rows() const { return m_lhs.rows(); }
    266   EIGEN_DEVICE_FUNC inline Index cols() const { return m_rhs.cols(); }
    267 
    268   template<typename Dest> EIGEN_DEVICE_FUNC void evalTo(Dest& dst) const
    269   {
    270     // FIXME investigate how to allow lazy evaluation of this product when possible
    271     dst = Block<const LhsMatrixTypeNested,
    272               LhsMatrixTypeNested::RowsAtCompileTime,
    273               LhsMatrixTypeNested::ColsAtCompileTime==Dynamic?Dynamic:LhsMatrixTypeNested::ColsAtCompileTime-1>
    274             (m_lhs,0,0,m_lhs.rows(),m_lhs.cols()-1) * m_rhs;
    275     dst += m_lhs.col(m_lhs.cols()-1).rowwise()
    276             .template replicate<MatrixType::ColsAtCompileTime>(m_rhs.cols());
    277   }
    278 
    279   typename LhsMatrixTypeCleaned::Nested m_lhs;
    280   typename MatrixType::Nested m_rhs;
    281 };
    282 
    283 template<typename MatrixType,typename Rhs>
    284 struct traits<homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs> >
    285 {
    286   typedef typename make_proper_matrix_type<typename traits<MatrixType>::Scalar,
    287                  MatrixType::RowsAtCompileTime,
    288                  Rhs::ColsAtCompileTime,
    289                  MatrixType::PlainObject::Options,
    290                  MatrixType::MaxRowsAtCompileTime,
    291                  Rhs::MaxColsAtCompileTime>::type ReturnType;
    292 };
    293 
    294 template<typename MatrixType,typename Rhs>
    295 struct homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs>
    296   : public ReturnByValue<homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs> >
    297 {
    298   typedef typename remove_all<typename Rhs::Nested>::type RhsNested;
    299   EIGEN_DEVICE_FUNC homogeneous_right_product_impl(const MatrixType& lhs, const Rhs& rhs)
    300     : m_lhs(lhs), m_rhs(rhs)
    301   {}
    302 
    303   EIGEN_DEVICE_FUNC inline Index rows() const { return m_lhs.rows(); }
    304   EIGEN_DEVICE_FUNC inline Index cols() const { return m_rhs.cols(); }
    305 
    306   template<typename Dest> EIGEN_DEVICE_FUNC void evalTo(Dest& dst) const
    307   {
    308     // FIXME investigate how to allow lazy evaluation of this product when possible
    309     dst = m_lhs * Block<const RhsNested,
    310                         RhsNested::RowsAtCompileTime==Dynamic?Dynamic:RhsNested::RowsAtCompileTime-1,
    311                         RhsNested::ColsAtCompileTime>
    312             (m_rhs,0,0,m_rhs.rows()-1,m_rhs.cols());
    313     dst += m_rhs.row(m_rhs.rows()-1).colwise()
    314             .template replicate<MatrixType::RowsAtCompileTime>(m_lhs.rows());
    315   }
    316 
    317   typename MatrixType::Nested m_lhs;
    318   typename Rhs::Nested m_rhs;
    319 };
    320 
    321 template<typename ArgType,int Direction>
    322 struct evaluator_traits<Homogeneous<ArgType,Direction> >
    323 {
    324   typedef typename storage_kind_to_evaluator_kind<typename ArgType::StorageKind>::Kind Kind;
    325   typedef HomogeneousShape Shape;
    326 };
    327 
    328 template<> struct AssignmentKind<DenseShape,HomogeneousShape> { typedef Dense2Dense Kind; };
    329 
    330 
    331 template<typename ArgType,int Direction>
    332 struct unary_evaluator<Homogeneous<ArgType,Direction>, IndexBased>
    333   : evaluator<typename Homogeneous<ArgType,Direction>::PlainObject >
    334 {
    335   typedef Homogeneous<ArgType,Direction> XprType;
    336   typedef typename XprType::PlainObject PlainObject;
    337   typedef evaluator<PlainObject> Base;
    338 
    339   EIGEN_DEVICE_FUNC explicit unary_evaluator(const XprType& op)
    340     : Base(), m_temp(op)
    341   {
    342     ::new (static_cast<Base*>(this)) Base(m_temp);
    343   }
    344 
    345 protected:
    346   PlainObject m_temp;
    347 };
    348 
    349 // dense = homogeneous
    350 template< typename DstXprType, typename ArgType, typename Scalar>
    351 struct Assignment<DstXprType, Homogeneous<ArgType,Vertical>, internal::assign_op<Scalar,typename ArgType::Scalar>, Dense2Dense>
    352 {
    353   typedef Homogeneous<ArgType,Vertical> SrcXprType;
    354   EIGEN_DEVICE_FUNC static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,typename ArgType::Scalar> &)
    355   {
    356     Index dstRows = src.rows();
    357     Index dstCols = src.cols();
    358     if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
    359       dst.resize(dstRows, dstCols);
    360 
    361     dst.template topRows<ArgType::RowsAtCompileTime>(src.nestedExpression().rows()) = src.nestedExpression();
    362     dst.row(dst.rows()-1).setOnes();
    363   }
    364 };
    365 
    366 // dense = homogeneous
    367 template< typename DstXprType, typename ArgType, typename Scalar>
    368 struct Assignment<DstXprType, Homogeneous<ArgType,Horizontal>, internal::assign_op<Scalar,typename ArgType::Scalar>, Dense2Dense>
    369 {
    370   typedef Homogeneous<ArgType,Horizontal> SrcXprType;
    371   EIGEN_DEVICE_FUNC static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,typename ArgType::Scalar> &)
    372   {
    373     Index dstRows = src.rows();
    374     Index dstCols = src.cols();
    375     if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
    376       dst.resize(dstRows, dstCols);
    377 
    378     dst.template leftCols<ArgType::ColsAtCompileTime>(src.nestedExpression().cols()) = src.nestedExpression();
    379     dst.col(dst.cols()-1).setOnes();
    380   }
    381 };
    382 
    383 template<typename LhsArg, typename Rhs, int ProductTag>
    384 struct generic_product_impl<Homogeneous<LhsArg,Horizontal>, Rhs, HomogeneousShape, DenseShape, ProductTag>
    385 {
    386   template<typename Dest>
    387   EIGEN_DEVICE_FUNC static void evalTo(Dest& dst, const Homogeneous<LhsArg,Horizontal>& lhs, const Rhs& rhs)
    388   {
    389     homogeneous_right_product_impl<Homogeneous<LhsArg,Horizontal>, Rhs>(lhs.nestedExpression(), rhs).evalTo(dst);
    390   }
    391 };
    392 
    393 template<typename Lhs,typename Rhs>
    394 struct homogeneous_right_product_refactoring_helper
    395 {
    396   enum {
    397     Dim  = Lhs::ColsAtCompileTime,
    398     Rows = Lhs::RowsAtCompileTime
    399   };
    400   typedef typename Rhs::template ConstNRowsBlockXpr<Dim>::Type          LinearBlockConst;
    401   typedef typename remove_const<LinearBlockConst>::type                 LinearBlock;
    402   typedef typename Rhs::ConstRowXpr                                     ConstantColumn;
    403   typedef Replicate<const ConstantColumn,Rows,1>                        ConstantBlock;
    404   typedef Product<Lhs,LinearBlock,LazyProduct>                          LinearProduct;
    405   typedef CwiseBinaryOp<internal::scalar_sum_op<typename Lhs::Scalar,typename Rhs::Scalar>, const LinearProduct, const ConstantBlock> Xpr;
    406 };
    407 
    408 template<typename Lhs, typename Rhs, int ProductTag>
    409 struct product_evaluator<Product<Lhs, Rhs, LazyProduct>, ProductTag, HomogeneousShape, DenseShape>
    410  : public evaluator<typename homogeneous_right_product_refactoring_helper<typename Lhs::NestedExpression,Rhs>::Xpr>
    411 {
    412   typedef Product<Lhs, Rhs, LazyProduct> XprType;
    413   typedef homogeneous_right_product_refactoring_helper<typename Lhs::NestedExpression,Rhs> helper;
    414   typedef typename helper::ConstantBlock ConstantBlock;
    415   typedef typename helper::Xpr RefactoredXpr;
    416   typedef evaluator<RefactoredXpr> Base;
    417 
    418   EIGEN_DEVICE_FUNC explicit product_evaluator(const XprType& xpr)
    419     : Base(  xpr.lhs().nestedExpression() .lazyProduct(  xpr.rhs().template topRows<helper::Dim>(xpr.lhs().nestedExpression().cols()) )
    420             + ConstantBlock(xpr.rhs().row(xpr.rhs().rows()-1),xpr.lhs().rows(), 1) )
    421   {}
    422 };
    423 
    424 template<typename Lhs, typename RhsArg, int ProductTag>
    425 struct generic_product_impl<Lhs, Homogeneous<RhsArg,Vertical>, DenseShape, HomogeneousShape, ProductTag>
    426 {
    427   template<typename Dest>
    428   EIGEN_DEVICE_FUNC static void evalTo(Dest& dst, const Lhs& lhs, const Homogeneous<RhsArg,Vertical>& rhs)
    429   {
    430     homogeneous_left_product_impl<Homogeneous<RhsArg,Vertical>, Lhs>(lhs, rhs.nestedExpression()).evalTo(dst);
    431   }
    432 };
    433 
    434 // TODO: the following specialization is to address a regression from 3.2 to 3.3
    435 // In the future, this path should be optimized.
    436 template<typename Lhs, typename RhsArg, int ProductTag>
    437 struct generic_product_impl<Lhs, Homogeneous<RhsArg,Vertical>, TriangularShape, HomogeneousShape, ProductTag>
    438 {
    439   template<typename Dest>
    440   static void evalTo(Dest& dst, const Lhs& lhs, const Homogeneous<RhsArg,Vertical>& rhs)
    441   {
    442     dst.noalias() = lhs * rhs.eval();
    443   }
    444 };
    445 
    446 template<typename Lhs,typename Rhs>
    447 struct homogeneous_left_product_refactoring_helper
    448 {
    449   enum {
    450     Dim = Rhs::RowsAtCompileTime,
    451     Cols = Rhs::ColsAtCompileTime
    452   };
    453   typedef typename Lhs::template ConstNColsBlockXpr<Dim>::Type          LinearBlockConst;
    454   typedef typename remove_const<LinearBlockConst>::type                 LinearBlock;
    455   typedef typename Lhs::ConstColXpr                                     ConstantColumn;
    456   typedef Replicate<const ConstantColumn,1,Cols>                        ConstantBlock;
    457   typedef Product<LinearBlock,Rhs,LazyProduct>                          LinearProduct;
    458   typedef CwiseBinaryOp<internal::scalar_sum_op<typename Lhs::Scalar,typename Rhs::Scalar>, const LinearProduct, const ConstantBlock> Xpr;
    459 };
    460 
    461 template<typename Lhs, typename Rhs, int ProductTag>
    462 struct product_evaluator<Product<Lhs, Rhs, LazyProduct>, ProductTag, DenseShape, HomogeneousShape>
    463  : public evaluator<typename homogeneous_left_product_refactoring_helper<Lhs,typename Rhs::NestedExpression>::Xpr>
    464 {
    465   typedef Product<Lhs, Rhs, LazyProduct> XprType;
    466   typedef homogeneous_left_product_refactoring_helper<Lhs,typename Rhs::NestedExpression> helper;
    467   typedef typename helper::ConstantBlock ConstantBlock;
    468   typedef typename helper::Xpr RefactoredXpr;
    469   typedef evaluator<RefactoredXpr> Base;
    470 
    471   EIGEN_DEVICE_FUNC explicit product_evaluator(const XprType& xpr)
    472     : Base(   xpr.lhs().template leftCols<helper::Dim>(xpr.rhs().nestedExpression().rows()) .lazyProduct( xpr.rhs().nestedExpression() )
    473             + ConstantBlock(xpr.lhs().col(xpr.lhs().cols()-1),1,xpr.rhs().cols()) )
    474   {}
    475 };
    476 
    477 template<typename Scalar, int Dim, int Mode,int Options, typename RhsArg, int ProductTag>
    478 struct generic_product_impl<Transform<Scalar,Dim,Mode,Options>, Homogeneous<RhsArg,Vertical>, DenseShape, HomogeneousShape, ProductTag>
    479 {
    480   typedef Transform<Scalar,Dim,Mode,Options> TransformType;
    481   template<typename Dest>
    482   EIGEN_DEVICE_FUNC static void evalTo(Dest& dst, const TransformType& lhs, const Homogeneous<RhsArg,Vertical>& rhs)
    483   {
    484     homogeneous_left_product_impl<Homogeneous<RhsArg,Vertical>, TransformType>(lhs, rhs.nestedExpression()).evalTo(dst);
    485   }
    486 };
    487 
    488 template<typename ExpressionType, int Side, bool Transposed>
    489 struct permutation_matrix_product<ExpressionType, Side, Transposed, HomogeneousShape>
    490   : public permutation_matrix_product<ExpressionType, Side, Transposed, DenseShape>
    491 {};
    492 
    493 } // end namespace internal
    494 
    495 } // end namespace Eigen
    496 
    497 #endif // EIGEN_HOMOGENEOUS_H
    498