Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Transform.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <[email protected]>
5 // Copyright (C) 2009 Benoit Jacob <[email protected]>
6 // Copyright (C) 2010 Hauke Heibel <[email protected]>
7 //
8 // This Source Code Form is subject to the terms of the Mozilla
9 // Public License v. 2.0. If a copy of the MPL was not distributed
10 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11 
12 #ifndef EIGEN_TRANSFORM_H
13 #define EIGEN_TRANSFORM_H
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 template<typename Transform>
20 struct transform_traits
21 {
22  enum
23  {
24  Dim = Transform::Dim,
25  HDim = Transform::HDim,
26  Mode = Transform::Mode,
27  IsProjective = (int(Mode)==int(Projective))
28  };
29 };
30 
31 template< typename TransformType,
32  typename MatrixType,
33  int Case = transform_traits<TransformType>::IsProjective ? 0
34  : int(MatrixType::RowsAtCompileTime) == int(transform_traits<TransformType>::HDim) ? 1
35  : 2>
36 struct transform_right_product_impl;
37 
38 template< typename Other,
39  int Mode,
40  int Options,
41  int Dim,
42  int HDim,
43  int OtherRows=Other::RowsAtCompileTime,
44  int OtherCols=Other::ColsAtCompileTime>
45 struct transform_left_product_impl;
46 
47 template< typename Lhs,
48  typename Rhs,
49  bool AnyProjective =
50  transform_traits<Lhs>::IsProjective ||
51  transform_traits<Rhs>::IsProjective>
52 struct transform_transform_product_impl;
53 
54 template< typename Other,
55  int Mode,
56  int Options,
57  int Dim,
58  int HDim,
59  int OtherRows=Other::RowsAtCompileTime,
60  int OtherCols=Other::ColsAtCompileTime>
61 struct transform_construct_from_matrix;
62 
63 template<typename TransformType> struct transform_take_affine_part;
64 
65 template<int Mode> struct transform_make_affine;
66 
67 } // end namespace internal
68 
183 template<typename _Scalar, int _Dim, int _Mode, int _Options>
184 class Transform
185 {
186 public:
188  enum {
189  Mode = _Mode,
190  Options = _Options,
191  Dim = _Dim,
192  HDim = _Dim+1,
193  Rows = int(Mode)==(AffineCompact) ? Dim : HDim
194  };
196  typedef _Scalar Scalar;
197  typedef DenseIndex Index;
201  typedef const MatrixType ConstMatrixType;
209  typedef typename internal::conditional<int(Mode)==int(AffineCompact),
210  MatrixType&,
213  typedef typename internal::conditional<int(Mode)==int(AffineCompact),
214  const MatrixType&,
224 
225  // this intermediate enum is needed to avoid an ICE with gcc 3.4 and 4.0
226  enum { TransformTimeDiagonalMode = ((Mode==int(Isometry))?Affine:int(Mode)) };
229 
230 protected:
231 
232  MatrixType m_matrix;
233 
234 public:
235 
238  inline Transform()
239  {
240  check_template_params();
241  internal::transform_make_affine<(int(Mode)==Affine) ? Affine : AffineCompact>::run(m_matrix);
242  }
243 
244  inline Transform(const Transform& other)
245  {
246  check_template_params();
247  m_matrix = other.m_matrix;
248  }
249 
250  inline explicit Transform(const TranslationType& t)
251  {
252  check_template_params();
253  *this = t;
254  }
255  inline explicit Transform(const UniformScaling<Scalar>& s)
256  {
257  check_template_params();
258  *this = s;
259  }
260  template<typename Derived>
261  inline explicit Transform(const RotationBase<Derived, Dim>& r)
262  {
263  check_template_params();
264  *this = r;
265  }
266 
267  inline Transform& operator=(const Transform& other)
268  { m_matrix = other.m_matrix; return *this; }
269 
270  typedef internal::transform_take_affine_part<Transform> take_affine_part;
271 
273  template<typename OtherDerived>
274  inline explicit Transform(const EigenBase<OtherDerived>& other)
275  {
276  EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
277  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
278 
279  check_template_params();
280  internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
281  }
282 
284  template<typename OtherDerived>
286  {
287  EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
288  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
289 
290  internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
291  return *this;
292  }
293 
294  template<int OtherOptions>
296  {
297  check_template_params();
298  // only the options change, we can directly copy the matrices
299  m_matrix = other.matrix();
300  }
301 
302  template<int OtherMode,int OtherOptions>
303  inline Transform(const Transform<Scalar,Dim,OtherMode,OtherOptions>& other)
304  {
305  check_template_params();
306  // prevent conversions as:
307  // Affine | AffineCompact | Isometry = Projective
308  EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Projective), Mode==int(Projective)),
309  YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
310 
311  // prevent conversions as:
312  // Isometry = Affine | AffineCompact
313  EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Affine)||OtherMode==int(AffineCompact), Mode!=int(Isometry)),
314  YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
315 
316  enum { ModeIsAffineCompact = Mode == int(AffineCompact),
317  OtherModeIsAffineCompact = OtherMode == int(AffineCompact)
318  };
319 
320  if(ModeIsAffineCompact == OtherModeIsAffineCompact)
321  {
322  // We need the block expression because the code is compiled for all
323  // combinations of transformations and will trigger a compile time error
324  // if one tries to assign the matrices directly
325  m_matrix.template block<Dim,Dim+1>(0,0) = other.matrix().template block<Dim,Dim+1>(0,0);
326  makeAffine();
327  }
328  else if(OtherModeIsAffineCompact)
329  {
330  typedef typename Transform<Scalar,Dim,OtherMode,OtherOptions>::MatrixType OtherMatrixType;
331  internal::transform_construct_from_matrix<OtherMatrixType,Mode,Options,Dim,HDim>::run(this, other.matrix());
332  }
333  else
334  {
335  // here we know that Mode == AffineCompact and OtherMode != AffineCompact.
336  // if OtherMode were Projective, the static assert above would already have caught it.
337  // So the only possibility is that OtherMode == Affine
338  linear() = other.linear();
339  translation() = other.translation();
340  }
341  }
342 
343  template<typename OtherDerived>
344  Transform(const ReturnByValue<OtherDerived>& other)
345  {
346  check_template_params();
347  other.evalTo(*this);
348  }
349 
350  template<typename OtherDerived>
351  Transform& operator=(const ReturnByValue<OtherDerived>& other)
352  {
353  other.evalTo(*this);
354  return *this;
355  }
356 
357  #ifdef EIGEN_QT_SUPPORT
358  inline Transform(const QMatrix& other);
359  inline Transform& operator=(const QMatrix& other);
360  inline QMatrix toQMatrix(void) const;
361  inline Transform(const QTransform& other);
362  inline Transform& operator=(const QTransform& other);
363  inline QTransform toQTransform(void) const;
364  #endif
365 
368  inline Scalar operator() (Index row, Index col) const { return m_matrix(row,col); }
371  inline Scalar& operator() (Index row, Index col) { return m_matrix(row,col); }
372 
374  inline const MatrixType& matrix() const { return m_matrix; }
376  inline MatrixType& matrix() { return m_matrix; }
377 
379  inline ConstLinearPart linear() const { return ConstLinearPart(m_matrix,0,0); }
381  inline LinearPart linear() { return LinearPart(m_matrix,0,0); }
382 
384  inline ConstAffinePart affine() const { return take_affine_part::run(m_matrix); }
386  inline AffinePart affine() { return take_affine_part::run(m_matrix); }
387 
389  inline ConstTranslationPart translation() const { return ConstTranslationPart(m_matrix,0,Dim); }
391  inline TranslationPart translation() { return TranslationPart(m_matrix,0,Dim); }
392 
417  // note: this function is defined here because some compilers cannot find the respective declaration
418  template<typename OtherDerived>
419  EIGEN_STRONG_INLINE const typename OtherDerived::PlainObject
421  { return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this,other.derived()); }
422 
430  template<typename OtherDerived> friend
431  inline const typename internal::transform_left_product_impl<OtherDerived,Mode,Options,_Dim,_Dim+1>::ResultType
433  { return internal::transform_left_product_impl<OtherDerived,Mode,Options,Dim,HDim>::run(a.derived(),b); }
434 
441  template<typename DiagonalDerived>
442  inline const TransformTimeDiagonalReturnType
443  operator * (const DiagonalBase<DiagonalDerived> &b) const
444  {
445  TransformTimeDiagonalReturnType res(*this);
446  res.linear() *= b;
447  return res;
448  }
449 
456  template<typename DiagonalDerived>
457  friend inline TransformTimeDiagonalReturnType
458  operator * (const DiagonalBase<DiagonalDerived> &a, const Transform &b)
459  {
460  TransformTimeDiagonalReturnType res;
461  res.linear().noalias() = a*b.linear();
462  res.translation().noalias() = a*b.translation();
463  if (Mode!=int(AffineCompact))
464  res.matrix().row(Dim) = b.matrix().row(Dim);
465  return res;
466  }
467 
468  template<typename OtherDerived>
469  inline Transform& operator*=(const EigenBase<OtherDerived>& other) { return *this = *this * other; }
470 
472  inline const Transform operator * (const Transform& other) const
473  {
474  return internal::transform_transform_product_impl<Transform,Transform>::run(*this,other);
475  }
476 
477  #ifdef __INTEL_COMPILER
478 private:
479  // this intermediate structure permits to workaround a bug in ICC 11:
480  // error: template instantiation resulted in unexpected function type of "Eigen::Transform<double, 3, 32, 0>
481  // (const Eigen::Transform<double, 3, 2, 0> &) const"
482  // (the meaning of a name may have changed since the template declaration -- the type of the template is:
483  // "Eigen::internal::transform_transform_product_impl<Eigen::Transform<double, 3, 32, 0>,
484  // Eigen::Transform<double, 3, Mode, Options>, <expression>>::ResultType (const Eigen::Transform<double, 3, Mode, Options> &) const")
485  //
486  template<int OtherMode,int OtherOptions> struct icc_11_workaround
487  {
488  typedef internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> > ProductType;
489  typedef typename ProductType::ResultType ResultType;
490  };
491 
492 public:
494  template<int OtherMode,int OtherOptions>
495  inline typename icc_11_workaround<OtherMode,OtherOptions>::ResultType
496  operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions>& other) const
497  {
498  typedef typename icc_11_workaround<OtherMode,OtherOptions>::ProductType ProductType;
499  return ProductType::run(*this,other);
500  }
501  #else
502 
503  template<int OtherMode,int OtherOptions>
504  inline typename internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::ResultType
506  {
507  return internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::run(*this,other);
508  }
509  #endif
510 
512  void setIdentity() { m_matrix.setIdentity(); }
513 
518  static const Transform Identity()
519  {
520  return Transform(MatrixType::Identity());
521  }
522 
523  template<typename OtherDerived>
524  inline Transform& scale(const MatrixBase<OtherDerived> &other);
525 
526  template<typename OtherDerived>
527  inline Transform& prescale(const MatrixBase<OtherDerived> &other);
528 
529  inline Transform& scale(const Scalar& s);
530  inline Transform& prescale(const Scalar& s);
531 
532  template<typename OtherDerived>
533  inline Transform& translate(const MatrixBase<OtherDerived> &other);
534 
535  template<typename OtherDerived>
536  inline Transform& pretranslate(const MatrixBase<OtherDerived> &other);
537 
538  template<typename RotationType>
539  inline Transform& rotate(const RotationType& rotation);
540 
541  template<typename RotationType>
542  inline Transform& prerotate(const RotationType& rotation);
543 
544  Transform& shear(const Scalar& sx, const Scalar& sy);
545  Transform& preshear(const Scalar& sx, const Scalar& sy);
546 
547  inline Transform& operator=(const TranslationType& t);
548  inline Transform& operator*=(const TranslationType& t) { return translate(t.vector()); }
549  inline Transform operator*(const TranslationType& t) const;
550 
551  inline Transform& operator=(const UniformScaling<Scalar>& t);
552  inline Transform& operator*=(const UniformScaling<Scalar>& s) { return scale(s.factor()); }
553  inline Transform<Scalar,Dim,(int(Mode)==int(Isometry)?int(Affine):int(Mode))> operator*(const UniformScaling<Scalar>& s) const
554  {
555  Transform<Scalar,Dim,(int(Mode)==int(Isometry)?int(Affine):int(Mode)),Options> res = *this;
556  res.scale(s.factor());
557  return res;
558  }
559 
560  inline Transform& operator*=(const DiagonalMatrix<Scalar,Dim>& s) { linear() *= s; return *this; }
561 
562  template<typename Derived>
563  inline Transform& operator=(const RotationBase<Derived,Dim>& r);
564  template<typename Derived>
565  inline Transform& operator*=(const RotationBase<Derived,Dim>& r) { return rotate(r.toRotationMatrix()); }
566  template<typename Derived>
567  inline Transform operator*(const RotationBase<Derived,Dim>& r) const;
568 
569  const LinearMatrixType rotation() const;
570  template<typename RotationMatrixType, typename ScalingMatrixType>
571  void computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const;
572  template<typename ScalingMatrixType, typename RotationMatrixType>
573  void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const;
574 
575  template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
576  Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived> &position,
577  const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale);
578 
579  inline Transform inverse(TransformTraits traits = (TransformTraits)Mode) const;
580 
582  const Scalar* data() const { return m_matrix.data(); }
584  Scalar* data() { return m_matrix.data(); }
585 
591  template<typename NewScalarType>
592  inline typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type cast() const
593  { return typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type(*this); }
594 
596  template<typename OtherScalarType>
598  {
599  check_template_params();
600  m_matrix = other.matrix().template cast<Scalar>();
601  }
602 
607  bool isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
608  { return m_matrix.isApprox(other.m_matrix, prec); }
609 
612  void makeAffine()
613  {
614  internal::transform_make_affine<int(Mode)>::run(m_matrix);
615  }
616 
622  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
627  inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt() const
628  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
629 
634  inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt()
635  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }
640  inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt() const
641  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }
642 
643 
644  #ifdef EIGEN_TRANSFORM_PLUGIN
645  #include EIGEN_TRANSFORM_PLUGIN
646  #endif
647 
648 protected:
649  #ifndef EIGEN_PARSED_BY_DOXYGEN
650  static EIGEN_STRONG_INLINE void check_template_params()
651  {
652  EIGEN_STATIC_ASSERT((Options & (DontAlign|RowMajor)) == Options, INVALID_MATRIX_TEMPLATE_PARAMETERS)
653  }
654  #endif
655 
656 };
657 
666 
675 
684 
693 
694 /**************************
695 *** Optional QT support ***
696 **************************/
697 
698 #ifdef EIGEN_QT_SUPPORT
699 
703 template<typename Scalar, int Dim, int Mode,int Options>
705 {
706  check_template_params();
707  *this = other;
708 }
709 
714 template<typename Scalar, int Dim, int Mode,int Options>
716 {
717  EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
718  m_matrix << other.m11(), other.m21(), other.dx(),
719  other.m12(), other.m22(), other.dy(),
720  0, 0, 1;
721  return *this;
722 }
723 
730 template<typename Scalar, int Dim, int Mode, int Options>
732 {
733  check_template_params();
734  EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
735  return QMatrix(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
736  m_matrix.coeff(0,1), m_matrix.coeff(1,1),
737  m_matrix.coeff(0,2), m_matrix.coeff(1,2));
738 }
739 
744 template<typename Scalar, int Dim, int Mode,int Options>
746 {
747  check_template_params();
748  *this = other;
749 }
750 
755 template<typename Scalar, int Dim, int Mode, int Options>
757 {
758  check_template_params();
759  EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
760  if (Mode == int(AffineCompact))
761  m_matrix << other.m11(), other.m21(), other.dx(),
762  other.m12(), other.m22(), other.dy();
763  else
764  m_matrix << other.m11(), other.m21(), other.dx(),
765  other.m12(), other.m22(), other.dy(),
766  other.m13(), other.m23(), other.m33();
767  return *this;
768 }
769 
774 template<typename Scalar, int Dim, int Mode, int Options>
776 {
777  EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
778  if (Mode == int(AffineCompact))
779  return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
780  m_matrix.coeff(0,1), m_matrix.coeff(1,1),
781  m_matrix.coeff(0,2), m_matrix.coeff(1,2));
782  else
783  return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0), m_matrix.coeff(2,0),
784  m_matrix.coeff(0,1), m_matrix.coeff(1,1), m_matrix.coeff(2,1),
785  m_matrix.coeff(0,2), m_matrix.coeff(1,2), m_matrix.coeff(2,2));
786 }
787 #endif
788 
789 /*********************
790 *** Procedural API ***
791 *********************/
792 
797 template<typename Scalar, int Dim, int Mode, int Options>
798 template<typename OtherDerived>
801 {
802  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
803  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
804  linearExt().noalias() = (linearExt() * other.asDiagonal());
805  return *this;
806 }
807 
812 template<typename Scalar, int Dim, int Mode, int Options>
814 {
815  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
816  linearExt() *= s;
817  return *this;
818 }
819 
824 template<typename Scalar, int Dim, int Mode, int Options>
825 template<typename OtherDerived>
828 {
829  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
830  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
831  m_matrix.template block<Dim,HDim>(0,0).noalias() = (other.asDiagonal() * m_matrix.template block<Dim,HDim>(0,0));
832  return *this;
833 }
834 
839 template<typename Scalar, int Dim, int Mode, int Options>
841 {
842  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
843  m_matrix.template topRows<Dim>() *= s;
844  return *this;
845 }
846 
851 template<typename Scalar, int Dim, int Mode, int Options>
852 template<typename OtherDerived>
855 {
856  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
857  translationExt() += linearExt() * other;
858  return *this;
859 }
860 
865 template<typename Scalar, int Dim, int Mode, int Options>
866 template<typename OtherDerived>
869 {
870  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
871  if(int(Mode)==int(Projective))
872  affine() += other * m_matrix.row(Dim);
873  else
874  translation() += other;
875  return *this;
876 }
877 
895 template<typename Scalar, int Dim, int Mode, int Options>
896 template<typename RotationType>
898 Transform<Scalar,Dim,Mode,Options>::rotate(const RotationType& rotation)
899 {
900  linearExt() *= internal::toRotationMatrix<Scalar,Dim>(rotation);
901  return *this;
902 }
903 
911 template<typename Scalar, int Dim, int Mode, int Options>
912 template<typename RotationType>
915 {
916  m_matrix.template block<Dim,HDim>(0,0) = internal::toRotationMatrix<Scalar,Dim>(rotation)
917  * m_matrix.template block<Dim,HDim>(0,0);
918  return *this;
919 }
920 
926 template<typename Scalar, int Dim, int Mode, int Options>
928 Transform<Scalar,Dim,Mode,Options>::shear(const Scalar& sx, const Scalar& sy)
929 {
930  EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
931  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
932  VectorType tmp = linear().col(0)*sy + linear().col(1);
933  linear() << linear().col(0) + linear().col(1)*sx, tmp;
934  return *this;
935 }
936 
942 template<typename Scalar, int Dim, int Mode, int Options>
944 Transform<Scalar,Dim,Mode,Options>::preshear(const Scalar& sx, const Scalar& sy)
945 {
946  EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
947  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
948  m_matrix.template block<Dim,HDim>(0,0) = LinearMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0);
949  return *this;
950 }
951 
952 /******************************************************
953 *** Scaling, Translation and Rotation compatibility ***
954 ******************************************************/
955 
956 template<typename Scalar, int Dim, int Mode, int Options>
958 {
959  linear().setIdentity();
960  translation() = t.vector();
961  makeAffine();
962  return *this;
963 }
964 
965 template<typename Scalar, int Dim, int Mode, int Options>
966 inline Transform<Scalar,Dim,Mode,Options> Transform<Scalar,Dim,Mode,Options>::operator*(const TranslationType& t) const
967 {
968  Transform res = *this;
969  res.translate(t.vector());
970  return res;
971 }
972 
973 template<typename Scalar, int Dim, int Mode, int Options>
974 inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const UniformScaling<Scalar>& s)
975 {
976  m_matrix.setZero();
977  linear().diagonal().fill(s.factor());
978  makeAffine();
979  return *this;
980 }
981 
982 template<typename Scalar, int Dim, int Mode, int Options>
983 template<typename Derived>
984 inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const RotationBase<Derived,Dim>& r)
985 {
986  linear() = internal::toRotationMatrix<Scalar,Dim>(r);
987  translation().setZero();
988  makeAffine();
989  return *this;
990 }
991 
992 template<typename Scalar, int Dim, int Mode, int Options>
993 template<typename Derived>
994 inline Transform<Scalar,Dim,Mode,Options> Transform<Scalar,Dim,Mode,Options>::operator*(const RotationBase<Derived,Dim>& r) const
995 {
996  Transform res = *this;
997  res.rotate(r.derived());
998  return res;
999 }
1000 
1001 /************************
1002 *** Special functions ***
1003 ************************/
1004 
1012 template<typename Scalar, int Dim, int Mode, int Options>
1013 const typename Transform<Scalar,Dim,Mode,Options>::LinearMatrixType
1015 {
1016  LinearMatrixType result;
1017  computeRotationScaling(&result, (LinearMatrixType*)0);
1018  return result;
1019 }
1020 
1021 
1033 template<typename Scalar, int Dim, int Mode, int Options>
1034 template<typename RotationMatrixType, typename ScalingMatrixType>
1035 void Transform<Scalar,Dim,Mode,Options>::computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const
1036 {
1038 
1039  Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant(); // so x has absolute value 1
1040  VectorType sv(svd.singularValues());
1041  sv.coeffRef(0) *= x;
1042  if(scaling) scaling->lazyAssign(svd.matrixV() * sv.asDiagonal() * svd.matrixV().adjoint());
1043  if(rotation)
1044  {
1045  LinearMatrixType m(svd.matrixU());
1046  m.col(0) /= x;
1047  rotation->lazyAssign(m * svd.matrixV().adjoint());
1048  }
1049 }
1050 
1062 template<typename Scalar, int Dim, int Mode, int Options>
1063 template<typename ScalingMatrixType, typename RotationMatrixType>
1064 void Transform<Scalar,Dim,Mode,Options>::computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const
1065 {
1067 
1068  Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant(); // so x has absolute value 1
1069  VectorType sv(svd.singularValues());
1070  sv.coeffRef(0) *= x;
1071  if(scaling) scaling->lazyAssign(svd.matrixU() * sv.asDiagonal() * svd.matrixU().adjoint());
1072  if(rotation)
1073  {
1074  LinearMatrixType m(svd.matrixU());
1075  m.col(0) /= x;
1076  rotation->lazyAssign(m * svd.matrixV().adjoint());
1077  }
1078 }
1079 
1083 template<typename Scalar, int Dim, int Mode, int Options>
1084 template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
1087  const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale)
1088 {
1089  linear() = internal::toRotationMatrix<Scalar,Dim>(orientation);
1090  linear() *= scale.asDiagonal();
1091  translation() = position;
1092  makeAffine();
1093  return *this;
1094 }
1095 
1096 namespace internal {
1097 
1098 template<int Mode>
1099 struct transform_make_affine
1100 {
1101  template<typename MatrixType>
1102  static void run(MatrixType &mat)
1103  {
1104  static const int Dim = MatrixType::ColsAtCompileTime-1;
1105  mat.template block<1,Dim>(Dim,0).setZero();
1106  mat.coeffRef(Dim,Dim) = typename MatrixType::Scalar(1);
1107  }
1108 };
1109 
1110 template<>
1111 struct transform_make_affine<AffineCompact>
1112 {
1113  template<typename MatrixType> static void run(MatrixType &) { }
1114 };
1115 
1116 // selector needed to avoid taking the inverse of a 3x4 matrix
1117 template<typename TransformType, int Mode=TransformType::Mode>
1118 struct projective_transform_inverse
1119 {
1120  static inline void run(const TransformType&, TransformType&)
1121  {}
1122 };
1123 
1124 template<typename TransformType>
1125 struct projective_transform_inverse<TransformType, Projective>
1126 {
1127  static inline void run(const TransformType& m, TransformType& res)
1128  {
1129  res.matrix() = m.matrix().inverse();
1130  }
1131 };
1132 
1133 } // end namespace internal
1134 
1135 
1156 template<typename Scalar, int Dim, int Mode, int Options>
1157 Transform<Scalar,Dim,Mode,Options>
1159 {
1160  Transform res;
1161  if (hint == Projective)
1162  {
1163  internal::projective_transform_inverse<Transform>::run(*this, res);
1164  }
1165  else
1166  {
1167  if (hint == Isometry)
1168  {
1169  res.matrix().template topLeftCorner<Dim,Dim>() = linear().transpose();
1170  }
1171  else if(hint&Affine)
1172  {
1173  res.matrix().template topLeftCorner<Dim,Dim>() = linear().inverse();
1174  }
1175  else
1176  {
1177  eigen_assert(false && "Invalid transform traits in Transform::Inverse");
1178  }
1179  // translation and remaining parts
1180  res.matrix().template topRightCorner<Dim,1>()
1181  = - res.matrix().template topLeftCorner<Dim,Dim>() * translation();
1182  res.makeAffine(); // we do need this, because in the beginning res is uninitialized
1183  }
1184  return res;
1185 }
1186 
1187 namespace internal {
1188 
1189 /*****************************************************
1190 *** Specializations of take affine part ***
1191 *****************************************************/
1192 
1193 template<typename TransformType> struct transform_take_affine_part {
1194  typedef typename TransformType::MatrixType MatrixType;
1195  typedef typename TransformType::AffinePart AffinePart;
1196  typedef typename TransformType::ConstAffinePart ConstAffinePart;
1197  static inline AffinePart run(MatrixType& m)
1198  { return m.template block<TransformType::Dim,TransformType::HDim>(0,0); }
1199  static inline ConstAffinePart run(const MatrixType& m)
1200  { return m.template block<TransformType::Dim,TransformType::HDim>(0,0); }
1201 };
1202 
1203 template<typename Scalar, int Dim, int Options>
1204 struct transform_take_affine_part<Transform<Scalar,Dim,AffineCompact, Options> > {
1205  typedef typename Transform<Scalar,Dim,AffineCompact,Options>::MatrixType MatrixType;
1206  static inline MatrixType& run(MatrixType& m) { return m; }
1207  static inline const MatrixType& run(const MatrixType& m) { return m; }
1208 };
1209 
1210 /*****************************************************
1211 *** Specializations of construct from matrix ***
1212 *****************************************************/
1213 
1214 template<typename Other, int Mode, int Options, int Dim, int HDim>
1215 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, Dim,Dim>
1216 {
1217  static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other)
1218  {
1219  transform->linear() = other;
1220  transform->translation().setZero();
1221  transform->makeAffine();
1222  }
1223 };
1224 
1225 template<typename Other, int Mode, int Options, int Dim, int HDim>
1226 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, Dim,HDim>
1227 {
1228  static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other)
1229  {
1230  transform->affine() = other;
1231  transform->makeAffine();
1232  }
1233 };
1234 
1235 template<typename Other, int Mode, int Options, int Dim, int HDim>
1236 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, HDim,HDim>
1237 {
1238  static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other)
1239  { transform->matrix() = other; }
1240 };
1241 
1242 template<typename Other, int Options, int Dim, int HDim>
1243 struct transform_construct_from_matrix<Other, AffineCompact,Options,Dim,HDim, HDim,HDim>
1244 {
1245  static inline void run(Transform<typename Other::Scalar,Dim,AffineCompact,Options> *transform, const Other& other)
1246  { transform->matrix() = other.template block<Dim,HDim>(0,0); }
1247 };
1248 
1249 /**********************************************************
1250 *** Specializations of operator* with rhs EigenBase ***
1251 **********************************************************/
1252 
1253 template<int LhsMode,int RhsMode>
1254 struct transform_product_result
1255 {
1256  enum
1257  {
1258  Mode =
1259  (LhsMode == (int)Projective || RhsMode == (int)Projective ) ? Projective :
1260  (LhsMode == (int)Affine || RhsMode == (int)Affine ) ? Affine :
1261  (LhsMode == (int)AffineCompact || RhsMode == (int)AffineCompact ) ? AffineCompact :
1262  (LhsMode == (int)Isometry || RhsMode == (int)Isometry ) ? Isometry : Projective
1263  };
1264 };
1265 
1266 template< typename TransformType, typename MatrixType >
1267 struct transform_right_product_impl< TransformType, MatrixType, 0 >
1268 {
1269  typedef typename MatrixType::PlainObject ResultType;
1270 
1271  static EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1272  {
1273  return T.matrix() * other;
1274  }
1275 };
1276 
1277 template< typename TransformType, typename MatrixType >
1278 struct transform_right_product_impl< TransformType, MatrixType, 1 >
1279 {
1280  enum {
1281  Dim = TransformType::Dim,
1282  HDim = TransformType::HDim,
1283  OtherRows = MatrixType::RowsAtCompileTime,
1284  OtherCols = MatrixType::ColsAtCompileTime
1285  };
1286 
1287  typedef typename MatrixType::PlainObject ResultType;
1288 
1289  static EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1290  {
1291  EIGEN_STATIC_ASSERT(OtherRows==HDim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1292 
1293  typedef Block<ResultType, Dim, OtherCols, int(MatrixType::RowsAtCompileTime)==Dim> TopLeftLhs;
1294 
1295  ResultType res(other.rows(),other.cols());
1296  TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() = T.affine() * other;
1297  res.row(OtherRows-1) = other.row(OtherRows-1);
1298 
1299  return res;
1300  }
1301 };
1302 
1303 template< typename TransformType, typename MatrixType >
1304 struct transform_right_product_impl< TransformType, MatrixType, 2 >
1305 {
1306  enum {
1307  Dim = TransformType::Dim,
1308  HDim = TransformType::HDim,
1309  OtherRows = MatrixType::RowsAtCompileTime,
1310  OtherCols = MatrixType::ColsAtCompileTime
1311  };
1312 
1313  typedef typename MatrixType::PlainObject ResultType;
1314 
1315  static EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other)
1316  {
1317  EIGEN_STATIC_ASSERT(OtherRows==Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1318 
1319  typedef Block<ResultType, Dim, OtherCols, true> TopLeftLhs;
1320  ResultType res(Replicate<typename TransformType::ConstTranslationPart, 1, OtherCols>(T.translation(),1,other.cols()));
1321  TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() += T.linear() * other;
1322 
1323  return res;
1324  }
1325 };
1326 
1327 /**********************************************************
1328 *** Specializations of operator* with lhs EigenBase ***
1329 **********************************************************/
1330 
1331 // generic HDim x HDim matrix * T => Projective
1332 template<typename Other,int Mode, int Options, int Dim, int HDim>
1333 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, HDim,HDim>
1334 {
1335  typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType;
1336  typedef typename TransformType::MatrixType MatrixType;
1337  typedef Transform<typename Other::Scalar,Dim,Projective,Options> ResultType;
1338  static ResultType run(const Other& other,const TransformType& tr)
1339  { return ResultType(other * tr.matrix()); }
1340 };
1341 
1342 // generic HDim x HDim matrix * AffineCompact => Projective
1343 template<typename Other, int Options, int Dim, int HDim>
1344 struct transform_left_product_impl<Other,AffineCompact,Options,Dim,HDim, HDim,HDim>
1345 {
1346  typedef Transform<typename Other::Scalar,Dim,AffineCompact,Options> TransformType;
1347  typedef typename TransformType::MatrixType MatrixType;
1348  typedef Transform<typename Other::Scalar,Dim,Projective,Options> ResultType;
1349  static ResultType run(const Other& other,const TransformType& tr)
1350  {
1351  ResultType res;
1352  res.matrix().noalias() = other.template block<HDim,Dim>(0,0) * tr.matrix();
1353  res.matrix().col(Dim) += other.col(Dim);
1354  return res;
1355  }
1356 };
1357 
1358 // affine matrix * T
1359 template<typename Other,int Mode, int Options, int Dim, int HDim>
1360 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, Dim,HDim>
1361 {
1362  typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType;
1363  typedef typename TransformType::MatrixType MatrixType;
1364  typedef TransformType ResultType;
1365  static ResultType run(const Other& other,const TransformType& tr)
1366  {
1367  ResultType res;
1368  res.affine().noalias() = other * tr.matrix();
1369  res.matrix().row(Dim) = tr.matrix().row(Dim);
1370  return res;
1371  }
1372 };
1373 
1374 // affine matrix * AffineCompact
1375 template<typename Other, int Options, int Dim, int HDim>
1376 struct transform_left_product_impl<Other,AffineCompact,Options,Dim,HDim, Dim,HDim>
1377 {
1378  typedef Transform<typename Other::Scalar,Dim,AffineCompact,Options> TransformType;
1379  typedef typename TransformType::MatrixType MatrixType;
1380  typedef TransformType ResultType;
1381  static ResultType run(const Other& other,const TransformType& tr)
1382  {
1383  ResultType res;
1384  res.matrix().noalias() = other.template block<Dim,Dim>(0,0) * tr.matrix();
1385  res.translation() += other.col(Dim);
1386  return res;
1387  }
1388 };
1389 
1390 // linear matrix * T
1391 template<typename Other,int Mode, int Options, int Dim, int HDim>
1392 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, Dim,Dim>
1393 {
1394  typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType;
1395  typedef typename TransformType::MatrixType MatrixType;
1396  typedef TransformType ResultType;
1397  static ResultType run(const Other& other, const TransformType& tr)
1398  {
1399  TransformType res;
1400  if(Mode!=int(AffineCompact))
1401  res.matrix().row(Dim) = tr.matrix().row(Dim);
1402  res.matrix().template topRows<Dim>().noalias()
1403  = other * tr.matrix().template topRows<Dim>();
1404  return res;
1405  }
1406 };
1407 
1408 /**********************************************************
1409 *** Specializations of operator* with another Transform ***
1410 **********************************************************/
1411 
1412 template<typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1413 struct transform_transform_product_impl<Transform<Scalar,Dim,LhsMode,LhsOptions>,Transform<Scalar,Dim,RhsMode,RhsOptions>,false >
1414 {
1415  enum { ResultMode = transform_product_result<LhsMode,RhsMode>::Mode };
1416  typedef Transform<Scalar,Dim,LhsMode,LhsOptions> Lhs;
1417  typedef Transform<Scalar,Dim,RhsMode,RhsOptions> Rhs;
1418  typedef Transform<Scalar,Dim,ResultMode,LhsOptions> ResultType;
1419  static ResultType run(const Lhs& lhs, const Rhs& rhs)
1420  {
1421  ResultType res;
1422  res.linear() = lhs.linear() * rhs.linear();
1423  res.translation() = lhs.linear() * rhs.translation() + lhs.translation();
1424  res.makeAffine();
1425  return res;
1426  }
1427 };
1428 
1429 template<typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1430 struct transform_transform_product_impl<Transform<Scalar,Dim,LhsMode,LhsOptions>,Transform<Scalar,Dim,RhsMode,RhsOptions>,true >
1431 {
1432  typedef Transform<Scalar,Dim,LhsMode,LhsOptions> Lhs;
1433  typedef Transform<Scalar,Dim,RhsMode,RhsOptions> Rhs;
1434  typedef Transform<Scalar,Dim,Projective> ResultType;
1435  static ResultType run(const Lhs& lhs, const Rhs& rhs)
1436  {
1437  return ResultType( lhs.matrix() * rhs.matrix() );
1438  }
1439 };
1440 
1441 template<typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1442 struct transform_transform_product_impl<Transform<Scalar,Dim,AffineCompact,LhsOptions>,Transform<Scalar,Dim,Projective,RhsOptions>,true >
1443 {
1444  typedef Transform<Scalar,Dim,AffineCompact,LhsOptions> Lhs;
1445  typedef Transform<Scalar,Dim,Projective,RhsOptions> Rhs;
1446  typedef Transform<Scalar,Dim,Projective> ResultType;
1447  static ResultType run(const Lhs& lhs, const Rhs& rhs)
1448  {
1449  ResultType res;
1450  res.matrix().template topRows<Dim>() = lhs.matrix() * rhs.matrix();
1451  res.matrix().row(Dim) = rhs.matrix().row(Dim);
1452  return res;
1453  }
1454 };
1455 
1456 template<typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1457 struct transform_transform_product_impl<Transform<Scalar,Dim,Projective,LhsOptions>,Transform<Scalar,Dim,AffineCompact,RhsOptions>,true >
1458 {
1459  typedef Transform<Scalar,Dim,Projective,LhsOptions> Lhs;
1460  typedef Transform<Scalar,Dim,AffineCompact,RhsOptions> Rhs;
1461  typedef Transform<Scalar,Dim,Projective> ResultType;
1462  static ResultType run(const Lhs& lhs, const Rhs& rhs)
1463  {
1464  ResultType res(lhs.matrix().template leftCols<Dim>() * rhs.matrix());
1465  res.matrix().col(Dim) += lhs.matrix().col(Dim);
1466  return res;
1467  }
1468 };
1469 
1470 } // end namespace internal
1471 
1472 } // end namespace Eigen
1473 
1474 #endif // EIGEN_TRANSFORM_H
Definition: Constants.h:270
QMatrix toQMatrix(void) const
Definition: Transform.h:731
QTransform toQTransform(void) const
Definition: Transform.h:775
Definition: Constants.h:398
ConstLinearPart linear() const
Definition: Transform.h:379
Block< MatrixType, Dim, Dim, int(Mode)==(AffineCompact)&&(Options &RowMajor)==0 > LinearPart
Definition: Transform.h:205
const MatrixType ConstMatrixType
Definition: Transform.h:201
RowXpr row(Index i)
Definition: DenseBase.h:750
const Scalar * data() const
Definition: PlainObjectBase.h:212
Transform()
Definition: Transform.h:238
Transform & shear(const Scalar &sx, const Scalar &sy)
Definition: Transform.h:928
Matrix< Scalar, Dim, 1 > VectorType
Definition: Transform.h:217
void setIdentity()
Definition: Transform.h:512
bool isApprox(const DenseBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Fuzzy.h:98
Transform< float, 2, Projective > Projective2f
Definition: Transform.h:686
Transform< float, 2, Isometry > Isometry2f
Definition: Transform.h:659
const OtherDerived::PlainObject operator*(const EigenBase< OtherDerived > &other) const
Definition: Transform.h:420
Transform< float, 2, AffineCompact > AffineCompact2f
Definition: Transform.h:677
Transform< Scalar, Dim, TransformTimeDiagonalMode > TransformTimeDiagonalReturnType
Definition: Transform.h:228
const LinearMatrixType rotation() const
Definition: Transform.h:1014
Definition: Constants.h:394
internal::conditional< int(Mode)==int(AffineCompact), const MatrixType &, const Block< const MatrixType, Dim, HDim > >::type ConstAffinePart
Definition: Transform.h:215
Block< MatrixType, Dim, 1, int(Mode)==(AffineCompact)> TranslationPart
Definition: Transform.h:219
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const internal::permut_matrix_product_retval< PermutationDerived, Derived, OnTheRight > operator*(const MatrixBase< Derived > &matrix, const PermutationBase< PermutationDerived > &permutation)
Definition: PermutationMatrix.h:539
const int Dynamic
Definition: Constants.h:21
Scalar * data()
Definition: Transform.h:584
Definition: Constants.h:331
Transform< double, 2, Affine > Affine2d
Definition: Transform.h:672
const SingularValuesType & singularValues() const
Definition: JacobiSVD.h:641
_Scalar Scalar
Definition: Transform.h:194
internal::cast_return_type< Transform, Transform< NewScalarType, Dim, Mode, Options > >::type cast() const
Definition: Transform.h:592
Definition: Constants.h:391
Transform< double, 3, AffineCompact > AffineCompact3d
Definition: Transform.h:683
Definition: EigenBase.h:26
Represents a translation transformation.
Definition: ForwardDeclarations.h:236
bool isApprox(const Transform &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Transform.h:607
Transform< double, 3, Affine > Affine3d
Definition: Transform.h:674
void makeAffine()
Definition: Transform.h:612
TransformTraits
Definition: Constants.h:389
void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const
Definition: Transform.h:1064
Translation< Scalar, Dim > TranslationType
Definition: Transform.h:223
Transform< float, 3, Affine > Affine3f
Definition: Transform.h:670
ConstAffinePart affine() const
Definition: Transform.h:384
Transform< float, 3, Projective > Projective3f
Definition: Transform.h:688
LinearPart linear()
Definition: Transform.h:381
Derived & derived()
Definition: EigenBase.h:34
Definition: Constants.h:396
const MatrixUType & matrixU() const
Definition: JacobiSVD.h:613
static const Transform Identity()
Returns an identity transformation.
Definition: Transform.h:518
const MatrixType & matrix() const
Definition: Transform.h:374
Transform< double, 2, Isometry > Isometry2d
Definition: Transform.h:663
AffinePart affine()
Definition: Transform.h:386
Transform(const EigenBase< OtherDerived > &other)
Definition: Transform.h:274
Transform inverse(TransformTraits traits=(TransformTraits) Mode) const
Definition: Transform.h:1158
internal::make_proper_matrix_type< Scalar, Rows, HDim, Options >::type MatrixType
Definition: Transform.h:199
Transform< double, 3, Projective > Projective3d
Definition: Transform.h:692
Transform< float, 3, AffineCompact > AffineCompact3f
Definition: Transform.h:679
Transform(const Transform< OtherScalarType, Dim, Mode, Options > &other)
Definition: Transform.h:597
const Block< ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact)&&(Options &RowMajor)==0 > ConstLinearPart
Definition: Transform.h:207
Transform< double, 2, AffineCompact > AffineCompact2d
Definition: Transform.h:681
const Scalar * data() const
Definition: Transform.h:582
const MatrixVType & matrixV() const
Definition: JacobiSVD.h:629
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
internal::conditional< int(Mode)==int(AffineCompact), MatrixType &, Block< MatrixType, Dim, HDim > >::type AffinePart
Definition: Transform.h:211
Transform< float, 3, Isometry > Isometry3f
Definition: Transform.h:661
Transform< float, 2, Affine > Affine2f
Definition: Transform.h:668
Transform & operator=(const EigenBase< OtherDerived > &other)
Definition: Transform.h:285
const Block< ConstMatrixType, Dim, 1, int(Mode)==(AffineCompact)> ConstTranslationPart
Definition: Transform.h:221
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:224
Derived & setIdentity()
Definition: CwiseNullaryOp.h:772
MatrixType & matrix()
Definition: Transform.h:376
Definition: Constants.h:266
Transform & preshear(const Scalar &sx, const Scalar &sy)
Definition: Transform.h:944
TranslationPart translation()
Definition: Transform.h:391
Transform< double, 2, Projective > Projective2d
Definition: Transform.h:690
void computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const
Definition: Transform.h:1035
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
Definition: Constants.h:327
const DiagonalWrapper< const Derived > asDiagonal() const
Definition: DiagonalMatrix.h:278
ConstTranslationPart translation() const
Definition: Transform.h:389
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar, _Dim==Dynamic?Dynamic:(_Dim+1)*(_Dim+1)) enum
Definition: Transform.h:187
Represents an homogeneous transformation in a N dimensional space.
Definition: ForwardDeclarations.h:264
Transform< double, 3, Isometry > Isometry3d
Definition: Transform.h:665
Matrix< Scalar, Dim, Dim, Options > LinearMatrixType
Definition: Transform.h:203
Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:515