[prev in list] [next in list] [prev in thread] [next in thread] 

List:       kde-commits
Subject:    kdesupport/eigen2
From:       Gael Guennebaud <g.gael () free ! fr>
Date:       2008-10-25 23:10:22
Message-ID: 1224976222.527470.314.nullmailer () svn ! kde ! org
[Download RAW message or body]

SVN commit 875906 by ggael:

Add isApprox in Geometry module's classes.
Complete unit tests wrt previous commits.


 M  +7 -0      Eigen/src/Geometry/AngleAxis.h  
 M  +7 -0      Eigen/src/Geometry/Hyperplane.h  
 M  +7 -0      Eigen/src/Geometry/ParametrizedLine.h  
 M  +8 -3      Eigen/src/Geometry/Quaternion.h  
 M  +7 -0      Eigen/src/Geometry/Rotation2D.h  
 M  +7 -0      Eigen/src/Geometry/Scaling.h  
 M  +7 -0      Eigen/src/Geometry/Transform.h  
 M  +7 -0      Eigen/src/Geometry/Translation.h  
 M  +64 -30    test/geometry.cpp  
 M  +11 -3     test/hyperplane.cpp  
 M  +9 -4      test/main.h  
 M  +10 -2     test/parametrizedline.cpp  


--- trunk/kdesupport/eigen2/Eigen/src/Geometry/AngleAxis.h #875905:875906
@@ -149,6 +149,13 @@
     m_axis = other.axis().template cast<OtherScalarType>();
     m_angle = other.angle();
   }
+
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const AngleAxis& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return m_axis.isApprox(other.m_axis, prec) && \
ei_isApprox(m_angle,other.m_angle, prec); }  };
 
 /** \ingroup GeometryModule
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/Hyperplane.h #875905:875906
@@ -256,6 +256,13 @@
   inline explicit Hyperplane(const \
Hyperplane<OtherScalarType,AmbientDimAtCompileTime>& other)  { m_coeffs = \
other.coeffs().template cast<OtherScalarType>(); }  
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const Hyperplane& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return m_coeffs.isApprox(other.m_coeffs, prec); }
+
 protected:
 
   Coefficients m_coeffs;
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/ParametrizedLine.h #875905:875906
@@ -122,6 +122,13 @@
     m_direction = other.direction().template cast<OtherScalarType>();
   }
 
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const ParametrizedLine& other, typename NumTraits<Scalar>::Real prec \
= precision<Scalar>()) const +  { return m_origin.isApprox(other.m_origin, prec) && \
m_direction.isApprox(other.m_direction, prec); } +
 protected:
 
   VectorType m_origin, m_direction;
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/Quaternion.h #875905:875906
@@ -207,10 +207,15 @@
   /** Copy constructor with scalar type conversion */
   template<typename OtherScalarType>
   inline explicit Quaternion(const Quaternion<OtherScalarType>& other)
-  {
-    m_coeffs = other.coeffs().template cast<OtherScalarType>();
-  }
+  { m_coeffs = other.coeffs().template cast<OtherScalarType>(); }
 
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const Quaternion& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return m_coeffs.isApprox(other.m_coeffs, prec); }
+
 };
 
 /** \ingroup GeometryModule
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/Rotation2D.h #875905:875906
@@ -116,6 +116,13 @@
   {
     m_angle = other.angle();
   }
+
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const Rotation2D& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return ei_isApprox(m_angle,other.m_angle, prec); }
 };
 
 /** \ingroup GeometryModule
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/Scaling.h #875905:875906
@@ -142,6 +142,13 @@
   inline explicit Scaling(const Scaling<OtherScalarType,Dim>& other)
   { m_coeffs = other.coeffs().template cast<OtherScalarType>(); }
 
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const Scaling& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return m_coeffs.isApprox(other.m_coeffs, prec); }
+
 };
 
 /** \addtogroup GeometryModule */
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/Transform.h #875905:875906
@@ -260,6 +260,13 @@
   inline explicit Transform(const Transform<OtherScalarType,Dim>& other)
   { m_matrix = other.matrix().template cast<OtherScalarType>(); }
 
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const Transform& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return m_matrix.isApprox(other.m_matrix, prec); }
+
 protected:
 
 };
--- trunk/kdesupport/eigen2/Eigen/src/Geometry/Translation.h #875905:875906
@@ -145,6 +145,13 @@
   inline explicit Translation(const Translation<OtherScalarType,Dim>& other)
   { m_coeffs = other.vector().template cast<OtherScalarType>(); }
 
