Re: [eigen] Enhanced AlignedBox

[ Thread Index | Date Index | More lists.tuxfamily.org/eigen Archives ]


Hi !
following my previous thoughts, I made a minor modification in the
file MathFunctions.h to add highest and lowest functions in the eigen
way.
For complex types, I choose to provide also a lowest and a highest
(the natural ones from a 2D bounding box point of view).

So here are the consequent patches for MathFunctions.h, AlignedBox.h
and geo_alignedbox.cpp.


P.S. I will open a topic on NumTraits and a possible link to boost as
soon as I have a clear design in head for what a class NumTraits
should look like (I have to look at other libraries that face the same
kind of problem to try to see the whole picture). In the mean time the
minimal change I propose seems ok for me.
# HG changeset patch
# User Manuel Yguel <manuel.yguel@xxxxxxxxx>
# Date 1262609869 -3600
# Node ID 7503d16a092beccddf662f922426c98e666a1ff3
# Parent  d96ed163a2749b8d3faa997b11ffd377cc8434e4
Add lowest and highest functions that return the minimum and maximum values of a certain numeric type. It is a workaround for the fact that with std::numeric_limits min stands for minimum positive value for floating numbers.

diff -r d96ed163a274 -r 7503d16a092b Eigen/src/Core/MathFunctions.h
--- a/Eigen/src/Core/MathFunctions.h	Sat Jan 02 18:01:39 2010 -0500
+++ b/Eigen/src/Core/MathFunctions.h	Mon Jan 04 13:57:49 2010 +0100
@@ -51,9 +51,21 @@
   return p * ei_sqrt(T(1) + qp*qp);
 }
 
+
+template<typename T> inline T highest();
+template<typename T> inline T lowest();
+
+template<typename T> inline T ei_default_float_highest(){ return std::numeric_limits<T>::max(); }
+template<typename T> inline T ei_default_float_lowest(){ return static_cast<T>(-std::numeric_limits<T>::max()); }
+template<typename T> inline T ei_default_integral_highest(){ return std::numeric_limits<T>::max(); }
+template<typename T> inline T ei_default_integral_lowest(){ return std::numeric_limits<T>::min(); }
+
 /**************
 ***   int   ***
 **************/
+
+template<> inline int highest<int>() { return ei_default_integral_highest<int>(); }
+template<> inline int lowest<int>() { return ei_default_integral_lowest<int>(); }
 
 template<> inline int dummy_precision<int>() { return 0; }
 inline int ei_real(int x)  { return x; }
@@ -108,6 +120,8 @@
 /**************
 *** float   ***
 **************/
+template<> inline float highest<float>() { return ei_default_float_highest<float>(); }
+template<> inline float lowest<float>() { return ei_default_float_lowest<float>(); }
 
 template<> inline float dummy_precision<float>() { return 1e-5f; }
 inline float ei_real(float x)  { return x; }
@@ -156,6 +170,8 @@
 /**************
 *** double  ***
 **************/
+template<> inline double highest<double>() { return ei_default_float_highest<double>(); }
+template<> inline double lowest<double>() { return ei_default_float_lowest<double>(); }
 
 template<> inline double dummy_precision<double>() { return 1e-12; }
 
@@ -205,6 +221,10 @@
 /*********************
 *** complex<float> ***
 *********************/
+template<> inline std::complex<float> highest<std::complex<float> >() {
+  return std::complex<float>( ei_default_float_highest<float>(), ei_default_float_highest<float>() ); }
+template<> inline std::complex<float> lowest<std::complex<float> >() {
+  return std::complex<float>( ei_default_float_lowest<float>(),ei_default_float_lowest<float>() ); }
 
 template<> inline float dummy_precision<std::complex<float> >() { return dummy_precision<float>(); }
 inline float ei_real(const std::complex<float>& x) { return std::real(x); }
@@ -242,6 +262,10 @@
 /**********************
 *** complex<double> ***
 **********************/
