blob: d352f1f2b8abcbf23c5a184b1e09f22fa16af89b [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// 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//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_SCALING_H
11#define EIGEN_SCALING_H
12
13namespace Eigen {
14
15/** \geometry_module \ingroup Geometry_Module
16 *
Austin Schuhc55b0172022-02-20 17:52:35 -080017 * \class UniformScaling
Brian Silverman72890c22015-09-19 14:37:37 -040018 *
19 * \brief Represents a generic uniform scaling transformation
20 *
Austin Schuh189376f2018-12-20 22:11:15 +110021 * \tparam _Scalar the scalar type, i.e., the type of the coefficients.
Brian Silverman72890c22015-09-19 14:37:37 -040022 *
23 * This class represent a uniform scaling transformation. It is the return
24 * type of Scaling(Scalar), and most of the time this is the only way it
25 * is used. In particular, this class is not aimed to be used to store a scaling transformation,
26 * but rather to make easier the constructions and updates of Transform objects.
27 *
28 * To represent an axis aligned scaling, use the DiagonalMatrix class.
29 *
30 * \sa Scaling(), class DiagonalMatrix, MatrixBase::asDiagonal(), class Translation, class Transform
31 */
Austin Schuhc55b0172022-02-20 17:52:35 -080032
33namespace internal
34{
35 // This helper helps nvcc+MSVC to properly parse this file.
36 // See bug 1412.
37 template <typename Scalar, int Dim, int Mode>
38 struct uniformscaling_times_affine_returntype
39 {
40 enum
41 {
42 NewMode = int(Mode) == int(Isometry) ? Affine : Mode
43 };
44 typedef Transform <Scalar, Dim, NewMode> type;
45 };
46}
47
Brian Silverman72890c22015-09-19 14:37:37 -040048template<typename _Scalar>
49class UniformScaling
50{
51public:
52 /** the scalar type of the coefficients */
53 typedef _Scalar Scalar;
54
55protected:
56
57 Scalar m_factor;
58
59public:
60
61 /** Default constructor without initialization. */
62 UniformScaling() {}
63 /** Constructs and initialize a uniform scaling transformation */
64 explicit inline UniformScaling(const Scalar& s) : m_factor(s) {}
65
66 inline const Scalar& factor() const { return m_factor; }
67 inline Scalar& factor() { return m_factor; }
68
69 /** Concatenates two uniform scaling */
70 inline UniformScaling operator* (const UniformScaling& other) const
71 { return UniformScaling(m_factor * other.factor()); }
72
73 /** Concatenates a uniform scaling and a translation */
74 template<int Dim>
75 inline Transform<Scalar,Dim,Affine> operator* (const Translation<Scalar,Dim>& t) const;
76
77 /** Concatenates a uniform scaling and an affine transformation */
78 template<int Dim, int Mode, int Options>
Austin Schuhc55b0172022-02-20 17:52:35 -080079 inline typename
80 internal::uniformscaling_times_affine_returntype<Scalar,Dim,Mode>::type
81 operator* (const Transform<Scalar, Dim, Mode, Options>& t) const
Brian Silverman72890c22015-09-19 14:37:37 -040082 {
Austin Schuhc55b0172022-02-20 17:52:35 -080083 typename internal::uniformscaling_times_affine_returntype<Scalar,Dim,Mode>::type res = t;
Austin Schuh189376f2018-12-20 22:11:15 +110084 res.prescale(factor());
85 return res;
86 }
Brian Silverman72890c22015-09-19 14:37:37 -040087
88 /** Concatenates a uniform scaling and a linear transformation matrix */
89 // TODO returns an expression
90 template<typename Derived>
Austin Schuhc55b0172022-02-20 17:52:35 -080091 inline typename Eigen::internal::plain_matrix_type<Derived>::type operator* (const MatrixBase<Derived>& other) const
Brian Silverman72890c22015-09-19 14:37:37 -040092 { return other * m_factor; }
93
94 template<typename Derived,int Dim>
95 inline Matrix<Scalar,Dim,Dim> operator*(const RotationBase<Derived,Dim>& r) const
96 { return r.toRotationMatrix() * m_factor; }
97
98 /** \returns the inverse scaling */
99 inline UniformScaling inverse() const
100 { return UniformScaling(Scalar(1)/m_factor); }
101
102 /** \returns \c *this with scalar type casted to \a NewScalarType
103 *
104 * Note that if \a NewScalarType is equal to the current scalar type of \c *this
105 * then this function smartly returns a const reference to \c *this.
106 */
107 template<typename NewScalarType>
108 inline UniformScaling<NewScalarType> cast() const
109 { return UniformScaling<NewScalarType>(NewScalarType(m_factor)); }
110
111 /** Copy constructor with scalar type conversion */
112 template<typename OtherScalarType>
113 inline explicit UniformScaling(const UniformScaling<OtherScalarType>& other)
114 { m_factor = Scalar(other.factor()); }
115
116 /** \returns \c true if \c *this is approximately equal to \a other, within the precision
117 * determined by \a prec.
118 *
119 * \sa MatrixBase::isApprox() */
120 bool isApprox(const UniformScaling& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
121 { return internal::isApprox(m_factor, other.factor(), prec); }
122
123};
124
Austin Schuh189376f2018-12-20 22:11:15 +1100125/** \addtogroup Geometry_Module */
126//@{
127
128/** Concatenates a linear transformation matrix and a uniform scaling
129 * \relates UniformScaling
130 */
Austin Schuhc55b0172022-02-20 17:52:35 -0800131// NOTE this operator is defined in MatrixBase and not as a friend function
Brian Silverman72890c22015-09-19 14:37:37 -0400132// of UniformScaling to fix an internal crash of Intel's ICC
Austin Schuh189376f2018-12-20 22:11:15 +1100133template<typename Derived,typename Scalar>
134EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Derived,Scalar,product)
135operator*(const MatrixBase<Derived>& matrix, const UniformScaling<Scalar>& s)
136{ return matrix.derived() * s.factor(); }
Brian Silverman72890c22015-09-19 14:37:37 -0400137
138/** Constructs a uniform scaling from scale factor \a s */
Austin Schuh189376f2018-12-20 22:11:15 +1100139inline UniformScaling<float> Scaling(float s) { return UniformScaling<float>(s); }
Brian Silverman72890c22015-09-19 14:37:37 -0400140/** Constructs a uniform scaling from scale factor \a s */
Austin Schuh189376f2018-12-20 22:11:15 +1100141inline UniformScaling<double> Scaling(double s) { return UniformScaling<double>(s); }
Brian Silverman72890c22015-09-19 14:37:37 -0400142/** Constructs a uniform scaling from scale factor \a s */
143template<typename RealScalar>
Austin Schuh189376f2018-12-20 22:11:15 +1100144inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s)
Brian Silverman72890c22015-09-19 14:37:37 -0400145{ return UniformScaling<std::complex<RealScalar> >(s); }
146
147/** Constructs a 2D axis aligned scaling */
148template<typename Scalar>
Austin Schuh189376f2018-12-20 22:11:15 +1100149inline DiagonalMatrix<Scalar,2> Scaling(const Scalar& sx, const Scalar& sy)
Brian Silverman72890c22015-09-19 14:37:37 -0400150{ return DiagonalMatrix<Scalar,2>(sx, sy); }
151/** Constructs a 3D axis aligned scaling */
152template<typename Scalar>
Austin Schuh189376f2018-12-20 22:11:15 +1100153inline DiagonalMatrix<Scalar,3> Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz)
Brian Silverman72890c22015-09-19 14:37:37 -0400154{ return DiagonalMatrix<Scalar,3>(sx, sy, sz); }
155
156/** Constructs an axis aligned scaling expression from vector expression \a coeffs
157 * This is an alias for coeffs.asDiagonal()
158 */
159template<typename Derived>
Austin Schuh189376f2018-12-20 22:11:15 +1100160inline const DiagonalWrapper<const Derived> Scaling(const MatrixBase<Derived>& coeffs)
Brian Silverman72890c22015-09-19 14:37:37 -0400161{ return coeffs.asDiagonal(); }
162
Brian Silverman72890c22015-09-19 14:37:37 -0400163/** \deprecated */
164typedef DiagonalMatrix<float, 2> AlignedScaling2f;
165/** \deprecated */
166typedef DiagonalMatrix<double,2> AlignedScaling2d;
167/** \deprecated */
168typedef DiagonalMatrix<float, 3> AlignedScaling3f;
169/** \deprecated */
170typedef DiagonalMatrix<double,3> AlignedScaling3d;
171//@}
172
173template<typename Scalar>
174template<int Dim>
175inline Transform<Scalar,Dim,Affine>
176UniformScaling<Scalar>::operator* (const Translation<Scalar,Dim>& t) const
177{
178 Transform<Scalar,Dim,Affine> res;
179 res.matrix().setZero();
180 res.linear().diagonal().fill(factor());
181 res.translation() = factor() * t.vector();
182 res(Dim,Dim) = Scalar(1);
183 return res;
184}
185
186} // end namespace Eigen
187
188#endif // EIGEN_SCALING_H