+  /** \returns \c true if \c *this is approximately equal to \a other, within the \
precision +    * determined by \a prec.
+    *
+    * \sa MatrixBase::isApprox() */
+  bool isApprox(const Translation& other, typename NumTraits<Scalar>::Real prec = \
precision<Scalar>()) const +  { return m_coeffs.isApprox(other.m_coeffs, prec); }
+
 };
 
 /** \addtogroup GeometryModule */
--- trunk/kdesupport/eigen2/test/geometry.cpp #875905:875906
@@ -39,8 +39,8 @@
   typedef Matrix<Scalar,2,1> Vector2;
   typedef Matrix<Scalar,3,1> Vector3;
   typedef Matrix<Scalar,4,1> Vector4;
-  typedef Quaternion<Scalar> Quaternion;
-  typedef AngleAxis<Scalar> AngleAxis;
+  typedef Quaternion<Scalar> Quaternionx;
+  typedef AngleAxis<Scalar> AngleAxisx;
   typedef Transform<Scalar,2> Transform2;
   typedef Transform<Scalar,3> Transform3;
   typedef Scaling<Scalar,2> Scaling2;
@@ -52,7 +52,7 @@
   if (ei_is_same_type<Scalar,float>::ret)
     largeEps = 1e-3f;
 
-  Quaternion q1, q2;
+  Quaternionx q1, q2;
   Vector3 v0 = Vector3::Random(),
     v1 = Vector3::Random(),
     v2 = Vector3::Random();
@@ -76,18 +76,18 @@
   VERIFY_IS_APPROX(v0.unitOrthogonal().norm(), Scalar(1));
 
 
-  VERIFY_IS_APPROX(v0, AngleAxis(a, v0.normalized()) * v0);
-  VERIFY_IS_APPROX(-v0, AngleAxis(M_PI, 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);
+  VERIFY_IS_APPROX(v0, AngleAxisx(a, v0.normalized()) * v0);
+  VERIFY_IS_APPROX(-v0, AngleAxisx(M_PI, v0.unitOrthogonal()) * v0);
+  VERIFY_IS_APPROX(ei_cos(a)*v0.norm2(), v0.dot(AngleAxisx(a, v0.unitOrthogonal()) * \
v0)); +  m = AngleAxisx(a, v0.normalized()).toRotationMatrix().adjoint();
+  VERIFY_IS_APPROX(Matrix3::Identity(), m * AngleAxisx(a, v0.normalized()));
+  VERIFY_IS_APPROX(Matrix3::Identity(), AngleAxisx(a, v0.normalized()) * m);
 
-  q1 = AngleAxis(a, v0.normalized());
-  q2 = AngleAxis(a, v1.normalized());
+  q1 = AngleAxisx(a, v0.normalized());
+  q2 = AngleAxisx(a, v1.normalized());
 
   // angular distance
-  Scalar refangle = ei_abs(AngleAxis(q1.inverse()*q2).angle());
+  Scalar refangle = ei_abs(AngleAxisx(q1.inverse()*q2).angle());
   if (refangle>M_PI)
     refangle = 2.*M_PI - refangle;
   VERIFY(ei_isApprox(q1.angularDistance(q2), refangle, largeEps));
@@ -101,18 +101,18 @@
   q2 = q1.toRotationMatrix();
   VERIFY_IS_APPROX(q1*v1,q2*v1);
 
-  matrot1 = AngleAxis(0.1, Vector3::UnitX())
-          * AngleAxis(0.2, Vector3::UnitY())
-          * AngleAxis(0.3, Vector3::UnitZ());
+  matrot1 = AngleAxisx(0.1, Vector3::UnitX())
+          * AngleAxisx(0.2, Vector3::UnitY())
+          * AngleAxisx(0.3, Vector3::UnitZ());
   VERIFY_IS_APPROX(matrot1 * v1,
-       AngleAxis(0.1, Vector3(1,0,0)).toRotationMatrix()
-    * (AngleAxis(0.2, Vector3(0,1,0)).toRotationMatrix()
-    * (AngleAxis(0.3, Vector3(0,0,1)).toRotationMatrix() * v1)));
+       AngleAxisx(0.1, Vector3(1,0,0)).toRotationMatrix()
+    * (AngleAxisx(0.2, Vector3(0,1,0)).toRotationMatrix()
+    * (AngleAxisx(0.3, Vector3(0,0,1)).toRotationMatrix() * v1)));
 
   // angle-axis conversion