+template<> inline std::complex<double> highest<std::complex<double> >() {
+  return std::complex<double>( ei_default_float_highest<double>(), ei_default_float_highest<double>() ); }
+template<> inline std::complex<double> lowest<std::complex<double> >() {
+  return std::complex<double>( ei_default_float_lowest<double>(),ei_default_float_lowest<double>() ); }
 
 template<> inline double dummy_precision<std::complex<double> >() { return dummy_precision<double>(); }
 inline double ei_real(const std::complex<double>& x) { return std::real(x); }
@@ -280,6 +304,8 @@
 /******************
 *** long double ***
 ******************/
+template<> inline long double highest<long double>() { return ei_default_float_highest<long double>(); }
+template<> inline long double lowest<long double>() { return ei_default_float_lowest<long double>(); }
 
 template<> inline long double dummy_precision<long double>() { return dummy_precision<double>(); }
 inline long double ei_real(long double x)  { return x; }
@@ -320,6 +346,8 @@
 /**************
 ***  bool  ***
 **************/
+template<> inline bool highest<bool>() { return ei_default_integral_highest<bool>(); }
+template<> inline bool lowest<bool>() { return ei_default_integral_lowest<bool>(); }
 
 template<> inline bool dummy_precision<bool>() { return 0; }
 inline bool ei_real(bool x)  { return x; }
# HG changeset patch
# User Manuel Yguel <manuel.yguel@xxxxxxxxx>
# Date 1262610056 -3600
# Node ID e75eb64d0ca189f7beb0a1f13e34b7545f2da2c3
# Parent  7503d16a092beccddf662f922426c98e666a1ff3
Enhance AlignedBox class by providing some simple functions (diagonal,volume,sample,sizes) and better handling integral types.

diff -r 7503d16a092b -r e75eb64d0ca1 Eigen/src/Geometry/AlignedBox.h
--- a/Eigen/src/Geometry/AlignedBox.h	Mon Jan 04 13:57:49 2010 +0100
+++ b/Eigen/src/Geometry/AlignedBox.h	Mon Jan 04 14:00:56 2010 +0100
@@ -43,23 +43,49 @@
 public:
 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
   enum { AmbientDimAtCompileTime = _AmbientDim };
-  typedef _Scalar Scalar;
-  typedef typename NumTraits<Scalar>::Real RealScalar;
-  typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
+  typedef _Scalar                                   Scalar;
+  typedef typename NumTraits<Scalar>::Real          RealScalar;
+  typedef typename NumTraits<Scalar>::FloatingPoint FloatingPoint;
+  typedef Matrix<Scalar,AmbientDimAtCompileTime,1>  VectorType;
+
+  /** Define constants to name the corners of a 1D, 2D or 3D axis aligned bounding box */
+  typedef enum
+  {
+  /** 1D names */
+  Min=0, Max=1,
+
+  /** Added names for 2D */
+  BottomLeft=0, BottomRight=1,
+  TopLeft=2, TopRight=3,
+
+  /** Added names for 3D */
+  BottomLeftFloor=0, BottomRightFloor=1,
+  TopLeftFloor=2, TopRightFloor=3,
+  BottomLeftCeil=4, BottomRightCeil=5,
+  TopLeftCeil=6, TopRightCeil=7
+  } CornerType;
+
 
   /** Default constructor initializing a null box. */
   inline explicit AlignedBox()
-  { if (AmbientDimAtCompileTime!=Dynamic) setNull(); }
+  { if (AmbientDimAtCompileTime!=Dynamic) setEmpty(); }
 
   /** Constructs a null box with \a _dim the dimension of the ambient space. */
   inline explicit AlignedBox(int _dim) : m_min(_dim), m_max(_dim)
-  { setNull(); }
+  { setEmpty(); }
 
   /** Constructs a box with extremities \a _min and \a _max. */
-  inline AlignedBox(const VectorType& _min, const VectorType& _max) : m_min(_min), m_max(_max) {}
+  template<typename OtherVectorType1, typename OtherVectorType2>
+  inline AlignedBox(const OtherVectorType1& _min, const OtherVectorType2& _max) : m_min(_min), m_max(_max) {}
 
   /** Constructs a box containing a single point \a p. */
