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 namespace Eigen {
14 
32 template<typename _Scalar>
33 class UniformScaling
34 {
35 public:
37  typedef _Scalar Scalar;
38 
39 protected:
40 
42 
43 public:
44 
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>
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>
74  { return other * m_factor; }
75 
76  template<typename Derived,int Dim>
78  { return r.toRotationMatrix() * m_factor; }
79 
81  inline UniformScaling inverse() const
82  { return UniformScaling(Scalar(1)/m_factor); }
83 
89  template<typename NewScalarType>
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>
117 operator*(const MatrixBase<Derived>& matrix, const UniformScaling<Scalar>& s)
118 { return matrix.derived() * s.factor(); }
119 
125 template<typename RealScalar>
126 inline UniformScaling<std::complex<RealScalar> > Scaling(const 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>
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::UniformScaling::Scalar
_Scalar Scalar
Definition: Eigen/src/Geometry/Scaling.h:37
Eigen
Definition: common.h:73
Eigen::Transform::linear
EIGEN_DEVICE_FUNC ConstLinearPart linear() const
Definition: Transform.h:400
Eigen::DiagonalMatrix
Represents a diagonal matrix with its storage.
Definition: DiagonalMatrix.h:116
Eigen::AlignedScaling2f
DiagonalMatrix< float, 2 > AlignedScaling2f
Definition: Eigen/src/Geometry/Scaling.h:146
Eigen::Transform
Represents an homogeneous transformation in a N dimensional space.
Definition: ForwardDeclarations.h:274
s
RealScalar s
Definition: level1_cplx_impl.h:104
Eigen::Affine
@ Affine
Definition: Constants.h:450
Eigen::UniformScaling::UniformScaling
UniformScaling(const Scalar &s)
Definition: Eigen/src/Geometry/Scaling.h:48
Eigen::Scaling
UniformScaling< float > Scaling(float s)
Definition: Eigen/src/Geometry/Scaling.h:121
Eigen::UniformScaling::factor
Scalar & factor()
Definition: Eigen/src/Geometry/Scaling.h:51
Eigen::internal::isApprox
EIGEN_DEVICE_FUNC bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition: Eigen/src/Core/MathFunctions.h:1361
Eigen::UniformScaling::cast
UniformScaling< NewScalarType > cast() const
Definition: Eigen/src/Geometry/Scaling.h:90
Scalar
SCALAR Scalar
Definition: common.h:84
Eigen::UniformScaling::UniformScaling
UniformScaling()
Definition: Eigen/src/Geometry/Scaling.h:46
Eigen::AlignedScaling3d
DiagonalMatrix< double, 3 > AlignedScaling3d
Definition: Eigen/src/Geometry/Scaling.h:152
Eigen::UniformScaling::operator*
UniformScaling operator*(const UniformScaling &other) const
Definition: Eigen/src/Geometry/Scaling.h:54
Eigen::internal::plain_matrix_type
Definition: XprHelper.h:275
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:102
Eigen::AlignedScaling3f
DiagonalMatrix< float, 3 > AlignedScaling3f
Definition: Eigen/src/Geometry/Scaling.h:150
Eigen::RotationBase::toRotationMatrix
EIGEN_DEVICE_FUNC RotationMatrixType toRotationMatrix() const
Definition: RotationBase.h:45
Eigen::UniformScaling::factor
const Scalar & factor() const
Definition: Eigen/src/Geometry/Scaling.h:50
Eigen::Isometry
@ Isometry
Definition: Constants.h:447
Eigen::UniformScaling::UniformScaling
UniformScaling(const UniformScaling< OtherScalarType > &other)
Definition: Eigen/src/Geometry/Scaling.h:95
Eigen::RotationBase
Common base class for compact rotation representations.
Definition: ForwardDeclarations.h:266
Eigen::UniformScaling::isApprox
bool isApprox(const UniformScaling &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Eigen/src/Geometry/Scaling.h:102
Eigen::MatrixBase::asDiagonal
const EIGEN_DEVICE_FUNC DiagonalWrapper< const Derived > asDiagonal() const
Definition: DiagonalMatrix.h:277
Eigen::Translation::vector
const EIGEN_DEVICE_FUNC VectorType & vector() const
Definition: Translation.h:87
Eigen::Transform::prescale
EIGEN_DEVICE_FUNC Transform & prescale(const MatrixBase< OtherDerived > &other)
Eigen::UniformScaling::operator*
Matrix< Scalar, Dim, Dim > operator*(const RotationBase< Derived, Dim > &r) const
Definition: Eigen/src/Geometry/Scaling.h:77
Eigen::DiagonalWrapper
Expression of a diagonal matrix.
Definition: DiagonalMatrix.h:245
int
return int(ret)+1
Eigen::UniformScaling::inverse
UniformScaling inverse() const
Definition: Eigen/src/Geometry/Scaling.h:81
Eigen::Transform::matrix
const EIGEN_DEVICE_FUNC MatrixType & matrix() const
Definition: Transform.h:395
Eigen::Transform::translation
EIGEN_DEVICE_FUNC ConstTranslationPart translation() const
Definition: Transform.h:410
Eigen::UniformScaling
Definition: ForwardDeclarations.h:277
Eigen::UniformScaling::m_factor
Scalar m_factor
Definition: Eigen/src/Geometry/Scaling.h:41
Eigen::AlignedScaling2d
DiagonalMatrix< double, 2 > AlignedScaling2d
Definition: Eigen/src/Geometry/Scaling.h:148
Eigen::Translation
Represents a translation transformation.
Definition: ForwardDeclarations.h:271
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:178
Eigen::PlainObjectBase::setZero
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:515
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Eigen::EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE
const EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Derived, typename Derived::Scalar, pow) pow(const Eigen
Definition: GlobalFunctions.h:113
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:150


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Wed Mar 2 2022 00:06:10