diff --git a/Eigen/Regression b/Eigen/Regression index 5504f922c..de9e0f04a 100644 --- a/Eigen/Regression +++ b/Eigen/Regression @@ -3,6 +3,7 @@ #include "LU" #include "QR" +#include "Geometry" namespace Eigen { diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h index 015f4e5b7..fca47bad4 100644 --- a/Eigen/src/Core/util/ForwardDeclarations.h +++ b/Eigen/src/Core/util/ForwardDeclarations.h @@ -106,6 +106,8 @@ template class Quaternion; template class Rotation2D; template class AngleAxis; template class Transform; +template class ParametrizedLine; +template class Hyperplane; template class Translation; template class Scaling; diff --git a/Eigen/src/Geometry/Hyperplane.h b/Eigen/src/Geometry/Hyperplane.h index f7049883e..65700cc32 100644 --- a/Eigen/src/Geometry/Hyperplane.h +++ b/Eigen/src/Geometry/Hyperplane.h @@ -2,6 +2,7 @@ // for linear algebra. Eigen itself is part of the KDE project. // // Copyright (C) 2008 Gael Guennebaud +// Copyright (C) 2008 Benoit Jacob // // Eigen is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public @@ -27,36 +28,76 @@ /** \geometry_module \ingroup GeometryModule * - * \class Hyperplane + * \class ParametrizedLine * - * \brief Represents an hyper plane in any dimensions + * \brief A parametrized line * * \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 - * - * This class represents an hyper-plane as the zero set of the implicit equation - * \f$ n \cdot x + d = 0 \f$ where \f$ n \f$ is the normal of the plane (linear part) - * and \f$ d \f$ is the distance (offset) to the origin. - * + * \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 -class Hyperplane +template +class ParametrizedLine { - public: - enum { DimAtCompileTime = _Dim }; + enum { AmbientDimAtCompileTime = _AmbientDim }; typedef _Scalar Scalar; typedef typename NumTraits::Real RealScalar; - typedef Matrix VectorType; - typedef Matrix VectorType; + + ParametrizedLine(const VectorType& origin, const VectorType& direction) + : m_origin(origin), m_direction(direction) {} + 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 + * + * \brief A hyperplane + * + * A hyperplane is an affine subspace of dimension n-1 in a space of dimension n. + * For example, a hyperplane in a plane is a line; a hyperplane in 3-space is a plane. + * + * \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. + * + * This class represents an hyperplane as the zero set of the implicit equation + * \f$ n \cdot x + d = 0 \f$ where \f$ n \f$ is a unit normal vector of the plane (linear part) + * and \f$ d \f$ is the distance (offset) to the origin. + */ +template +class Hyperplane +{ + public: + + enum { AmbientDimAtCompileTime = _AmbientDim }; + typedef _Scalar Scalar; + typedef typename NumTraits::Real RealScalar; + typedef Matrix VectorType; + typedef Matrix Coefficients; - typedef Block NormalReturnType; + : AmbientDimAtCompileTime+1,1> Coefficients; + typedef Block NormalReturnType; /** Default constructor without initialization */ - inline Hyperplane(int _dim = DimAtCompileTime) : m_coeffs(_dim+1) {} + inline 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. @@ -64,8 +105,8 @@ class Hyperplane inline Hyperplane(const VectorType& n, const VectorType e) : m_coeffs(n.size()+1) { - _normal() = n; - _offset() = -e.dot(n); + normal() = n; + offset() = -e.dot(n); } /** Constructs a plane from its normal \a n and distance to the origin \a d. @@ -74,85 +115,152 @@ class Hyperplane inline Hyperplane(const VectorType& n, Scalar d) : m_coeffs(n.size()+1) { - _normal() = n; - _offset() = d; + normal() = n; + offset() = d; + } + + /** Constructs a hyperplane passing through the two points. If the dimension of the ambient space + * is greater than 2, then there isn't uniqueness, so an arbitrary choice is made. + */ + static inline Hyperplane Through(const VectorType& p0, const VectorType& p1) + { + Hyperplane result(p0.size()); + result.normal() = (p1 - p0).unitOrthogonal(); + result.offset() = -result.normal().dot(p0); + return result; + } + + /** Constructs a hyperplane passing through the three points. The dimension of the ambient space + * is required to be exactly 3. + */ + static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2) + { + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3); + Hyperplane result(p0.size()); + result.normal() = (p2 - p0).cross(p1 - p0).normalized(); + result.offset() = -result.normal().dot(p0); + return result; + } + + Hyperplane(const ParametrizedLine& parametrized) + { + normal() = parametrized.direction().unitOrthogonal(); + offset() = -normal().dot(parametrized.origin()); } ~Hyperplane() {} /** \returns the dimension in which the plane holds */ - inline int dim() const { return DimAtCompileTime==Dynamic ? m_coeffs.size()-1 : DimAtCompileTime; } - - void normalize(void); + inline int dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : AmbientDimAtCompileTime; } + + /** normalizes \c *this */ + void normalize(void) + { + m_coeffs /= normal().norm(); + } /** \returns the signed distance between the plane \c *this and a point \a p. */ - inline Scalar distanceTo(const VectorType& p) const { return p.dot(normal()) + offset(); } + inline Scalar signedDistance(const VectorType& p) const { return p.dot(normal()) + offset(); } + + /** \returns the absolute distance between the plane \c *this and a point \a p. + */ + inline Scalar absDistance(const VectorType& p) const { return ei_abs(signedDistance(p)); } /** \returns the projection of a point \a p onto the plane \c *this. */ - inline VectorType project(const VectorType& p) const { return p - distanceTo(p) * normal(); } - /** \returns the normal of the plane, which corresponds to the linear part of the implicit equation. */ + 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 + * to the linear part of the implicit equation. + */ inline const NormalReturnType normal() const { return NormalReturnType(m_coeffs,0,0,dim(),1); } - /** \returns the distance to the origin, which is also the constant part + /** \returns a non-constant reference to the unit normal vector of the plane, which corresponds + * to the linear part of the implicit equation. + */ + inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); } + + /** \returns the distance to the origin, which is also the "constant term" of the implicit equation + * \warning the vector normal is assumed to be normalized. + */ + inline const Scalar& offset() const { return m_coeffs.coeff(dim()); } + + /** \returns a non-constant reference to the distance to the origin, which is also the constant part * of the implicit equation */ - inline Scalar offset() const { return m_coeffs(dim()); } + inline Scalar& offset() { return m_coeffs(dim()); } - /** Set the normal of the plane. - * \warning the vector normal is assumed to be normalized. */ - inline void setNormal(const VectorType& normal) { _normal() = normal; } - - /** Set the distance to origin */ - inline void setOffset(Scalar d) { _offset() = d; } - - /** \returns the coefficients c_i of the plane equation: + /** \returns a constant reference to the coefficients c_i of the plane equation: * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$ */ - // FIXME name: equation vs coeffs vs coefficients ??? - inline Coefficients equation(void) const { return m_coeffs; } - - /** \brief Plane/ray intersection. - Returns the parameter value of the intersection between the plane \a *this - and the parametric ray of origin \a rayOrigin and axis \a rayDir - */ - inline Scalar rayIntersection(const VectorType& rayOrigin, const VectorType& rayDir) + inline const Coefficients& coeffs() const { return m_coeffs; } + + /** \returns a non-constant reference to the coefficients c_i of the plane equation: + * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$ + */ + inline Coefficients& coeffs() { return m_coeffs; } + + /** \returns the intersection of *this with \a other. + * + * \warning The ambient space must be a plane, i.e. have dimension 2, so that *this and \a other are lines. + * + * \note If \a other is approximately parallel to *this, this method will return any point on *this. + */ + VectorType intersection(const Hyperplane& other) { - return -(_offset()+rayOrigin.dot(_normal()))/(rayDir.dot(_normal())); + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2); + Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0); + // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests + // whether the two lines are approximately parallel. + if(ei_isMuchSmallerThan(det, Scalar(1))) + { // special case where the two lines are approximately parallel. Pick any point on the first line. + if(ei_abs(coeffs().coeff(1))>ei_abs(coeffs().coeff(0))) + return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0)); + else + return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0)); + } + else + { // general case + Scalar invdet = Scalar(1) / det; + return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)), + invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2))); + } } +#if 0 template inline Hyperplane operator* (const MatrixBase& mat) const - { return Hyperplane(mat.inverse().transpose() * _normal(), _offset()); } + { return Hyperplane(mat.inverse().transpose() * normal(), offset()); } template inline Hyperplane& operator*= (const MatrixBase& mat) const - { _normal() = mat.inverse().transpose() * _normal(); return *this; } - - // TODO some convenient functions to fit a 3D plane on 3 points etc... -// void makePassBy(const VectorType& p0, const VectorType& p1, const VectorType& p2) -// { -// EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(3); -// m_normal = (p2 - p0).cross(p1 - p0).normalized(); -// m_offset = -m_normal.dot(p0); -// } -// -// void makePassBy(const VectorType& p0, const VectorType& p1) -// { -// EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(2); -// m_normal = (p2 - p0).cross(p1 - p0).normalized(); -// m_offset = -m_normal.dot(p0); -// } + { normal() = mat.inverse().transpose() * normal(); return *this; } +#endif protected: - inline NormalReturnType _normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); } - inline Scalar& _offset() { return m_coeffs(dim()); } - Coefficients m_coeffs; }; +template +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 +inline _Scalar ParametrizedLine<_Scalar, _AmbientDim>::intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane) +{ + return -(hyperplane.offset()+origin().dot(hyperplane.normal())) + /(direction().dot(hyperplane.normal())); +} + +#if 0 /** \addtogroup GeometryModule */ //@{ typedef Hyperplane Hyperplane2f; @@ -166,14 +274,6 @@ typedef Hyperplane Planed; typedef Hyperplane HyperplaneXf; typedef Hyperplane HyperplaneXd; //@} - -/** normalizes \c *this */ -template -void Hyperplane<_Scalar,_Dim>::normalize(void) -{ - RealScalar l = Scalar(1)/_normal().norm(); - _normal() *= l; - _offset() *= l; -} +#endif #endif // EIGEN_Hyperplane_H diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h index 7475b75df..a857eff53 100644 --- a/Eigen/src/Geometry/Transform.h +++ b/Eigen/src/Geometry/Transform.h @@ -539,8 +539,9 @@ Transform::extractRotation(TransformTraits traits) const } else if (traits == NoScaling) // though that's stupid let's handle it ! return linear(); - else + else { ei_assert("invalid traits value in Transform::inverse()"); + } return LinearMatrixType(); } diff --git a/Eigen/src/Regression/Regression.h b/Eigen/src/Regression/Regression.h index b03799f49..7f189a77c 100644 --- a/Eigen/src/Regression/Regression.h +++ b/Eigen/src/Regression/Regression.h @@ -155,21 +155,20 @@ void linearRegression(int numPoints, * * \sa linearRegression() */ -template +template void fitHyperplane(int numPoints, VectorType **points, - BigVectorType *result, + HyperplaneType *result, typename NumTraits::Real* soundness = 0) { typedef typename VectorType::Scalar Scalar; typedef Matrix CovMatrixType; EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType) - EIGEN_STATIC_ASSERT_VECTOR_ONLY(BigVectorType) ei_assert(numPoints >= 1); int size = points[0]->size(); - ei_assert(size+1 == result->size()); + ei_assert(size+1 == result->coeffs().size()); - // compue the mean of the data + // compute the mean of the data VectorType mean = VectorType::Zero(size); for(int i = 0; i < numPoints; i++) mean += *(points[i]); @@ -186,13 +185,13 @@ void fitHyperplane(int numPoints, // now we just have to pick the eigen vector with smallest eigen value SelfAdjointEigenSolver eig(covMat); - result->start(size) = eig.eigenvectors().col(0); + result->normal() = eig.eigenvectors().col(0); if (soundness) *soundness = eig.eigenvalues().coeff(0)/eig.eigenvalues().coeff(1); // let's compute the constant coefficient such that the // plane pass trough the mean point: - result->coeffRef(size) = - (result->start(size).cwise()* mean).sum(); + result->offset() = - (result->normal().cwise()* mean).sum(); } diff --git a/test/hyperplane.cpp b/test/hyperplane.cpp index af98652a5..3d001c7d9 100644 --- a/test/hyperplane.cpp +++ b/test/hyperplane.cpp @@ -26,16 +26,16 @@ #include #include -template void hyperplane(const PlaneType& _plane) +template void hyperplane(const HyperplaneType& _plane) { /* this test covers the following files: Hyperplane.h */ const int dim = _plane.dim(); - typedef typename PlaneType::Scalar Scalar; + typedef typename HyperplaneType::Scalar Scalar; typedef typename NumTraits::Real RealScalar; - typedef Matrix VectorType; + typedef Matrix VectorType; VectorType p0 = VectorType::Random(dim); VectorType p1 = VectorType::Random(dim); @@ -43,8 +43,8 @@ template void hyperplane(const PlaneType& _plane) VectorType n0 = VectorType::Random(dim).normalized(); VectorType n1 = VectorType::Random(dim).normalized(); - PlaneType pl0(n0, p0); - PlaneType pl1(n1, p1); + HyperplaneType pl0(n0, p0); + HyperplaneType pl1(n1, p1); Scalar s0 = ei_random(); Scalar s1 = ei_random(); @@ -52,11 +52,43 @@ template void hyperplane(const PlaneType& _plane) VERIFY_IS_APPROX( n1.dot(n1), Scalar(1) ); VERIFY_IS_APPROX( n1.dot(n1), Scalar(1) ); - VERIFY_IS_MUCH_SMALLER_THAN( pl0.distanceTo(p0), Scalar(1) ); - VERIFY_IS_APPROX( pl1.distanceTo(p1 + n1 * s0), s0 ); - VERIFY_IS_MUCH_SMALLER_THAN( pl1.distanceTo(pl1.project(p0)), Scalar(1) ); - VERIFY_IS_MUCH_SMALLER_THAN( pl1.distanceTo(p1 + pl1.normal().unitOrthogonal() * s1), Scalar(1) ); + VERIFY_IS_MUCH_SMALLER_THAN( pl0.absDistance(p0), Scalar(1) ); + VERIFY_IS_APPROX( pl1.signedDistance(p1 + n1 * s0), s0 ); + VERIFY_IS_MUCH_SMALLER_THAN( pl1.signedDistance(pl1.projection(p0)), Scalar(1) ); + VERIFY_IS_MUCH_SMALLER_THAN( pl1.absDistance(p1 + pl1.normal().unitOrthogonal() * s1), Scalar(1) ); +} +template void lines() +{ + typedef Hyperplane HLine; + typedef ParametrizedLine PLine; + typedef Matrix Vector; + typedef Matrix CoeffsType; + + for(int i = 0; i < 10; i++) + { + Vector center = Vector::Random(); + Vector u = Vector::Random(); + Vector v = Vector::Random(); + Scalar a = ei_random(); + + HLine line_u = HLine::Through(center + u, center + a*u); + HLine line_v = HLine::Through(center + v, center + a*v); + + // the line equations should be normalized so that a^2+b^2=1 + VERIFY_IS_APPROX(line_u.normal().norm(), Scalar(1)); + VERIFY_IS_APPROX(line_v.normal().norm(), Scalar(1)); + + Vector result = line_u.intersection(line_v); + + // the lines should intersect at the point we called "center" + VERIFY_IS_APPROX(result, center); + + // check conversions between two types of lines + CoeffsType converted_coeffs(HLine(PLine(line_u)).coeffs()); + converted_coeffs *= line_u.coeffs()(0)/converted_coeffs(0); + VERIFY(line_u.coeffs().isApprox(converted_coeffs)); + } } void test_hyperplane() @@ -66,6 +98,7 @@ void test_hyperplane() CALL_SUBTEST( hyperplane(Hyperplane()) ); CALL_SUBTEST( hyperplane(Hyperplane()) ); CALL_SUBTEST( hyperplane(Hyperplane,5>()) ); - CALL_SUBTEST( hyperplane(Hyperplane(13)) ); + CALL_SUBTEST( lines() ); + CALL_SUBTEST( lines() ); } } diff --git a/test/regression.cpp b/test/regression.cpp index 98d58255e..8bbe0816f 100644 --- a/test/regression.cpp +++ b/test/regression.cpp @@ -26,21 +26,21 @@ #include template + typename HyperplaneType> void makeNoisyCohyperplanarPoints(int numPoints, VectorType **points, - BigVecType *coeffs, + HyperplaneType *hyperplane, typename VectorType::Scalar noiseAmplitude ) { typedef typename VectorType::Scalar Scalar; const int size = points[0]->size(); // pick a random hyperplane, store the coefficients of its equation - coeffs->resize(size + 1); + hyperplane->coeffs().resize(size + 1); for(int j = 0; j < size + 1; j++) { do { - coeffs->coeffRef(j) = ei_random(); - } while(ei_abs(coeffs->coeffRef(j)) < 0.5); + hyperplane->coeffs().coeffRef(j) = ei_random(); + } while(ei_abs(hyperplane->coeffs().coeff(j)) < 0.5); } // now pick numPoints random points on this hyperplane @@ -51,8 +51,8 @@ void makeNoisyCohyperplanarPoints(int numPoints, { cur_point = VectorType::Random(size)/*.normalized()*/; // project cur_point onto the hyperplane - Scalar x = - (coeffs->start(size).cwise()*cur_point).sum(); - cur_point *= coeffs->coeff(size) / x; + Scalar x = - (hyperplane->coeffs().start(size).cwise()*cur_point).sum(); + cur_point *= hyperplane->coeffs().coeff(size) / x; } while( ei_abs(cur_point.norm()) < 0.5 || ei_abs(cur_point.norm()) > 2.0 ); } @@ -63,18 +63,17 @@ void makeNoisyCohyperplanarPoints(int numPoints, } template + typename HyperplaneType> void check_fitHyperplane(int numPoints, VectorType **points, - BigVecType *coeffs, + const HyperplaneType& original, typename VectorType::Scalar tolerance) { int size = points[0]->size(); - BigVecType result(size + 1); + HyperplaneType result(size); fitHyperplane(numPoints, points, &result); - result /= result.coeff(size); - result *= coeffs->coeff(size); - typename VectorType::Scalar error = (result - *coeffs).norm() / coeffs->norm(); + result.coeffs() *= original.coeffs().coeff(size)/result.coeffs().coeff(size); + typename VectorType::Scalar error = (result.coeffs() - original.coeffs()).norm() / original.coeffs().norm(); VERIFY(ei_abs(error) < ei_abs(tolerance)); } @@ -86,31 +85,33 @@ void test_regression() Vector2f points2f [1000]; Vector2f *points2f_ptrs [1000]; for(int i = 0; i < 1000; i++) points2f_ptrs[i] = &(points2f[i]); - Vector3f coeffs3f; + Hyperplane coeffs3f; makeNoisyCohyperplanarPoints(1000, points2f_ptrs, &coeffs3f, 0.01f); - CALL_SUBTEST(check_fitHyperplane(10, points2f_ptrs, &coeffs3f, 0.05f)); - CALL_SUBTEST(check_fitHyperplane(100, points2f_ptrs, &coeffs3f, 0.01f)); - CALL_SUBTEST(check_fitHyperplane(1000, points2f_ptrs, &coeffs3f, 0.002f)); + CALL_SUBTEST(check_fitHyperplane(10, points2f_ptrs, coeffs3f, 0.05f)); + CALL_SUBTEST(check_fitHyperplane(100, points2f_ptrs, coeffs3f, 0.01f)); + CALL_SUBTEST(check_fitHyperplane(1000, points2f_ptrs, coeffs3f, 0.002f)); } { Vector4d points4d [1000]; Vector4d *points4d_ptrs [1000]; for(int i = 0; i < 1000; i++) points4d_ptrs[i] = &(points4d[i]); - Matrix coeffs5d; + Hyperplane coeffs5d; makeNoisyCohyperplanarPoints(1000, points4d_ptrs, &coeffs5d, 0.01); - CALL_SUBTEST(check_fitHyperplane(10, points4d_ptrs, &coeffs5d, 0.05)); - CALL_SUBTEST(check_fitHyperplane(100, points4d_ptrs, &coeffs5d, 0.01)); - CALL_SUBTEST(check_fitHyperplane(1000, points4d_ptrs, &coeffs5d, 0.002)); + CALL_SUBTEST(check_fitHyperplane(10, points4d_ptrs, coeffs5d, 0.05)); + CALL_SUBTEST(check_fitHyperplane(100, points4d_ptrs, coeffs5d, 0.01)); + CALL_SUBTEST(check_fitHyperplane(1000, points4d_ptrs, coeffs5d, 0.002)); } { VectorXcd *points11cd_ptrs[1000]; for(int i = 0; i < 1000; i++) points11cd_ptrs[i] = new VectorXcd(11); - VectorXcd *coeffs12cd = new VectorXcd(12); + Hyperplane,Dynamic> *coeffs12cd = new Hyperplane,Dynamic>(11); makeNoisyCohyperplanarPoints(1000, points11cd_ptrs, coeffs12cd, 0.01); - CALL_SUBTEST(check_fitHyperplane(100, points11cd_ptrs, coeffs12cd, 0.025)); - CALL_SUBTEST(check_fitHyperplane(1000, points11cd_ptrs, coeffs12cd, 0.006)); + CALL_SUBTEST(check_fitHyperplane(100, points11cd_ptrs, *coeffs12cd, 0.025)); + CALL_SUBTEST(check_fitHyperplane(1000, points11cd_ptrs, *coeffs12cd, 0.006)); + delete coeffs12cd; + for(int i = 0; i < 1000; i++) delete points11cd_ptrs[i]; } } } diff --git a/test/submatrices.cpp b/test/submatrices.cpp index cc00ae810..74bf8c7af 100644 --- a/test/submatrices.cpp +++ b/test/submatrices.cpp @@ -128,7 +128,6 @@ template void submatrices(const MatrixType& m) VERIFY(ones.row(r1).sum() == Scalar(cols)); VERIFY(ones.col(c1).dot(ones.col(c2)) == Scalar(rows)); - std::cerr << ones.row(r1).dot(ones.row(r2)) << " == " << cols << "\n"; VERIFY(ones.row(r1).dot(ones.row(r2)) == Scalar(cols)); }