Scaling.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//
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
32template<typename _Scalar>
33class UniformScaling
34{
35public:
37 typedef _Scalar Scalar;
38
39protected:
40
41 Scalar m_factor;
42
43public:
44
46 UniformScaling() {}
48 explicit inline UniformScaling(const Scalar& s) : m_factor(s) {}
49
50 inline const Scalar& factor() const { return m_factor; }
51 inline Scalar& factor() { return m_factor; }
52
54 inline UniformScaling operator* (const UniformScaling& other) const
55 { return UniformScaling(m_factor * other.factor()); }
56
58 template<int Dim>
59 inline Transform<Scalar,Dim,Affine> operator* (const Translation<Scalar,Dim>& t) const;
60
62 template<int Dim, int Mode, int Options>
63 inline Transform<Scalar,Dim,(int(Mode)==int(Isometry)?Affine:Mode)> operator* (const Transform<Scalar,Dim, Mode, Options>& t) const
64 {
65 Transform<Scalar,Dim,(int(Mode)==int(Isometry)?Affine:Mode)> res = t;
66 res.prescale(factor());
67 return res;
68}
69
71 // TODO returns an expression
72 template<typename Derived>
73 inline typename internal::plain_matrix_type<Derived>::type operator* (const MatrixBase<Derived>& other) const
74 { return other * m_factor; }
75
76 template<typename Derived,int Dim>
77 inline Matrix<Scalar,Dim,Dim> operator*(const RotationBase<Derived,Dim>& r) const
78 { return r.toRotationMatrix() * m_factor; }
79
81 inline UniformScaling inverse() const
82 { return UniformScaling(Scalar(1)/m_factor); }
83
89 template<typename NewScalarType>
90 inline UniformScaling<NewScalarType> cast() const
91 { return UniformScaling<NewScalarType>(NewScalarType(m_factor)); }
92
94 template<typename OtherScalarType>
95 inline explicit UniformScaling(const UniformScaling<OtherScalarType>& other)
96 { m_factor = Scalar(other.factor()); }
97
102 bool isApprox(const UniformScaling& other, typename NumTraits<Scalar>::Real prec = NumTraits<Scalar>::dummy_precision()) const
103 { return internal::isApprox(m_factor, other.factor(), prec); }
104
105};
106
108// NOTE this operator is defiend in MatrixBase and not as a friend function
109// of UniformScaling to fix an internal crash of Intel's ICC
110template<typename Derived> typename MatrixBase<Derived>::ScalarMultipleReturnType
111MatrixBase<Derived>::operator*(const UniformScaling<Scalar>& s) const
112{ return derived() * s.factor(); }
113
115static inline UniformScaling<float> Scaling(float s) { return UniformScaling<float>(s); }
117static inline UniformScaling<double> Scaling(double s) { return UniformScaling<double>(s); }
119template<typename RealScalar>
120static inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s)
121{ return UniformScaling<std::complex<RealScalar> >(s); }
122
124template<typename Scalar>
125static inline DiagonalMatrix<Scalar,2> Scaling(Scalar sx, Scalar sy)
126{ return DiagonalMatrix<Scalar,2>(sx, sy); }
128template<typename Scalar>
129static inline DiagonalMatrix<Scalar,3> Scaling(Scalar sx, Scalar sy, Scalar sz)
130{ return DiagonalMatrix<Scalar,3>(sx, sy, sz); }
131
135template<typename Derived>
136static inline const DiagonalWrapper<const Derived> Scaling(const MatrixBase<Derived>& coeffs)
137{ return coeffs.asDiagonal(); }
138
141
142typedef DiagonalMatrix<float, 2> AlignedScaling2f;
144typedef DiagonalMatrix<double,2> AlignedScaling2d;
146typedef DiagonalMatrix<float, 3> AlignedScaling3f;
148typedef DiagonalMatrix<double,3> AlignedScaling3d;
150
151template<typename Scalar>
152template<int Dim>
154UniformScaling<Scalar>::operator* (const Translation<Scalar,Dim>& t) const
155{
157 res.matrix().setZero();
158 res.linear().diagonal().fill(factor());
159 res.translation() = factor() * t.vector();
160 res(Dim,Dim) = Scalar(1);
161 return res;
162}
163
164} // end namespace Eigen
165
166#endif // EIGEN_SCALING_H
Represents a diagonal matrix with its storage.
Definition DiagonalMatrix.h:132
Expression of a diagonal matrix.
Definition DiagonalMatrix.h:245
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
const ScalarMultipleReturnType operator*(const Scalar &scalar) const
Definition MatrixBase.h:50
Represents an homogeneous transformation in a N dimensional space.
Definition Transform.h:177
const MatrixType & matrix() const
Definition Transform.h:367
ConstLinearPart linear() const
Definition Transform.h:372
Represents a translation transformation.
Definition Translation.h:31
Represents a generic uniform scaling transformation.
@ Affine
Definition Constants.h:387
@ Isometry
Definition Constants.h:384