-  inline explicit AlignedBox(const VectorType& p) : m_min(p), m_max(p) {}
+  template<typename Derived>
+  inline explicit AlignedBox(const MatrixBase<Derived>& a_p)
+  {
+    const typename ei_nested<Derived,2>::type p(a_p.derived());
+    m_min = p;
+    m_max = p;
+  }
 
   ~AlignedBox() {}
 
@@ -67,14 +93,21 @@
   inline int dim() const { return AmbientDimAtCompileTime==Dynamic ? m_min.size()-1 : AmbientDimAtCompileTime; }
 
   /** \returns true if the box is null, i.e, empty. */
-  inline bool isNull() const { return (m_min.cwise() > m_max).any(); }
+  inline bool isNull() const { return isEmpty(); }
 
   /** Makes \c *this a null/empty box. */
-  inline void setNull()
+  inline void setNull() { setEmpty(); }
+
+  /** \returns true if the box is empty. */
+  inline bool isEmpty() const { return (m_min.cwise() > m_max).any(); }
+
+  /** Makes \c *this an empty box. */
+  inline void setEmpty()
   {
-    m_min.setConstant( std::numeric_limits<Scalar>::max());
-    m_max.setConstant(-std::numeric_limits<Scalar>::max());
+    m_min.setConstant( highest<Scalar>() );
+    m_max.setConstant( lowest<Scalar>() );
   }
+
 
   /** \returns the minimal corner */
   inline const VectorType& min() const { return m_min; }
@@ -86,65 +119,162 @@
   inline VectorType& max() { return m_max; }
 
   /** \returns the center of the box */
-  inline VectorType center() const { return (m_min + m_max) / 2; }
+  inline
+    const CwiseUnaryOp<ei_scalar_quotient1_op<Scalar>,
+          CwiseBinaryOp<
+            ei_scalar_sum_op<Scalar>, VectorType, VectorType> > center() const
+  {
+    return (m_min+m_max)/2;
+  }
+
+  /** \returns the lengths of the sides of the bounding box.
+   * Note that this function does not get the same
+   * result for integral or floating scalar types: see
+   */
+  inline
+  const CwiseBinaryOp< ei_scalar_difference_op<Scalar>, VectorType, VectorType> sizes() const
+  {
+    return m_max - m_min;
+  }
+
+  /** \returns the volume of the bounding box */
+  inline Scalar volume() const {
+    return sizes().prod(); }
+
+  /** \returns an expression for the bounding box diagonal vector
+   * if the length of the diagonal is needed: diagonal().norm()
+   * will provide it.
+   * */
+  inline CwiseBinaryOp< ei_scalar_difference_op<Scalar>, VectorType, VectorType> diagonal() const {
+    return sizes(); }
+
+  /** \returns the vertex of the bounding box at the corner defined by
+   * the corner-id corner. It works only for a 1D, 2D or 3D bounding box.
+   * For 1D bounding boxes corners are named by 2 enum constants:
+   * BottomLeft and BottomRight.
+   * For 2D bounding boxes, corners are named by 4 enum constants:
+   * BottomLeft, BottomRight, TopLeft, TopRight.
+   * For 3D bounding boxes, the following names are added:
+   * BottomLeftCeil, BottomRightCeil, TopLeftCeil, TopRightCeil
+   * */
+  inline
+  VectorType corner( CornerType corner ) const
+  {
+    EIGEN_STATIC_ASSERT( _AmbientDim <= 3,
+        THIS_METHOD_IS_ONLY_FOR_VECTORS_OF_A_SPECIFIC_SIZE );
+
+    VectorType res;
+
+    int mult=1;
+    for( int d=0; d<dim(); ++d )
+    {
+      if( mult & corner ){
+        res[d] = m_max[d]; }
+      else{
+        res[d] = m_min[d]; }
+      mult*=2;
+    }
+    return res;
+  }
+
+  /** \returns a random point inside the bounding box sampled with
+   * a uniform distribution */
+  inline
+  VectorType sample() const
+  {
+    VectorType r;
+    for( int d=0; d<dim(); ++d )
+    {
+      if( NumTraits<Scalar>::HasFloatingPoint )
+      {
+        r[d] = m_min[d] + (m_max[d]-m_min[d])
+            * ( ei_random<Scalar>() + ei_random_amplitude<Scalar>() )
+            / (2*ei_random_amplitude<Scalar>() );
+      }
+      else{
+        r[d] = ei_random( m_min[d], m_max[d] ); }
+    }
+    return r;
+  }
 
   /** \returns true if the point \a p is inside the box \c *this. */
