* Extend a bit ParametrizedLine and move it to a separate file,

add unit-tests for it.
* remove "using namespace std" in test/main.h such that the compilation
  bug found today in SVD won't happen again.
This commit is contained in:
Gael Guennebaud 2008-09-03 22:35:45 +00:00
parent c29c7b0ea9
commit 52406aecad
8 changed files with 206 additions and 78 deletions

View File

@ -12,10 +12,11 @@ namespace Eigen {
/** \defgroup GeometryModule Geometry module
* This module provides support for:
* - fixed-size homogeneous transformations
* - 2D and 3D rotations
* - translation, scaling, 2D and 3D rotations
* - quaternions
* - \ref MatrixBase::cross() "cross product"
* - \ref MatrixBase::unitOrthogonal() "orthognal vector generation"
* - some linear components: parametrized-lines and hyperplanes
*
* \code
* #include <Eigen/Geometry>
@ -36,6 +37,7 @@ namespace Eigen {
#include "src/Geometry/Translation.h"
#include "src/Geometry/Scaling.h"
#include "src/Geometry/Hyperplane.h"
#include "src/Geometry/ParametrizedLine.h"
} // namespace Eigen

View File

@ -26,48 +26,6 @@
#ifndef EIGEN_HYPERPLANE_H
#define EIGEN_HYPERPLANE_H
/** \geometry_module \ingroup GeometryModule
*
* \class ParametrizedLine
*
* \brief A parametrized line
*
* \param _Scalar the scalar type, i.e., the type of the coefficients
* \param _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
* Notice that the dimension of the hyperplane is _AmbientDim-1.
*/
template <typename _Scalar, int _AmbientDim>
class ParametrizedLine
#ifdef EIGEN_VECTORIZE
: public ei_with_aligned_operator_new<_Scalar,_AmbientDim>
#endif
{
public:
enum { AmbientDimAtCompileTime = _AmbientDim };
typedef _Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
ParametrizedLine(const VectorType& origin, const VectorType& direction)
: m_origin(origin), m_direction(direction) {}
explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane);
~ParametrizedLine() {}
const VectorType& origin() const { return m_origin; }
VectorType& origin() { return m_origin; }
const VectorType& direction() const { return m_direction; }
VectorType& direction() { return m_direction; }
Scalar intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane);
protected:
VectorType m_origin, m_direction;
};
/** \geometry_module \ingroup GeometryModule
*
* \class Hyperplane
@ -103,7 +61,7 @@ class Hyperplane
typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
/** Default constructor without initialization */
inline Hyperplane(int _dim = AmbientDimAtCompileTime) : m_coeffs(_dim+1) {}
inline explicit Hyperplane(int _dim = AmbientDimAtCompileTime) : m_coeffs(_dim+1) {}
/** Construct a plane from its normal \a n and a point \a e onto the plane.
* \warning the vector normal is assumed to be normalized.
@ -175,7 +133,6 @@ class Hyperplane
/** \returns the projection of a point \a p onto the plane \c *this.
*/
inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
/** \returns a constant reference to the unit normal vector of the plane, which corresponds
@ -261,25 +218,4 @@ protected:
Coefficients m_coeffs;
};
/** Construct a parametrized line from a 2D hyperplane
*
* \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line
*/
template <typename _Scalar, int _AmbientDim>
ParametrizedLine<_Scalar, _AmbientDim>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2);
direction() = hyperplane.normal().unitOrthogonal();
origin() = -hyperplane.normal()*hyperplane.offset();
}
/** \returns the parameter value of the intersection between *this and the given hyperplane
*/
template <typename _Scalar, int _AmbientDim>
inline _Scalar ParametrizedLine<_Scalar, _AmbientDim>::intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
{
return -(hyperplane.offset()+origin().dot(hyperplane.normal()))
/(direction().dot(hyperplane.normal()));
}
#endif // EIGEN_HYPERPLANE_H

View File

