xref: /aosp_15_r20/external/eigen/Eigen/src/Geometry/Scaling.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
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 //
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 
13 namespace Eigen {
14 
15 /** \geometry_module \ingroup Geometry_Module
16   *
17   * \class UniformScaling
18   *
19   * \brief Represents a generic uniform scaling transformation
20   *
21   * \tparam _Scalar the scalar type, i.e., the type of the coefficients.
22   *
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   */
32 
33 namespace 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 
48 template<typename _Scalar>
49 class UniformScaling
50 {
51 public:
52   /** the scalar type of the coefficients */
53   typedef _Scalar Scalar;
54 
55 protected:
56 
57   Scalar m_factor;
58 
59 public:
60 
61   /** Default constructor without initialization. */
UniformScaling()62   UniformScaling() {}
63   /** Constructs and initialize a uniform scaling transformation */
UniformScaling(const Scalar & s)64   explicit inline UniformScaling(const Scalar& s) : m_factor(s) {}
65 
factor()66   inline const Scalar& factor() const { return m_factor; }
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>
79   inline typename
80 	internal::uniformscaling_times_affine_returntype<Scalar,Dim,Mode>::type
81 	operator* (const Transform<Scalar, Dim, Mode, Options>& t) const
82   {
83     typename internal::uniformscaling_times_affine_returntype<Scalar,Dim,Mode>::type res = t;
84     res.prescale(factor());
85     return res;
86   }
87 
88   /** Concatenates a uniform scaling and a linear transformation matrix */
89   // TODO returns an expression
90   template<typename Derived>
91   inline typename Eigen::internal::plain_matrix_type<Derived>::type operator* (const MatrixBase<Derived>& other) const
92   { 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 */
inverse()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>
cast()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>
UniformScaling(const UniformScaling<OtherScalarType> & other)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 
125 /** \addtogroup Geometry_Module */
126 //@{
127 
128 /** Concatenates a linear transformation matrix and a uniform scaling
129   * \relates UniformScaling
130   */
131 // NOTE this operator is defined in MatrixBase and not as a friend function
132 // of UniformScaling to fix an internal crash of Intel's ICC
133 template<typename Derived,typename Scalar>
EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Derived,Scalar,product)134 EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Derived,Scalar,product)
135 operator*(const MatrixBase<Derived>& matrix, const UniformScaling<Scalar>& s)
136 { return matrix.derived() * s.factor(); }
137 
138 /** Constructs a uniform scaling from scale factor \a s */
Scaling(float s)139 inline UniformScaling<float> Scaling(float s) { return UniformScaling<float>(s); }
140 /** Constructs a uniform scaling from scale factor \a s */
Scaling(double s)141 inline UniformScaling<double> Scaling(double s) { return UniformScaling<double>(s); }
142 /** Constructs a uniform scaling from scale factor \a s */
143 template<typename RealScalar>
Scaling(const std::complex<RealScalar> & s)144 inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s)
145 { return UniformScaling<std::complex<RealScalar> >(s); }
146 
147 /** Constructs a 2D axis aligned scaling */
148 template<typename Scalar>
Scaling(const Scalar & sx,const Scalar & sy)149 inline DiagonalMatrix<Scalar,2> Scaling(const Scalar& sx, const Scalar& sy)
150 { return DiagonalMatrix<Scalar,2>(sx, sy); }
151 /** Constructs a 3D axis aligned scaling */
152 template<typename Scalar>
Scaling(const Scalar & sx,const Scalar & sy,const Scalar & sz)153 inline DiagonalMatrix<Scalar,3> Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz)
154 { 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   */
159 template<typename Derived>
Scaling(const MatrixBase<Derived> & coeffs)160 inline const DiagonalWrapper<const Derived> Scaling(const MatrixBase<Derived>& coeffs)
161 { return coeffs.asDiagonal(); }
162 
163 /** \deprecated */
164 typedef DiagonalMatrix<float, 2> AlignedScaling2f;
165 /** \deprecated */
166 typedef DiagonalMatrix<double,2> AlignedScaling2d;
167 /** \deprecated */
168 typedef DiagonalMatrix<float, 3> AlignedScaling3f;
169 /** \deprecated */
170 typedef DiagonalMatrix<double,3> AlignedScaling3d;
171 //@}
172 
173 template<typename Scalar>
174 template<int Dim>
175 inline Transform<Scalar,Dim,Affine>
176 UniformScaling<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
189