-  inline bool contains(const VectorType& p) const
-  { return (m_min.cwise()<=p).all() && (p.cwise()<=m_max).all(); }
+  template<typename Derived>
+  inline bool contains(const MatrixBase<Derived>& a_p) const
+  {
+    const typename ei_nested<Derived,2>::type p(a_p.derived());
+    return (m_min.cwise()<=p).all() && (p.cwise()<=m_max).all();
+  }
 
   /** \returns true if the box \a b is entirely inside the box \c *this. */
-  inline bool contains(const AlignedBox& b) const
+  template<typename OtherScalar>
+  inline bool contains(const AlignedBox<OtherScalar,_AmbientDim>& b) const
   { return (m_min.cwise()<=b.min()).all() && (b.max().cwise()<=m_max).all(); }
 
   /** Extends \c *this such that it contains the point \a p and returns a reference to \c *this. */
-  inline AlignedBox& extend(const VectorType& p)
-  { m_min = m_min.cwise().min(p); m_max = m_max.cwise().max(p); return *this; }
+  template<typename Derived>
+  inline AlignedBox& extend(const MatrixBase<Derived>& a_p)
+  {
+    const typename ei_nested<Derived,2>::type p(a_p.derived());
+    m_min = m_min.cwise().min(p); m_max = m_max.cwise().max(p); return *this;
+  }
 
   /** Extends \c *this such that it contains the box \a b and returns a reference to \c *this. */
-  inline AlignedBox& extend(const AlignedBox& b)
+  template<typename OtherScalar>
+  inline AlignedBox& extend(const AlignedBox<OtherScalar,_AmbientDim>& b)
   { m_min = m_min.cwise().min(b.m_min); m_max = m_max.cwise().max(b.m_max); return *this; }
 
   /** Clamps \c *this by the box \a b and returns a reference to \c *this. */
-  inline AlignedBox& clamp(const AlignedBox& b)
+  template<typename OtherScalar>
+  inline AlignedBox& clamp(const AlignedBox<OtherScalar,_AmbientDim>& b)
   { m_min = m_min.cwise().max(b.m_min); m_max = m_max.cwise().min(b.m_max); return *this; }
 
   /** Returns an AlignedBox that is the intersection of \a b and \c *this */
-  inline AlignedBox intersection(const AlignedBox &b) const
+  template<typename OtherScalar>
+  inline AlignedBox intersection(const AlignedBox<OtherScalar,_AmbientDim> &b) const
   { return AlignedBox(m_min.cwise().max(b.m_min), m_max.cwise().min(b.m_max)); }
 
   /** Returns an AlignedBox that is the union of \a b and \c *this */
-  inline AlignedBox merged(const AlignedBox &b) const
+  template<typename OtherScalar>
+  inline AlignedBox merged(const AlignedBox<OtherScalar,_AmbientDim> &b) const
   { return AlignedBox(m_min.cwise().min(b.m_min), m_max.cwise().max(b.m_max)); }
 
   /** Translate \c *this by the vector \a t and returns a reference to \c *this. */
-  inline AlignedBox& translate(const VectorType& t)
-  { m_min += t; m_max += t; return *this; }
+  template<typename Derived>
+  inline AlignedBox& translate(const MatrixBase<Derived>& a_t)
+  {
+    const typename ei_nested<Derived,2>::type t(a_t.derived());
+    m_min += t; m_max += t; return *this;
+  }
 
   /** \returns the squared distance between the point \a p and the box \c *this,
     * and zero if \a p is inside the box.
     * \sa exteriorDistance()
     */
