Please, help us to better know about our user community by answering the following short survey: https://forms.gle/wpyrxWi18ox9Z5ae9
Eigen  3.4.0
 
Loading...
Searching...
No Matches
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
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
48template<typename _Scalar>
50{
51public:
53 typedef _Scalar Scalar;
54
55protected:
56
57 Scalar m_factor;
58
59public:
60
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
70 inline UniformScaling operator* (const UniformScaling& other) const
71 { return UniformScaling(m_factor * other.factor()); }
72
74 template<int Dim>
76
78 template<int Dim, int Mode, int Options>
79 inline typename
82 {
84 res.prescale(factor());
85 return res;
86 }
87
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>
96 { return r.toRotationMatrix() * m_factor; }
97
99 inline UniformScaling inverse() const
100 { return UniformScaling(Scalar(1)/m_factor); }
101
107 template<typename NewScalarType>
109 { return UniformScaling<NewScalarType>(NewScalarType(m_factor)); }
110
112 template<typename OtherScalarType>
114 { m_factor = Scalar(other.factor()); }
115
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
127
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
133template<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(); }
137
143template<typename RealScalar>
144inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s)
146
148template<typename Scalar>
149inline DiagonalMatrix<Scalar,2> Scaling(const Scalar& sx, const Scalar& sy)
150{ return DiagonalMatrix<Scalar,2>(sx, sy); }
152template<typename Scalar>
153inline DiagonalMatrix<Scalar,3> Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz)
154{ return DiagonalMatrix<Scalar,3>(sx, sy, sz); }
155
159template<typename Derived>
161{ return coeffs.asDiagonal(); }
162
172
173template<typename Scalar>
174template<int Dim>
177{
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
Represents a diagonal matrix with its storage.
Definition: DiagonalMatrix.h:142
Expression of a diagonal matrix.
Definition: DiagonalMatrix.h:295
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
const DiagonalWrapper< const Derived > asDiagonal() const
Definition: DiagonalMatrix.h:325
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:562
Common base class for compact rotation representations.
Definition: RotationBase.h:30
RotationMatrixType toRotationMatrix() const
Definition: RotationBase.h:45
Represents an homogeneous transformation in a N dimensional space.
Definition: Transform.h:205
const MatrixType & matrix() const
Definition: Transform.h:389
ConstTranslationPart translation() const
Definition: Transform.h:404
ConstLinearPart linear() const
Definition: Transform.h:394
Represents a translation transformation.
Definition: Translation.h:31
Represents a generic uniform scaling transformation.
Definition: Scaling.h:50
_Scalar Scalar
Definition: Scaling.h:53
UniformScaling< NewScalarType > cast() const
Definition: Scaling.h:108
UniformScaling(const Scalar &s)
Definition: Scaling.h:64
UniformScaling inverse() const
Definition: Scaling.h:99
bool isApprox(const UniformScaling &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Scaling.h:120
UniformScaling(const UniformScaling< OtherScalarType > &other)
Definition: Scaling.h:113
UniformScaling()
Definition: Scaling.h:62
UniformScaling operator*(const UniformScaling &other) const
Definition: Scaling.h:70
@ Affine
Definition: Constants.h:460
@ Isometry
Definition: Constants.h:457
Namespace containing all symbols from the Eigen library.
Definition: Core:141
UniformScaling< float > Scaling(float s)
Definition: Scaling.h:139
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:233