-  AngleAxis aa = q1;
-  VERIFY_IS_APPROX(q1 * v1, Quaternion(aa) * v1);
-  VERIFY_IS_NOT_APPROX(q1 * v1, Quaternion(AngleAxis(aa.angle()*2,aa.axis())) * v1);
+  AngleAxisx aa = q1;
+  VERIFY_IS_APPROX(q1 * v1, Quaternionx(aa) * v1);
+  VERIFY_IS_NOT_APPROX(q1 * v1, Quaternionx(AngleAxisx(aa.angle()*2,aa.axis())) * \
v1);  
   // from two vector creation
   VERIFY_IS_APPROX(v2.normalized(),(q2.setFromTwoVectors(v1,v2)*v1).normalized());
@@ -123,21 +123,21 @@
   VERIFY_IS_APPROX(q1 * (q1.conjugate() * v1), v1);
 
   // AngleAxis
-  VERIFY_IS_APPROX(AngleAxis(a,v1.normalized()).toRotationMatrix(),
-    Quaternion(AngleAxis(a,v1.normalized())).toRotationMatrix());
+  VERIFY_IS_APPROX(AngleAxisx(a,v1.normalized()).toRotationMatrix(),
+    Quaternionx(AngleAxisx(a,v1.normalized())).toRotationMatrix());
 
-  AngleAxis aa1;
+  AngleAxisx aa1;
   m = q1.toRotationMatrix();
   aa1 = m;
-  VERIFY_IS_APPROX(AngleAxis(m).toRotationMatrix(),
-    Quaternion(m).toRotationMatrix());
+  VERIFY_IS_APPROX(AngleAxisx(m).toRotationMatrix(),
+    Quaternionx(m).toRotationMatrix());
 
   // Transform
   // TODO complete the tests !
   a = 0;
   while (ei_abs(a)<0.1)
     a = ei_random<Scalar>(-0.4*M_PI, 0.4*M_PI);
-  q1 = AngleAxis(a, v0.normalized());
+  q1 = AngleAxisx(a, v0.normalized());
   Transform3 t0, t1, t2;
   t0.setIdentity();
   t0.linear() = q1.toRotationMatrix();
@@ -171,7 +171,7 @@
   t1.setIdentity(); t1.scale(v0).rotate(q1);
   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
 
-  t0.setIdentity(); t0.scale(v0).rotate(AngleAxis(q1));
+  t0.setIdentity(); t0.scale(v0).rotate(AngleAxisx(q1));
   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
 
   VERIFY_IS_APPROX(t0.scale(a).matrix(), t1.scale(Vector3::Constant(a)).matrix());
@@ -212,7 +212,7 @@
   // 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
@@ -226,7 +226,7 @@
   t0.translate(v0);
   t1 = t1 * Translation3(v0);
   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
-  // translation * transformation 
+  // translation * transformation
   t0.pretranslate(v0);
   t1 = Translation3(v0) * t1;
   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
@@ -284,6 +284,40 @@
   t0.setIdentity();
   t0.translate(v0).rotate(q1).scale(v1);
   VERIFY_IS_APPROX(t0.extractRotation(Affine) * v1, Matrix3(q1) * v1);
+
+  // test casting
+  Transform<float,3> t1f = t1.template cast<float>();
+  VERIFY_IS_APPROX(t1f.template cast<Scalar>(),t1);
+  Transform<double,3> t1d = t1.template cast<double>();
+  VERIFY_IS_APPROX(t1d.template cast<Scalar>(),t1);
+
+  Translation3 tr1(v0);
+  Translation<float,3> tr1f = tr1.template cast<float>();
+  VERIFY_IS_APPROX(tr1f.template cast<Scalar>(),tr1);
+  Translation<double,3> tr1d = tr1.template cast<double>();
+  VERIFY_IS_APPROX(tr1d.template cast<Scalar>(),tr1);
+
+  Scaling3 sc1(v0);
+  Scaling<float,3> sc1f = sc1.template cast<float>();
+  VERIFY_IS_APPROX(sc1f.template cast<Scalar>(),sc1);
+  Scaling<double,3> sc1d = sc1.template cast<double>();
+  VERIFY_IS_APPROX(sc1d.template cast<Scalar>(),sc1);
+
+  Quaternion<float> q1f = q1.template cast<float>();
+  VERIFY_IS_APPROX(q1f.template cast<Scalar>(),q1);
+  Quaternion<double> q1d = q1.template cast<double>();
+  VERIFY_IS_APPROX(q1d.template cast<Scalar>(),q1);
+
+  AngleAxis<float> aa1f = aa1.template cast<float>();
+  VERIFY_IS_APPROX(aa1f.template cast<Scalar>(),aa1);
+  AngleAxis<double> aa1d = aa1.template cast<double>();
+  VERIFY_IS_APPROX(aa1d.template cast<Scalar>(),aa1);
+
+  Rotation2D<Scalar> r2d1(ei_random<Scalar>());
+  Rotation2D<float> r2d1f = r2d1.template cast<float>();
+  VERIFY_IS_APPROX(r2d1f.template cast<Scalar>(),r2d1);
+  Rotation2D<double> r2d1d = r2d1.template cast<double>();
+  VERIFY_IS_APPROX(r2d1d.template cast<Scalar>(),r2d1);
 }
 
 void test_geometry()
