Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/Geometry/Transform.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
0006 // Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
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_TRANSFORM_H
0013 #define EIGEN_TRANSFORM_H
0014 
0015 namespace Eigen {
0016 
0017 namespace internal {
0018 
0019 template<typename Transform>
0020 struct transform_traits
0021 {
0022   enum
0023   {
0024     Dim = Transform::Dim,
0025     HDim = Transform::HDim,
0026     Mode = Transform::Mode,
0027     IsProjective = (int(Mode)==int(Projective))
0028   };
0029 };
0030 
0031 template< typename TransformType,
0032           typename MatrixType,
0033           int Case = transform_traits<TransformType>::IsProjective ? 0
0034                    : int(MatrixType::RowsAtCompileTime) == int(transform_traits<TransformType>::HDim) ? 1
0035                    : 2,
0036           int RhsCols = MatrixType::ColsAtCompileTime>
0037 struct transform_right_product_impl;
0038 
0039 template< typename Other,
0040           int Mode,
0041           int Options,
0042           int Dim,
0043           int HDim,
0044           int OtherRows=Other::RowsAtCompileTime,
0045           int OtherCols=Other::ColsAtCompileTime>
0046 struct transform_left_product_impl;
0047 
0048 template< typename Lhs,
0049           typename Rhs,
0050           bool AnyProjective =
0051             transform_traits<Lhs>::IsProjective ||
0052             transform_traits<Rhs>::IsProjective>
0053 struct transform_transform_product_impl;
0054 
0055 template< typename Other,
0056           int Mode,
0057           int Options,
0058           int Dim,
0059           int HDim,
0060           int OtherRows=Other::RowsAtCompileTime,
0061           int OtherCols=Other::ColsAtCompileTime>
0062 struct transform_construct_from_matrix;
0063 
0064 template<typename TransformType> struct transform_take_affine_part;
0065 
0066 template<typename _Scalar, int _Dim, int _Mode, int _Options>
0067 struct traits<Transform<_Scalar,_Dim,_Mode,_Options> >
0068 {
0069   typedef _Scalar Scalar;
0070   typedef Eigen::Index StorageIndex;
0071   typedef Dense StorageKind;
0072   enum {
0073     Dim1 = _Dim==Dynamic ? _Dim : _Dim + 1,
0074     RowsAtCompileTime = _Mode==Projective ? Dim1 : _Dim,
0075     ColsAtCompileTime = Dim1,
0076     MaxRowsAtCompileTime = RowsAtCompileTime,
0077     MaxColsAtCompileTime = ColsAtCompileTime,
0078     Flags = 0
0079   };
0080 };
0081 
0082 template<int Mode> struct transform_make_affine;
0083 
0084 } // end namespace internal
0085 
0086 /** \geometry_module \ingroup Geometry_Module
0087   *
0088   * \class Transform
0089   *
0090   * \brief Represents an homogeneous transformation in a N dimensional space
0091   *
0092   * \tparam _Scalar the scalar type, i.e., the type of the coefficients
0093   * \tparam _Dim the dimension of the space
0094   * \tparam _Mode the type of the transformation. Can be:
0095   *              - #Affine: the transformation is stored as a (Dim+1)^2 matrix,
0096   *                         where the last row is assumed to be [0 ... 0 1].
0097   *              - #AffineCompact: the transformation is stored as a (Dim)x(Dim+1) matrix.
0098   *              - #Projective: the transformation is stored as a (Dim+1)^2 matrix
0099   *                             without any assumption.
0100   *              - #Isometry: same as #Affine with the additional assumption that
0101   *                           the linear part represents a rotation. This assumption is exploited
0102   *                           to speed up some functions such as inverse() and rotation().
0103   * \tparam _Options has the same meaning as in class Matrix. It allows to specify DontAlign and/or RowMajor.
0104   *                  These Options are passed directly to the underlying matrix type.
0105   *
0106   * The homography is internally represented and stored by a matrix which
0107   * is available through the matrix() method. To understand the behavior of
0108   * this class you have to think a Transform object as its internal
0109   * matrix representation. The chosen convention is right multiply:
0110   *
0111   * \code v' = T * v \endcode
0112   *
0113   * Therefore, an affine transformation matrix M is shaped like this:
0114   *
0115   * \f$ \left( \begin{array}{cc}
0116   * linear & translation\\
0117   * 0 ... 0 & 1
0118   * \end{array} \right) \f$
0119   *
0120   * Note that for a projective transformation the last row can be anything,
0121   * and then the interpretation of different parts might be slightly different.
0122   *
0123   * However, unlike a plain matrix, the Transform class provides many features
0124   * simplifying both its assembly and usage. In particular, it can be composed
0125   * with any other transformations (Transform,Translation,RotationBase,DiagonalMatrix)
0126   * and can be directly used to transform implicit homogeneous vectors. All these
0127   * operations are handled via the operator*. For the composition of transformations,
0128   * its principle consists to first convert the right/left hand sides of the product
0129   * to a compatible (Dim+1)^2 matrix and then perform a pure matrix product.
0130   * Of course, internally, operator* tries to perform the minimal number of operations
0131   * according to the nature of each terms. Likewise, when applying the transform
0132   * to points, the latters are automatically promoted to homogeneous vectors
0133   * before doing the matrix product. The conventions to homogeneous representations
0134   * are performed as follow:
0135   *
0136   * \b Translation t (Dim)x(1):
0137   * \f$ \left( \begin{array}{cc}
0138   * I & t \\
0139   * 0\,...\,0 & 1
0140   * \end{array} \right) \f$
0141   *
0142   * \b Rotation R (Dim)x(Dim):
0143   * \f$ \left( \begin{array}{cc}
0144   * R & 0\\
0145   * 0\,...\,0 & 1
0146   * \end{array} \right) \f$
0147   *<!--
0148   * \b Linear \b Matrix L (Dim)x(Dim):
0149   * \f$ \left( \begin{array}{cc}
0150   * L & 0\\
0151   * 0\,...\,0 & 1
0152   * \end{array} \right) \f$
0153   *
0154   * \b Affine \b Matrix A (Dim)x(Dim+1):
0155   * \f$ \left( \begin{array}{c}
0156   * A\\
0157   * 0\,...\,0\,1
0158   * \end{array} \right) \f$
0159   *-->
0160   * \b Scaling \b DiagonalMatrix S (Dim)x(Dim):
0161   * \f$ \left( \begin{array}{cc}
0162   * S & 0\\
0163   * 0\,...\,0 & 1
0164   * \end{array} \right) \f$
0165   *
0166   * \b Column \b point v (Dim)x(1):
0167   * \f$ \left( \begin{array}{c}
0168   * v\\
0169   * 1
0170   * \end{array} \right) \f$
0171   *
0172   * \b Set \b of \b column \b points V1...Vn (Dim)x(n):
0173   * \f$ \left( \begin{array}{ccc}
0174   * v_1 & ... & v_n\\
0175   * 1 & ... & 1
0176   * \end{array} \right) \f$
0177   *
0178   * The concatenation of a Transform object with any kind of other transformation
0179   * always returns a Transform object.
0180   *
0181   * A little exception to the "as pure matrix product" rule is the case of the
0182   * transformation of non homogeneous vectors by an affine transformation. In
0183   * that case the last matrix row can be ignored, and the product returns non
0184   * homogeneous vectors.
0185   *
0186   * Since, for instance, a Dim x Dim matrix is interpreted as a linear transformation,
0187   * it is not possible to directly transform Dim vectors stored in a Dim x Dim matrix.
0188   * The solution is either to use a Dim x Dynamic matrix or explicitly request a
0189   * vector transformation by making the vector homogeneous:
0190   * \code
0191   * m' = T * m.colwise().homogeneous();
0192   * \endcode
0193   * Note that there is zero overhead.
0194   *
0195   * Conversion methods from/to Qt's QMatrix and QTransform are available if the
0196   * preprocessor token EIGEN_QT_SUPPORT is defined.
0197   *
0198   * This class can be extended with the help of the plugin mechanism described on the page
0199   * \ref TopicCustomizing_Plugins by defining the preprocessor symbol \c EIGEN_TRANSFORM_PLUGIN.
0200   *
0201   * \sa class Matrix, class Quaternion
0202   */
0203 template<typename _Scalar, int _Dim, int _Mode, int _Options>
0204 class Transform
0205 {
0206 public:
0207   EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_Dim==Dynamic ? Dynamic : (_Dim+1)*(_Dim+1))
0208   enum {
0209     Mode = _Mode,
0210     Options = _Options,
0211     Dim = _Dim,     ///< space dimension in which the transformation holds
0212     HDim = _Dim+1,  ///< size of a respective homogeneous vector
0213     Rows = int(Mode)==(AffineCompact) ? Dim : HDim
0214   };
0215   /** the scalar type of the coefficients */
0216   typedef _Scalar Scalar;
0217   typedef Eigen::Index StorageIndex;
0218   typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
0219   /** type of the matrix used to represent the transformation */
0220   typedef typename internal::make_proper_matrix_type<Scalar,Rows,HDim,Options>::type MatrixType;
0221   /** constified MatrixType */
0222   typedef const MatrixType ConstMatrixType;
0223   /** type of the matrix used to represent the linear part of the transformation */
0224   typedef Matrix<Scalar,Dim,Dim,Options> LinearMatrixType;
0225   /** type of read/write reference to the linear part of the transformation */
0226   typedef Block<MatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (int(Options)&RowMajor)==0> LinearPart;
0227   /** type of read reference to the linear part of the transformation */
0228   typedef const Block<ConstMatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (int(Options)&RowMajor)==0> ConstLinearPart;
0229   /** type of read/write reference to the affine part of the transformation */
0230   typedef typename internal::conditional<int(Mode)==int(AffineCompact),
0231                               MatrixType&,
0232                               Block<MatrixType,Dim,HDim> >::type AffinePart;
0233   /** type of read reference to the affine part of the transformation */
0234   typedef typename internal::conditional<int(Mode)==int(AffineCompact),
0235                               const MatrixType&,
0236                               const Block<const MatrixType,Dim,HDim> >::type ConstAffinePart;
0237   /** type of a vector */
0238   typedef Matrix<Scalar,Dim,1> VectorType;
0239   /** type of a read/write reference to the translation part of the rotation */
0240   typedef Block<MatrixType,Dim,1,!(internal::traits<MatrixType>::Flags & RowMajorBit)> TranslationPart;
0241   /** type of a read reference to the translation part of the rotation */
0242   typedef const Block<ConstMatrixType,Dim,1,!(internal::traits<MatrixType>::Flags & RowMajorBit)> ConstTranslationPart;
0243   /** corresponding translation type */
0244   typedef Translation<Scalar,Dim> TranslationType;
0245 
0246   // this intermediate enum is needed to avoid an ICE with gcc 3.4 and 4.0
0247   enum { TransformTimeDiagonalMode = ((Mode==int(Isometry))?Affine:int(Mode)) };
0248   /** The return type of the product between a diagonal matrix and a transform */
0249   typedef Transform<Scalar,Dim,TransformTimeDiagonalMode> TransformTimeDiagonalReturnType;
0250 
0251 protected:
0252 
0253   MatrixType m_matrix;
0254 
0255 public:
0256 
0257   /** Default constructor without initialization of the meaningful coefficients.
0258     * If Mode==Affine or Mode==Isometry, then the last row is set to [0 ... 0 1] */
0259   EIGEN_DEVICE_FUNC inline Transform()
0260   {
0261     check_template_params();
0262     internal::transform_make_affine<(int(Mode)==Affine || int(Mode)==Isometry) ? Affine : AffineCompact>::run(m_matrix);
0263   }
0264 
0265   EIGEN_DEVICE_FUNC inline explicit Transform(const TranslationType& t)
0266   {
0267     check_template_params();
0268     *this = t;
0269   }
0270   EIGEN_DEVICE_FUNC inline explicit Transform(const UniformScaling<Scalar>& s)
0271   {
0272     check_template_params();
0273     *this = s;
0274   }
0275   template<typename Derived>
0276   EIGEN_DEVICE_FUNC inline explicit Transform(const RotationBase<Derived, Dim>& r)
0277   {
0278     check_template_params();
0279     *this = r;
0280   }
0281 
0282   typedef internal::transform_take_affine_part<Transform> take_affine_part;
0283 
0284   /** Constructs and initializes a transformation from a Dim^2 or a (Dim+1)^2 matrix. */
0285   template<typename OtherDerived>
0286   EIGEN_DEVICE_FUNC inline explicit Transform(const EigenBase<OtherDerived>& other)
0287   {
0288     EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
0289       YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
0290 
0291     check_template_params();
0292     internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
0293   }
0294 
0295   /** Set \c *this from a Dim^2 or (Dim+1)^2 matrix. */
0296   template<typename OtherDerived>
0297   EIGEN_DEVICE_FUNC inline Transform& operator=(const EigenBase<OtherDerived>& other)
0298   {
0299     EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
0300       YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
0301 
0302     internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
0303     return *this;
0304   }
0305 
0306   template<int OtherOptions>
0307   EIGEN_DEVICE_FUNC inline Transform(const Transform<Scalar,Dim,Mode,OtherOptions>& other)
0308   {
0309     check_template_params();
0310     // only the options change, we can directly copy the matrices
0311     m_matrix = other.matrix();
0312   }
0313 
0314   template<int OtherMode,int OtherOptions>
0315   EIGEN_DEVICE_FUNC inline Transform(const Transform<Scalar,Dim,OtherMode,OtherOptions>& other)
0316   {
0317     check_template_params();
0318     // prevent conversions as:
0319     // Affine | AffineCompact | Isometry = Projective
0320     EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Projective), Mode==int(Projective)),
0321                         YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
0322 
0323     // prevent conversions as:
0324     // Isometry = Affine | AffineCompact
0325     EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Affine)||OtherMode==int(AffineCompact), Mode!=int(Isometry)),
0326                         YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
0327 
0328     enum { ModeIsAffineCompact = Mode == int(AffineCompact),
0329            OtherModeIsAffineCompact = OtherMode == int(AffineCompact)
0330     };
0331 
0332     if(EIGEN_CONST_CONDITIONAL(ModeIsAffineCompact == OtherModeIsAffineCompact))
0333     {
0334       // We need the block expression because the code is compiled for all
0335       // combinations of transformations and will trigger a compile time error
0336       // if one tries to assign the matrices directly
0337       m_matrix.template block<Dim,Dim+1>(0,0) = other.matrix().template block<Dim,Dim+1>(0,0);
0338       makeAffine();
0339     }
0340     else if(EIGEN_CONST_CONDITIONAL(OtherModeIsAffineCompact))
0341     {
0342       typedef typename Transform<Scalar,Dim,OtherMode,OtherOptions>::MatrixType OtherMatrixType;
0343       internal::transform_construct_from_matrix<OtherMatrixType,Mode,Options,Dim,HDim>::run(this, other.matrix());
0344     }
0345     else
0346     {
0347       // here we know that Mode == AffineCompact and OtherMode != AffineCompact.
0348       // if OtherMode were Projective, the static assert above would already have caught it.
0349       // So the only possibility is that OtherMode == Affine
0350       linear() = other.linear();
0351       translation() = other.translation();
0352     }
0353   }
0354 
0355   template<typename OtherDerived>
0356   EIGEN_DEVICE_FUNC Transform(const ReturnByValue<OtherDerived>& other)
0357   {
0358     check_template_params();
0359     other.evalTo(*this);
0360   }
0361 
0362   template<typename OtherDerived>
0363   EIGEN_DEVICE_FUNC Transform& operator=(const ReturnByValue<OtherDerived>& other)
0364   {
0365     other.evalTo(*this);
0366     return *this;
0367   }
0368 
0369   #ifdef EIGEN_QT_SUPPORT
0370   inline Transform(const QMatrix& other);
0371   inline Transform& operator=(const QMatrix& other);
0372   inline QMatrix toQMatrix(void) const;
0373   inline Transform(const QTransform& other);
0374   inline Transform& operator=(const QTransform& other);
0375   inline QTransform toQTransform(void) const;
0376   #endif
0377 
0378   EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return int(Mode)==int(Projective) ? m_matrix.cols() : (m_matrix.cols()-1); }
0379   EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
0380 
0381   /** shortcut for m_matrix(row,col);
0382     * \sa MatrixBase::operator(Index,Index) const */
0383   EIGEN_DEVICE_FUNC inline Scalar operator() (Index row, Index col) const { return m_matrix(row,col); }
0384   /** shortcut for m_matrix(row,col);
0385     * \sa MatrixBase::operator(Index,Index) */
0386   EIGEN_DEVICE_FUNC inline Scalar& operator() (Index row, Index col) { return m_matrix(row,col); }
0387 
0388   /** \returns a read-only expression of the transformation matrix */
0389   EIGEN_DEVICE_FUNC inline const MatrixType& matrix() const { return m_matrix; }
0390   /** \returns a writable expression of the transformation matrix */
0391   EIGEN_DEVICE_FUNC inline MatrixType& matrix() { return m_matrix; }
0392 
0393   /** \returns a read-only expression of the linear part of the transformation */
0394   EIGEN_DEVICE_FUNC inline ConstLinearPart linear() const { return ConstLinearPart(m_matrix,0,0); }
0395   /** \returns a writable expression of the linear part of the transformation */
0396   EIGEN_DEVICE_FUNC inline LinearPart linear() { return LinearPart(m_matrix,0,0); }
0397 
0398   /** \returns a read-only expression of the Dim x HDim affine part of the transformation */
0399   EIGEN_DEVICE_FUNC inline ConstAffinePart affine() const { return take_affine_part::run(m_matrix); }
0400   /** \returns a writable expression of the Dim x HDim affine part of the transformation */
0401   EIGEN_DEVICE_FUNC inline AffinePart affine() { return take_affine_part::run(m_matrix); }
0402 
0403   /** \returns a read-only expression of the translation vector of the transformation */
0404   EIGEN_DEVICE_FUNC inline ConstTranslationPart translation() const { return ConstTranslationPart(m_matrix,0,Dim); }
0405   /** \returns a writable expression of the translation vector of the transformation */
0406   EIGEN_DEVICE_FUNC inline TranslationPart translation() { return TranslationPart(m_matrix,0,Dim); }
0407 
0408   /** \returns an expression of the product between the transform \c *this and a matrix expression \a other.
0409     *
0410     * The right-hand-side \a other can be either:
0411     * \li an homogeneous vector of size Dim+1,
0412     * \li a set of homogeneous vectors of size Dim+1 x N,
0413     * \li a transformation matrix of size Dim+1 x Dim+1.
0414     *
0415     * Moreover, if \c *this represents an affine transformation (i.e., Mode!=Projective), then \a other can also be:
0416     * \li a point of size Dim (computes: \code this->linear() * other + this->translation()\endcode),
0417     * \li a set of N points as a Dim x N matrix (computes: \code (this->linear() * other).colwise() + this->translation()\endcode),
0418     *
0419     * In all cases, the return type is a matrix or vector of same sizes as the right-hand-side \a other.
0420     *
0421     * If you want to interpret \a other as a linear or affine transformation, then first convert it to a Transform<> type,
0422     * or do your own cooking.
0423     *
0424     * Finally, if you want to apply Affine transformations to vectors, then explicitly apply the linear part only:
0425     * \code
0426     * Affine3f A;
0427     * Vector3f v1, v2;
0428     * v2 = A.linear() * v1;
0429     * \endcode
0430     *
0431     */
0432   // note: this function is defined here because some compilers cannot find the respective declaration
0433   template<typename OtherDerived>
0434   EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const typename internal::transform_right_product_impl<Transform, OtherDerived>::ResultType
0435   operator * (const EigenBase<OtherDerived> &other) const
0436   { return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this,other.derived()); }
0437 
0438   /** \returns the product expression of a transformation matrix \a a times a transform \a b
0439     *
0440     * The left hand side \a other can be either:
0441     * \li a linear transformation matrix of size Dim x Dim,
0442     * \li an affine transformation matrix of size Dim x Dim+1,
0443     * \li a general transformation matrix of size Dim+1 x Dim+1.
0444     */
0445   template<typename OtherDerived> friend
0446   EIGEN_DEVICE_FUNC inline const typename internal::transform_left_product_impl<OtherDerived,Mode,Options,_Dim,_Dim+1>::ResultType
0447     operator * (const EigenBase<OtherDerived> &a, const Transform &b)
0448   { return internal::transform_left_product_impl<OtherDerived,Mode,Options,Dim,HDim>::run(a.derived(),b); }
0449 
0450   /** \returns The product expression of a transform \a a times a diagonal matrix \a b
0451     *
0452     * The rhs diagonal matrix is interpreted as an affine scaling transformation. The
0453     * product results in a Transform of the same type (mode) as the lhs only if the lhs
0454     * mode is no isometry. In that case, the returned transform is an affinity.
0455     */
0456   template<typename DiagonalDerived>
0457   EIGEN_DEVICE_FUNC inline const TransformTimeDiagonalReturnType
0458     operator * (const DiagonalBase<DiagonalDerived> &b) const
0459   {
0460     TransformTimeDiagonalReturnType res(*this);
0461     res.linearExt() *= b;
0462     return res;
0463   }
0464 
0465   /** \returns The product expression of a diagonal matrix \a a times a transform \a b
0466     *
0467     * The lhs diagonal matrix is interpreted as an affine scaling transformation. The
0468     * product results in a Transform of the same type (mode) as the lhs only if the lhs
0469     * mode is no isometry. In that case, the returned transform is an affinity.
0470     */
0471   template<typename DiagonalDerived>
0472   EIGEN_DEVICE_FUNC friend inline TransformTimeDiagonalReturnType
0473     operator * (const DiagonalBase<DiagonalDerived> &a, const Transform &b)
0474   {
0475     TransformTimeDiagonalReturnType res;
0476     res.linear().noalias() = a*b.linear();
0477     res.translation().noalias() = a*b.translation();
0478     if (EIGEN_CONST_CONDITIONAL(Mode!=int(AffineCompact)))
0479       res.matrix().row(Dim) = b.matrix().row(Dim);
0480     return res;
0481   }
0482 
0483   template<typename OtherDerived>
0484   EIGEN_DEVICE_FUNC inline Transform& operator*=(const EigenBase<OtherDerived>& other) { return *this = *this * other; }
0485 
0486   /** Concatenates two transformations */
0487   EIGEN_DEVICE_FUNC inline const Transform operator * (const Transform& other) const
0488   {
0489     return internal::transform_transform_product_impl<Transform,Transform>::run(*this,other);
0490   }
0491 
0492   #if EIGEN_COMP_ICC
0493 private:
0494   // this intermediate structure permits to workaround a bug in ICC 11:
0495   //   error: template instantiation resulted in unexpected function type of "Eigen::Transform<double, 3, 32, 0>
0496   //             (const Eigen::Transform<double, 3, 2, 0> &) const"
0497   //  (the meaning of a name may have changed since the template declaration -- the type of the template is:
0498   // "Eigen::internal::transform_transform_product_impl<Eigen::Transform<double, 3, 32, 0>,
0499   //     Eigen::Transform<double, 3, Mode, Options>, <expression>>::ResultType (const Eigen::Transform<double, 3, Mode, Options> &) const")
0500   //
0501   template<int OtherMode,int OtherOptions> struct icc_11_workaround
0502   {
0503     typedef internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> > ProductType;
0504     typedef typename ProductType::ResultType ResultType;
0505   };
0506 
0507 public:
0508   /** Concatenates two different transformations */
0509   template<int OtherMode,int OtherOptions>
0510   inline typename icc_11_workaround<OtherMode,OtherOptions>::ResultType
0511     operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions>& other) const
0512   {
0513     typedef typename icc_11_workaround<OtherMode,OtherOptions>::ProductType ProductType;
0514     return ProductType::run(*this,other);
0515   }
0516   #else
0517   /** Concatenates two different transformations */
0518   template<int OtherMode,int OtherOptions>
0519   EIGEN_DEVICE_FUNC inline typename internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::ResultType
0520     operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions>& other) const
0521   {
0522     return internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::run(*this,other);
0523   }
0524   #endif
0525 
0526   /** \sa MatrixBase::setIdentity() */
0527   EIGEN_DEVICE_FUNC void setIdentity() { m_matrix.setIdentity(); }
0528 
0529   /**
0530    * \brief Returns an identity transformation.
0531    * \todo In the future this function should be returning a Transform expression.
0532    */
0533   EIGEN_DEVICE_FUNC static const Transform Identity()
0534   {
0535     return Transform(MatrixType::Identity());
0536   }
0537 
0538   template<typename OtherDerived>
0539   EIGEN_DEVICE_FUNC
0540   inline Transform& scale(const MatrixBase<OtherDerived> &other);
0541 
0542   template<typename OtherDerived>
0543   EIGEN_DEVICE_FUNC
0544   inline Transform& prescale(const MatrixBase<OtherDerived> &other);
0545 
0546   EIGEN_DEVICE_FUNC inline Transform& scale(const Scalar& s);
0547   EIGEN_DEVICE_FUNC inline Transform& prescale(const Scalar& s);
0548 
0549   template<typename OtherDerived>
0550   EIGEN_DEVICE_FUNC
0551   inline Transform& translate(const MatrixBase<OtherDerived> &other);
0552 
0553   template<typename OtherDerived>
0554   EIGEN_DEVICE_FUNC
0555   inline Transform& pretranslate(const MatrixBase<OtherDerived> &other);
0556 
0557   template<typename RotationType>
0558   EIGEN_DEVICE_FUNC
0559   inline Transform& rotate(const RotationType& rotation);
0560 
0561   template<typename RotationType>
0562   EIGEN_DEVICE_FUNC
0563   inline Transform& prerotate(const RotationType& rotation);
0564 
0565   EIGEN_DEVICE_FUNC Transform& shear(const Scalar& sx, const Scalar& sy);
0566   EIGEN_DEVICE_FUNC Transform& preshear(const Scalar& sx, const Scalar& sy);
0567 
0568   EIGEN_DEVICE_FUNC inline Transform& operator=(const TranslationType& t);
0569 
0570   EIGEN_DEVICE_FUNC
0571   inline Transform& operator*=(const TranslationType& t) { return translate(t.vector()); }
0572 
0573   EIGEN_DEVICE_FUNC inline Transform operator*(const TranslationType& t) const;
0574 
0575   EIGEN_DEVICE_FUNC
0576   inline Transform& operator=(const UniformScaling<Scalar>& t);
0577 
0578   EIGEN_DEVICE_FUNC
0579   inline Transform& operator*=(const UniformScaling<Scalar>& s) { return scale(s.factor()); }
0580 
0581   EIGEN_DEVICE_FUNC
0582   inline TransformTimeDiagonalReturnType operator*(const UniformScaling<Scalar>& s) const
0583   {
0584     TransformTimeDiagonalReturnType res = *this;
0585     res.scale(s.factor());
0586     return res;
0587   }
0588 
0589   EIGEN_DEVICE_FUNC
0590   inline Transform& operator*=(const DiagonalMatrix<Scalar,Dim>& s) { linearExt() *= s; return *this; }
0591 
0592   template<typename Derived>
0593   EIGEN_DEVICE_FUNC inline Transform& operator=(const RotationBase<Derived,Dim>& r);
0594   template<typename Derived>
0595   EIGEN_DEVICE_FUNC inline Transform& operator*=(const RotationBase<Derived,Dim>& r) { return rotate(r.toRotationMatrix()); }
0596   template<typename Derived>
0597   EIGEN_DEVICE_FUNC inline Transform operator*(const RotationBase<Derived,Dim>& r) const;
0598 
0599   typedef typename internal::conditional<int(Mode)==Isometry,ConstLinearPart,const LinearMatrixType>::type RotationReturnType;
0600   EIGEN_DEVICE_FUNC RotationReturnType rotation() const;
0601 
0602   template<typename RotationMatrixType, typename ScalingMatrixType>
0603   EIGEN_DEVICE_FUNC
0604   void computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const;
0605   template<typename ScalingMatrixType, typename RotationMatrixType>
0606   EIGEN_DEVICE_FUNC
0607   void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const;
0608 
0609   template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
0610   EIGEN_DEVICE_FUNC
0611   Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived> &position,
0612     const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale);
0613 
0614   EIGEN_DEVICE_FUNC
0615   inline Transform inverse(TransformTraits traits = (TransformTraits)Mode) const;
0616 
0617   /** \returns a const pointer to the column major internal matrix */
0618   EIGEN_DEVICE_FUNC const Scalar* data() const { return m_matrix.data(); }
0619   /** \returns a non-const pointer to the column major internal matrix */
0620   EIGEN_DEVICE_FUNC Scalar* data() { return m_matrix.data(); }
0621 
0622   /** \returns \c *this with scalar type casted to \a NewScalarType
0623     *
0624     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
0625     * then this function smartly returns a const reference to \c *this.
0626     */
0627   template<typename NewScalarType>
0628   EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type cast() const
0629   { return typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type(*this); }
0630 
0631   /** Copy constructor with scalar type conversion */
0632   template<typename OtherScalarType>
0633   EIGEN_DEVICE_FUNC inline explicit Transform(const Transform<OtherScalarType,Dim,Mode,Options>& other)
0634   {
0635     check_template_params();
0636     m_matrix = other.matrix().template cast<Scalar>();
0637   }
0638 
0639   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
0640     * determined by \a prec.
0641     *
0642     * \sa MatrixBase::isApprox() */
0643   EIGEN_DEVICE_FUNC bool isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
0644   { return m_matrix.isApprox(other.m_matrix, prec); }
0645 
0646   /** Sets the last row to [0 ... 0 1]
0647     */
0648   EIGEN_DEVICE_FUNC void makeAffine()
0649   {
0650     internal::transform_make_affine<int(Mode)>::run(m_matrix);
0651   }
0652 
0653   /** \internal
0654     * \returns the Dim x Dim linear part if the transformation is affine,
0655     *          and the HDim x Dim part for projective transformations.
0656     */
0657   EIGEN_DEVICE_FUNC inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt()
0658   { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
0659   /** \internal
0660     * \returns the Dim x Dim linear part if the transformation is affine,
0661     *          and the HDim x Dim part for projective transformations.
0662     */
0663   EIGEN_DEVICE_FUNC inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt() const
0664   { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
0665 
0666   /** \internal
0667     * \returns the translation part if the transformation is affine,
0668     *          and the last column for projective transformations.
0669     */
0670   EIGEN_DEVICE_FUNC inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt()
0671   { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }
0672   /** \internal
0673     * \returns the translation part if the transformation is affine,
0674     *          and the last column for projective transformations.
0675     */
0676   EIGEN_DEVICE_FUNC inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt() const
0677   { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }
0678 
0679 
0680   #ifdef EIGEN_TRANSFORM_PLUGIN
0681   #include EIGEN_TRANSFORM_PLUGIN
0682   #endif
0683 
0684 protected:
0685   #ifndef EIGEN_PARSED_BY_DOXYGEN
0686     EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE void check_template_params()
0687     {
0688       EIGEN_STATIC_ASSERT((Options & (DontAlign|RowMajor)) == Options, INVALID_MATRIX_TEMPLATE_PARAMETERS)
0689     }
0690   #endif
0691 
0692 };
0693 
0694 /** \ingroup Geometry_Module */
0695 typedef Transform<float,2,Isometry> Isometry2f;
0696 /** \ingroup Geometry_Module */
0697 typedef Transform<float,3,Isometry> Isometry3f;
0698 /** \ingroup Geometry_Module */
0699 typedef Transform<double,2,Isometry> Isometry2d;
0700 /** \ingroup Geometry_Module */
0701 typedef Transform<double,3,Isometry> Isometry3d;
0702 
0703 /** \ingroup Geometry_Module */
0704 typedef Transform<float,2,Affine> Affine2f;
0705 /** \ingroup Geometry_Module */
0706 typedef Transform<float,3,Affine> Affine3f;
0707 /** \ingroup Geometry_Module */
0708 typedef Transform<double,2,Affine> Affine2d;
0709 /** \ingroup Geometry_Module */
0710 typedef Transform<double,3,Affine> Affine3d;
0711 
0712 /** \ingroup Geometry_Module */
0713 typedef Transform<float,2,AffineCompact> AffineCompact2f;
0714 /** \ingroup Geometry_Module */
0715 typedef Transform<float,3,AffineCompact> AffineCompact3f;
0716 /** \ingroup Geometry_Module */
0717 typedef Transform<double,2,AffineCompact> AffineCompact2d;
0718 /** \ingroup Geometry_Module */
0719 typedef Transform<double,3,AffineCompact> AffineCompact3d;
0720 
0721 /** \ingroup Geometry_Module */
0722 typedef Transform<float,2,Projective> Projective2f;
0723 /** \ingroup Geometry_Module */
0724 typedef Transform<float,3,Projective> Projective3f;
0725 /** \ingroup Geometry_Module */
0726 typedef Transform<double,2,Projective> Projective2d;
0727 /** \ingroup Geometry_Module */
0728 typedef Transform<double,3,Projective> Projective3d;
0729 
0730 /**************************
0731 *** Optional QT support ***
0732 **************************/
0733 
0734 #ifdef EIGEN_QT_SUPPORT
0735 /** Initializes \c *this from a QMatrix assuming the dimension is 2.
0736   *
0737   * This function is available only if the token EIGEN_QT_SUPPORT is defined.
0738   */
0739 template<typename Scalar, int Dim, int Mode,int Options>
0740 Transform<Scalar,Dim,Mode,Options>::Transform(const QMatrix& other)
0741 {
0742   check_template_params();
0743   *this = other;
0744 }
0745 
0746 /** Set \c *this from a QMatrix assuming the dimension is 2.
0747   *
0748   * This function is available only if the token EIGEN_QT_SUPPORT is defined.
0749   */
0750 template<typename Scalar, int Dim, int Mode,int Options>
0751 Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const QMatrix& other)
0752 {
0753   EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
0754   if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
0755     m_matrix << other.m11(), other.m21(), other.dx(),
0756                 other.m12(), other.m22(), other.dy();
0757   else
0758     m_matrix << other.m11(), other.m21(), other.dx(),
0759                 other.m12(), other.m22(), other.dy(),
0760                 0, 0, 1;
0761   return *this;
0762 }
0763 
0764 /** \returns a QMatrix from \c *this assuming the dimension is 2.
0765   *
0766   * \warning this conversion might loss data if \c *this is not affine
0767   *
0768   * This function is available only if the token EIGEN_QT_SUPPORT is defined.
0769   */
0770 template<typename Scalar, int Dim, int Mode, int Options>
0771 QMatrix Transform<Scalar,Dim,Mode,Options>::toQMatrix(void) const
0772 {
0773   check_template_params();
0774   EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
0775   return QMatrix(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
0776                  m_matrix.coeff(0,1), m_matrix.coeff(1,1),
0777                  m_matrix.coeff(0,2), m_matrix.coeff(1,2));
0778 }
0779 
0780 /** Initializes \c *this from a QTransform assuming the dimension is 2.
0781   *
0782   * This function is available only if the token EIGEN_QT_SUPPORT is defined.
0783   */
0784 template<typename Scalar, int Dim, int Mode,int Options>
0785 Transform<Scalar,Dim,Mode,Options>::Transform(const QTransform& other)
0786 {
0787   check_template_params();
0788   *this = other;
0789 }
0790 
0791 /** Set \c *this from a QTransform assuming the dimension is 2.
0792   *
0793   * This function is available only if the token EIGEN_QT_SUPPORT is defined.
0794   */
0795 template<typename Scalar, int Dim, int Mode, int Options>
0796 Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const QTransform& other)
0797 {
0798   check_template_params();
0799   EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
0800   if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
0801     m_matrix << other.m11(), other.m21(), other.dx(),
0802                 other.m12(), other.m22(), other.dy();
0803   else
0804     m_matrix << other.m11(), other.m21(), other.dx(),
0805                 other.m12(), other.m22(), other.dy(),
0806                 other.m13(), other.m23(), other.m33();
0807   return *this;
0808 }
0809 
0810 /** \returns a QTransform from \c *this assuming the dimension is 2.
0811   *
0812   * This function is available only if the token EIGEN_QT_SUPPORT is defined.
0813   */
0814 template<typename Scalar, int Dim, int Mode, int Options>
0815 QTransform Transform<Scalar,Dim,Mode,Options>::toQTransform(void) const
0816 {
0817   EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
0818   if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
0819     return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
0820                       m_matrix.coeff(0,1), m_matrix.coeff(1,1),
0821                       m_matrix.coeff(0,2), m_matrix.coeff(1,2));
0822   else
0823     return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0), m_matrix.coeff(2,0),
0824                       m_matrix.coeff(0,1), m_matrix.coeff(1,1), m_matrix.coeff(2,1),
0825                       m_matrix.coeff(0,2), m_matrix.coeff(1,2), m_matrix.coeff(2,2));
0826 }
0827 #endif
0828 
0829 /*********************
0830 *** Procedural API ***
0831 *********************/
0832 
0833 /** Applies on the right the non uniform scale transformation represented
0834   * by the vector \a other to \c *this and returns a reference to \c *this.
0835   * \sa prescale()
0836   */
0837 template<typename Scalar, int Dim, int Mode, int Options>
0838 template<typename OtherDerived>
0839 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0840 Transform<Scalar,Dim,Mode,Options>::scale(const MatrixBase<OtherDerived> &other)
0841 {
0842   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
0843   EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
0844   linearExt().noalias() = (linearExt() * other.asDiagonal());
0845   return *this;
0846 }
0847 
0848 /** Applies on the right a uniform scale of a factor \a c to \c *this
0849   * and returns a reference to \c *this.
0850   * \sa prescale(Scalar)
0851   */
0852 template<typename Scalar, int Dim, int Mode, int Options>
0853 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::scale(const Scalar& s)
0854 {
0855   EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
0856   linearExt() *= s;
0857   return *this;
0858 }
0859 
0860 /** Applies on the left the non uniform scale transformation represented
0861   * by the vector \a other to \c *this and returns a reference to \c *this.
0862   * \sa scale()
0863   */
0864 template<typename Scalar, int Dim, int Mode, int Options>
0865 template<typename OtherDerived>
0866 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0867 Transform<Scalar,Dim,Mode,Options>::prescale(const MatrixBase<OtherDerived> &other)
0868 {
0869   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
0870   EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
0871   affine().noalias() = (other.asDiagonal() * affine());
0872   return *this;
0873 }
0874 
0875 /** Applies on the left a uniform scale of a factor \a c to \c *this
0876   * and returns a reference to \c *this.
0877   * \sa scale(Scalar)
0878   */
0879 template<typename Scalar, int Dim, int Mode, int Options>
0880 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::prescale(const Scalar& s)
0881 {
0882   EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
0883   m_matrix.template topRows<Dim>() *= s;
0884   return *this;
0885 }
0886 
0887 /** Applies on the right the translation matrix represented by the vector \a other
0888   * to \c *this and returns a reference to \c *this.
0889   * \sa pretranslate()
0890   */
0891 template<typename Scalar, int Dim, int Mode, int Options>
0892 template<typename OtherDerived>
0893 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0894 Transform<Scalar,Dim,Mode,Options>::translate(const MatrixBase<OtherDerived> &other)
0895 {
0896   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
0897   translationExt() += linearExt() * other;
0898   return *this;
0899 }
0900 
0901 /** Applies on the left the translation matrix represented by the vector \a other
0902   * to \c *this and returns a reference to \c *this.
0903   * \sa translate()
0904   */
0905 template<typename Scalar, int Dim, int Mode, int Options>
0906 template<typename OtherDerived>
0907 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0908 Transform<Scalar,Dim,Mode,Options>::pretranslate(const MatrixBase<OtherDerived> &other)
0909 {
0910   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
0911   if(EIGEN_CONST_CONDITIONAL(int(Mode)==int(Projective)))
0912     affine() += other * m_matrix.row(Dim);
0913   else
0914     translation() += other;
0915   return *this;
0916 }
0917 
0918 /** Applies on the right the rotation represented by the rotation \a rotation
0919   * to \c *this and returns a reference to \c *this.
0920   *
0921   * The template parameter \a RotationType is the type of the rotation which
0922   * must be known by internal::toRotationMatrix<>.
0923   *
0924   * Natively supported types includes:
0925   *   - any scalar (2D),
0926   *   - a Dim x Dim matrix expression,
0927   *   - a Quaternion (3D),
0928   *   - a AngleAxis (3D)
0929   *
0930   * This mechanism is easily extendable to support user types such as Euler angles,
0931   * or a pair of Quaternion for 4D rotations.
0932   *
0933   * \sa rotate(Scalar), class Quaternion, class AngleAxis, prerotate(RotationType)
0934   */
0935 template<typename Scalar, int Dim, int Mode, int Options>
0936 template<typename RotationType>
0937 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0938 Transform<Scalar,Dim,Mode,Options>::rotate(const RotationType& rotation)
0939 {
0940   linearExt() *= internal::toRotationMatrix<Scalar,Dim>(rotation);
0941   return *this;
0942 }
0943 
0944 /** Applies on the left the rotation represented by the rotation \a rotation
0945   * to \c *this and returns a reference to \c *this.
0946   *
0947   * See rotate() for further details.
0948   *
0949   * \sa rotate()
0950   */
0951 template<typename Scalar, int Dim, int Mode, int Options>
0952 template<typename RotationType>
0953 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0954 Transform<Scalar,Dim,Mode,Options>::prerotate(const RotationType& rotation)
0955 {
0956   m_matrix.template block<Dim,HDim>(0,0) = internal::toRotationMatrix<Scalar,Dim>(rotation)
0957                                          * m_matrix.template block<Dim,HDim>(0,0);
0958   return *this;
0959 }
0960 
0961 /** Applies on the right the shear transformation represented
0962   * by the vector \a other to \c *this and returns a reference to \c *this.
0963   * \warning 2D only.
0964   * \sa preshear()
0965   */
0966 template<typename Scalar, int Dim, int Mode, int Options>
0967 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0968 Transform<Scalar,Dim,Mode,Options>::shear(const Scalar& sx, const Scalar& sy)
0969 {
0970   EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
0971   EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
0972   VectorType tmp = linear().col(0)*sy + linear().col(1);
0973   linear() << linear().col(0) + linear().col(1)*sx, tmp;
0974   return *this;
0975 }
0976 
0977 /** Applies on the left the shear transformation represented
0978   * by the vector \a other to \c *this and returns a reference to \c *this.
0979   * \warning 2D only.
0980   * \sa shear()
0981   */
0982 template<typename Scalar, int Dim, int Mode, int Options>
0983 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
0984 Transform<Scalar,Dim,Mode,Options>::preshear(const Scalar& sx, const Scalar& sy)
0985 {
0986   EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
0987   EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
0988   m_matrix.template block<Dim,HDim>(0,0) = LinearMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0);
0989   return *this;
0990 }
0991 
0992 /******************************************************
0993 *** Scaling, Translation and Rotation compatibility ***
0994 ******************************************************/
0995 
0996 template<typename Scalar, int Dim, int Mode, int Options>
0997 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const TranslationType& t)
0998 {
0999   linear().setIdentity();
1000   translation() = t.vector();
1001   makeAffine();
1002   return *this;
1003 }
1004 
1005 template<typename Scalar, int Dim, int Mode, int Options>
1006 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options> Transform<Scalar,Dim,Mode,Options>::operator*(const TranslationType& t) const
1007 {
1008   Transform res = *this;
1009   res.translate(t.vector());
1010   return res;
1011 }
1012 
1013 template<typename Scalar, int Dim, int Mode, int Options>
1014 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const UniformScaling<Scalar>& s)
1015 {
1016   m_matrix.setZero();
1017   linear().diagonal().fill(s.factor());
1018   makeAffine();
1019   return *this;
1020 }
1021 
1022 template<typename Scalar, int Dim, int Mode, int Options>
1023 template<typename Derived>
1024 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const RotationBase<Derived,Dim>& r)
1025 {
1026   linear() = internal::toRotationMatrix<Scalar,Dim>(r);
1027   translation().setZero();
1028   makeAffine();
1029   return *this;
1030 }
1031 
1032 template<typename Scalar, int Dim, int Mode, int Options>
1033 template<typename Derived>
1034 EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options> Transform<Scalar,Dim,Mode,Options>::operator*(const RotationBase<Derived,Dim>& r) const
1035 {
1036   Transform res = *this;
1037   res.rotate(r.derived());
1038   return res;
1039 }
1040 
1041 /************************
1042 *** Special functions ***
1043 ************************/
1044 
1045 namespace internal {
1046 template<int Mode> struct transform_rotation_impl {
1047   template<typename TransformType>
1048   EIGEN_DEVICE_FUNC static inline
1049   const typename TransformType::LinearMatrixType run(const TransformType& t)
1050   {
1051     typedef typename TransformType::LinearMatrixType LinearMatrixType;
1052     LinearMatrixType result;
1053     t.computeRotationScaling(&result, (LinearMatrixType*)0);
1054     return result;
1055   }
1056 };
1057 template<> struct transform_rotation_impl<Isometry> {
1058   template<typename TransformType>
1059   EIGEN_DEVICE_FUNC static inline
1060   typename TransformType::ConstLinearPart run(const TransformType& t)
1061   {
1062     return t.linear();
1063   }
1064 };
1065 }
1066 /** \returns the rotation part of the transformation
1067   *
1068   * If Mode==Isometry, then this method is an alias for linear(),
1069   * otherwise it calls computeRotationScaling() to extract the rotation
1070   * through a SVD decomposition.
1071   *
1072   * \svd_module
1073   *
1074   * \sa computeRotationScaling(), computeScalingRotation(), class SVD
1075   */
1076 template<typename Scalar, int Dim, int Mode, int Options>
1077 EIGEN_DEVICE_FUNC
1078 typename Transform<Scalar,Dim,Mode,Options>::RotationReturnType
1079 Transform<Scalar,Dim,Mode,Options>::rotation() const
1080 {
1081   return internal::transform_rotation_impl<Mode>::run(*this);
1082 }
1083 
1084 
1085 /** decomposes the linear part of the transformation as a product rotation x scaling, the scaling being
1086   * not necessarily positive.
1087   *
1088   * If either pointer is zero, the corresponding computation is skipped.
1089   *
1090   *
1091   *
1092   * \svd_module
1093   *
1094   * \sa computeScalingRotation(), rotation(), class SVD
1095   */
1096 template<typename Scalar, int Dim, int Mode, int Options>
1097 template<typename RotationMatrixType, typename ScalingMatrixType>
1098 EIGEN_DEVICE_FUNC void Transform<Scalar,Dim,Mode,Options>::computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const
1099 {
1100   // Note that JacobiSVD is faster than BDCSVD for small matrices.
1101   JacobiSVD<LinearMatrixType> svd(linear(), ComputeFullU | ComputeFullV);
1102 
1103   Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0) ? Scalar(-1) : Scalar(1); // so x has absolute value 1
1104   VectorType sv(svd.singularValues());
1105   sv.coeffRef(Dim-1) *= x;
1106   if(scaling) *scaling = svd.matrixV() * sv.asDiagonal() * svd.matrixV().adjoint();
1107   if(rotation)
1108   {
1109     LinearMatrixType m(svd.matrixU());
1110     m.col(Dim-1) *= x;
1111     *rotation = m * svd.matrixV().adjoint();
1112   }
1113 }
1114 
1115 /** decomposes the linear part of the transformation as a product scaling x rotation, the scaling being
1116   * not necessarily positive.
1117   *
1118   * If either pointer is zero, the corresponding computation is skipped.
1119   *
1120   *
1121   *
1122   * \svd_module
1123   *
1124   * \sa computeRotationScaling(), rotation(), class SVD
1125   */
1126 template<typename Scalar, int Dim, int Mode, int Options>
1127 template<typename ScalingMatrixType, typename RotationMatrixType>
1128 EIGEN_DEVICE_FUNC void Transform<Scalar,Dim,Mode,Options>::computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const
1129 {
1130   // Note that JacobiSVD is faster than BDCSVD for small matrices.
1131   JacobiSVD<LinearMatrixType> svd(linear(), ComputeFullU | ComputeFullV);
1132 
1133   Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0) ? Scalar(-1) : Scalar(1); // so x has absolute value 1
1134   VectorType sv(svd.singularValues());
1135   sv.coeffRef(Dim-1) *= x;
1136   if(scaling) *scaling = svd.matrixU() * sv.asDiagonal() * svd.matrixU().adjoint();
1137   if(rotation)
1138   {
1139     LinearMatrixType m(svd.matrixU());
1140     m.col(Dim-1) *= x;
1141     *rotation = m * svd.matrixV().adjoint();
1142   }
1143 }
1144 
1145 /** Convenient method to set \c *this from a position, orientation and scale
1146   * of a 3D object.
1147   */
1148 template<typename Scalar, int Dim, int Mode, int Options>
1149 template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
1150 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>&
1151 Transform<Scalar,Dim,Mode,Options>::fromPositionOrientationScale(const MatrixBase<PositionDerived> &position,
1152   const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale)
1153 {
1154   linear() = internal::toRotationMatrix<Scalar,Dim>(orientation);
1155   linear() *= scale.asDiagonal();
1156   translation() = position;
1157   makeAffine();
1158   return *this;
1159 }
1160 
1161 namespace internal {
1162 
1163 template<int Mode>
1164 struct transform_make_affine
1165 {
1166   template<typename MatrixType>
1167   EIGEN_DEVICE_FUNC static void run(MatrixType &mat)
1168   {
1169     static const int Dim = MatrixType::ColsAtCompileTime-1;
1170     mat.template block<1,Dim>(Dim,0).setZero();
1171     mat.coeffRef(Dim,Dim) = typename MatrixType::Scalar(1);
1172   }
1173 };
1174 
1175 template<>
1176 struct transform_make_affine<AffineCompact>
1177 {
1178   template<typename MatrixType> EIGEN_DEVICE_FUNC static void run(MatrixType &) { }
1179 };
1180 
1181 // selector needed to avoid taking the inverse of a 3x4 matrix
1182 template<typename TransformType, int Mode=TransformType::Mode>
1183 struct projective_transform_inverse
1184 {
1185   EIGEN_DEVICE_FUNC static inline void run(const TransformType&, TransformType&)
1186   {}
1187 };
1188 
1189 template<typename TransformType>
1190 struct projective_transform_inverse<TransformType, Projective>
1191 {
1192   EIGEN_DEVICE_FUNC static inline void run(const TransformType& m, TransformType& res)
1193   {
1194     res.matrix() = m.matrix().inverse();
1195   }
1196 };
1197 
1198 } // end namespace internal
1199 
1200 
1201 /**
1202   *
1203   * \returns the inverse transformation according to some given knowledge
1204   * on \c *this.
1205   *
1206   * \param hint allows to optimize the inversion process when the transformation
1207   * is known to be not a general transformation (optional). The possible values are:
1208   *  - #Projective if the transformation is not necessarily affine, i.e., if the
1209   *    last row is not guaranteed to be [0 ... 0 1]
1210   *  - #Affine if the last row can be assumed to be [0 ... 0 1]
1211   *  - #Isometry if the transformation is only a concatenations of translations
1212   *    and rotations.
1213   *  The default is the template class parameter \c Mode.
1214   *
1215   * \warning unless \a traits is always set to NoShear or NoScaling, this function
1216   * requires the generic inverse method of MatrixBase defined in the LU module. If
1217   * you forget to include this module, then you will get hard to debug linking errors.
1218   *
1219   * \sa MatrixBase::inverse()
1220   */
1221 template<typename Scalar, int Dim, int Mode, int Options>
1222 EIGEN_DEVICE_FUNC Transform<Scalar,Dim,Mode,Options>
1223 Transform<Scalar,Dim,Mode,Options>::inverse(TransformTraits hint) const
1224 {
1225   Transform res;
1226   if (hint == Projective)
1227   {
1228     internal::projective_transform_inverse<Transform>::run(*this, res);
1229   }
1230   else
1231   {
1232     if (hint == Isometry)
1233     {
1234       res.matrix().template topLeftCorner<Dim,Dim>() = linear().transpose();
1235     }
1236     else if(hint&Affine)
1237     {
1238       res.matrix().template topLeftCorner<Dim,Dim>() = linear().inverse();
1239     }
1240     else
1241     {
1242       eigen_assert(false && "Invalid transform traits in Transform::Inverse");
1243     }
1244     // translation and remaining parts
1245     res.matrix().template topRightCorner<Dim,1>()
1246       = - res.matrix().template topLeftCorner<Dim,Dim>() * translation();
1247     res.makeAffine(); // we do need this, because in the beginning res is uninitialized
1248   }
1249   return res;
1250 }
1251 
1252 namespace internal {
1253 
1254 /*****************************************************
1255 *** Specializations of take affine part            ***
1256 *****************************************************/
1257 
1258 template<typename TransformType> struct transform_take_affine_part {
1259   typedef typename TransformType::MatrixType MatrixType;
1260   typedef typename TransformType::AffinePart AffinePart;
1261   typedef typename TransformType::ConstAffinePart ConstAffinePart;
1262   static inline AffinePart run(MatrixType& m)
1263   { return m.template block<TransformType::Dim,TransformType::HDim>(0,0); }
1264   static inline ConstAffinePart run(const MatrixType& m)
1265   { return m.template block<TransformType::Dim,TransformType::HDim>(0,0); }
1266 };
1267 
1268 template<typename Scalar, int Dim, int Options>
1269 struct transform_take_affine_part<Transform<Scalar,Dim,AffineCompact, Options> > {
1270   typedef typename Transform<Scalar,Dim,AffineCompact,Options>::MatrixType MatrixType;
1271   static inline MatrixType& run(MatrixType& m) { return m; }
1272   static inline const MatrixType& run(const MatrixType& m) { return m; }
1273 };
1274 
1275 /*****************************************************
1276 *** Specializations of construct from matrix       ***
1277 *****************************************************/
1278 
1279 template<typename Other, int Mode, int Options, int Dim, int HDim>
1280 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, Dim,Dim>
1281 {
1282   static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other)
1283   {
1284     transform->linear() = other;
1285     transform->translation().setZero();
1286     transform->makeAffine();
1287   }
1288 };
1289 
1290 template<typename Other, int Mode, int Options, int Dim, int HDim>
1291 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, Dim,HDim>
1292 {
1293   static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other)
1294   {
1295     transform->affine() = other;
1296     transform->makeAffine();
1297   }
1298 };
1299 
1300 template<typename Other, int Mode, int Options, int Dim, int HDim>
1301 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, HDim,HDim>
1302 {
1303   static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other)
1304   { transform->matrix() = other; }
1305 };
1306 
1307 template<typename Other, int Options, int Dim, int HDim>
1308 struct transform_construct_from_matrix<Other, AffineCompact,Options,Dim,HDim, HDim,HDim>
1309 {
1310   static inline void run(Transform<typename Other::Scalar,Dim,AffineCompact,Options> *transform, const Other& other)
1311   { transform->matrix() = other.template block<Dim,HDim>(0,0); }
1312 };
1313 
1314 /**********************************************************
1315 ***   Specializations of operator* with rhs EigenBase   ***
1316 **********************************************************/
1317 
1318 template<int LhsMode,int RhsMode>
1319 struct transform_product_result
1320 {
1321   enum
1322   {
1323     Mode =
1324       (LhsMode == (int)Projective    || RhsMode == (int)Projective    ) ? Projective :
1325       (LhsMode == (int)Affine        || RhsMode == (int)Affine        ) ? Affine :
1326       (LhsMode == (int)AffineCompact || RhsMode == (int)AffineCompact ) ? AffineCompact :
1327       (LhsMode == (int)Isometry      || RhsMode == (int)Isometry      ) ? Isometry : Projective
1328   };
1329 };
1330 
1331 template< typename TransformType, typename MatrixType, int RhsCols>
1332 struct transform_right_product_impl< TransformType, MatrixType, 0, RhsCols>
1333 {
1334   typedef typename MatrixType::PlainObject ResultType;
1335 
1336   static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1337   {
1338     return T.matrix() * other;
1339   }
1340 };
1341 
1342 template< typename TransformType, typename MatrixType, int RhsCols>
1343 struct transform_right_product_impl< TransformType, MatrixType, 1, RhsCols>
1344 {
1345   enum {
1346     Dim = TransformType::Dim,
1347     HDim = TransformType::HDim,
1348     OtherRows = MatrixType::RowsAtCompileTime,
1349     OtherCols = MatrixType::ColsAtCompileTime
1350   };
1351 
1352   typedef typename MatrixType::PlainObject ResultType;
1353 
1354   static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1355   {
1356     EIGEN_STATIC_ASSERT(OtherRows==HDim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1357 
1358     typedef Block<ResultType, Dim, OtherCols, int(MatrixType::RowsAtCompileTime)==Dim> TopLeftLhs;
1359 
1360     ResultType res(other.rows(),other.cols());
1361     TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() = T.affine() * other;
1362     res.row(OtherRows-1) = other.row(OtherRows-1);
1363 
1364     return res;
1365   }
1366 };
1367 
1368 template< typename TransformType, typename MatrixType, int RhsCols>
1369 struct transform_right_product_impl< TransformType, MatrixType, 2, RhsCols>
1370 {
1371   enum {
1372     Dim = TransformType::Dim,
1373     HDim = TransformType::HDim,
1374     OtherRows = MatrixType::RowsAtCompileTime,
1375     OtherCols = MatrixType::ColsAtCompileTime
1376   };
1377 
1378   typedef typename MatrixType::PlainObject ResultType;
1379 
1380   static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1381   {
1382     EIGEN_STATIC_ASSERT(OtherRows==Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1383 
1384     typedef Block<ResultType, Dim, OtherCols, true> TopLeftLhs;
1385     ResultType res(Replicate<typename TransformType::ConstTranslationPart, 1, OtherCols>(T.translation(),1,other.cols()));
1386     TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() += T.linear() * other;
1387 
1388     return res;
1389   }
1390 };
1391 
1392 template< typename TransformType, typename MatrixType >
1393 struct transform_right_product_impl< TransformType, MatrixType, 2, 1> // rhs is a vector of size Dim
1394 {
1395   typedef typename TransformType::MatrixType TransformMatrix;
1396   enum {
1397     Dim = TransformType::Dim,
1398     HDim = TransformType::HDim,
1399     OtherRows = MatrixType::RowsAtCompileTime,
1400     WorkingRows = EIGEN_PLAIN_ENUM_MIN(TransformMatrix::RowsAtCompileTime,HDim)
1401   };
1402 
1403   typedef typename MatrixType::PlainObject ResultType;
1404 
1405   static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1406   {
1407     EIGEN_STATIC_ASSERT(OtherRows==Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1408 
1409     Matrix<typename ResultType::Scalar, Dim+1, 1> rhs;
1410     rhs.template head<Dim>() = other; rhs[Dim] = typename ResultType::Scalar(1);
1411     Matrix<typename ResultType::Scalar, WorkingRows, 1> res(T.matrix() * rhs);
1412     return res.template head<Dim>();
1413   }
1414 };
1415 
1416 /**********************************************************
1417 ***   Specializations of operator* with lhs EigenBase   ***
1418 **********************************************************/
1419 
1420 // generic HDim x HDim matrix * T => Projective
1421 template<typename Other,int Mode, int Options, int Dim, int HDim>
1422 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, HDim,HDim>
1423 {
1424   typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType;
1425   typedef typename TransformType::MatrixType MatrixType;
1426   typedef Transform<typename Other::Scalar,Dim,Projective,Options> ResultType;
1427   static ResultType run(const Other& other,const TransformType& tr)
1428   { return ResultType(other * tr.matrix()); }
1429 };
1430 
1431 // generic HDim x HDim matrix * AffineCompact => Projective
1432 template<typename Other, int Options, int Dim, int HDim>
1433 struct transform_left_product_impl<Other,AffineCompact,Options,Dim,HDim, HDim,HDim>
1434 {
1435   typedef Transform<typename Other::Scalar,Dim,AffineCompact,Options> TransformType;
1436   typedef typename TransformType::MatrixType MatrixType;
1437   typedef Transform<typename Other::Scalar,Dim,Projective,Options> ResultType;
1438   static ResultType run(const Other& other,const TransformType& tr)
1439   {
1440     ResultType res;
1441     res.matrix().noalias() = other.template block<HDim,Dim>(0,0) * tr.matrix();
1442     res.matrix().col(Dim) += other.col(Dim);
1443     return res;
1444   }
1445 };
1446 
1447 // affine matrix * T
1448 template<typename Other,int Mode, int Options, int Dim, int HDim>
1449 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, Dim,HDim>
1450 {
1451   typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType;
1452   typedef typename TransformType::MatrixType MatrixType;
1453   typedef TransformType ResultType;
1454   static ResultType run(const Other& other,const TransformType& tr)
1455   {
1456     ResultType res;
1457     res.affine().noalias() = other * tr.matrix();
1458     res.matrix().row(Dim) = tr.matrix().row(Dim);
1459     return res;
1460   }
1461 };
1462 
1463 // affine matrix * AffineCompact
1464 template<typename Other, int Options, int Dim, int HDim>
1465 struct transform_left_product_impl<Other,AffineCompact,Options,Dim,HDim, Dim,HDim>
1466 {
1467   typedef Transform<typename Other::Scalar,Dim,AffineCompact,Options> TransformType;
1468   typedef typename TransformType::MatrixType MatrixType;
1469   typedef TransformType ResultType;
1470   static ResultType run(const Other& other,const TransformType& tr)
1471   {
1472     ResultType res;
1473     res.matrix().noalias() = other.template block<Dim,Dim>(0,0) * tr.matrix();
1474     res.translation() += other.col(Dim);
1475     return res;
1476   }
1477 };
1478 
1479 // linear matrix * T
1480 template<typename Other,int Mode, int Options, int Dim, int HDim>
1481 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, Dim,Dim>
1482 {
1483   typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType;
1484   typedef typename TransformType::MatrixType MatrixType;
1485   typedef TransformType ResultType;
1486   static ResultType run(const Other& other, const TransformType& tr)
1487   {
1488     TransformType res;
1489     if(Mode!=int(AffineCompact))
1490       res.matrix().row(Dim) = tr.matrix().row(Dim);
1491     res.matrix().template topRows<Dim>().noalias()
1492       = other * tr.matrix().template topRows<Dim>();
1493     return res;
1494   }
1495 };
1496 
1497 /**********************************************************
1498 *** Specializations of operator* with another Transform ***
1499 **********************************************************/
1500 
1501 template<typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1502 struct transform_transform_product_impl<Transform<Scalar,Dim,LhsMode,LhsOptions>,Transform<Scalar,Dim,RhsMode,RhsOptions>,false >
1503 {
1504   enum { ResultMode = transform_product_result<LhsMode,RhsMode>::Mode };
1505   typedef Transform<Scalar,Dim,LhsMode,LhsOptions> Lhs;
1506   typedef Transform<Scalar,Dim,RhsMode,RhsOptions> Rhs;
1507   typedef Transform<Scalar,Dim,ResultMode,LhsOptions> ResultType;
1508   static ResultType run(const Lhs& lhs, const Rhs& rhs)
1509   {
1510     ResultType res;
1511     res.linear() = lhs.linear() * rhs.linear();
1512     res.translation() = lhs.linear() * rhs.translation() + lhs.translation();
1513     res.makeAffine();
1514     return res;
1515   }
1516 };
1517 
1518 template<typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1519 struct transform_transform_product_impl<Transform<Scalar,Dim,LhsMode,LhsOptions>,Transform<Scalar,Dim,RhsMode,RhsOptions>,true >
1520 {
1521   typedef Transform<Scalar,Dim,LhsMode,LhsOptions> Lhs;
1522   typedef Transform<Scalar,Dim,RhsMode,RhsOptions> Rhs;
1523   typedef Transform<Scalar,Dim,Projective> ResultType;
1524   static ResultType run(const Lhs& lhs, const Rhs& rhs)
1525   {
1526     return ResultType( lhs.matrix() * rhs.matrix() );
1527   }
1528 };
1529 
1530 template<typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1531 struct transform_transform_product_impl<Transform<Scalar,Dim,AffineCompact,LhsOptions>,Transform<Scalar,Dim,Projective,RhsOptions>,true >
1532 {
1533   typedef Transform<Scalar,Dim,AffineCompact,LhsOptions> Lhs;
1534   typedef Transform<Scalar,Dim,Projective,RhsOptions> Rhs;
1535   typedef Transform<Scalar,Dim,Projective> ResultType;
1536   static ResultType run(const Lhs& lhs, const Rhs& rhs)
1537   {
1538     ResultType res;
1539     res.matrix().template topRows<Dim>() = lhs.matrix() * rhs.matrix();
1540     res.matrix().row(Dim) = rhs.matrix().row(Dim);
1541     return res;
1542   }
1543 };
1544 
1545 template<typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1546 struct transform_transform_product_impl<Transform<Scalar,Dim,Projective,LhsOptions>,Transform<Scalar,Dim,AffineCompact,RhsOptions>,true >
1547 {
1548   typedef Transform<Scalar,Dim,Projective,LhsOptions> Lhs;
1549   typedef Transform<Scalar,Dim,AffineCompact,RhsOptions> Rhs;
1550   typedef Transform<Scalar,Dim,Projective> ResultType;
1551   static ResultType run(const Lhs& lhs, const Rhs& rhs)
1552   {
1553     ResultType res(lhs.matrix().template leftCols<Dim>() * rhs.matrix());
1554     res.matrix().col(Dim) += lhs.matrix().col(Dim);
1555     return res;
1556   }
1557 };
1558 
1559 } // end namespace internal
1560 
1561 } // end namespace Eigen
1562 
1563 #endif // EIGEN_TRANSFORM_H