Add Scaling and Translation class as discussed on ML, still missing:

* handling Quaternion, AngleAxis and Rotation2D, 2 options here:
 1- make all of them inheriting a common base class Rotation such that we can
    have a single version of operator* for all the rotation type (they all get converted to a matrix)
 2- write a version for all type (so 3 rotations types * 3 for Transform,Translation and Scaling)
* real documentation
This commit is contained in:
Gael Guennebaud 2008-08-30 00:08:23 +00:00
parent 13a9d93bc0
commit 9e7a9cde14
6 changed files with 429 additions and 16 deletions

View File

@ -31,6 +31,8 @@ namespace Eigen {
#include "src/Geometry/AngleAxis.h"
#include "src/Geometry/Rotation.h"
#include "src/Geometry/Transform.h"
#include "src/Geometry/Translation.h"
#include "src/Geometry/Scaling.h"
#include "src/Geometry/Hyperplane.h"
} // namespace Eigen

View File

@ -105,5 +105,7 @@ template<typename Scalar> class Quaternion;
template<typename Scalar> class Rotation2D;
template<typename Scalar> class AngleAxis;
template<typename Scalar,int Dim> class Transform;
template<typename Scalar,int Dim> class Translation;
template<typename Scalar,int Dim> class Scaling;
#endif // EIGEN_FORWARDDECLARATIONS_H

View File

@ -0,0 +1,148 @@
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#ifndef EIGEN_SCALING_H
#define EIGEN_SCALING_H
/** \geometry_module \ingroup GeometryModule
*
* \class Scaling
*
* \brief Represents a possibly non uniform scaling transformation
*
* \param _Scalar the scalar type, i.e., the type of the coefficients.
* \param _Dim the dimension of the space, can be a compile time value or Dynamic
*
*
* \sa class Translate, class Transform
*/
template<typename _Scalar, int _Dim>
class Scaling
{
public:
enum { Dim = _Dim };
/** the scalar type of the coefficients */
typedef _Scalar Scalar;
typedef Matrix<Scalar,Dim,Dim> LinearMatrixType;
typedef Matrix<Scalar,Dim,1> VectorType;
typedef Translation<Scalar,Dim> TranslationType;
typedef Transform<Scalar,Dim> TransformType;
protected:
VectorType m_coeffs;
public:
/** Default constructor without initialization. */
Scaling() {}
/** Constructs and initialize a uniform scaling transformation */
explicit inline Scaling(const Scalar& s) { m_coeffs.setConstant(s); }
/** 2D only */
inline Scaling(const Scalar& sx, const Scalar& sy)
{
ei_assert(Dim==2);
m_coeffs.x() = sx;
m_coeffs.y() = sy;
}
/** 3D only */
inline Scaling(const Scalar& sx, const Scalar& sy, const Scalar& sz)
{
ei_assert(Dim==3);
m_coeffs.x() = sx;
m_coeffs.y() = sy;
m_coeffs.z() = sz;
}
/** Constructs and initialize the scaling transformation from a vector of scaling coefficients */
explicit inline Scaling(const VectorType& coeffs) : m_coeffs(coeffs) {}
const VectorType& coeffs() const { return m_coeffs; }
VectorType& coeffs() { return m_coeffs; }
/** Concatenates two scaling */
inline Scaling operator* (const Scaling& other) const
{ return Scaling(coeffs().cwise() * other.coeffs()); }
/** Concatenates a scaling and a translation */
inline TransformType operator* (const TranslationType& t) const;
/** Concatenates a scaling and an affine transformation */
inline TransformType operator* (const TransformType& t) const;
/** Concatenates a scaling and a linear transformation matrix */
// TODO returns an expression
inline LinearMatrixType operator* (const LinearMatrixType& other) const
{ return coeffs().asDiagonal() * other; }
/** Concatenates a linear transformation matrix and a scaling */
// TODO returns an expression
friend inline LinearMatrixType operator* (const LinearMatrixType& other, const Scaling& s)
{ return other * s.coeffs().asDiagonal(); }
/** Applies scaling to vector */
inline VectorType operator* (const VectorType& other) const
{ return coeffs().asDiagonal() * other; }
/** \returns the inverse scaling */
inline Scaling inverse() const
{ return Scaling(coeffs.cwise().inverse()); }
inline Scaling& operator=(const Scaling& other)
{
m_coeffs = other.m_coeffs;
return *this;
}
};
/** \addtogroup GeometryModule */
//@{
typedef Scaling<float, 2> Scaling2f;
typedef Scaling<double,2> Scaling2d;
typedef Scaling<float, 3> Scaling3f;
typedef Scaling<double,3> Scaling3d;
//@}
template<typename Scalar, int Dim>
inline typename Scaling<Scalar,Dim>::TransformType
Scaling<Scalar,Dim>::operator* (const TranslationType& t) const
{
TransformType res;
res.matrix().setZero();
res.linear().diagonal() = coeffs();
res.translation() = m_coeffs.cwise() * t.vector();
res(Dim,Dim) = Scalar(1);
return res;
}
template<typename Scalar, int Dim>
inline typename Scaling<Scalar,Dim>::TransformType
Scaling<Scalar,Dim>::operator* (const TransformType& t) const
{
TransformType res = t;
res.prescale(m_coeffs);
return res;
}
#endif // EIGEN_SCALING_H

