Eigen  5.0.1-dev+7c7d8473
 
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 <gael.guennebaud@inria.fr>
5// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
6// Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
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 };
203
204 typedef Scalar_ Scalar;
205 typedef Eigen::Index StorageIndex;
208 typedef typename internal::make_proper_matrix_type<Scalar, Rows, HDim, Options>::type MatrixType;
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 constexpr Index rows() const noexcept {
357 return int(Mode) == int(Projective) ? m_matrix.cols() : (m_matrix.cols() - 1);
358 }
359 EIGEN_DEVICE_FUNC constexpr Index cols() const 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
368#ifdef EIGEN_MULTIDIMENSIONAL_SUBSCRIPT
371 EIGEN_DEVICE_FUNC inline Scalar operator[](Index row, Index col) const { return m_matrix[row, col]; }
374 EIGEN_DEVICE_FUNC inline Scalar& operator[](Index row, Index col) { return m_matrix[row, col]; }
375#endif
376
378 EIGEN_DEVICE_FUNC inline const MatrixType& matrix() const { return m_matrix; }
380 EIGEN_DEVICE_FUNC inline MatrixType& matrix() { return m_matrix; }
381
383 EIGEN_DEVICE_FUNC inline ConstLinearPart linear() const { return ConstLinearPart(m_matrix, 0, 0); }
385 EIGEN_DEVICE_FUNC inline LinearPart linear() { return LinearPart(m_matrix, 0, 0); }
386
388 EIGEN_DEVICE_FUNC inline ConstAffinePart affine() const { return take_affine_part::run(m_matrix); }
390 EIGEN_DEVICE_FUNC inline AffinePart affine() { return take_affine_part::run(m_matrix); }
391
393 EIGEN_DEVICE_FUNC inline ConstTranslationPart translation() const { return ConstTranslationPart(m_matrix, 0, Dim); }
395 EIGEN_DEVICE_FUNC inline TranslationPart translation() { return TranslationPart(m_matrix, 0, Dim); }
396
422 // note: this function is defined here because some compilers cannot find the respective declaration
423 template <typename OtherDerived>
424 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const typename internal::transform_right_product_impl<Transform,
425 OtherDerived>::ResultType
426 operator*(const EigenBase<OtherDerived>& other) const {
427 return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this, other.derived());
428 }
429
437 template <typename OtherDerived>
438 friend EIGEN_DEVICE_FUNC inline const typename internal::transform_left_product_impl<OtherDerived, Mode, Options,
439 Dim_, Dim_ + 1>::ResultType
441 return internal::transform_left_product_impl<OtherDerived, Mode, Options, Dim, HDim>::run(a.derived(), b);
442 }
443
450 template <typename DiagonalDerived>
451 EIGEN_DEVICE_FUNC inline const TransformTimeDiagonalReturnType operator*(
452 const DiagonalBase<DiagonalDerived>& b) const {
454 res.linearExt() *= b;
455 return res;
456 }
457
464 template <typename DiagonalDerived>
466 const Transform& b) {
468 res.linear().noalias() = a * b.linear();
469 res.translation().noalias() = a * b.translation();
470 if (EIGEN_CONST_CONDITIONAL(Mode != int(AffineCompact))) res.matrix().row(Dim) = b.matrix().row(Dim);
471 return res;
472 }
473
474 template <typename OtherDerived>
475 EIGEN_DEVICE_FUNC inline Transform& operator*=(const EigenBase<OtherDerived>& other) {
476 return *this = *this * other;
477 }
478
480 EIGEN_DEVICE_FUNC inline const Transform operator*(const Transform& other) const {
481 return internal::transform_transform_product_impl<Transform, Transform>::run(*this, other);
482 }
483
484#if EIGEN_COMP_ICC
485 private:
486 // this intermediate structure permits to workaround a bug in ICC 11:
487 // error: template instantiation resulted in unexpected function type of "Eigen::Transform<double, 3, 32, 0>
488 // (const Eigen::Transform<double, 3, 2, 0> &) const"
489 // (the meaning of a name may have changed since the template declaration -- the type of the template is:
490 // "Eigen::internal::transform_transform_product_impl<Eigen::Transform<double, 3, 32, 0>,
491 // Eigen::Transform<double, 3, Mode, Options>, <expression>>::ResultType (const Eigen::Transform<double, 3, Mode,
492 // Options> &) const")
493 //
494 template <int OtherMode, int OtherOptions>
495 struct icc_11_workaround {
496 typedef internal::transform_transform_product_impl<Transform, Transform<Scalar, Dim, OtherMode, OtherOptions> >
497 ProductType;
498 typedef typename ProductType::ResultType ResultType;
499 };
500
501 public:
503 template <int OtherMode, int OtherOptions>
504 inline typename icc_11_workaround<OtherMode, OtherOptions>::ResultType operator*(
506 typedef typename icc_11_workaround<OtherMode, OtherOptions>::ProductType ProductType;
507 return ProductType::run(*this, other);
508 }
509#else
511 template <int OtherMode, int OtherOptions>
512 EIGEN_DEVICE_FUNC inline
513 typename internal::transform_transform_product_impl<Transform,
516 return internal::transform_transform_product_impl<Transform, Transform<Scalar, Dim, OtherMode, OtherOptions> >::run(
517 *this, other);
518 }
519#endif
520
522 EIGEN_DEVICE_FUNC void setIdentity() { m_matrix.setIdentity(); }
523
528 EIGEN_DEVICE_FUNC static const Transform Identity() { return Transform(MatrixType::Identity()); }
529
530 template <typename OtherDerived>
531 EIGEN_DEVICE_FUNC inline Transform& scale(const MatrixBase<OtherDerived>& other);
532
533 template <typename OtherDerived>
534 EIGEN_DEVICE_FUNC inline Transform& prescale(const MatrixBase<OtherDerived>& other);
535
536 EIGEN_DEVICE_FUNC inline Transform& scale(const Scalar& s);
537 EIGEN_DEVICE_FUNC inline Transform& prescale(const Scalar& s);
538
539 template <typename OtherDerived>
540 EIGEN_DEVICE_FUNC inline Transform& translate(const MatrixBase<OtherDerived>& other);
541
542 template <typename OtherDerived>
543 EIGEN_DEVICE_FUNC inline Transform& pretranslate(const MatrixBase<OtherDerived>& other);
544
545 template <typename RotationType>
546 EIGEN_DEVICE_FUNC inline Transform& rotate(const RotationType& rotation);
547
548 template <typename RotationType>
549 EIGEN_DEVICE_FUNC inline Transform& prerotate(const RotationType& rotation);
550
551 EIGEN_DEVICE_FUNC Transform& shear(const Scalar& sx, const Scalar& sy);
552 EIGEN_DEVICE_FUNC Transform& preshear(const Scalar& sx, const Scalar& sy);
553
554 EIGEN_DEVICE_FUNC inline Transform& operator=(const TranslationType& t);
555
556 EIGEN_DEVICE_FUNC inline Transform& operator*=(const TranslationType& t) { return translate(t.vector()); }
557
558 EIGEN_DEVICE_FUNC inline Transform operator*(const TranslationType& t) const;
559
560 EIGEN_DEVICE_FUNC inline Transform& operator=(const UniformScaling<Scalar>& t);
561
562 EIGEN_DEVICE_FUNC inline Transform& operator*=(const UniformScaling<Scalar>& s) { return scale(s.factor()); }
563
564 EIGEN_DEVICE_FUNC inline TransformTimeDiagonalReturnType operator*(const UniformScaling<Scalar>& s) const {
566 res.scale(s.factor());
567 return res;
568 }
569
570 EIGEN_DEVICE_FUNC inline Transform& operator*=(const DiagonalMatrix<Scalar, Dim>& s) {
571 linearExt() *= s;
572 return *this;
573 }
574
575 template <typename Derived>
576 EIGEN_DEVICE_FUNC inline Transform& operator=(const RotationBase<Derived, Dim>& r);
577 template <typename Derived>
578 EIGEN_DEVICE_FUNC inline Transform& operator*=(const RotationBase<Derived, Dim>& r) {
579 return rotate(r.toRotationMatrix());
580 }
581 template <typename Derived>
582 EIGEN_DEVICE_FUNC inline Transform operator*(const RotationBase<Derived, Dim>& r) const;
583
584 typedef std::conditional_t<int(Mode) == Isometry, ConstLinearPart, const LinearMatrixType> RotationReturnType;
585 EIGEN_DEVICE_FUNC RotationReturnType rotation() const;
586
587 template <typename RotationMatrixType, typename ScalingMatrixType>
588 EIGEN_DEVICE_FUNC void computeRotationScaling(RotationMatrixType* rotation, ScalingMatrixType* scaling) const;
589 template <typename ScalingMatrixType, typename RotationMatrixType>
590 EIGEN_DEVICE_FUNC void computeScalingRotation(ScalingMatrixType* scaling, RotationMatrixType* rotation) const;
591
592 template <typename PositionDerived, typename OrientationType, typename ScaleDerived>
593 EIGEN_DEVICE_FUNC Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived>& position,
594 const OrientationType& orientation,
595 const MatrixBase<ScaleDerived>& scale);
596
597 EIGEN_DEVICE_FUNC inline Transform inverse(TransformTraits traits = (TransformTraits)Mode) const;
598
600 EIGEN_DEVICE_FUNC constexpr const Scalar* data() const { return m_matrix.data(); }
602 EIGEN_DEVICE_FUNC constexpr Scalar* data() { return m_matrix.data(); }
603
609 template <typename NewScalarType>
610 EIGEN_DEVICE_FUNC inline
611 typename internal::cast_return_type<Transform, Transform<NewScalarType, Dim, Mode, Options> >::type
612 cast() const {
613 return typename internal::cast_return_type<Transform, Transform<NewScalarType, Dim, Mode, Options> >::type(*this);
614 }
615
617 template <typename OtherScalarType>
618 EIGEN_DEVICE_FUNC inline explicit Transform(const Transform<OtherScalarType, Dim, Mode, Options>& other) {
619 check_template_params();
620 m_matrix = other.matrix().template cast<Scalar>();
621 }
622
627 EIGEN_DEVICE_FUNC bool isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec =
628 NumTraits<Scalar>::dummy_precision()) const {
629 return m_matrix.isApprox(other.m_matrix, prec);
630 }
631
634 EIGEN_DEVICE_FUNC void makeAffine() { internal::transform_make_affine<int(Mode)>::run(m_matrix); }
635
640 EIGEN_DEVICE_FUNC inline Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, Dim> linearExt() {
641 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, Dim > (0, 0);
642 }
647 EIGEN_DEVICE_FUNC inline const Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, Dim> linearExt() const {
648 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, Dim > (0, 0);
649 }
650
655 EIGEN_DEVICE_FUNC inline Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, 1> translationExt() {
656 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, 1 > (0, Dim);
657 }
662 EIGEN_DEVICE_FUNC inline const Block<MatrixType, int(Mode) == int(Projective) ? HDim : Dim, 1> translationExt()
663 const {
664 return m_matrix.template block < int(Mode) == int(Projective) ? HDim : Dim, 1 > (0, Dim);
665 }
666
667#ifdef EIGEN_TRANSFORM_PLUGIN
668#include EIGEN_TRANSFORM_PLUGIN
669#endif
670
671 protected:
672#ifndef EIGEN_PARSED_BY_DOXYGEN
673 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE void check_template_params() {
674 EIGEN_STATIC_ASSERT((Options & (DontAlign | RowMajor)) == Options, INVALID_MATRIX_TEMPLATE_PARAMETERS)
675 }
676#endif
677};
678
680typedef Transform<float, 2, Isometry> Isometry2f;
682typedef Transform<float, 3, Isometry> Isometry3f;
684typedef Transform<double, 2, Isometry> Isometry2d;
686typedef Transform<double, 3, Isometry> Isometry3d;
687
689typedef Transform<float, 2, Affine> Affine2f;
691typedef Transform<float, 3, Affine> Affine3f;
693typedef Transform<double, 2, Affine> Affine2d;
695typedef Transform<double, 3, Affine> Affine3d;
696
698typedef Transform<float, 2, AffineCompact> AffineCompact2f;
700typedef Transform<float, 3, AffineCompact> AffineCompact3f;
702typedef Transform<double, 2, AffineCompact> AffineCompact2d;
704typedef Transform<double, 3, AffineCompact> AffineCompact3d;
705
707typedef Transform<float, 2, Projective> Projective2f;
709typedef Transform<float, 3, Projective> Projective3f;
711typedef Transform<double, 2, Projective> Projective2d;
713typedef Transform<double, 3, Projective> Projective3d;
714
715/**************************
716*** Optional QT support ***
717**************************/
718
719#ifdef EIGEN_QT_SUPPORT
720
721#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
726template <typename Scalar, int Dim, int Mode, int Options>
728 check_template_params();
729 *this = other;
730}
731
736template <typename Scalar, int Dim, int Mode, int Options>
738 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
739 if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
740 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy();
741 else
742 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy(), 0, 0, 1;
743 return *this;
744}
745
752template <typename Scalar, int Dim, int Mode, int Options>
754 check_template_params();
755 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
756 return QMatrix(m_matrix.coeff(0, 0), m_matrix.coeff(1, 0), m_matrix.coeff(0, 1), m_matrix.coeff(1, 1),
757 m_matrix.coeff(0, 2), m_matrix.coeff(1, 2));
758}
759#endif
760
765template <typename Scalar, int Dim, int Mode, int Options>
767 check_template_params();
768 *this = other;
769}
770
775template <typename Scalar, int Dim, int Mode, int Options>
777 check_template_params();
778 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
779 if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
780 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy();
781 else
782 m_matrix << other.m11(), other.m21(), other.dx(), other.m12(), other.m22(), other.dy(), other.m13(), other.m23(),
783 other.m33();
784 return *this;
785}
786
791template <typename Scalar, int Dim, int Mode, int Options>
793 EIGEN_STATIC_ASSERT(Dim == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
794 if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
795 return QTransform(m_matrix.coeff(0, 0), m_matrix.coeff(1, 0), m_matrix.coeff(0, 1), m_matrix.coeff(1, 1),
796 m_matrix.coeff(0, 2), m_matrix.coeff(1, 2));
797 else
798 return QTransform(m_matrix.coeff(0, 0), m_matrix.coeff(1, 0), m_matrix.coeff(2, 0), m_matrix.coeff(0, 1),
799 m_matrix.coeff(1, 1), m_matrix.coeff(2, 1), m_matrix.coeff(0, 2), m_matrix.coeff(1, 2),
800 m_matrix.coeff(2, 2));
801}
802#endif
803
804/*********************
805*** Procedural API ***
806*********************/
807
812template <typename Scalar, int Dim, int Mode, int Options>
813template <typename OtherDerived>
814EIGEN_DEVICE_FUNC Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::scale(
815 const MatrixBase<OtherDerived>& other) {
816 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
817 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
818 linearExt().noalias() = (linearExt() * other.asDiagonal());
819 return *this;
820}
821
826template <typename Scalar, int Dim, int Mode, int Options>
827EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::scale(
828 const Scalar& s) {
829 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
830 linearExt() *= s;
831 return *this;
832}
833
838template <typename Scalar, int Dim, int Mode, int Options>
839template <typename OtherDerived>
840EIGEN_DEVICE_FUNC Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::prescale(
841 const MatrixBase<OtherDerived>& other) {
842 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
843 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
844 affine().noalias() = (other.asDiagonal() * affine());
845 return *this;
846}
847
852template <typename Scalar, int Dim, int Mode, int Options>
853EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::prescale(
854 const Scalar& s) {
855 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
856 m_matrix.template topRows<Dim>() *= s;
857 return *this;
858}
859
864template <typename Scalar, int Dim, int Mode, int Options>
865template <typename OtherDerived>
866EIGEN_DEVICE_FUNC Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::translate(
867 const MatrixBase<OtherDerived>& other) {
868 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
869 translationExt() += linearExt() * other;
870 return *this;
871}
872
877template <typename Scalar, int Dim, int Mode, int Options>
878template <typename OtherDerived>
879EIGEN_DEVICE_FUNC Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::pretranslate(
880 const MatrixBase<OtherDerived>& other) {
881 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, int(Dim))
882 if (EIGEN_CONST_CONDITIONAL(int(Mode) == int(Projective)))
883 affine() += other * m_matrix.row(Dim);
884 else
885 translation() += other;
886 return *this;
887}
888
906template <typename Scalar, int Dim, int Mode, int Options>
907template <typename RotationType>
908EIGEN_DEVICE_FUNC Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::rotate(
909 const RotationType& rotation) {
910 linearExt() *= internal::toRotationMatrix<Scalar, Dim>(rotation);
911 return *this;
912}
913
921template <typename Scalar, int Dim, int Mode, int Options>
922template <typename RotationType>
923EIGEN_DEVICE_FUNC Transform<Scalar, Dim, Mode, Options>& Transform<Scalar, Dim, Mode, Options>::prerotate(
924 const RotationType& rotation) {
925 m_matrix.template block<Dim, HDim>(0, 0) =
926 internal::toRotationMatrix<Scalar, Dim>(rotation) * m_matrix.template block<Dim, HDim>(0, 0);
927 return *this;
928}
929
935template <typename Scalar, int Dim, int Mode, int Options>
937 const Scalar& sx, const Scalar& sy) {
938 EIGEN_STATIC_ASSERT(int(Dim) == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
939 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
940 VectorType tmp = linear().col(0) * sy + linear().col(1);
941 linear() << linear().col(0) + linear().col(1) * sx, tmp;
942 return *this;
943}
944
950template <typename Scalar, int Dim, int Mode, int Options>
952 const Scalar& sx, const Scalar& sy) {
953 EIGEN_STATIC_ASSERT(int(Dim) == 2, YOU_MADE_A_PROGRAMMING_MISTAKE)
954 EIGEN_STATIC_ASSERT(Mode != int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
955 m_matrix.template block<Dim, HDim>(0, 0) =
956 LinearMatrixType({{1, sy}, {sx, 1}}) * m_matrix.template block<Dim, HDim>(0, 0);
957 return *this;
958}
959
960/******************************************************
961*** Scaling, Translation and Rotation compatibility ***
962******************************************************/
963
964template <typename Scalar, int Dim, int Mode, int Options>
966 const TranslationType& t) {
967 linear().setIdentity();
968 translation() = t.vector();
969 makeAffine();
970 return *this;
971}
972
973template <typename Scalar, int Dim, int Mode, int Options>
974EIGEN_DEVICE_FUNC inline Transform<Scalar, Dim, Mode, Options> Transform<Scalar, Dim, Mode, Options>::operator*(
975 const TranslationType& t) const {
976 Transform res = *this;
977 res.translate(t.vector());
978 return res;
979}
980
981template <typename Scalar, int Dim, int Mode, int Options>
983 const UniformScaling<Scalar>& s) {
984 m_matrix.setZero();
985 linear().diagonal().fill(s.factor());
986 makeAffine();
987 return *this;
988}
989
990template <typename Scalar, int Dim, int Mode, int Options>
991template <typename Derived>
994 linear() = internal::toRotationMatrix<Scalar, Dim>(r);
995 translation().setZero();
996 makeAffine();
997 return *this;
998}
999
1000template <typename Scalar, int Dim, int Mode, int Options>
1001template <typename Derived>
1003 const RotationBase<Derived, Dim>& r) const {
1004 Transform res = *this;
1005 res.rotate(r.derived());
1006 return res;
1007}
1008
1009/************************
1010*** Special functions ***
1011************************/
1012
1013namespace internal {
1014template <int Mode>
1015struct transform_rotation_impl {
1016 template <typename TransformType>
1017 EIGEN_DEVICE_FUNC static inline const typename TransformType::LinearMatrixType run(const TransformType& t) {
1018 typedef typename TransformType::LinearMatrixType LinearMatrixType;
1019 LinearMatrixType result;
1020 t.computeRotationScaling(&result, (LinearMatrixType*)0);
1021 return result;
1022 }
1023};
1024template <>
1025struct transform_rotation_impl<Isometry> {
1026 template <typename TransformType>
1027 EIGEN_DEVICE_FUNC static inline typename TransformType::ConstLinearPart run(const TransformType& t) {
1028 return t.linear();
1029 }
1030};
1031} // namespace internal
1042template <typename Scalar, int Dim, int Mode, int Options>
1043EIGEN_DEVICE_FUNC typename Transform<Scalar, Dim, Mode, Options>::RotationReturnType
1045 return internal::transform_rotation_impl<Mode>::run(*this);
1046}
1047
1059template <typename Scalar, int Dim, int Mode, int Options>
1060template <typename RotationMatrixType, typename ScalingMatrixType>
1062 ScalingMatrixType* scaling) const {
1063 // Note that JacobiSVD is faster than BDCSVD for small matrices.
1065
1066 Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0)
1067 ? Scalar(-1)
1068 : Scalar(1); // so x has absolute value 1
1069 VectorType sv(svd.singularValues());
1070 sv.coeffRef(Dim - 1) *= x;
1071 if (scaling) (*scaling).noalias() = svd.matrixV() * sv.asDiagonal() * svd.matrixV().adjoint();
1072 if (rotation) {
1073 LinearMatrixType m(svd.matrixU());
1074 m.col(Dim - 1) *= x;
1075 (*rotation).noalias() = m * svd.matrixV().adjoint();
1076 }
1077}
1078
1090template <typename Scalar, int Dim, int Mode, int Options>
1091template <typename ScalingMatrixType, typename RotationMatrixType>
1093 ScalingMatrixType* scaling, RotationMatrixType* rotation) const {
1094 // Note that JacobiSVD is faster than BDCSVD for small matrices.
1096
1097 Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0)
1098 ? Scalar(-1)
1099 : Scalar(1); // so x has absolute value 1
1100 VectorType sv(svd.singularValues());
1101 sv.coeffRef(Dim - 1) *= x;
1102 if (scaling) *scaling = svd.matrixU() * sv.asDiagonal() * svd.matrixU().adjoint();
1103 if (rotation) {
1104 LinearMatrixType m(svd.matrixU());
1105 m.col(Dim - 1) *= x;
1106 *rotation = m * svd.matrixV().adjoint();
1107 }
1108}
1109
1113template <typename Scalar, int Dim, int Mode, int Options>
1114template <typename PositionDerived, typename OrientationType, typename ScaleDerived>
1116Transform<Scalar, Dim, Mode, Options>::fromPositionOrientationScale(const MatrixBase<PositionDerived>& position,
1117 const OrientationType& orientation,
1118 const MatrixBase<ScaleDerived>& scale) {
1119 linear() = internal::toRotationMatrix<Scalar, Dim>(orientation);
1120 linear() *= scale.asDiagonal();
1121 translation() = position;
1122 makeAffine();
1123 return *this;
1124}
1125
1126namespace internal {
1127
1128template <int Mode>
1129struct transform_make_affine {
1130 template <typename MatrixType>
1131 EIGEN_DEVICE_FUNC static void run(MatrixType& mat) {
1132 static const int Dim = MatrixType::ColsAtCompileTime - 1;
1133 mat.template block<1, Dim>(Dim, 0).setZero();
1134 mat.coeffRef(Dim, Dim) = typename MatrixType::Scalar(1);
1135 }
1136};
1137
1138template <>
1139struct transform_make_affine<AffineCompact> {
1140 template <typename MatrixType>
1141 EIGEN_DEVICE_FUNC static void run(MatrixType&) {}
1142};
1143
1144// selector needed to avoid taking the inverse of a 3x4 matrix
1145template <typename TransformType, int Mode = TransformType::Mode>
1146struct projective_transform_inverse {
1147 EIGEN_DEVICE_FUNC static inline void run(const TransformType&, TransformType&) {}
1148};
1149
1150template <typename TransformType>
1151struct projective_transform_inverse<TransformType, Projective> {
1152 EIGEN_DEVICE_FUNC static inline void run(const TransformType& m, TransformType& res) {
1153 res.matrix() = m.matrix().inverse();
1154 }
1155};
1156
1157} // end namespace internal
1158
1179template <typename Scalar, int Dim, int Mode, int Options>
1181 TransformTraits hint) const {
1182 Transform res;
1183 if (hint == Projective) {
1184 internal::projective_transform_inverse<Transform>::run(*this, res);
1185 } else {
1186 if (hint == Isometry) {
1187 res.matrix().template topLeftCorner<Dim, Dim>() = linear().transpose();
1188 } else if (hint & Affine) {
1189 res.matrix().template topLeftCorner<Dim, Dim>() = linear().inverse();
1190 } else {
1191 eigen_assert(false && "Invalid transform traits in Transform::Inverse");
1192 }
1193 // translation and remaining parts
1194 res.matrix().template topRightCorner<Dim, 1>().noalias() =
1195 -res.matrix().template topLeftCorner<Dim, Dim>() * translation();
1196 res.makeAffine(); // we do need this, because in the beginning res is uninitialized
1197 }
1198 return res;
1199}
1200
1201namespace internal {
1202
1203/*****************************************************
1204*** Specializations of take affine part ***
1205*****************************************************/
1206
1207template <typename TransformType>
1208struct transform_take_affine_part {
1209 typedef typename TransformType::MatrixType MatrixType;
1210 typedef typename TransformType::AffinePart AffinePart;
1211 typedef typename TransformType::ConstAffinePart ConstAffinePart;
1212 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE AffinePart run(MatrixType& m) {
1213 return m.template block<TransformType::Dim, TransformType::HDim>(0, 0);
1214 }
1215 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ConstAffinePart run(const MatrixType& m) {
1216 return m.template block<TransformType::Dim, TransformType::HDim>(0, 0);
1217 }
1218};
1219
1220template <typename Scalar, int Dim, int Options>
1221struct transform_take_affine_part<Transform<Scalar, Dim, AffineCompact, Options> > {
1223 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE MatrixType& run(MatrixType& m) { return m; }
1224 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const MatrixType& run(const MatrixType& m) { return m; }
1225};
1226
1227/*****************************************************
1228*** Specializations of construct from matrix ***
1229*****************************************************/
1230
1231template <typename Other, int Mode, int Options, int Dim, int HDim>
1232struct transform_construct_from_matrix<Other, Mode, Options, Dim, HDim, Dim, Dim> {
1233 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1234 Transform<typename Other::Scalar, Dim, Mode, Options>* transform, const Other& other) {
1235 transform->linear() = other;
1236 transform->translation().setZero();
1237 transform->makeAffine();
1238 }
1239};
1240
1241template <typename Other, int Mode, int Options, int Dim, int HDim>
1242struct transform_construct_from_matrix<Other, Mode, Options, Dim, HDim, Dim, HDim> {
1243 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1244 Transform<typename Other::Scalar, Dim, Mode, Options>* transform, const Other& other) {
1245 transform->affine() = other;
1246 transform->makeAffine();
1247 }
1248};
1249
1250template <typename Other, int Mode, int Options, int Dim, int HDim>
1251struct transform_construct_from_matrix<Other, Mode, Options, Dim, HDim, HDim, HDim> {
1252 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1253 Transform<typename Other::Scalar, Dim, Mode, Options>* transform, const Other& other) {
1254 transform->matrix() = other;
1255 }
1256};
1257
1258template <typename Other, int Options, int Dim, int HDim>
1259struct transform_construct_from_matrix<Other, AffineCompact, Options, Dim, HDim, HDim, HDim> {
1260 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void run(
1261 Transform<typename Other::Scalar, Dim, AffineCompact, Options>* transform, const Other& other) {
1262 transform->matrix() = other.template block<Dim, HDim>(0, 0);
1263 }
1264};
1265
1266/**********************************************************
1267*** Specializations of operator* with rhs EigenBase ***
1268**********************************************************/
1269
1270template <int LhsMode, int RhsMode>
1271struct transform_product_result {
1272 enum {
1273 Mode = (LhsMode == (int)Projective || RhsMode == (int)Projective) ? Projective
1274 : (LhsMode == (int)Affine || RhsMode == (int)Affine) ? Affine
1275 : (LhsMode == (int)AffineCompact || RhsMode == (int)AffineCompact) ? AffineCompact
1276 : (LhsMode == (int)Isometry || RhsMode == (int)Isometry) ? Isometry
1277 : Projective
1278 };
1279};
1280
1281template <typename TransformType, typename MatrixType, int RhsCols>
1282struct transform_right_product_impl<TransformType, MatrixType, 0, RhsCols> {
1283 typedef typename MatrixType::PlainObject ResultType;
1284
1285 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1286 return T.matrix() * other;
1287 }
1288};
1289
1290template <typename TransformType, typename MatrixType, int RhsCols>
1291struct transform_right_product_impl<TransformType, MatrixType, 1, RhsCols> {
1292 enum {
1293 Dim = TransformType::Dim,
1294 HDim = TransformType::HDim,
1295 OtherRows = MatrixType::RowsAtCompileTime,
1296 OtherCols = MatrixType::ColsAtCompileTime
1297 };
1298
1299 typedef typename MatrixType::PlainObject ResultType;
1300
1301 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1302 EIGEN_STATIC_ASSERT(OtherRows == HDim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1303
1304 typedef Block<ResultType, Dim, OtherCols, int(MatrixType::RowsAtCompileTime) == Dim> TopLeftLhs;
1305
1306 ResultType res(other.rows(), other.cols());
1307 TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() = T.affine() * other;
1308 res.row(OtherRows - 1) = other.row(OtherRows - 1);
1309
1310 return res;
1311 }
1312};
1313
1314template <typename TransformType, typename MatrixType, int RhsCols>
1315struct transform_right_product_impl<TransformType, MatrixType, 2, RhsCols> {
1316 enum {
1317 Dim = TransformType::Dim,
1318 HDim = TransformType::HDim,
1319 OtherRows = MatrixType::RowsAtCompileTime,
1320 OtherCols = MatrixType::ColsAtCompileTime
1321 };
1322
1323 typedef typename MatrixType::PlainObject ResultType;
1324
1325 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1326 EIGEN_STATIC_ASSERT(OtherRows == Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1327
1328 typedef Block<ResultType, Dim, OtherCols, true> TopLeftLhs;
1329 ResultType res(
1330 Replicate<typename TransformType::ConstTranslationPart, 1, OtherCols>(T.translation(), 1, other.cols()));
1331 TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() += T.linear() * other;
1332
1333 return res;
1334 }
1335};
1336
1337template <typename TransformType, typename MatrixType>
1338struct transform_right_product_impl<TransformType, MatrixType, 2, 1> // rhs is a vector of size Dim
1339{
1340 typedef typename TransformType::MatrixType TransformMatrix;
1341 enum {
1342 Dim = TransformType::Dim,
1343 HDim = TransformType::HDim,
1344 OtherRows = MatrixType::RowsAtCompileTime,
1345 WorkingRows = plain_enum_min(TransformMatrix::RowsAtCompileTime, HDim)
1346 };
1347
1348 typedef typename MatrixType::PlainObject ResultType;
1349
1350 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) {
1351 EIGEN_STATIC_ASSERT(OtherRows == Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES);
1352
1353 Matrix<typename ResultType::Scalar, Dim + 1, 1> rhs;
1354 rhs.template head<Dim>() = other;
1355 rhs[Dim] = typename ResultType::Scalar(1);
1356 Matrix<typename ResultType::Scalar, WorkingRows, 1> res(T.matrix() * rhs);
1357 return res.template head<Dim>();
1358 }
1359};
1360
1361/**********************************************************
1362*** Specializations of operator* with lhs EigenBase ***
1363**********************************************************/
1364
1365// generic HDim x HDim matrix * T => Projective
1366template <typename Other, int Mode, int Options, int Dim, int HDim>
1367struct transform_left_product_impl<Other, Mode, Options, Dim, HDim, HDim, HDim> {
1368 typedef Transform<typename Other::Scalar, Dim, Mode, Options> TransformType;
1369 typedef typename TransformType::MatrixType MatrixType;
1370 typedef Transform<typename Other::Scalar, Dim, Projective, Options> ResultType;
1371 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1372 return ResultType(other * tr.matrix());
1373 }
1374};
1375
1376// generic HDim x HDim matrix * AffineCompact => Projective
1377template <typename Other, int Options, int Dim, int HDim>
1378struct transform_left_product_impl<Other, AffineCompact, Options, Dim, HDim, HDim, HDim> {
1379 typedef Transform<typename Other::Scalar, Dim, AffineCompact, Options> TransformType;
1380 typedef typename TransformType::MatrixType MatrixType;
1381 typedef Transform<typename Other::Scalar, Dim, Projective, Options> ResultType;
1382 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1383 ResultType res;
1384 res.matrix().noalias() = other.template block<HDim, Dim>(0, 0) * tr.matrix();
1385 res.matrix().col(Dim) += other.col(Dim);
1386 return res;
1387 }
1388};
1389
1390// affine matrix * T
1391template <typename Other, int Mode, int Options, int Dim, int HDim>
1392struct transform_left_product_impl<Other, Mode, Options, Dim, HDim, Dim, HDim> {
1393 typedef Transform<typename Other::Scalar, Dim, Mode, Options> TransformType;
1394 typedef typename TransformType::MatrixType MatrixType;
1395 typedef TransformType ResultType;
1396 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1397 ResultType res;
1398 res.affine().noalias() = other * tr.matrix();
1399 res.matrix().row(Dim) = tr.matrix().row(Dim);
1400 return res;
1401 }
1402};
1403
1404// affine matrix * AffineCompact
1405template <typename Other, int Options, int Dim, int HDim>
1406struct transform_left_product_impl<Other, AffineCompact, Options, Dim, HDim, Dim, HDim> {
1407 typedef Transform<typename Other::Scalar, Dim, AffineCompact, Options> TransformType;
1408 typedef typename TransformType::MatrixType MatrixType;
1409 typedef TransformType ResultType;
1410 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1411 ResultType res;
1412 res.matrix().noalias() = other.template block<Dim, Dim>(0, 0) * tr.matrix();
1413 res.translation() += other.col(Dim);
1414 return res;
1415 }
1416};
1417
1418// linear matrix * T
1419template <typename Other, int Mode, int Options, int Dim, int HDim>
1420struct transform_left_product_impl<Other, Mode, Options, Dim, HDim, Dim, Dim> {
1421 typedef Transform<typename Other::Scalar, Dim, Mode, Options> TransformType;
1422 typedef typename TransformType::MatrixType MatrixType;
1423 typedef TransformType ResultType;
1424 static EIGEN_DEVICE_FUNC ResultType run(const Other& other, const TransformType& tr) {
1425 TransformType res;
1426 if (Mode != int(AffineCompact)) res.matrix().row(Dim) = tr.matrix().row(Dim);
1427 res.matrix().template topRows<Dim>().noalias() = other * tr.matrix().template topRows<Dim>();
1428 return res;
1429 }
1430};
1431
1432/**********************************************************
1433*** Specializations of operator* with another Transform ***
1434**********************************************************/
1435
1436template <typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1437struct transform_transform_product_impl<Transform<Scalar, Dim, LhsMode, LhsOptions>,
1438 Transform<Scalar, Dim, RhsMode, RhsOptions>, false> {
1439 enum { ResultMode = transform_product_result<LhsMode, RhsMode>::Mode };
1440 typedef Transform<Scalar, Dim, LhsMode, LhsOptions> Lhs;
1441 typedef Transform<Scalar, Dim, RhsMode, RhsOptions> Rhs;
1442 typedef Transform<Scalar, Dim, ResultMode, LhsOptions> ResultType;
1443 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1444 ResultType res;
1445 res.linear().noalias() = lhs.linear() * rhs.linear();
1446 res.translation() = lhs.linear() * rhs.translation() + lhs.translation();
1447 res.makeAffine();
1448 return res;
1449 }
1450};
1451
1452template <typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions>
1453struct transform_transform_product_impl<Transform<Scalar, Dim, LhsMode, LhsOptions>,
1454 Transform<Scalar, Dim, RhsMode, RhsOptions>, true> {
1455 typedef Transform<Scalar, Dim, LhsMode, LhsOptions> Lhs;
1456 typedef Transform<Scalar, Dim, RhsMode, RhsOptions> Rhs;
1457 typedef Transform<Scalar, Dim, Projective> ResultType;
1458 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1459 return ResultType(lhs.matrix() * rhs.matrix());
1460 }
1461};
1462
1463template <typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1464struct transform_transform_product_impl<Transform<Scalar, Dim, AffineCompact, LhsOptions>,
1465 Transform<Scalar, Dim, Projective, RhsOptions>, true> {
1466 typedef Transform<Scalar, Dim, AffineCompact, LhsOptions> Lhs;
1467 typedef Transform<Scalar, Dim, Projective, RhsOptions> Rhs;
1468 typedef Transform<Scalar, Dim, Projective> ResultType;
1469 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1470 ResultType res;
1471 res.matrix().template topRows<Dim>() = lhs.matrix() * rhs.matrix();
1472 res.matrix().row(Dim) = rhs.matrix().row(Dim);
1473 return res;
1474 }
1475};
1476
1477template <typename Scalar, int Dim, int LhsOptions, int RhsOptions>
1478struct transform_transform_product_impl<Transform<Scalar, Dim, Projective, LhsOptions>,
1479 Transform<Scalar, Dim, AffineCompact, RhsOptions>, true> {
1480 typedef Transform<Scalar, Dim, Projective, LhsOptions> Lhs;
1481 typedef Transform<Scalar, Dim, AffineCompact, RhsOptions> Rhs;
1482 typedef Transform<Scalar, Dim, Projective> ResultType;
1483 static EIGEN_DEVICE_FUNC ResultType run(const Lhs& lhs, const Rhs& rhs) {
1484 ResultType res(lhs.matrix().template leftCols<Dim>() * rhs.matrix());
1485 res.matrix().col(Dim) += lhs.matrix().col(Dim);
1486 return res;
1487 }
1488};
1489
1490} // end namespace internal
1491
1492} // end namespace Eigen
1493
1494#endif // EIGEN_TRANSFORM_H
Expression of a fixed-size or dynamic-size block.
Definition Block.h:110
RowXpr row(Index i)
Definition DenseBase.h:1085
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:347
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:191
constexpr const Scalar & coeff(Index rowId, Index colId) const
Definition PlainObjectBase.h:172
Derived & setZero(Index size)
Definition CwiseNullaryOp.h:567
Common base class for compact rotation representations.
Definition RotationBase.h:32
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:388
ConstLinearPart linear() const
Definition Transform.h:383
const internal::transform_right_product_impl< Transform, OtherDerived >::ResultType operator*(const EigenBase< OtherDerived > &other) const
Definition Transform.h:426
static const Transform Identity()
Returns an identity transformation.
Definition Transform.h:528
Transform & preshear(const Scalar &sx, const Scalar &sy)
Definition Transform.h:951
Block< MatrixType, Dim, 1, !(internal::traits< MatrixType >::Flags &RowMajorBit)> TranslationPart
Definition Transform.h:227
void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const
Definition Transform.h:1092
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
const Block< ConstMatrixType, Dim, 1, !(internal::traits< MatrixType >::Flags &RowMajorBit)> ConstTranslationPart
Definition Transform.h:230
ConstTranslationPart translation() const
Definition Transform.h:393
Transform & operator=(const QTransform &other)
Definition Transform.h:776
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:627
AffinePart affine()
Definition Transform.h:390
Transform(const QTransform &other)
Definition Transform.h:766
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:1061
Transform & scale(const Scalar &s)
Definition Transform.h:827
void setIdentity()
Definition Transform.h:522
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:440
constexpr Scalar * data()
Definition Transform.h:602
const Transform operator*(const Transform &other) const
Definition Transform.h:480
friend TransformTimeDiagonalReturnType operator*(const DiagonalBase< DiagonalDerived > &a, const Transform &b)
Definition Transform.h:465
LinearPart linear()
Definition Transform.h:385
const MatrixType & matrix() const
Definition Transform.h:378
const TransformTimeDiagonalReturnType operator*(const DiagonalBase< DiagonalDerived > &b) const
Definition Transform.h:451
Transform & prescale(const Scalar &s)
Definition Transform.h:853
Transform & operator=(const EigenBase< OtherDerived > &other)
Definition Transform.h:280
Transform(const Transform< OtherScalarType, Dim, Mode, Options > &other)
Definition Transform.h:618
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:612
Transform inverse(TransformTraits traits=(TransformTraits) Mode) const
Definition Transform.h:1180
Matrix< Scalar, Dim, Dim, Options > LinearMatrixType
Definition Transform.h:212
const MatrixType ConstMatrixType
Definition Transform.h:210
void makeAffine()
Definition Transform.h:634
MatrixType & matrix()
Definition Transform.h:380
internal::transform_transform_product_impl< Transform, Transform< Scalar, Dim, OtherMode, OtherOptions > >::ResultType operator*(const Transform< Scalar, Dim, OtherMode, OtherOptions > &other) const
Definition Transform.h:515
Eigen::Index Index
Definition Transform.h:206
Transform()
Definition Transform.h:245
constexpr const Scalar * data() const
Definition Transform.h:600
Transform & shear(const Scalar &sx, const Scalar &sy)
Definition Transform.h:936
QTransform toQTransform(void) const
Definition Transform.h:792
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:395
Represents a translation transformation.
Definition Translation.h:33
Represents a generic uniform scaling transformation.
Definition Scaling.h:47
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 B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const int Dynamic
Definition Constants.h:25
Definition EigenBase.h:33
constexpr Derived & derived()
Definition EigenBase.h:49