Eigen/src/Geometry/Scaling.h
Go to the documentation of this file.
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 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
36 namespace internal {
37 // This helper helps nvcc+MSVC to properly parse this file.
38 // See bug 1412.
39 template <typename Scalar, int Dim, int Mode>
41  enum { NewMode = int(Mode) == int(Isometry) ? Affine : Mode };
43 };
44 } // namespace internal
45 
46 template <typename Scalar_>
48  public:
50  typedef Scalar_ Scalar;
51 
52  protected:
54 
55  public:
59  explicit inline UniformScaling(const Scalar& s) : m_factor(s) {}
60 
61  inline const Scalar& factor() const { return m_factor; }
62  inline Scalar& factor() { return m_factor; }
63 
65  inline UniformScaling operator*(const UniformScaling& other) const {
66  return UniformScaling(m_factor * other.factor());
67  }
68 
70  template <int Dim>
72 
74  template <int Dim, int Mode, int Options>
78  res.prescale(factor());
79  return res;
80  }
81 
83  // TODO returns an expression
84  template <typename Derived>
86  return other * m_factor;
87  }
88 
89  template <typename Derived, int Dim>
91  return r.toRotationMatrix() * m_factor;
92  }
93 
95  inline UniformScaling inverse() const { return UniformScaling(Scalar(1) / m_factor); }
96 
102  template <typename NewScalarType>
104  return UniformScaling<NewScalarType>(NewScalarType(m_factor));
105  }
106 
108  template <typename OtherScalarType>
109  inline explicit UniformScaling(const UniformScaling<OtherScalarType>& other) {
110  m_factor = Scalar(other.factor());
111  }
112 
117  bool isApprox(const UniformScaling& other,
118  const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const {
119  return internal::isApprox(m_factor, other.factor(), prec);
120  }
121 };
122 
125 
129 // NOTE this operator is defined in MatrixBase and not as a friend function
130 // of UniformScaling to fix an internal crash of Intel's ICC
131 template <typename Derived, typename Scalar>
133 operator*(const MatrixBase<Derived>& matrix, const UniformScaling<Scalar>& s) {
134  return matrix.derived() * s.factor();
135 }
136 
142 template <typename RealScalar>
143 inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s) {
145 }
146 
148 template <typename Scalar>
149 inline DiagonalMatrix<Scalar, 2> Scaling(const Scalar& sx, const Scalar& sy) {
150  return DiagonalMatrix<Scalar, 2>(sx, sy);
151 }
153 template <typename Scalar>
154 inline DiagonalMatrix<Scalar, 3> Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz) {
155  return DiagonalMatrix<Scalar, 3>(sx, sy, sz);
156 }
157 
161 template <typename Derived>
163  return coeffs.asDiagonal();
164 }
165 
167 template <typename Derived>
169  return typename DiagonalWrapper<const Derived>::PlainObject(std::move(coeffs.derived()));
170 }
171 
181 
182 template <typename Scalar>
183 template <int Dim>
186  res.matrix().setZero();
187  res.linear().diagonal().fill(factor());
188  res.translation() = factor() * t.vector();
189  res(Dim, Dim) = Scalar(1);
190  return res;
191 }
192 
193 } // end namespace Eigen
194 
195 #endif // EIGEN_SCALING_H
#define EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(EXPR, SCALAR, OPNAME)
Definition: Macros.h:1200
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
SCALAR Scalar
Definition: bench_gemm.cpp:45
DiagonalMatrix< Scalar, DiagonalVectorType::SizeAtCompileTime, DiagonalVectorType::MaxSizeAtCompileTime > PlainObject
Definition: DiagonalMatrix.h:54
Represents a diagonal matrix with its storage.
Definition: DiagonalMatrix.h:172
Expression of a diagonal matrix.
Definition: DiagonalMatrix.h:320
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
EIGEN_DEVICE_FUNC const DiagonalWrapper< const Derived > asDiagonal() const
Definition: DiagonalMatrix.h:347
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
Common base class for compact rotation representations.
Definition: RotationBase.h:32
Represents an homogeneous transformation in a N dimensional space.
Definition: Transform.h:192
Represents a translation transformation.
Definition: Translation.h:33
Represents a generic uniform scaling transformation.
Definition: Eigen/src/Geometry/Scaling.h:47
Scalar & factor()
Definition: Eigen/src/Geometry/Scaling.h:62
Eigen::internal::plain_matrix_type< Derived >::type operator*(const MatrixBase< Derived > &other) const
Definition: Eigen/src/Geometry/Scaling.h:85
UniformScaling inverse() const
Definition: Eigen/src/Geometry/Scaling.h:95
Scalar_ Scalar
Definition: Eigen/src/Geometry/Scaling.h:50
bool isApprox(const UniformScaling &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Eigen/src/Geometry/Scaling.h:117
Scalar m_factor
Definition: Eigen/src/Geometry/Scaling.h:53
internal::uniformscaling_times_affine_returntype< Scalar, Dim, Mode >::type operator*(const Transform< Scalar, Dim, Mode, Options > &t) const
Definition: Eigen/src/Geometry/Scaling.h:75
UniformScaling< NewScalarType > cast() const
Definition: Eigen/src/Geometry/Scaling.h:103
UniformScaling operator*(const UniformScaling &other) const
Definition: Eigen/src/Geometry/Scaling.h:65
Matrix< Scalar, Dim, Dim > operator*(const RotationBase< Derived, Dim > &r) const
Definition: Eigen/src/Geometry/Scaling.h:90
UniformScaling()
Definition: Eigen/src/Geometry/Scaling.h:57
UniformScaling(const Scalar &s)
Definition: Eigen/src/Geometry/Scaling.h:59
const Scalar & factor() const
Definition: Eigen/src/Geometry/Scaling.h:61
UniformScaling(const UniformScaling< OtherScalarType > &other)
Definition: Eigen/src/Geometry/Scaling.h:109
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
@ Affine
Definition: Constants.h:458
@ Isometry
Definition: Constants.h:455
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
EIGEN_DEVICE_FUNC bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition: MathFunctions.h:1923
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
UniformScaling< float > Scaling(float s)
Definition: Eigen/src/Geometry/Scaling.h:138
DiagonalMatrix< float, 3 > AlignedScaling3f
Definition: Eigen/src/Geometry/Scaling.h:177
DiagonalMatrix< double, 3 > AlignedScaling3d
Definition: Eigen/src/Geometry/Scaling.h:179
DiagonalMatrix< double, 2 > AlignedScaling2d
Definition: Eigen/src/Geometry/Scaling.h:175
DiagonalMatrix< float, 2 > AlignedScaling2f
Definition: Eigen/src/Geometry/Scaling.h:173
static const unsigned Dim
Problem dimension.
Definition: two_d_tilted_square.cc:62
r
Definition: UniformPSDSelfTest.py:20
Definition: Eigen_Colamd.h:49
t
Definition: plotPSD.py:36
void product(const MatrixType &m)
Definition: product.h:42
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: XprHelper.h:389
Definition: Eigen/src/Geometry/Scaling.h:40
Transform< Scalar, Dim, NewMode > type
Definition: Eigen/src/Geometry/Scaling.h:42
@ NewMode
Definition: Eigen/src/Geometry/Scaling.h:41