Eigne Matrix Class Library

Dependents:   Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more

Embed: (wiki syntax)

« Back to documentation index

Show/hide line numbers Scaling.h Source File

Scaling.h

00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_SCALING_H
00011 #define EIGEN_SCALING_H
00012 
00013 namespace Eigen { 
00014 
00015 /** \geometry_module \ingroup Geometry_Module
00016   *
00017   * \class Scaling
00018   *
00019   * \brief Represents a generic uniform scaling transformation
00020   *
00021   * \param _Scalar the scalar type, i.e., the type of the coefficients.
00022   *
00023   * This class represent a uniform scaling transformation. It is the return
00024   * type of Scaling(Scalar), and most of the time this is the only way it
00025   * is used. In particular, this class is not aimed to be used to store a scaling transformation,
00026   * but rather to make easier the constructions and updates of Transform objects.
00027   *
00028   * To represent an axis aligned scaling, use the DiagonalMatrix class.
00029   *
00030   * \sa Scaling(), class DiagonalMatrix, MatrixBase::asDiagonal(), class Translation, class Transform
00031   */
00032 template<typename _Scalar>
00033 class UniformScaling
00034 {
00035 public:
00036   /** the scalar type of the coefficients */
00037   typedef _Scalar Scalar;
00038 
00039 protected:
00040 
00041   Scalar m_factor;
00042 
00043 public:
00044 
00045   /** Default constructor without initialization. */
00046   UniformScaling() {}
00047   /** Constructs and initialize a uniform scaling transformation */
00048   explicit inline UniformScaling(const Scalar& s) : m_factor(s) {}
00049 
00050   inline const Scalar& factor() const { return m_factor; }
00051   inline Scalar& factor() { return m_factor; }
00052 
00053   /** Concatenates two uniform scaling */
00054   inline UniformScaling operator*  (const UniformScaling& other) const
00055   { return UniformScaling(m_factor * other.factor()); }
00056 
00057   /** Concatenates a uniform scaling and a translation */
00058   template<int Dim>
00059   inline Transform<Scalar,Dim,Affine> operator*  (const Translation<Scalar,Dim>& t) const;
00060 
00061   /** Concatenates a uniform scaling and an affine transformation */
00062   template<int Dim, int Mode, int Options>
00063   inline Transform<Scalar,Dim,(int(Mode)==int(Isometry)?Affine:Mode)> operator*  (const Transform<Scalar,Dim, Mode, Options>& t) const
00064   {
00065    Transform<Scalar,Dim,(int(Mode)==int(Isometry)?Affine:Mode)> res = t;
00066    res.prescale(factor());
00067    return res;
00068 }
00069 
00070   /** Concatenates a uniform scaling and a linear transformation matrix */
00071   // TODO returns an expression
00072   template<typename Derived>
00073   inline typename internal::plain_matrix_type<Derived>::type operator*  (const MatrixBase<Derived>& other) const
00074   { return other * m_factor; }
00075 
00076   template<typename Derived,int Dim>
00077   inline Matrix<Scalar,Dim,Dim> operator* (const RotationBase<Derived,Dim>& r) const
00078   { return r.toRotationMatrix() * m_factor; }
00079 
00080   /** \returns the inverse scaling */
00081   inline UniformScaling inverse() const
00082   { return UniformScaling(Scalar(1)/m_factor); }
00083 
00084   /** \returns \c *this with scalar type casted to \a NewScalarType
00085     *
00086     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
00087     * then this function smartly returns a const reference to \c *this.
00088     */
00089   template<typename NewScalarType>
00090   inline UniformScaling<NewScalarType> cast() const
00091   { return UniformScaling<NewScalarType>(NewScalarType(m_factor)); }
00092 
00093   /** Copy constructor with scalar type conversion */
00094   template<typename OtherScalarType>
00095   inline explicit UniformScaling(const UniformScaling<OtherScalarType>& other)
00096   { m_factor = Scalar(other.factor()); }
00097 
00098   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
00099     * determined by \a prec.
00100     *
00101     * \sa MatrixBase::isApprox() */
00102   bool isApprox(const UniformScaling& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
00103   { return internal::isApprox(m_factor, other.factor(), prec); }
00104 
00105 };
00106 
00107 /** Concatenates a linear transformation matrix and a uniform scaling */
00108 // NOTE this operator is defiend in MatrixBase and not as a friend function
00109 // of UniformScaling to fix an internal crash of Intel's ICC
00110 template<typename Derived> typename MatrixBase<Derived>::ScalarMultipleReturnType
00111 MatrixBase<Derived>::operator* (const UniformScaling<Scalar>& s) const
00112 { return derived() * s.factor(); }
00113 
00114 /** Constructs a uniform scaling from scale factor \a s */
00115 static inline UniformScaling<float> Scaling(float s) { return UniformScaling<float>(s); }
00116 /** Constructs a uniform scaling from scale factor \a s */
00117 static inline UniformScaling<double> Scaling(double s) { return UniformScaling<double>(s); }
00118 /** Constructs a uniform scaling from scale factor \a s */
00119 template<typename RealScalar>
00120 static inline UniformScaling<std::complex<RealScalar> > Scaling(const std::complex<RealScalar>& s)
00121 { return UniformScaling<std::complex<RealScalar> >(s); }
00122 
00123 /** Constructs a 2D axis aligned scaling */
00124 template<typename Scalar>
00125 static inline DiagonalMatrix<Scalar,2> Scaling(const Scalar& sx, const Scalar& sy)
00126 { return DiagonalMatrix<Scalar,2>(sx, sy); }
00127 /** Constructs a 3D axis aligned scaling */
00128 template<typename Scalar>
00129 static inline DiagonalMatrix<Scalar,3> Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz)
00130 { return DiagonalMatrix<Scalar,3>(sx, sy, sz); }
00131 
00132 /** Constructs an axis aligned scaling expression from vector expression \a coeffs
00133   * This is an alias for coeffs.asDiagonal()
00134   */
00135 template<typename Derived>
00136 static inline const DiagonalWrapper<const Derived> Scaling(const MatrixBase<Derived>& coeffs)
00137 { return coeffs.asDiagonal (); }
00138 
00139 /** \addtogroup Geometry_Module */
00140 //@{
00141 /** \deprecated */
00142 typedef DiagonalMatrix<float, 2> AlignedScaling2f ;
00143 /** \deprecated */
00144 typedef DiagonalMatrix<double,2> AlignedScaling2d ;
00145 /** \deprecated */
00146 typedef DiagonalMatrix<float, 3> AlignedScaling3f ;
00147 /** \deprecated */
00148 typedef DiagonalMatrix<double,3> AlignedScaling3d ;
00149 //@}
00150 
00151 template<typename Scalar>
00152 template<int Dim>
00153 inline Transform<Scalar,Dim,Affine> 
00154 UniformScaling<Scalar>::operator*  (const Translation<Scalar,Dim> & t) const
00155 {
00156   Transform<Scalar,Dim,Affine>  res;
00157   res.matrix ().setZero();
00158   res.linear ().diagonal().fill(factor());
00159   res.translation () = factor() * t.vector();
00160   res(Dim,Dim) = Scalar(1);
00161   return res;
00162 }
00163 
00164 } // end namespace Eigen
00165 
00166 #endif // EIGEN_SCALING_H