--- trunk/kdesupport/eigen2/test/hyperplane.cpp #875905:875906
@@ -46,7 +46,7 @@
 
   VectorType n0 = VectorType::Random(dim).normalized();
   VectorType n1 = VectorType::Random(dim).normalized();
-  
+
   HyperplaneType pl0(n0, p0);
   HyperplaneType pl1(n1, p1);
   HyperplaneType pl2 = pl1;
@@ -55,7 +55,7 @@
   Scalar s1 = ei_random<Scalar>();
 
   VERIFY_IS_APPROX( n1.dot(n1), 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) );
@@ -67,7 +67,7 @@
     MatrixType rot = MatrixType::Random(dim,dim).qr().matrixQ();
     Scaling<Scalar,HyperplaneType::AmbientDimAtCompileTime> \
                scaling(VectorType::Random());
     Translation<Scalar,HyperplaneType::AmbientDimAtCompileTime> \
                translation(VectorType::Random());
-    
+
     pl2 = pl1;
     VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot).absDistance(rot * p1), Scalar(1) \
);  pl2 = pl1;
@@ -81,6 +81,14 @@
     VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot*translation,Isometry)
                                  .absDistance((rot*translation) * p1), Scalar(1) );
   }
+
+  // casting
+  const int Dim = HyperplaneType::AmbientDimAtCompileTime;
+  typedef typename GetDifferentType<Scalar>::type OtherScalar;
+  Hyperplane<OtherScalar,Dim> hp1f = pl1.template cast<OtherScalar>();
+  VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),pl1);
+  Hyperplane<Scalar,Dim> hp1d = pl1.template cast<Scalar>();
+  VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),pl1);
 }
 
 template<typename Scalar> void lines()
--- trunk/kdesupport/eigen2/test/main.h #875905:875906
@@ -209,11 +209,10 @@
 inline bool test_ei_isApproxOrLessThan(const long double& a, const long double& b)
 { return ei_isApproxOrLessThan(a, b, test_precision<long double>()); }
 
-template<typename Derived1, typename Derived2>
-inline bool test_ei_isApprox(const MatrixBase<Derived1>& m1,
-                   const MatrixBase<Derived2>& m2)
+template<typename Type1, typename Type2>
+inline bool test_ei_isApprox(const Type1& a, const Type2& b)
 {
-  return m1.isApprox(m2, test_precision<typename ei_traits<Derived1>::Scalar>());
+  return a.isApprox(b, test_precision<typename Type1::Scalar>());
 }
 
 template<typename Derived1, typename Derived2>
@@ -232,7 +231,13 @@
 
 } // end namespace Eigen
 
+template<typename T> struct GetDifferentType;
 
+template<> struct GetDifferentType<float> { typedef double type; };
+template<> struct GetDifferentType<double> { typedef float type; };
+template<typename T> struct GetDifferentType<std::complex<T> >
+{ typedef std::complex<typename GetDifferentType<T>::type> type; };
+
 // forward declaration of the main test function
 void EI_PP_CAT(test_,EIGEN_TEST_FUNC)();
 
--- trunk/kdesupport/eigen2/test/parametrizedline.cpp #875905:875906
@@ -45,7 +45,7 @@
   VectorType p1 = VectorType::Random(dim);
 
   VectorType d0 = VectorType::Random(dim).normalized();
-  
+
   LineType l0(p0, d0);
 
   Scalar s0 = ei_random<Scalar>();
@@ -55,7 +55,15 @@
   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 );
+  VERIFY_IS_APPROX( Scalar(l0.distance((p0+s0*d0) + d0.unitOrthogonal() * s1)), s1 \
); +
+  // casting
+  const int Dim = LineType::AmbientDimAtCompileTime;
+  typedef typename GetDifferentType<Scalar>::type OtherScalar;
+  ParametrizedLine<OtherScalar,Dim> hp1f = l0.template cast<OtherScalar>();
+  VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),l0);
+  ParametrizedLine<Scalar,Dim> hp1d = l0.template cast<Scalar>();
+  VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),l0);
 }
 
 void test_parametrizedline()


[prev in list] [next in list] [prev in thread] [next in thread] 

Configure | About | News | Add a list | Sponsored by KoreLogic