Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
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// IWYU pragma: private
16#include "./InternalHeaderCheck.h"
17
18namespace Eigen {
19
20namespace internal {
21
22template <typename Transform>
23struct transform_traits {
24 enum {
25 Dim = Transform::Dim,
26 HDim = Transform::HDim,
27 Mode = Transform::Mode,
28 IsProjective = (int(Mode) == int(Projective))
29 };
30};
31
32template <typename TransformType, typename MatrixType,
33 int Case = transform_traits<TransformType>::IsProjective ? 0
34 : int(MatrixType::RowsAtCompileTime) == int(transform_traits<TransformType>::HDim) ? 1
35 : 2,
36 int RhsCols = MatrixType::ColsAtCompileTime>
37struct transform_right_product_impl;
38
39template <typename Other, int Mode, int Options, int Dim, int HDim, int OtherRows = Other::RowsAtCompileTime,
40 int OtherCols = Other::ColsAtCompileTime>
41struct transform_left_product_impl;
42
43template <typename Lhs, typename Rhs,
44 bool AnyProjective = transform_traits<Lhs>::IsProjective || transform_traits<Rhs>::IsProjective>
45struct transform_transform_product_impl;
46
47template <typename Other, int Mode, int Options, int Dim, int HDim, int OtherRows = Other::RowsAtCompileTime,
48 int OtherCols = Other::ColsAtCompileTime>
49struct transform_construct_from_matrix;
50
51template <typename TransformType>
52struct transform_take_affine_part;
53
54template <typename Scalar_, int Dim_, int Mode_, int Options_>
55struct traits<Transform<Scalar_, Dim_, Mode_, Options_> > {
56 typedef Scalar_ Scalar;
57 typedef Eigen::Index StorageIndex;
58 typedef Dense StorageKind;
59 enum {
60 Dim1 = Dim_ == Dynamic ? Dim_ : Dim_ + 1,
61 RowsAtCompileTime = Mode_ == Projective ? Dim1 : Dim_,
62 ColsAtCompileTime = Dim1,
63 MaxRowsAtCompileTime = RowsAtCompileTime,
64 MaxColsAtCompileTime = ColsAtCompileTime,
65 Flags = 0
66 };
67};
68
69template <int Mode>
70struct transform_make_affine;
71
72} // end namespace internal
73
191template <typename Scalar_, int Dim_, int Mode_, int Options_>
193 public:
195 Dim_ == Dynamic ? Dynamic : (Dim_ + 1) * (Dim_ + 1))
196 enum {
197 Mode = Mode_,
198 Options = Options_,
199 Dim = Dim_,
200 HDim = Dim_ + 1,
201 Rows = int(Mode) == (AffineCompact) ? Dim : HDim
202 };
204 typedef Scalar_ Scalar;
205 typedef Eigen::Index StorageIndex;
214 typedef Block<MatrixType, Dim, Dim, int(Mode) == (AffineCompact) && (int(Options) & RowMajor) == 0> LinearPart;
216 typedef const Block<ConstMatrixType, Dim, Dim, int(Mode) == (AffineCompact) && (int(Options) & RowMajor) == 0>
219 typedef std::conditional_t<int(Mode) == int(AffineCompact), MatrixType&, Block<MatrixType, Dim, HDim> > AffinePart;
221 typedef std::conditional_t<int(Mode) == int(AffineCompact), const MatrixType&,
227 typedef Block<MatrixType, Dim, 1, !(internal::traits<MatrixType>::Flags & RowMajorBit)> TranslationPart;
229 typedef const Block<ConstMatrixType, Dim, 1, !(internal::traits<MatrixType>::Flags & RowMajorBit)>
233
234 // this intermediate enum is needed to avoid an ICE with gcc 3.4 and 4.0
235 enum { TransformTimeDiagonalMode = ((Mode == int(Isometry)) ? Affine : int(Mode)) };
238
239 protected:
240 MatrixType m_matrix;
241
242 public:
245 EIGEN_DEVICE_FUNC inline Transform() {
246 check_template_params();
247 internal::transform_make_affine<(int(Mode) == Affine || int(Mode) == Isometry) ? Affine : AffineCompact>::run(
248 m_matrix);
249 }
250
251 EIGEN_DEVICE_FUNC inline explicit Transform(const TranslationType& t) {
252 check_template_params();
253 *this = t;
254 }
255 EIGEN_DEVICE_FUNC inline explicit Transform(const UniformScaling<Scalar>& s) {
256 check_template_params();
257 *this = s;
258 }
259 template <typename Derived>
260 EIGEN_DEVICE_FUNC inline explicit Transform(const RotationBase<Derived, Dim>& r) {
261 check_template_params();
262 *this = r;
263 }
264
265 typedef internal::transform_take_affine_part<Transform> take_affine_part;
266
268 template <typename OtherDerived>
269 EIGEN_DEVICE_FUNC inline explicit Transform(const EigenBase<OtherDerived>& other) {
270 EIGEN_STATIC_ASSERT(
271 (internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
272 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
273
274 check_template_params();
275 internal::transform_construct_from_matrix<OtherDerived, Mode, Options, Dim, HDim>::run(this, other.derived());
276 }
277
279 template <typename OtherDerived>
280 EIGEN_DEVICE_FUNC inline Transform& operator=(const EigenBase<OtherDerived>& other) {
281 EIGEN_STATIC_ASSERT(
282 (internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
283 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
284
285 internal::transform_construct_from_matrix<OtherDerived, Mode, Options, Dim, HDim>::run(this, other.derived());
286 return *this;
287 }
288
289 template <int OtherOptions>
290 EIGEN_DEVICE_FUNC inline Transform(const Transform<Scalar, Dim, Mode, OtherOptions>& other) {
291 check_template_params();
292 // only the options change, we can directly copy the matrices
293 m_matrix = other.matrix();
294 }
295
296 template <int OtherMode, int OtherOptions>
297 EIGEN_DEVICE_FUNC inline Transform(const Transform<Scalar, Dim, OtherMode, OtherOptions>& other) {
298 check_template_params();
299 // prevent conversions as:
300 // Affine | AffineCompact | Isometry = Projective
301 EIGEN_STATIC_ASSERT(internal::check_implication(OtherMode == int(Projective), Mode == int(Projective)),
302 YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
303
304 // prevent conversions as:
305 // Isometry = Affine | AffineCompact
306 EIGEN_STATIC_ASSERT(
307 internal::check_implication(OtherMode == int(Affine) || OtherMode == int(AffineCompact), Mode != int(Isometry)),
308 YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
309
310 enum {
311 ModeIsAffineCompact = Mode == int(AffineCompact),
312 OtherModeIsAffineCompact = OtherMode == int(AffineCompact)
313 };
314
315 if (EIGEN_CONST_CONDITIONAL(ModeIsAffineCompact == OtherModeIsAffineCompact)) {
316 // We need the block expression because the code is compiled for all
317 // combinations of transformations and will trigger a compile time error
318 // if one tries to assign the matrices directly
319 m_matrix.template block<Dim, Dim + 1>(0, 0) = other.matrix().template block<Dim, Dim + 1>(0, 0);
320 makeAffine();
321 } else if (EIGEN_CONST_CONDITIONAL(OtherModeIsAffineCompact)) {
323 internal::transform_construct_from_matrix<OtherMatrixType, Mode, Options, Dim, HDim>::run(this, other.matrix());
324 } else {
325 // here we know that Mode == AffineCompact and OtherMode != AffineCompact.
326 // if OtherMode were Projective, the static assert above would already have caught it.
327 // So the only possibility is that OtherMode == Affine
328 linear() = other.linear();
329 translation() = other.translation();
330 }
331 }
332
333 template <typename OtherDerived>
334 EIGEN_DEVICE_FUNC Transform(const ReturnByValue<OtherDerived>& other) {
335 check_template_params();
336 other.evalTo(*this);
337 }
338
339 template <typename OtherDerived>
340 EIGEN_DEVICE_FUNC Transform& operator=(const ReturnByValue<OtherDerived>& other) {
341 other.evalTo(*this);
342 return *this;
343 }
344
345#ifdef EIGEN_QT_SUPPORT
346#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
347 inline Transform(const QMatrix& other);
348 inline Transform& operator=(const QMatrix& other);
349 inline QMatrix toQMatrix(void) const;
350#endif
351 inline Transform(const QTransform& other);
352 inline Transform& operator=(const QTransform& other);
353 inline QTransform toQTransform(void) const;
354#endif
355
356 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT {
357 return int(Mode) == int(Projective) ? m_matrix.cols() : (m_matrix.cols() - 1);
358 }
359 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
360
363 EIGEN_DEVICE_FUNC inline Scalar operator()(Index row, Index col) const { return m_matrix(row, col); }
366 EIGEN_DEVICE_FUNC inline Scalar& operator()(Index row, Index col) { return m_matrix(row, col); }
367
369 EIGEN_DEVICE_FUNC inline const MatrixType& matrix() const { return m_matrix; }
371 EIGEN_DEVICE_FUNC inline MatrixType& matrix() { return m_matrix; }
372
374 EIGEN_DEVICE_FUNC inline ConstLinearPart linear() const { return ConstLinearPart(m_matrix, 0, 0); }
376 EIGEN_DEVICE_FUNC inline LinearPart linear() { return LinearPart(m_matrix, 0, 0); }
377
379 EIGEN_DEVICE_FUNC inline ConstAffinePart affine() const { return take_affine_part::run(m_matrix); }
381 EIGEN_DEVICE_FUNC inline AffinePart affine() { return take_affine_part::run(m_matrix); }
382
384 EIGEN_DEVICE_FUNC inline ConstTranslationPart translation() const { return ConstTranslationPart(m_matrix, 0, Dim); }
386 EIGEN_DEVICE_FUNC inline TranslationPart translation() { return TranslationPart(m_matrix, 0, Dim); }
387
413 // note: this function is defined here because some compilers cannot find the respective declaration
414 template <typename OtherDerived>
415 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const typename internal::transform_right_product_impl<Transform,
416 OtherDerived>::ResultType
417 operator*(const EigenBase<OtherDerived>& other) const {
418 return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this, other.derived());
419 }
420
428 template <typename OtherDerived>
429 friend EIGEN_DEVICE_FUNC inline const typename internal::transform_left_product_impl<OtherDerived, Mode, Options,
430 Dim_, Dim_ + 1>::ResultType
432 return internal::transform_left_product_impl<OtherDerived, Mode, Options, Dim, HDim>::run(a.derived(), b);
433 }
434
441 template <typename DiagonalDerived>
442 EIGEN_DEVICE_FUNC inline const TransformTimeDiagonalReturnType operator*(
443 const DiagonalBase<DiagonalDerived>& b) const {
445 res.linearExt() *= b;
446 return res;
447 }
448
455 template <typename DiagonalDerived>
457 const Transform& b) {
459 res.linear().noalias() = a * b.linear();
460 res.translation().noalias() = a * b.translation();
461 if (EIGEN_CONST_CONDITIONAL(Mode != int(AffineCompact))) res.matrix().row(Dim) = b.matrix().row(Dim);
462 return res;
463 }
464
465 template <typename OtherDerived>
466 EIGEN_DEVICE_FUNC inline Transform& operator*=(const EigenBase<OtherDerived>& other) {
467 return *this = *this * other;
468 }
469
471 EIGEN_DEVICE_FUNC inline const Transform operator*(const Transform& other) const {
472 return internal::transform_transform_product_impl<Transform, Transform>::run(*this, other);
473 }
474
475#if EIGEN_COMP_ICC
476 private:
477 // this intermediate structure permits to workaround a bug in ICC 11:
478 // error: template instantiation resulted in unexpected function type of "Eigen::Transform<double, 3, 32, 0>
479 // (const Eigen::Transform<double, 3, 2, 0> &) const"
480 // (the meaning of a name may have changed since the template declaration -- the type of the template is:
481 // "Eigen::internal::transform_transform_product_impl<Eigen::Transform<double, 3, 32, 0>,
482 // Eigen::Transform<double, 3, Mode, Options>, <expression>>::ResultType (const Eigen::Transform<double, 3, Mode,
483 // Options> &) const")
484 //
485 template <int OtherMode, int OtherOptions>
486 struct icc_11_workaround {
487 typedef internal::transform_transform_product_impl<Transform, Transform<Scalar, Dim, OtherMode, OtherOptions> >
488 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 operator*(
497 typedef typename icc_11_workaround<OtherMode, OtherOptions>::ProductType ProductType;
498 return ProductType::run(*this, other);
499 }
500#else
502 template <int OtherMode, int OtherOptions>
503 EIGEN_DEVICE_FUNC inline
504 typename internal::transform_transform_product_impl<Transform,
507 return internal::transform_transform_product_impl<Transform, Transform<Scalar, Dim, OtherMode, OtherOptions> >::run(
508 *this, other);
509 }
510#endif
511
513 EIGEN_DEVICE_FUNC void setIdentity() { m_matrix.setIdentity(); }
514
519 EIGEN_DEVICE_FUNC static const Transform Identity() { return Transform(MatrixType::Identity()); }
520
521 template <typename OtherDerived>
522 EIGEN_DEVICE_FUNC inline Transform& scale(const MatrixBase<OtherDerived>& other);
523
524 template <typename OtherDerived>
525 EIGEN_DEVICE_FUNC inline Transform& prescale(const MatrixBase<OtherDerived>& other);
526
527 EIGEN_DEVICE_FUNC inline Transform& scale(const Scalar& s);
528 EIGEN_DEVICE_FUNC inline Transform& prescale(const Scalar& s);
529
530 template <typename OtherDerived>
531 EIGEN_DEVICE_FUNC inline Transform& translate(const MatrixBase<OtherDerived>& other);
532
533 template <typename OtherDerived>
534 EIGEN_DEVICE_FUNC inline Transform& pretranslate(const MatrixBase<OtherDerived>& other);
535
536 template <typename RotationType>
537 EIGEN_DEVICE_FUNC inline Transform& rotate(const RotationType& rotation);
538
539 template <typename RotationType>
540 EIGEN_DEVICE_FUNC inline Transform& prerotate(const RotationType& rotation);
541
542 EIGEN_DEVICE_FUNC Transform& shear(const Scalar& sx, const Scalar& sy);
543 EIGEN_DEVICE_FUNC Transform& preshear(const Scalar& sx, const Scalar& sy);
544
545 EIGEN_DEVICE_FUNC inline Transform& operator=(const TranslationType& t);
546
547 EIGEN_DEVICE_FUNC inline Transform& operator*=(const TranslationType& t) { return translate(t.vector()); }
548
549 EIGEN_DEVICE_FUNC inline Transform operator*(const TranslationType& t) const;
550
551 EIGEN_DEVICE_FUNC inline Transform& operator=(const UniformScaling<Scalar>& t);
552
553 EIGEN_DEVICE_FUNC inline Transform& operator*=(const UniformScaling<Scalar>& s) { return scale(s.factor()); }
554
555 EIGEN_DEVICE_FUNC inline TransformTimeDiagonalReturnType operator*(const UniformScaling<Scalar>& s) const {
557 res.scale(s.factor());
558 return res;
559 }
560
561 EIGEN_DEVICE_FUNC inline Transform& operator*=(const DiagonalMatrix<Scalar, Dim>& s) {
562 linearExt() *= s;
563 return *this;
564 }
565
566 template <typename Derived>
567 EIGEN_DEVICE_FUNC inline Transform& operator=(const RotationBase<Derived, Dim>& r);
568 template <typename Derived>
569 EIGEN_DEVICE_FUNC inline Transform& operator*=(const RotationBase<Derived, Dim>& r) {
570 return rotate(r.toRotationMatrix());
571 }
572 template <typename Derived>
573 EIGEN_DEVICE_FUNC inline Transform operator*(const RotationBase<Derived, Dim>& r) const;
574
575 typedef std::conditional_t<int(Mode) == Isometry, ConstLinearPart, const LinearMatrixType> RotationReturnType;
576 EIGEN_DEVICE_FUNC RotationReturnType rotation() const;
577
578 template <typename RotationMatrixType, typename ScalingMatrixType>
579 EIGEN_DEVICE_FUNC void computeRotationScaling(RotationMatrixType* rotation, ScalingMatrixType* scaling) const;
580 template <typename ScalingMatrixType, typename RotationMatrixType>
581 EIGEN_DEVICE_FUNC void computeScalingRotation(ScalingMatrixType* scaling, RotationMatrixType* rotation) const;
582
583 template <typename PositionDerived, typename OrientationType, typename ScaleDerived>
584 EIGEN_DEVICE_FUNC Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived>& position,
585 const OrientationType& orientation,
586 const MatrixBase<ScaleDerived>& scale);
587
588 EIGEN_DEVICE_FUNC inline Transform inverse(TransformTraits traits = (TransformTraits)Mode) const;
589
591 EIGEN_DEVICE_FUNC const Scalar* data() const { return m_matrix.data(); }
593 EIGEN_DEVICE_FUNC Scalar* data() { return m_matrix.data(); }
594
600 template <typename NewScalarType>
601 EIGEN_DEVICE_FUNC inline
602 typename internal::cast_return_type<Transform, Transform<NewScalarType, Dim, Mode, Options> >::type
603 cast() const {
604 return typename internal::cast_return_type<Transform, Transform<NewScalarType, Dim, Mode, Options> >::type(*this);
605 }
606
608 template <typename OtherScalarType>
609 EIGEN_DEVICE_FUNC inline explicit Transform(const Transform<OtherScalarType, Dim, Mode, Options>& other) {
610 check_template_params();
611 m_matrix = other.matrix().template cast<Scalar>();
612 }
613
618 EIGEN_DEVICE_FUNC bool isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec =
620 return m_matrix.isApprox(other.m_matrix, prec);
621 }
622
625 EIGEN_DEVICE_FUNC void makeAffine() { internal::transform_make_affine<int(Mode)>::run(m_matrix); }
626
631 EIGEN_DEVICE_FUNC inline Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, Dim> linearExt() {
632 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, Dim > (0, 0);
633 }
638 EIGEN_DEVICE_FUNC inline const Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, Dim> linearExt() const {
639 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, Dim > (0, 0);
640 }
641
646 EIGEN_DEVICE_FUNC inline Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, 1> translationExt() {
647 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, 1 > (0, Dim);
648 }
653 EIGEN_DEVICE_FUNC inline const Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, 1> translationExt()
654 const {
655 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, 1 > (0, Dim);
656 }
657
658#ifdef EIGEN_TRANSFORM_PLUGIN
659#include EIGEN_TRANSFORM_PLUGIN
660#endif
661
662 protected:
663#ifndef EIGEN_PARSED_BY_DOXYGEN
664 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE void check_template_params() {
665 EIGEN_STATIC_ASSERT((Options & (DontAlign | RowMajor)) == Options, INVALID_MATRIX_TEMPLATE_PARAMETERS)
666 }
667#endif
668};
669
671typedef Transform<float, 2, Isometry> Isometry2f;
673typedef Transform<float, 3, Isometry> Isometry3f;
675typedef Transform<double, 2, Isometry> Isometry2d;
677typedef Transform<double, 3, Isometry> Isometry3d;
678
680typedef Transform<float, 2, Affine> Affine2f;
682typedef Transform<float, 3, Affine> Affine3f;
684typedef Transform<double, 2, Affine> Affine2d;
686typedef Transform<double, 3, Affine> Affine3d;
687
689typedef Transform<float, 2, AffineCompact> AffineCompact2f;
691typedef Transform<float, 3, AffineCompact> AffineCompact3f;
693typedef Transform<double, 2, AffineCompact> AffineCompact2d;
695typedef Transform<double, 3, AffineCompact> AffineCompact3d;
696
698typedef Transform<float, 2, Projective> Projective2f;
700typedef Transform<float, 3, Projective> Projective3f;
702typedef Transform<double, 2, Projective> Projective2d;
704typedef Transform<double, 3, Projective> Projective3d;
705
706/**************************
707*** Optional QT support ***
708**************************/
709
710#ifdef EIGEN_QT_SUPPORT
711
712#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
717template <typename Scalar, int Dim, int Mode, int Options>
718Transform<Scalar, Dim, Mode, Options>::Transform(const QMatrix& other) {
719 check_template_params();
720 *this = other;
721}
722
727template <typename Scalar, int Dim, int Mode, int Options>
728Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::operator=(const QMatrix& other) {
729 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
730 if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
731 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy();
732 else
733 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy(), 0, 0, 1;
734 return *this;
735}
736
743template <typename Scalar, int Dim, int Mode, int Options>
744QMatrix Transform<Scalar, Dim, Mode, Options>::toQMatrix(void) const {
745 check_template_params();
746 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
747 return QMatrix(m_matrix.coeff(0, 0), m_matrix.coeff(1, 0), m_matrix.coeff(0, 1), m_matrix.coeff(1, 1),
748 m_matrix.coeff(0, 2), m_matrix.coeff(1, 2));
749}
750#endif
751
756template <typename Scalar, int Dim, int Mode, int Options>
758 check_template_params();
759 *this = other;
760}
761
766template <typename Scalar, int Dim, int Mode, int Options>
768 check_template_params();
769 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
770 if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
771 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy();
772 else
773 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy(), other.m13(), other.m23(),
774 other.m33();
775 return *this;
776}
777
782template <typename Scalar, int Dim, int Mode, int Options>
784 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
785 if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
786 return QTransform(m_matrix.coeff(0, 0), m_matrix.coeff(1, 0), m_matrix.coeff(0, 1), m_matrix.coeff(1, 1),
787 m_matrix.coeff(0, 2), m_matrix.coeff(1, 2));
788 else
789 return QTransform(m_matrix.coeff(0, 0), m_matrix.coeff(1, 0), m_matrix.coeff(2, 0), m_matrix.coeff(0, 1),
790 m_matrix.coeff(1, 1), m_matrix.coeff(2, 1), m_matrix.coeff(0, 2), m_matrix.coeff(1, 2),
791 m_matrix.coeff(2, 2));
792}
793#endif
794
795/*********************
796*** Procedural API ***
797*********************/
798
803template <typename Scalar, int Dim, int Mode, int Options>
804template <typename OtherDerived>
806 const MatrixBase<OtherDerived>& other) {
807 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
808 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
809 linearExt().noalias() = (linearExt() * other.asDiagonal());
810 return *this;
811}
812
817template <typename Scalar, int Dim, int Mode, int Options>
819 const Scalar& s) {
820 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
821 linearExt() *= s;
822 return *this;
823}
824
829template <typename Scalar, int Dim, int Mode, int Options>
830template <typename OtherDerived>
832 const MatrixBase<OtherDerived>& other) {
833 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
834 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
835 affine().noalias() = (other.asDiagonal() * affine());
836 return *this;
837}
838
843template <typename Scalar, int Dim, int Mode, int Options>
845 const Scalar& s) {
846 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
847 m_matrix.template topRows<Dim>() *= s;
848 return *this;
849}
850
855template <typename Scalar, int Dim, int Mode, int Options>
856template <typename OtherDerived>
858 const MatrixBase<OtherDerived>& other) {
859 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
860 translationExt() += linearExt() * other;
861 return *this;
862}
863
868template <typename Scalar, int Dim, int Mode, int Options>
869template <typename OtherDerived>
871 const MatrixBase<OtherDerived>& other) {
872 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
873 if (EIGEN_CONST_CONDITIONAL(int(Mode) == int(Projective)))
874 affine() += other * m_matrix.row(Dim);
875 else
876 translation() += other;
877 return *this;
878}
879
897template <typename Scalar, int Dim, int Mode, int Options>
898template <typename RotationType>
900 const RotationType& rotation) {
901 linearExt() *= internal::toRotationMatrix<Scalar, Dim>(rotation);
902 return *this;
903}
904
912template <typename Scalar, int Dim, int Mode, int Options>
913template <typename RotationType>
915 const RotationType& rotation) {
916 m_matrix.template block<Dim, HDim>(0, 0) =
917 internal::toRotationMatrix<Scalar, Dim>(rotation) * m_matrix.template block<Dim, HDim>(0, 0);
918 return *this;
919}
920
926template <typename Scalar, int Dim, int Mode, int Options>
928 const Scalar& sx, const Scalar& sy) {
929 EIGEN_STATIC_ASSERT(int(Dim) == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
930 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
931 VectorType tmp = linear().col(0) * sy + linear().col(1);
932 linear() << linear().col(0) + linear().col(1) * sx, tmp;
933 return *this;
934}
935
941template <typename Scalar, int Dim, int Mode, int Options>
943 const Scalar& sx, const Scalar& sy) {
944 EIGEN_STATIC_ASSERT(int(Dim) == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
945 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
946 m_matrix.template block<Dim, HDim>(0, 0) =
947 LinearMatrixType({{1, sy}, {sx, 1}}) * m_matrix.template block<Dim, HDim>(0, 0);
948 return *this;
949}
950
951/******************************************************
952*** Scaling, Translation and Rotation compatibility ***
953******************************************************/
954
955template <typename Scalar, int Dim, int Mode, int Options>
957 const TranslationType& t) {
958 linear().setIdentity();
959 translation() = t.vector();
960 makeAffine();
961 return *this;
962}
963
964template <typename Scalar, int Dim, int Mode, int Options>
965EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options> Transform<Scalar, Dim, Mode, Options>::operator*(
966 const TranslationType& t) const {
967 Transform res = *this;
968 res.translate(t.vector());
969 return res;
970}
971
972template <typename Scalar, int Dim, int Mode, int Options>
973EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::operator=(
974 const UniformScaling<Scalar>& s) {
975 m_matrix.setZero();
976 linear().diagonal().fill(s.factor());
977 makeAffine();
978 return *this;
979}
980
981template <typename Scalar, int Dim, int Mode, int Options>
982template <typename Derived>
983EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::operator=(
984 const RotationBase<Derived, Dim>& r) {
985 linear() = internal::toRotationMatrix<Scalar, Dim>(r);
986 translation().setZero();
987 makeAffine();
988 return *this;
989}
990
991template <typename Scalar, int Dim, int Mode, int Options>
992template <typename Derived>
993EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options> Transform<Scalar, Dim, Mode, Options>::operator*(
994 const RotationBase<Derived, Dim>& r) const {
995 Transform res = *this;
996 res.rotate(r.derived());
997 return res;
998}
999
1000/************************
1001*** Special functions ***
1002************************/
1003
1004namespace internal {
1005template <int Mode>
1006struct transform_rotation_impl {
1007 template <typename TransformType>
1008 EIGEN_DEVICE_FUNC static inline const typename TransformType::LinearMatrixType run(const TransformType& t) {
1009 typedef typename TransformType::LinearMatrixType LinearMatrixType;
1010 LinearMatrixType result;
1011 t.computeRotationScaling(&result, (LinearMatrixType*)0);
1012 return result;
1013 }
1014};
1015template <>
1016struct transform_rotation_impl<Isometry> {
1017 template <typename TransformType>
1018 EIGEN_DEVICE_FUNC static inline typename TransformType::ConstLinearPart run(const TransformType& t) {
1019 return t.linear();
1020 }
1021};
1022} // namespace internal
1033template <typename Scalar, int Dim, int Mode, int Options>
1034EIGEN_DEVICE_FUNC typename Transform<Scalar, Dim, Mode, Options>::RotationReturnType
1036 return internal::transform_rotation_impl<Mode>::run(*this);
1037}
1038
1050template <typename Scalar, int Dim, int Mode, int Options>
1051template <typename RotationMatrixType, typename ScalingMatrixType>
1052EIGEN_DEVICE_FUNC void Transform<Scalar, Dim, Mode, Options>::computeRotationScaling(RotationMatrixType* rotation,
1053 ScalingMatrixType* scaling) const {
1054 // Note that JacobiSVD is faster than BDCSVD for small matrices.
1056
1057 Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0)
1058 ? Scalar(-1)
1059 : Scalar(1); // so x has absolute value 1
1060 VectorType sv(svd.singularValues());
1061 sv.coeffRef(Dim - 1) *= x;
1062 if (scaling) *scaling = svd.matrixV() * sv.asDiagonal() * svd.matrixV().adjoint();
1063 if (rotation) {
1064 LinearMatrixType m(svd.matrixU());
1065 m.col(Dim - 1) *= x;
1066 *rotation = m * svd.matrixV().adjoint();
1067 }
1068}
1069
1081template <typename Scalar, int Dim, int Mode, int Options>
1082template <typename ScalingMatrixType, typename RotationMatrixType>
1084 ScalingMatrixType* scaling, RotationMatrixType* rotation) const {
1085 // Note that JacobiSVD is faster than BDCSVD for small matrices.
1087
1088 Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0)
1089 ? Scalar(-1)
1090 : Scalar(1); // so x has absolute value 1
1091 VectorType sv(svd.singularValues());
1092 sv.coeffRef(Dim - 1) *= x;
1093 if (scaling) *scaling = svd.matrixU() * sv.asDiagonal() * svd.matrixU().adjoint();
1094 if (rotation) {
1095 LinearMatrixType m(svd.matrixU());
1096 m.col(Dim - 1) *= x;
1097 *rotation = m * svd.matrixV().adjoint();
1098 }
1099}
1100
1104template <typename Scalar, int Dim, int Mode, int Options>
1105template <typename PositionDerived, typename OrientationType, typename ScaleDerived>
1108 const OrientationType& orientation,
1109 const MatrixBase<ScaleDerived>& scale) {
1110 linear() = internal::toRotationMatrix<Scalar, Dim>(orientation);
1111 linear() *= scale.asDiagonal();
1112 translation() = position;
1113 makeAffine();
1114 return *this;
1115}
1116
1117namespace internal {
1118
1119template <int Mode>
1120struct transform_make_affine {
1121 template <typename MatrixType>
1122 EIGEN_DEVICE_FUNC static void run(MatrixType& mat) {
1123 static const int Dim = MatrixType::ColsAtCompileTime - 1;
1124 mat.template block<1, Dim>(Dim, 0).setZero();
1125 mat.coeffRef(Dim, Dim) = typename MatrixType::Scalar(1);
1126 }
1127};
1128
1129template <>
1130struct transform_make_affine<AffineCompact> {
1131 template <typename MatrixType>
1132 EIGEN_DEVICE_FUNC static void run(MatrixType&) {}
1133};
1134
1135// selector needed to avoid taking the inverse of a 3x4 matrix
1136template <typename TransformType, int Mode = TransformType::Mode>
1137struct projective_transform_inverse {
1138 EIGEN_DEVICE_FUNC static inline void run(const TransformType&, TransformType&) {}
1139};
1140
1141template <typename TransformType>
1142struct projective_transform_inverse<TransformType, Projective> {
1143 EIGEN_DEVICE_FUNC static inline void run(const TransformType& m, TransformType& res) {
1144 res.matrix() = m.matrix().inverse();
1145 }
1146};
1147
1148} // end namespace internal
1149
1170template <typename Scalar, int Dim, int Mode, int Options>
1172 TransformTraits hint) const {
1173 Transform res;
1174 if (hint == Projective) {
1175 internal::projective_transform_inverse<Transform>::run(*this, res);
1176 } else {
1177 if (hint == Isometry) {
1178 res.matrix().template topLeftCorner<Dim, Dim>() = linear().transpose();
1179 } else if (hint & Affine) {
1180 res.matrix().template topLeftCorner<Dim, Dim>() = linear().inverse();
1181 } else {
1182 eigen_assert(false && "Invalid transform traits in Transform::Inverse");
1183 }
1184 // translation and remaining parts
1185 res.matrix().template topRightCorner<Dim, 1>() = -res.matrix().template topLeftCorner<Dim, Dim>() * translation();
1186 res.makeAffine(); // we do need this, because in the beginning res is uninitialized
1187 }
1188 return res;
1189}
1190
1191namespace internal {
1192
1193/*****************************************************
1194*** Specializations of take affine part ***
1195*****************************************************/
1196
1197template <typename TransformType>
1198struct transform_take_affine_part {
1199 typedef typename TransformType::MatrixType MatrixType;
1200 typedef typename TransformType::AffinePart AffinePart;
1201 typedef typename TransformType::ConstAffinePart ConstAffinePart;
1202 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE AffinePart run(MatrixType& m) {
1203 return m.template block<TransformType::Dim, TransformType::HDim>(0, 0);
1204 }
1205 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ConstAffinePart run(const MatrixType& m) {
1206 return m.template block<TransformType::Dim, TransformType::HDim>(0, 0);
1207 }
1208};
1209
1210template <typename Scalar, int Dim, int Options>
1211struct transform_take_affine_part<Transform<Scalar, Dim, AffineCompact, Options> > {
1213 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE MatrixType& run(MatrixType& m) { return m; }
1214 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const MatrixType& run(const MatrixType& m) { return m; }
1215};
1216
1217/*****************************************************
1218*** Specializations of construct from matrix ***
1219*****************************************************/
1220
1221template <typename Other, int Mode, int Options, int Dim, int HDim>
1222struct transform_construct_from_matrix<Other, Mode, Options, Dim, HDim, Dim, Dim> {
1223 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1224 Transform<typename Other::Scalar, Dim, Mode, Options>* transform, const Other& other) {
1225 transform->linear() = other;
1226 transform->translation().setZero();
1227 transform->makeAffine();
1228 }
1229};
1230
1231template <typename Other, int Mode, int Options, int Dim, int HDim>
1232struct transform_construct_from_matrix<Other, Mode, Options, Dim, HDim, Dim, HDim> {
1233 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1234 Transform<typename Other::Scalar, Dim, Mode, Options>* transform, const Other& other) {
1235 transform->affine() = other;
1236 transform->makeAffine();
1237 }
1238};
1239
1240template <typename Other, int Mode, int Options, int Dim, int HDim>
1241struct transform_construct_from_matrix<Other, Mode, Options, Dim, HDim, HDim, HDim> {
1242 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1243 Transform<typename Other::Scalar, Dim, Mode, Options>* transform, const Other& other) {
1244 transform->matrix() = other;
1245 }
1246};
1247
1248template <typename Other, int Options, int Dim, int HDim>
1249struct transform_construct_from_matrix<Other, AffineCompact, Options, Dim, HDim, HDim, HDim> {
1250 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1251 Transform<typename Other::Scalar, Dim, AffineCompact, Options>* transform, const Other& other) {
1252 transform->matrix() = other.template block<Dim, HDim>(0, 0);
1253 }
1254};
1255
1256/**********************************************************
1257*** Specializations of operator* with rhs EigenBase ***
1258**********************************************************/
1259
1260template <int LhsMode, int RhsMode>
1261struct transform_product_result {
1262 enum {
1263 Mode = (LhsMode == (int)Projective || RhsMode == (int)Projective) ? Projective
1264 : (LhsMode == (int)Affine || RhsMode == (int)Affine) ? Affine
1265 : (LhsMode == (int)AffineCompact || RhsMode == (int)AffineCompact) ? AffineCompact
1266 : (LhsMode == (int)Isometry || RhsMode == (int)Isometry) ? Isometry
1267 : Projective
1268 };
1269};
1270
1271template <typename TransformType, typename MatrixType, int RhsCols>
1272struct transform_right_product_impl<TransformType, MatrixType, 0, RhsCols> {
1273 typedef typename MatrixType::PlainObject ResultType;
1274
1275 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1276 return T.matrix() * other;
1277 }
1278};
1279
1280template <typename TransformType, typename MatrixType, int RhsCols>
1281struct transform_right_product_impl<TransformType, MatrixType, 1, RhsCols> {
1282 enum {
1283 Dim = TransformType::Dim,
1284 HDim = TransformType::HDim,
1285 OtherRows = MatrixType::RowsAtCompileTime,
1286 OtherCols = MatrixType::ColsAtCompileTime
1287 };
1288
1289 typedef typename MatrixType::PlainObject ResultType;
1290
1291 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1292 EIGEN_STATIC_ASSERT(OtherRows == HDim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1293
1294 typedef Block<ResultType, Dim, OtherCols, int(MatrixType::RowsAtCompileTime) == Dim> TopLeftLhs;
1295
1296 ResultType res(other.rows(), other.cols());
1297 TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() = T.affine() * other;
1298 res.row(OtherRows - 1) = other.row(OtherRows - 1);
1299
1300 return res;
1301 }
1302};
1303
1304template <typename TransformType, typename MatrixType, int RhsCols>
1305struct transform_right_product_impl<TransformType, MatrixType, 2, RhsCols> {
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_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1316 EIGEN_STATIC_ASSERT(OtherRows == Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1317
1318 typedef Block<ResultType, Dim, OtherCols, true> TopLeftLhs;
1319 ResultType res(
1320 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
1327template <typename TransformType, typename MatrixType>
1328struct transform_right_product_impl<TransformType, MatrixType, 2, 1> // rhs is a vector of size Dim
1329{
1330 typedef typename TransformType::MatrixType TransformMatrix;
1331 enum {
1332 Dim = TransformType::Dim,
1333 HDim = TransformType::HDim,
1334 OtherRows = MatrixType::RowsAtCompileTime,
1335 WorkingRows = plain_enum_min(TransformMatrix::RowsAtCompileTime, HDim)
1336 };
1337
1338 typedef typename MatrixType::PlainObject ResultType;
1339
1340 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1341 EIGEN_STATIC_ASSERT(OtherRows == Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1342
1343 Matrix<typename ResultType::Scalar, Dim + 1, 1> rhs;
1344 rhs.template head<Dim>() = other;
1345 rhs[Dim] = typename ResultType::Scalar(1);
1346 Matrix<typename ResultType::Scalar, WorkingRows, 1> res(T.matrix() * rhs);
1347 return res.template head<Dim>();
1348 }
1349};
1350
1351/**********************************************************
1352*** Specializations of operator* with lhs EigenBase ***
1353**********************************************************/
1354
1355// generic HDim x HDim matrix * T => Projective
1356template <typename Other, int Mode, int Options, int Dim, int HDim>
1357struct transform_left_product_impl<Other, Mode, Options, Dim, HDim, HDim, HDim> {
1358 typedef Transform<typename Other::Scalar, Dim, Mode, Options> TransformType;
1359 typedef typename TransformType::MatrixType MatrixType;
1360 typedef Transform<typename Other::Scalar, Dim, Projective, Options> ResultType;
1361 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1362 return ResultType(other * tr.matrix());
1363 }
1364};
1365
1366// generic HDim x HDim matrix * AffineCompact => Projective
1367template <typename Other, int Options, int Dim, int HDim>
1368struct transform_left_product_impl<Other, AffineCompact, Options, Dim, HDim, HDim, HDim> {
1369 typedef Transform<typename Other::Scalar, Dim, AffineCompact, Options> TransformType;
1370 typedef typename TransformType::MatrixType MatrixType;
1371 typedef Transform<typename Other::Scalar, Dim, Projective, Options> ResultType;
1372 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1373 ResultType res;
1374 res.matrix().noalias() = other.template block<HDim, Dim>(0, 0) * tr.matrix();
1375 res.matrix().col(Dim) += other.col(Dim);
1376 return res;
1377 }
1378};
1379
1380// affine matrix * T
1381template <typename Other, int Mode, int Options, int Dim, int HDim>
1382struct transform_left_product_impl<Other, Mode, Options, Dim, HDim, Dim, HDim> {
1383 typedef Transform<typename Other::Scalar, Dim, Mode, Options> TransformType;
1384 typedef typename TransformType::MatrixType MatrixType;
1385 typedef TransformType ResultType;
1386 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1387 ResultType res;
1388 res.affine().noalias() = other * tr.matrix();
1389 res.matrix().row(Dim) = tr.matrix().row(Dim);
1390 return res;
1391 }
1392};
1393
1394// affine matrix * AffineCompact
1395template <typename Other, int Options, int Dim, int HDim>
1396struct transform_left_product_impl<Other, AffineCompact, Options, Dim, HDim, Dim, HDim> {
1397 typedef Transform<typename Other::Scalar, Dim, AffineCompact, Options> TransformType;
1398 typedef typename TransformType::MatrixType MatrixType;
1399 typedef TransformType ResultType;
1400 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1401 ResultType res;
1402 res.matrix().noalias() = other.template block<Dim, Dim>(0, 0) * tr.matrix();
1403 res.translation() += other.col(Dim);
1404 return res;
1405 }
1406};
1407
1408// linear matrix * T
1409template <typename Other, int Mode, int Options, int Dim, int HDim>
1410struct transform_left_product_impl<Other, Mode, Options, Dim, HDim, Dim, Dim> {
1411 typedef Transform<typename Other::Scalar, Dim, Mode, Options> TransformType;
1412 typedef typename TransformType::MatrixType MatrixType;
1413 typedef TransformType ResultType;
1414 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1415 TransformType res;
1416 if (Mode != int(AffineCompact)) res.matrix().row(Dim) = tr.matrix().row(Dim);
1417 res.matrix().template topRows<Dim>().noalias() = other * tr.matrix().template topRows<Dim>();
1418 return res;
1419 }
1420};
1421
1422/**********************************************************
1423*** Specializations of operator* with another Transform ***
1424**********************************************************/
1425
1426template <typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1427struct transform_transform_product_impl<Transform<Scalar, Dim, LhsMode, LhsOptions>,
1428 Transform<Scalar, Dim, RhsMode, RhsOptions>, false> {
1429 enum { ResultMode = transform_product_result<LhsMode, RhsMode>::Mode };
1430 typedef Transform<Scalar, Dim, LhsMode, LhsOptions> Lhs;
1431 typedef Transform<Scalar, Dim, RhsMode, RhsOptions> Rhs;
1432 typedef Transform<Scalar, Dim, ResultMode, LhsOptions> ResultType;
1433 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1434 ResultType res;
1435 res.linear() = lhs.linear() * rhs.linear();
1436 res.translation() = lhs.linear() * rhs.translation() + lhs.translation();
1437 res.makeAffine();
1438 return res;
1439 }
1440};
1441
1442template <typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1443struct transform_transform_product_impl<Transform<Scalar, Dim, LhsMode, LhsOptions>,
1444 Transform<Scalar, Dim, RhsMode, RhsOptions>, true> {
1445 typedef Transform<Scalar, Dim, LhsMode, LhsOptions> Lhs;
1446 typedef Transform<Scalar, Dim, RhsMode, RhsOptions> Rhs;
1447 typedef Transform<Scalar, Dim, Projective> ResultType;
1448 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1449 return ResultType(lhs.matrix() * rhs.matrix());
1450 }
1451};
1452
1453template <typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1454struct transform_transform_product_impl<Transform<Scalar, Dim, AffineCompact, LhsOptions>,
1455 Transform<Scalar, Dim, Projective, RhsOptions>, true> {
1456 typedef Transform<Scalar, Dim, AffineCompact, LhsOptions> Lhs;
1457 typedef Transform<Scalar, Dim, Projective, RhsOptions> Rhs;
1458 typedef Transform<Scalar, Dim, Projective> ResultType;
1459 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1460 ResultType res;
1461 res.matrix().template topRows<Dim>() = lhs.matrix() * rhs.matrix();
1462 res.matrix().row(Dim) = rhs.matrix().row(Dim);
1463 return res;
1464 }
1465};
1466
1467template <typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1468struct transform_transform_product_impl<Transform<Scalar, Dim, Projective, LhsOptions>,
1469 Transform<Scalar, Dim, AffineCompact, RhsOptions>, true> {
1470 typedef Transform<Scalar, Dim, Projective, LhsOptions> Lhs;
1471 typedef Transform<Scalar, Dim, AffineCompact, RhsOptions> Rhs;
1472 typedef Transform<Scalar, Dim, Projective> ResultType;
1473 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1474 ResultType res(lhs.matrix().template leftCols<Dim>() * rhs.matrix());
1475 res.matrix().col(Dim) += lhs.matrix().col(Dim);
1476 return res;
1477 }
1478};
1479
1480} // end namespace internal
1481
1482} // end namespace Eigen
1483
1484#endif // EIGEN_TRANSFORM_H
Expression of a fixed-size or dynamic-size block.
Definition Block.h:110
TransposeReturnType transpose()
Definition Transpose.h:160
bool isApprox(const DenseBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition Fuzzy.h:89
Base class for diagonal matrices and expressions.
Definition DiagonalMatrix.h:33
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition JacobiSVD.h:500
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
const DiagonalWrapper< const Derived > asDiagonal() const
Definition DiagonalMatrix.h:344
Derived & setIdentity()
Definition CwiseNullaryOp.h:839
const Inverse< Derived > inverse() const
Definition InverseImpl.h:279
static const IdentityReturnType Identity()
Definition CwiseNullaryOp.h:781
const AdjointReturnType adjoint() const
Definition Transpose.h:195
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:217
constexpr const Scalar & coeff(Index rowId, Index colId) const
Definition PlainObjectBase.h:198
Derived & setZero(Index size)
Definition CwiseNullaryOp.h:563
const Scalar * data() const
Definition PlainObjectBase.h:273
const SingularValuesType & singularValues() const
Definition SVDBase.h:200
const MatrixUType & matrixU() const
Definition SVDBase.h:173
const MatrixVType & matrixV() const
Definition SVDBase.h:189
Represents an homogeneous transformation in a N dimensional space.
Definition Transform.h:192
Block< MatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(int(Options) &RowMajor)==0 > LinearPart
Definition Transform.h:214
std::conditional_t< int(Mode)==int(AffineCompact), MatrixType &, Block< MatrixType, Dim, HDim > > AffinePart
Definition Transform.h:219
ConstAffinePart affine() const
Definition Transform.h:379
ConstLinearPart linear() const
Definition Transform.h:374
const internal::transform_right_product_impl< Transform, OtherDerived >::ResultType operator*(const EigenBase< OtherDerived > &other) const
Definition Transform.h:417
static const Transform Identity()
Returns an identity transformation.
Definition Transform.h:519
Transform & preshear(const Scalar &sx, const Scalar &sy)
Definition Transform.h:942
const Scalar * data() const
Definition Transform.h:591
Block< MatrixType, Dim, 1, !(internal::traits< MatrixType >::Flags &RowMajorBit)> TranslationPart
Definition Transform.h:227
void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const
Definition Transform.h:1083
internal::make_proper_matrix_type< Scalar, Rows, HDim, Options >::type MatrixType
Definition Transform.h:208
Scalar operator()(Index row, Index col) const
Definition Transform.h:363
Matrix< Scalar, Dim, 1 > VectorType
Definition Transform.h:225
Scalar_ Scalar
Definition Transform.h:204
Scalar & operator()(Index row, Index col)
Definition Transform.h:366
Scalar * data()
Definition Transform.h:593
const Block< ConstMatrixType, Dim, 1, !(internal::traits< MatrixType >::Flags &RowMajorBit)> ConstTranslationPart
Definition Transform.h:230
ConstTranslationPart translation() const
Definition Transform.h:384
Transform< Scalar, Dim, TransformTimeDiagonalMode > TransformTimeDiagonalReturnType
Definition Transform.h:237
bool isApprox(const Transform &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition Transform.h:618
AffinePart affine()
Definition Transform.h:381
std::conditional_t< int(Mode)==int(AffineCompact), const MatrixType &, const Block< const MatrixType, Dim, HDim > > ConstAffinePart
Definition Transform.h:223
void computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const
Definition Transform.h:1052
void setIdentity()
Definition Transform.h:513
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar_, Dim_==Dynamic ? Dynamic :(Dim_+1) *(Dim_+1)) enum
Definition Transform.h:194
friend const internal::transform_left_product_impl< OtherDerived, Mode, Options, Dim_, Dim_+1 >::ResultType operator*(const EigenBase< OtherDerived > &a, const Transform &b)
Definition Transform.h:431
const Transform operator*(const Transform &other) const
Definition Transform.h:471
friend TransformTimeDiagonalReturnType operator*(const DiagonalBase< DiagonalDerived > &a, const Transform &b)
Definition Transform.h:456
LinearPart linear()
Definition Transform.h:376
const MatrixType & matrix() const
Definition Transform.h:369
const TransformTimeDiagonalReturnType operator*(const DiagonalBase< DiagonalDerived > &b) const
Definition Transform.h:442
Transform & operator=(const EigenBase< OtherDerived > &other)
Definition Transform.h:280
Transform(const Transform< OtherScalarType, Dim, Mode, Options > &other)
Definition Transform.h:609
Transform(const EigenBase< OtherDerived > &other)
Definition Transform.h:269
internal::cast_return_type< Transform, Transform< NewScalarType, Dim, Mode, Options > >::type cast() const
Definition Transform.h:603
Transform inverse(TransformTraits traits=(TransformTraits) Mode) const
Definition Transform.h:1171
Matrix< Scalar, Dim, Dim, Options > LinearMatrixType
Definition Transform.h:212
RotationReturnType rotation() const
Definition Transform.h:1035
const MatrixType ConstMatrixType
Definition Transform.h:210
void makeAffine()
Definition Transform.h:625
MatrixType & matrix()
Definition Transform.h:371
internal::transform_transform_product_impl< Transform, Transform< Scalar, Dim, OtherMode, OtherOptions > >::ResultType operator*(const Transform< Scalar, Dim, OtherMode, OtherOptions > &other) const
Definition Transform.h:506
Eigen::Index Index
Definition Transform.h:206
Transform()
Definition Transform.h:245
Transform & shear(const Scalar &sx, const Scalar &sy)
Definition Transform.h:927
QTransform toQTransform(void) const
Definition Transform.h:783
const Block< ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(int(Options) &RowMajor)==0 > ConstLinearPart
Definition Transform.h:217
Translation< Scalar, Dim > TranslationType
Definition Transform.h:232
TranslationPart translation()
Definition Transform.h:386
Represents a translation transformation.
Definition Translation.h:33
TransformTraits
Definition Constants.h:453
@ DontAlign
Definition Constants.h:324
@ RowMajor
Definition Constants.h:320
@ Affine
Definition Constants.h:458
@ Projective
Definition Constants.h:462
@ AffineCompact
Definition Constants.h:460
@ Isometry
Definition Constants.h:455
const unsigned int RowMajorBit
Definition Constants.h:70
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const int Dynamic
Definition Constants.h:25
Definition EigenBase.h:33
Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523