Please, help us to better know about our user community by answering the following short survey: https://forms.gle/wpyrxWi18ox9Z5ae9
Eigen  3.3.8
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 
13 namespace Eigen {
14 
32 template<typename _Scalar>
33 class UniformScaling
34 {
35 public:
37  typedef _Scalar Scalar;
38 
39 protected:
40 
41  Scalar m_factor;
42 
43 public:
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, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
103  { return internal::isApprox(m_factor, other.factor(), prec); }
104 
105 };
106 
109 
113 // NOTE this operator is defiend in MatrixBase and not as a friend function
114 // of UniformScaling to fix an internal crash of Intel's ICC
115 template<typename Derived,typename Scalar>
116 EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Derived,Scalar,product)
117 operator*(const MatrixBase<Derived>& matrix, const UniformScaling<Scalar>& s)
118 { return matrix.derived() * s.factor(); }
119 
121 inline UniformScaling<float> Scaling(float s) { return UniformScaling<float>(s); }
123 inline UniformScaling<double> Scaling(double s) { return UniformScaling<double>(s); }
125 template<typename RealScalar>
126 inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s)
127 { return UniformScaling<std::complex<RealScalar> >(s); }
128 
130 template<typename Scalar>
131 inline DiagonalMatrix<Scalar,2> Scaling(const Scalar& sx, const Scalar& sy)
132 { return DiagonalMatrix<Scalar,2>(sx, sy); }
134 template<typename Scalar>
135 inline DiagonalMatrix<Scalar,3> Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz)
136 { return DiagonalMatrix<Scalar,3>(sx, sy, sz); }
137 
141 template<typename Derived>
143 { return coeffs.asDiagonal(); }
144 
154 
155 template<typename Scalar>
156 template<int Dim>
158 UniformScaling<Scalar>::operator* (const Translation<Scalar,Dim>& t) const
159 {
161  res.matrix().setZero();
162  res.linear().diagonal().fill(factor());
163  res.translation() = factor() * t.vector();
164  res(Dim,Dim) = Scalar(1);
165  return res;
166 }
167 
168 } // end namespace Eigen
169 
170 #endif // EIGEN_SCALING_H
Eigen
Namespace containing all symbols from the Eigen library.
Definition: Core:309
Eigen::DiagonalMatrix
Represents a diagonal matrix with its storage.
Definition: DiagonalMatrix.h:116
Eigen::AlignedScaling2f
DiagonalMatrix< float, 2 > AlignedScaling2f
Definition: Scaling.h:146
Eigen::Transform
Represents an homogeneous transformation in a N dimensional space.
Definition: ForwardDeclarations.h:270
Eigen::Affine
@ Affine
Definition: Constants.h:450
Eigen::Scaling
UniformScaling< float > Scaling(float s)
Definition: Scaling.h:121
Eigen::AlignedScaling3d
DiagonalMatrix< double, 3 > AlignedScaling3d
Definition: Scaling.h:152
Eigen::Transform::linear
ConstLinearPart linear() const
Definition: Transform.h:400
Eigen::AlignedScaling3f
DiagonalMatrix< float, 3 > AlignedScaling3f
Definition: Scaling.h:150
Eigen::inverse
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_inverse_op< typename Derived::Scalar >, const Derived > inverse(const Eigen::ArrayBase< Derived > &x)
Eigen::Isometry
@ Isometry
Definition: Constants.h:447
Eigen::PlainObjectBase::setZero
Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:515
Eigen::MatrixBase::asDiagonal
const DiagonalWrapper< const Derived > asDiagonal() const
Definition: DiagonalMatrix.h:277
Eigen::DiagonalWrapper
Expression of a diagonal matrix.
Definition: DiagonalMatrix.h:245
Eigen::Transform::translation
ConstTranslationPart translation() const
Definition: Transform.h:410
Eigen::AlignedScaling2d
DiagonalMatrix< double, 2 > AlignedScaling2d
Definition: Scaling.h:148
Eigen::Translation
Represents a translation transformation.
Definition: ForwardDeclarations.h:267
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Eigen::operator*
const Product< MatrixDerived, PermutationDerived, AliasFreeProduct > operator*(const MatrixBase< MatrixDerived > &matrix, const PermutationBase< PermutationDerived > &permutation)
Definition: PermutationMatrix.h:543
Eigen::Transform::matrix
const MatrixType & matrix() const
Definition: Transform.h:395