View File

@ -73,6 +73,8 @@ public:
typedef Matrix<Scalar,Dim,1> VectorType;
/** type of a read/write reference to the translation part of the rotation */
typedef Block<MatrixType,Dim,1> TranslationPart;
typedef Translation<Scalar,Dim> TranslationType;
typedef Scaling<Scalar,Dim> ScalingType;
protected:
@ -81,7 +83,7 @@ protected:
public:
/** Default constructor without initialization of the coefficients. */
Transform() { }
inline Transform() { }
inline Transform(const Transform& other)
{ m_matrix = other.m_matrix; }
@ -129,10 +131,10 @@ public:
/** shortcut for m_matrix(row,col);
* \sa MatrixBase::operaror(int,int) const */
Scalar operator() (int row, int col) const { return m_matrix(row,col); }
inline Scalar operator() (int row, int col) const { return m_matrix(row,col); }
/** shortcut for m_matrix(row,col);
* \sa MatrixBase::operaror(int,int) */
Scalar& operator() (int row, int col) { return m_matrix(row,col); }
inline Scalar& operator() (int row, int col) { return m_matrix(row,col); }
/** \returns a read-only expression of the transformation matrix */
inline const MatrixType& matrix() const { return m_matrix; }
@ -158,12 +160,12 @@ public:
*/
// note: this function is defined here because some compilers cannot find the respective declaration
template<typename OtherDerived>
const typename ei_transform_product_impl<OtherDerived,_Dim,_Dim+1>::ResultType
inline const typename ei_transform_product_impl<OtherDerived,_Dim,_Dim+1>::ResultType
operator * (const MatrixBase<OtherDerived> &other) const
{ return ei_transform_product_impl<OtherDerived,Dim,HDim>::run(*this,other.derived()); }
/** Contatenates two transformations */
const typename ProductReturnType<MatrixType,MatrixType>::Type
inline const typename ProductReturnType<MatrixType,MatrixType>::Type
operator * (const Transform& other) const
{ return m_matrix * other.matrix(); }
@ -171,26 +173,38 @@ public:
void setIdentity() { m_matrix.setIdentity(); }
template<typename OtherDerived>
Transform& scale(const MatrixBase<OtherDerived> &other);
inline Transform& scale(const MatrixBase<OtherDerived> &other);
template<typename OtherDerived>
Transform& prescale(const MatrixBase<OtherDerived> &other);
inline Transform& prescale(const MatrixBase<OtherDerived> &other);
template<typename OtherDerived>
Transform& translate(const MatrixBase<OtherDerived> &other);
inline Transform& translate(const MatrixBase<OtherDerived> &other);
template<typename OtherDerived>
Transform& pretranslate(const MatrixBase<OtherDerived> &other);
inline Transform& pretranslate(const MatrixBase<OtherDerived> &other);
template<typename RotationType>
Transform& rotate(const RotationType& rotation);
inline Transform& rotate(const RotationType& rotation);
template<typename RotationType>
Transform& prerotate(const RotationType& rotation);
inline Transform& prerotate(const RotationType& rotation);
Transform& shear(Scalar sx, Scalar sy);
Transform& preshear(Scalar sx, Scalar sy);
inline Transform& operator*=(const TranslationType& t) { return translate(t.vector()); }
inline Transform operator*(const TranslationType& t) const;
inline Transform& operator*=(const ScalingType& s) { return scale(s.coeffs()); }
inline Transform operator*(const ScalingType& s) const;
friend inline Transform operator*(const LinearMatrixType& mat, const Transform& t)
{
Transform res = t;
res.matrix().row(Dim) = t.matrix().row(Dim);
res.matrix().template block<Dim,HDim>(0,0) = (mat * t.matrix().template block<Dim,HDim>(0,0)).lazy();
return res;
}
LinearMatrixType extractRotation() const;
LinearMatrixType extractRotationNoShear() const;
@ -385,6 +399,22 @@ Transform<Scalar,Dim>::preshear(Scalar sx, Scalar sy)
return *this;
}
template<typename Scalar, int Dim>
inline Transform<Scalar,Dim> Transform<Scalar,Dim>::operator*(const TranslationType& t) const
{
Transform res = *this;
res.translate(t.vector());
return res;
}
template<typename Scalar, int Dim>
inline Transform<Scalar,Dim> Transform<Scalar,Dim>::operator*(const ScalingType& s) const
{
Transform res = *this;
res.scale(s.coeffs());
return res;
}
/** \returns the rotation part of the transformation using a QR decomposition.
* \sa extractRotationNoShear(), class QR
*/
@ -438,6 +468,22 @@ struct ei_transform_product_impl<Other,Dim,HDim, HDim,HDim>
{ return tr.matrix() * other; }
};
template<typename Other, int Dim, int HDim>
struct ei_transform_product_impl<Other,Dim,HDim, Dim,Dim>
{
typedef Transform<typename Other::Scalar,Dim> TransformType;
typedef typename TransformType::MatrixType MatrixType;
typedef TransformType ResultType;
static ResultType run(const TransformType& tr, const Other& other)
{
TransformType res;
res.translation() = tr.translation();
res.matrix().row(Dim) = tr.matrix().row(Dim);
res.linear() = (tr.linear() * other).lazy();
return res;
}
};
template<typename Other, int Dim, int HDim>
struct ei_transform_product_impl<Other,Dim,HDim, HDim,1>
{

View File

@ -0,0 +1,165 @@
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#ifndef EIGEN_TRANSLATION_H
#define EIGEN_TRANSLATION_H
/** \geometry_module \ingroup GeometryModule
*
* \class Translation
*
* \brief Represents a translation transformation
*
* \param _Scalar the scalar type, i.e., the type of the coefficients.
* \param _Dim the dimension of the space, can be a compile time value or Dynamic
*
*
* \sa class Scaling, class Transform
*/
template<typename _Scalar, int _Dim>
class Translation
{
public:
enum { Dim = _Dim };
/** the scalar type of the coefficients */
typedef _Scalar Scalar;
typedef Matrix<Scalar,Dim,Dim> LinearMatrixType;
typedef Matrix<Scalar,Dim,1> VectorType;
typedef Scaling<Scalar,Dim> ScalingType;
typedef Transform<Scalar,Dim> TransformType;
protected:
VectorType m_coeffs;
public:
/** Default constructor without initialization. */
Translation() {}
/** */
inline Translation(const Scalar& sx, const Scalar& sy)
{
ei_assert(Dim==2);
m_coeffs.x() = sx;
m_coeffs.y() = sy;
}
/** */
inline Translation(const Scalar& sx, const Scalar& sy, const Scalar& sz)
{
ei_assert(Dim==3);
m_coeffs.x() = sx;
m_coeffs.y() = sy;
m_coeffs.z() = sz;
}
/** Constructs and initialize the scaling transformation from a vector of scaling coefficients */
explicit inline Translation(const VectorType& vector) : m_coeffs(vector) {}
const VectorType& vector() const { return m_coeffs; }
VectorType& vector() { return m_coeffs; }
/** Concatenates two translation */
inline Translation operator* (const Translation& other) const
{ return Translation(m_coeffs + other.m_coeffs); }
/** Concatenates a translation and a scaling */
inline TransformType operator* (const ScalingType& other) const;
/** Concatenates a translation and a linear transformation */
inline TransformType operator* (const LinearMatrixType& linear) const;
/** Concatenates a linear transformation and a translation */
// its a nightmare to define a templated friend function outside its declaration
friend inline TransformType operator* (const LinearMatrixType& linear, const Translation& t)
{
TransformType res;
res.matrix().setZero();
res.linear() = linear;
res.translation() = linear * t.m_coeffs;
res.matrix().row(Dim).setZero();
res(Dim,Dim) = Scalar(1);
return res;
}
/** Concatenates a translation and an affine transformation */
inline TransformType operator* (const TransformType& t) const;
/** Applies translation to vector */
inline VectorType operator* (const VectorType& other) const
{ return m_coeffs + other; }
/** \returns the inverse translation (opposite) */
Translation inverse() const { return Translation(-m_coeffs); }
Translation& operator=(const Translation& other)
{
m_coeffs = other.m_coeffs;
return *this;
}
};
/** \addtogroup GeometryModule */
//@{
typedef Translation<float, 2> Translation2f;
typedef Translation<double,2> Translation2d;
typedef Translation<float, 3> Translation3f;
typedef Translation<double,3> Translation3d;
//@}
template<typename Scalar, int Dim>
inline typename Translation<Scalar,Dim>::TransformType
Translation<Scalar,Dim>::operator* (const ScalingType& other) const
{
TransformType res;
res.matrix().setZero();
res.linear().diagonal() = other.coeffs();
res.translation() = m_coeffs;
res(Dim,Dim) = Scalar(1);
return res;
}
template<typename Scalar, int Dim>
inline typename Translation<Scalar,Dim>::TransformType
Translation<Scalar,Dim>::operator* (const LinearMatrixType& linear) const
{
TransformType res;
res.matrix().setZero();
res.linear() = linear;
res.translation() = m_coeffs;
res.matrix().row(Dim).setZero();
res(Dim,Dim) = Scalar(1);
return res;
}
template<typename Scalar, int Dim>
inline typename Translation<Scalar,Dim>::TransformType
Translation<Scalar,Dim>::operator* (const TransformType& t) const
{
TransformType res = t;
res.pretranslate(m_coeffs);
return res;
}
#endif // EIGEN_TRANSLATION_H

View File

@ -40,6 +40,12 @@ template<typename Scalar> void geometry(void)
typedef Matrix<Scalar,4,1> Vector4;
typedef Quaternion<Scalar> Quaternion;
typedef AngleAxis<Scalar> AngleAxis;
typedef Transform<Scalar,2> Transform2;
typedef Transform<Scalar,3> Transform3;
typedef Scaling<Scalar,2> Scaling2;
typedef Scaling<Scalar,3> Scaling3;
typedef Translation<Scalar,2> Translation2;
typedef Translation<Scalar,3> Translation3;
Quaternion q1, q2;
Vector3 v0 = test_random_matrix<Vector3>(),
@ -115,12 +121,8 @@ template<typename Scalar> void geometry(void)
VERIFY_IS_APPROX(AngleAxis(m).toRotationMatrix(),
Quaternion(m).toRotationMatrix());
// Transform
// TODO complete the tests !
typedef Transform<Scalar,2> Transform2;
typedef Transform<Scalar,3> Transform3;
a = 0;
while (ei_abs(a)<0.1)
a = ei_random<Scalar>(-0.4*M_PI, 0.4*M_PI);
@ -169,12 +171,60 @@ template<typename Scalar> void geometry(void)
t21.linear() = Rotation2D<Scalar>(-a).toRotationMatrix();
VERIFY( (t20.fromPositionOrientationScale(v20,a,v21)
* (t21.prescale(v21.cwise().inverse()).translate(-v20))).isIdentity(test_precision<Scalar>()) );
// Transform - new API
// 3D
t0.setIdentity();
t0.rotate(q1).scale(v0).translate(v0);
// mat * scaling and mat * translation
t1 = (Matrix3(q1) * Scaling3(v0)) * Translation3(v0);
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
// mat * transformation and scaling * translation
t1 = Matrix3(q1) * (Scaling3(v0) * Translation3(v0));
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
t0.setIdentity();
t0.prerotate(q1).prescale(v0).pretranslate(v0);
// translation * scaling and transformation * mat
t1 = (Translation3(v0) * Scaling3(v0)) * Matrix3(q1);
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
// scaling * mat and translation * mat
t1 = Translation3(v0) * (Scaling3(v0) * Matrix3(q1));
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
t0.setIdentity();
t0.scale(v0).translate(v0).rotate(q1);
// translation * mat and scaling * transformation
t1 = Scaling3(v0) * (Translation3(v0) * Matrix3(q1));
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
// transformation * scaling
t0.scale(v0);
t1 = t1 * Scaling3(v0);
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
// transformation * translation
t0.translate(v0);
t1 = t1 * Translation3(v0);
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
// translation * transformation
t0.pretranslate(v0);
t1 = Translation3(v0) * t1;
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
// translation * vector
t0.setIdentity();
t0.translate(v0);
VERIFY_IS_APPROX(t0 * v1, Translation3(v0) * v1);
// scaling * vector
t0.setIdentity();
t0.scale(v0);
VERIFY_IS_APPROX(t0 * v1, Scaling3(v0) * v1);
}
void test_geometry()
{
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( geometry<float>() );
CALL_SUBTEST( geometry<double>() );
// CALL_SUBTEST( geometry<double>() );
}
}