@ -0,0 +1,118 @@
// 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>
// Copyright (C) 2008 Benoit Jacob <jacob@math.jussieu.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_PARAMETRIZEDLINE_H
#define EIGEN_PARAMETRIZEDLINE_H
/** \geometry_module \ingroup GeometryModule
*
* \class ParametrizedLine
*
* \brief A parametrized line
*
* \param _Scalar the scalar type, i.e., the type of the coefficients
* \param _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
* Notice that the dimension of the hyperplane is _AmbientDim-1.
*/
template <typename _Scalar, int _AmbientDim>
class ParametrizedLine
#ifdef EIGEN_VECTORIZE
: public ei_with_aligned_operator_new<_Scalar,_AmbientDim>
#endif
{
public:
enum { AmbientDimAtCompileTime = _AmbientDim };
typedef _Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
/** Default constructor without initialization */
inline explicit ParametrizedLine(int _dim = AmbientDimAtCompileTime)
: m_origin(_dim), m_direction(_dim)
{}
ParametrizedLine(const VectorType& origin, const VectorType& direction)
: m_origin(origin), m_direction(direction) {}
explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane);
~ParametrizedLine() {}
/** \returns the dimension in which the line holds */
inline int dim() const { return m_direction.size(); }
const VectorType& origin() const { return m_origin; }
VectorType& origin() { return m_origin; }
const VectorType& direction() const { return m_direction; }
VectorType& direction() { return m_direction; }
/** \returns the squared distance of a point \a p to its projection onto the line \c *this.
* \sa distance()
*/
RealScalar squaredDistance(const VectorType& p) const
{
VectorType diff = p-origin();
return (diff - diff.dot(direction())* direction()).norm2();
}
/** \returns the distance of a point \a p to its projection onto the line \c *this.
* \sa squaredDistance()
*/
RealScalar distance(const VectorType& p) const { return ei_sqrt(squaredDistance(p)); }
/** \returns the projection of a point \a p onto the line \c *this.
*/
VectorType projection(const VectorType& p) const
{ return origin() + (p-origin()).dot(direction()) * direction(); }
Scalar intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane);
protected:
VectorType m_origin, m_direction;
};
/** Construct a parametrized line from a 2D hyperplane
*
* \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line
*/
template <typename _Scalar, int _AmbientDim>
inline ParametrizedLine<_Scalar, _AmbientDim>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2);
direction() = hyperplane.normal().unitOrthogonal();
origin() = -hyperplane.normal()*hyperplane.offset();
}
/** \returns the parameter value of the intersection between *this and the given hyperplane
*/
template <typename _Scalar, int _AmbientDim>
inline _Scalar ParametrizedLine<_Scalar, _AmbientDim>::intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
{
return -(hyperplane.offset()+origin().dot(hyperplane.normal()))
/(direction().dot(hyperplane.normal()));
}
#endif // EIGEN_PARAMETRIZEDLINE_H

View File

@ -120,6 +120,7 @@ EI_ADD_TEST(eigensolver)
EI_ADD_TEST(svd)
EI_ADD_TEST(geometry)
EI_ADD_TEST(hyperplane)
EI_ADD_TEST(parametrizedline)
EI_ADD_TEST(regression)
EI_ADD_TEST(sparse)

View File

@ -74,7 +74,7 @@ template<typename Scalar> void geometry(void)
VERIFY_IS_APPROX(v0, AngleAxis(a, v0.normalized()) * v0);
VERIFY_IS_APPROX(-v0, AngleAxis(M_PI, v0.unitOrthogonal()) * v0);
VERIFY_IS_APPROX(cos(a)*v0.norm2(), v0.dot(AngleAxis(a, v0.unitOrthogonal()) * v0));
VERIFY_IS_APPROX(ei_cos(a)*v0.norm2(), v0.dot(AngleAxis(a, v0.unitOrthogonal()) * v0));
m = AngleAxis(a, v0.normalized()).toRotationMatrix().adjoint();
VERIFY_IS_APPROX(Matrix3::Identity(), m * AngleAxis(a, v0.normalized()));
VERIFY_IS_APPROX(Matrix3::Identity(), AngleAxis(a, v0.normalized()) * m);

View File

@ -2,6 +2,7 @@
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
// Copyright (C) 2008 Benoit Jacob <jacob@math.jussieu.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
@ -53,7 +54,6 @@ template<typename HyperplaneType> void hyperplane(const HyperplaneType& _plane)
Scalar s0 = ei_random<Scalar>();
Scalar s1 = ei_random<Scalar>();
VERIFY_IS_APPROX( n1.dot(n1), Scalar(1) );
VERIFY_IS_APPROX( n1.dot(n1), Scalar(1) );
VERIFY_IS_MUCH_SMALLER_THAN( pl0.absDistance(p0), Scalar(1) );
@ -96,6 +96,9 @@ template<typename Scalar> void lines()
Vector u = Vector::Random();
Vector v = Vector::Random();
Scalar a = ei_random<Scalar>();
while (ei_abs(a-1) < 1e-4) a = ei_random<Scalar>();
while (u.norm() < 1e-4) u = Vector::Random();
while (v.norm() < 1e-4) v = Vector::Random();
HLine line_u = HLine::Through(center + u, center + a*u);
HLine line_v = HLine::Through(center + v, center + a*v);

View File