-  inline Scalar squaredExteriorDistance(const VectorType& p) const;
+  template<typename Derived>
+  inline Scalar squaredExteriorDistance(const MatrixBase<Derived>& a_p) const;
 
   /** \returns the squared distance between the boxes \a b and \c *this,
     * and zero if the boxes intersect.
     * \sa exteriorDistance()
     */
-  inline Scalar squaredExteriorDistance(const AlignedBox& b) const;
+  template<typename OtherScalar>
+  inline Scalar squaredExteriorDistance(const AlignedBox<OtherScalar,_AmbientDim>& b) const;
 
   /** \returns the distance between the point \a p and the box \c *this,
     * and zero if \a p is inside the box.
     * \sa squaredExteriorDistance()
     */
-  inline Scalar exteriorDistance(const VectorType& p) const
-  { return ei_sqrt(squaredExteriorDistance(p)); }
+  template<typename Derived>
+  inline FloatingPoint exteriorDistance(const MatrixBase<Derived>& p) const
+  { return ei_sqrt(FloatingPoint(squaredExteriorDistance(p))); }
 
   /** \returns the distance between the boxes \a b and \c *this,
     * and zero if the boxes intersect.
     * \sa squaredExteriorDistance()
     */
-  inline Scalar exteriorDistance(const AlignedBox& b) const
-  { return ei_sqrt(squaredExteriorDistance(b)); }
+  template<typename OtherScalar>
+  inline FloatingPoint exteriorDistance(const AlignedBox<OtherScalar,_AmbientDim>& b) const
+  { return ei_sqrt(FloatingPoint(squaredExteriorDistance(b))); }
 
   /** \returns \c *this with scalar type casted to \a NewScalarType
     *
@@ -171,7 +301,7 @@
     * determined by \a prec.
     *
     * \sa MatrixBase::isApprox() */
-  bool isApprox(const AlignedBox& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
+  bool isApprox(const AlignedBox& other, RealScalar prec = dummy_precision<Scalar>()) const
   { return m_min.isApprox(other.m_min, prec) && m_max.isApprox(other.m_max, prec); }
 
 protected:
@@ -179,32 +309,51 @@
   VectorType m_min, m_max;
 };
 
-template<typename Scalar,int AmbiantDim>
-inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const VectorType& p) const
+
+
+template<typename _Scalar,int _AmbientDim>
+template<typename Derived>
+inline _Scalar AlignedBox<_Scalar,_AmbientDim>::squaredExteriorDistance(
+    const MatrixBase<Derived>& a_p) const
+{
+  const typename ei_nested<Derived,2*_AmbientDim>::type p(a_p.derived());
+  Scalar dist2 = 0.;
+  Scalar aux;
+  for (int k=0; k<dim(); ++k)
+  {
+    if( m_min[k] > p[k] )
+    {
+      aux = m_min[k] - p[k];
+      dist2 += aux*aux;
+    }
+    else if( p[k] > m_max[k] )
+    {
+      aux = p[k] - m_max[k];
+      dist2 += aux*aux;
+    }
+  }
+  return dist2;
+}
+
+template<typename _Scalar,int _AmbientDim>
+template<typename OtherScalarType>
+inline _Scalar AlignedBox<_Scalar,_AmbientDim>::squaredExteriorDistance(
+    const AlignedBox<OtherScalarType,_AmbientDim>& b) const
 {
   Scalar dist2 = 0.;
   Scalar aux;
   for (int k=0; k<dim(); ++k)
   {
-    if ((aux = (p[k]-m_min[k]))<0.)
+    if( m_min[k] > b.m_max[k] )
+    {
+      aux = m_min[k] - b.m_max[k];
       dist2 += aux*aux;
-    else if ( (aux = (m_max[k]-p[k]))<0. )
+    }
+    else if( b.m_min[k] > m_max[k] )
+    {
+      aux = b.m_min[k] - m_max[k];
       dist2 += aux*aux;
-  }
-  return dist2;
-}
-
-template<typename Scalar,int AmbiantDim>
-inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const AlignedBox& b) const
-{
-  Scalar dist2 = 0.;
-  Scalar aux;
-  for (int k=0; k<dim(); ++k)
-  {
-    if ((aux = (b.m_min[k]-m_max[k]))>0.)
-      dist2 += aux*aux;
-    else if ( (aux = (m_min[k]-b.m_max[k]))>0. )
-      dist2 += aux*aux;
+    }
   }
   return dist2;
 }