@ -229,7 +229,6 @@ inline bool test_ei_isMuchSmallerThan(const MatrixBase<Derived>& m,
void EI_PP_CAT(test_,EIGEN_TEST_FUNC)();
using namespace Eigen;
using namespace std;
int main(int argc, char *argv[])
{
@ -245,14 +244,14 @@ int main(int argc, char *argv[])
{
if(has_set_repeat)
{
cout << "Argument " << argv[i] << " conflicting with a former argument" << endl;
std::cout << "Argument " << argv[i] << " conflicting with a former argument" << std::endl;
return 1;
}
repeat = atoi(argv[i]+1);
has_set_repeat = true;
if(repeat <= 0)
{
cout << "Invalid \'repeat\' value " << argv[i]+1 << endl;
std::cout << "Invalid \'repeat\' value " << argv[i]+1 << std::endl;
return 1;
}
}
@ -260,7 +259,7 @@ int main(int argc, char *argv[])
{
if(has_set_seed)
{
cout << "Argument " << argv[i] << " conflicting with a former argument" << endl;
std::cout << "Argument " << argv[i] << " conflicting with a former argument" << std::endl;
return 1;
}
seed = strtoul(argv[i]+1, 0, 10);
@ -268,7 +267,7 @@ int main(int argc, char *argv[])
bool ok = seed!=0;
if(!ok)
{
cout << "Invalid \'seed\' value " << argv[i]+1 << endl;
std::cout << "Invalid \'seed\' value " << argv[i]+1 << std::endl;
return 1;
}
}
@ -280,18 +279,18 @@ int main(int argc, char *argv[])
if(need_help)
{
cout << "This test application takes the following optional arguments:" << endl;
cout << " rN Repeat each test N times (default: " << DEFAULT_REPEAT << ")" << endl;
cout << " sN Use N as seed for random numbers (default: based on current time)" << endl;
std::cout << "This test application takes the following optional arguments:" << std::endl;
std::cout << " rN Repeat each test N times (default: " << DEFAULT_REPEAT << ")" << std::endl;
std::cout << " sN Use N as seed for random numbers (default: based on current time)" << std::endl;
return 1;
}
if(!has_set_seed) seed = (unsigned int) time(NULL);
if(!has_set_repeat) repeat = DEFAULT_REPEAT;
cout << "Initializing random number generator with seed " << seed << endl;
std::cout << "Initializing random number generator with seed " << seed << std::endl;
srand(seed);
cout << "Repeating each test " << repeat << " times" << endl;
std::cout << "Repeating each test " << repeat << " times" << std::endl;
Eigen::g_repeat = repeat;
Eigen::g_test_stack.push_back(EI_PP_MAKE_STRING(EIGEN_TEST_FUNC));

69
test/parametrizedline.cpp Normal file
View File

@ -0,0 +1,69 @@
// 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>
// Copyright (C) 2008 Benoit Jacob <jacob@math.jussieu.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/>.
#include "main.h"
#include <Eigen/Geometry>
#include <Eigen/LU>
#include <Eigen/QR>
template<typename LineType> void parametrizedline(const LineType& _line)
{
/* this test covers the following files:
ParametrizedLine.h
*/
const int dim = _line.dim();
typedef typename LineType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, LineType::AmbientDimAtCompileTime, 1> VectorType;
typedef Matrix<Scalar, LineType::AmbientDimAtCompileTime,
LineType::AmbientDimAtCompileTime> MatrixType;
VectorType p0 = VectorType::Random(dim);
VectorType p1 = VectorType::Random(dim);
VectorType d0 = VectorType::Random(dim).normalized();
LineType l0(p0, d0);
Scalar s0 = ei_random<Scalar>();
Scalar s1 = ei_abs(ei_random<Scalar>());
VERIFY_IS_MUCH_SMALLER_THAN( l0.distance(p0), RealScalar(1) );
VERIFY_IS_MUCH_SMALLER_THAN( l0.distance(p0+s0*d0), RealScalar(1) );
VERIFY_IS_APPROX( (l0.projection(p1)-p1).norm(), l0.distance(p1) );
VERIFY_IS_MUCH_SMALLER_THAN( l0.distance(l0.projection(p1)), RealScalar(1) );
VERIFY_IS_APPROX( l0.distance((p0+s0*d0) + d0.unitOrthogonal() * s1), s1 );
}
void test_parametrizedline()
{
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( parametrizedline(ParametrizedLine<float,2>()) );
CALL_SUBTEST( parametrizedline(ParametrizedLine<float,3>()) );
CALL_SUBTEST( parametrizedline(ParametrizedLine<double,4>()) );
CALL_SUBTEST( parametrizedline(ParametrizedLine<std::complex<double>,5>()) );
}
}