# HG changeset patch
# User Manuel Yguel <manuel.yguel@xxxxxxxxx>
# Date 1262610114 -3600
# Node ID f9210d03ff0399a37e9d3cde66d8ac2881329c30
# Parent  e75eb64d0ca189f7beb0a1f13e34b7545f2da2c3
Provide some tests for integral aligned bounding boxes.

diff -r e75eb64d0ca1 -r f9210d03ff03 test/geo_alignedbox.cpp
--- a/test/geo_alignedbox.cpp	Mon Jan 04 14:00:56 2010 +0100
+++ b/test/geo_alignedbox.cpp	Mon Jan 04 14:01:54 2010 +0100
@@ -27,6 +27,9 @@
 #include <Eigen/LU>
 #include <Eigen/QR>

+#include<iostream>
+using namespace std;
+
 template<typename BoxType> void alignedbox(const BoxType& _box)
 {
   /* this test covers the following files:
@@ -40,6 +43,8 @@

   VectorType p0 = VectorType::Random(dim);
   VectorType p1 = VectorType::Random(dim);
+  while( p1 == p0 ){
+      p1 =  VectorType::Random(dim); }
   RealScalar s1 = ei_random<RealScalar>(0,1);

   BoxType b0(dim);
@@ -49,20 +54,13 @@
   b0.extend(p0);
   b0.extend(p1);
   VERIFY(b0.contains(p0*s1+(Scalar(1)-s1)*p1));
-  VERIFY(!b0.contains(p0 + (1+s1)*(p1-p0)));
+  VERIFY(!b0.contains(p0 + (2+s1)*(p1-p0)));

   (b2 = b0).extend(b1);
   VERIFY(b2.contains(b0));
   VERIFY(b2.contains(b1));
   VERIFY_IS_APPROX(b2.clamp(b0), b0);

-  // casting
-  const int Dim = BoxType::AmbientDimAtCompileTime;
-  typedef typename GetDifferentType<Scalar>::type OtherScalar;
-  AlignedBox<OtherScalar,Dim> hp1f = b0.template cast<OtherScalar>();
-  VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),b0);
-  AlignedBox<Scalar,Dim> hp1d = b0.template cast<Scalar>();
-  VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),b0);

   // alignment -- make sure there is no memory alignment assertion
   BoxType *bp0 = new BoxType(dim);
@@ -70,13 +68,117 @@
   bp0->extend(*bp1);
   delete bp0;
   delete bp1;
+
+  // sampling
+  for( int i=0; i<10; ++i )
+  {
+      VectorType r = b0.sample();
+      VERIFY(b0.contains(r));
+  }
+
 }
+
+
+
+template<typename BoxType>
+void alignedboxCastTests(const BoxType& _box)
+{
+  // casting
+  const int dim = _box.dim();
+  typedef typename BoxType::Scalar Scalar;
+  typedef typename NumTraits<Scalar>::Real RealScalar;
+  typedef Matrix<Scalar, BoxType::AmbientDimAtCompileTime, 1> VectorType;
+
+  VectorType p0 = VectorType::Random(dim);
+  VectorType p1 = VectorType::Random(dim);
+
+  BoxType b0(dim);
+
+  b0.extend(p0);
+  b0.extend(p1);
+
+  const int Dim = BoxType::AmbientDimAtCompileTime;
+  typedef typename GetDifferentType<Scalar>::type OtherScalar;
+  AlignedBox<OtherScalar,Dim> hp1f = b0.template cast<OtherScalar>();
+  VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),b0);
+  AlignedBox<Scalar,Dim> hp1d = b0.template cast<Scalar>();
+  VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),b0);
+}
+
+
+void specificTest1()
+{
+    Vector2f m; m << -1.0f, -2.0f;
+    Vector2f M; M <<  1.0f,  5.0f;
+
+    typedef AlignedBox<float,2>  BoxType;
+    BoxType box( m, M );
+
+    Vector2f sides = M-m;
+    VERIFY_IS_APPROX(sides, box.sizes() );
+    VERIFY_IS_APPROX(sides[1], box.sizes()[1] );
+    VERIFY_IS_APPROX(sides[1], box.sizes().maxCoeff() );
+    VERIFY_IS_APPROX(sides[0], box.sizes().minCoeff() );
+
+    VERIFY_IS_APPROX( 14.0f, box.volume() );
+    VERIFY_IS_APPROX( 53.0f, box.diagonal().squaredNorm() );
+    VERIFY_IS_APPROX( ei_sqrt( 53.0f ), box.diagonal().norm() );
+
+    VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeft ) );
+    VERIFY_IS_APPROX( M, box.corner( BoxType::TopRight ) );
+    Vector2f bottomRight; bottomRight << M[0], m[1];
+    Vector2f topLeft; topLeft << m[0], M[1];
+    VERIFY_IS_APPROX( bottomRight, box.corner( BoxType::BottomRight ) );
+    VERIFY_IS_APPROX( topLeft, box.corner( BoxType::TopLeft ) );
+}
+
+
+void specificTest2()
+{
+    Vector3i m; m << -1, -2, 0;
+    Vector3i M; M <<  1,  5, 3;
+
+    typedef AlignedBox<int,3>  BoxType;
+    BoxType box( m, M );
+
+    Vector3i sides = M-m;
+    VERIFY_IS_APPROX(sides, box.sizes() );
+    VERIFY_IS_APPROX(sides[1], box.sizes()[1] );
+    VERIFY_IS_APPROX(sides[1], box.sizes().maxCoeff() );
+    VERIFY_IS_APPROX(sides[0], box.sizes().minCoeff() );
+
+    VERIFY_IS_APPROX( 42, box.volume() );
+    VERIFY_IS_APPROX( 62, box.diagonal().squaredNorm() );
+
+    VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeftFloor ) );
+    VERIFY_IS_APPROX( M, box.corner( BoxType::TopRightCeil ) );
+    Vector3i bottomRightFloor; bottomRightFloor << M[0], m[1], m[2];
+    Vector3i topLeftFloor; topLeftFloor << m[0], M[1], m[2];
+    VERIFY_IS_APPROX( bottomRightFloor, box.corner( BoxType::BottomRightFloor ) );
+    VERIFY_IS_APPROX( topLeftFloor, box.corner( BoxType::TopLeftFloor ) );
+}
+

 void test_geo_alignedbox()
 {
-  for(int i = 0; i < g_repeat; i++) {
+  for(int i = 0; i < g_repeat; i++)
+  {
     CALL_SUBTEST_1( alignedbox(AlignedBox<float,2>()) );
-    CALL_SUBTEST_2( alignedbox(AlignedBox<float,3>()) );
-    CALL_SUBTEST_3( alignedbox(AlignedBox<double,4>()) );
+    CALL_SUBTEST_2( alignedboxCastTests(AlignedBox<float,2>()) );
+
+    CALL_SUBTEST_3( alignedbox(AlignedBox<float,3>()) );
+    CALL_SUBTEST_4( alignedboxCastTests(AlignedBox<float,3>()) );
+
+    CALL_SUBTEST_5( alignedbox(AlignedBox<double,4>()) );
+    CALL_SUBTEST_6( alignedboxCastTests(AlignedBox<double,4>()) );
+
+    CALL_SUBTEST_7( alignedbox(AlignedBox<double,1>()) );
+    CALL_SUBTEST_8( alignedboxCastTests(AlignedBox<double,1>()) );
+
+    CALL_SUBTEST_9( alignedbox(AlignedBox<int,1>()) );
+    CALL_SUBTEST_10( alignedbox(AlignedBox<int,2>()) );
+    CALL_SUBTEST_11( alignedbox(AlignedBox<int,3>()) );
   }
+  CALL_SUBTEST_12( specificTest1() );
+  CALL_SUBTEST_13( specificTest2() );
 }


Mail converted by MHonArc 2.6.19+ http://listengine.tuxfamily.org/