[eigen] Enhanced AlignedBox |
[ Thread Index |
Date Index
| More lists.tuxfamily.org/eigen Archives
]
- To: eigen@xxxxxxxxxxxxxxxxxxx
- Subject: [eigen] Enhanced AlignedBox
- From: Manuel Yguel <manuel.yguel@xxxxxxxxx>
- Date: Sun, 20 Dec 2009 23:30:24 +0100
- Dkim-signature: v=1; a=rsa-sha256; c=relaxed/relaxed; d=gmail.com; s=gamma; h=domainkey-signature:mime-version:received:from:date:message-id :subject:to:content-type; bh=MjWnpIBf8BzevJ9fMoxjSNHd4vM+Mdu7AMGEU7pMFW4=; b=q+30eIIi55lgDoYJe3/uiTOe4HeSi34xW4EqVv0GobR0ufPcY92Mu8OoJ46+Ls43lX qoCZ8XZgQuDnl0Q2lOmJYxzoQpYsvwAcMtRfiZfidILxvIdK3uqPbgasCL5LGffC3Vjd EW6fEKJ6rXCYFxVsFOC1J8DSyHSyfIrhkg1do=
- Domainkey-signature: a=rsa-sha1; c=nofws; d=gmail.com; s=gamma; h=mime-version:from:date:message-id:subject:to:content-type; b=T1pxah2Vd73r+7LFXKiuKI0K9BLLKZmiwm/EPqfQyMj6LN2wq0de5D9VsWbqUoW2nc xsfRMkwbZbdwpKfepsJmqlwtiFj7DZWMh3RKUmr5X2q0GJxovr7fUWxDb/CwlHk8mYNO S58BLClyFiq/sP2/dYJosU0+DmdK1jY0FNKu8=
Hello,
I have written some small pieces of codes to enhance AlignedBox.
The following changes I propose are as follow:
1) Define some new functions:
Scalar size( size_t k ) --> returns the size of the k-th side
Vector sizes() --> returns a vector with all the sides
Vector corner( CornerType t ) --> defines some enum and return the
corresponding corner for 1D,2D and 3D boundingBox
Scalar squaredDiagonale()
Real diagonal()
Scalar volume()
Vector sample() --> sample a point in the boundingBox volume with a
uniform distribution
2) Add the function isEmpty (I do not like very much the isNull name,
it is a matter of taste, so if it does not fit you, let me know...)
3) Semantic: I think that it is usefull to have bounding boxes with
integers: in this case the min and the max are considered inside the
bounding box and a side of the bounding box is max - min + 1 instead
of max - min.
This is a complex semantic problem, but I think it is the more natural choice.
4) Following that line, I have changed the code for the distances
functions such that non-negative integers types are safely handled and
I trust the compiler to produce as optimized code as the one that Gael
has written.
I know that such types are not natively handled by Eigen, but I think
that it is a safe change.
What do you think about those proposals?
If you agree with that lines, I plan to write an OrientedBox class
with almost same interface.
Also, I would like to launch a topic about the way random numbers are
generated because when writing sample(), I found that we can enhance
the way it is done.
Is this already an active topic ?
- best regards,
Manuel
# HG changeset patch
# User Manuel Yguel <manuel.yguel@xxxxxxxxx>
# Date 1261346899 -3600
# Node ID 8d31bd6406a8f9302d977172fb5a3d7d0172110e
# Parent 550cdad8b2a2a6c0a9df4d7a45c8f6aea5ac33a0
Enhanced AlignedBox: takes into account integer AlignedBox and add methods for side, diagonal, volume and corners for 1D, 2D and 3D cases.
diff -r 550cdad8b2a2 -r 8d31bd6406a8 Eigen/src/Geometry/AlignedBox.h
--- a/Eigen/src/Geometry/AlignedBox.h Tue Dec 15 08:09:14 2009 +0100
+++ b/Eigen/src/Geometry/AlignedBox.h Sun Dec 20 23:08:19 2009 +0100
@@ -43,9 +43,25 @@
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 */
+ BottomLeft=0, BottomRight=1,
+
+ /** Added names for 2D */
+ TopLeft=2, TopRight=3,
+
+ /** Added names for 3D */
+ BottomLeftCeil=4, BottomRightCeil=5,
+ TopLeftCeil=6, TopRightCeil=7
+ } CornerType;
+
/** Default constructor initializing a null box. */
inline explicit AlignedBox()
@@ -53,13 +69,15 @@
/** 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 OtherVectorType>
+ inline explicit AlignedBox(const OtherVectorType& p) : m_min(p), m_max(p) {}
~AlignedBox() {}
@@ -76,6 +94,17 @@
m_max.setConstant(-std::numeric_limits<Scalar>::max());
}
+ /** \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());
+ }
+
+
/** \returns the minimal corner */
inline const VectorType& min() const { return m_min; }
/** \returns a non const reference to the minimal corner */
@@ -88,63 +117,163 @@
/** \returns the center of the box */
inline VectorType center() const { return (m_min + m_max) / 2; }
+ /** \returns the k-th coordinate of the center of the
+ * bounding box as (max()[k]+min()[k])/2 */
+ inline Scalar center( size_t k ) const { return (m_min[k]+m_max[k])/2; }
+
+ /** \returns the length of the k-th side of the bounding box
+ * Note that this function does not get the same
+ * result for integral or floating scalar types:
+ * - floating types: max()[k]-min()[k]
+ * - integral types: max()[k]-min()[k]+1
+ *
+ */
+ Scalar side( size_t k ) const;
+
+ /** \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
+ */
+ VectorType sides() const;
+
+ /** \returns the length of the largest side of the bounding box */
+ inline Scalar maxSide() const {
+ return sides().maxCoeff(); }
+
+ /** \returns the length of the smallest side of the bounding box */
+ inline Scalar minSide() const {
+ return sides().minCoeff(); }
+
+ /** \returns the volume of the bounding box */
+ inline Scalar volume() const {
+ return sides().prod(); }
+
+ /** \returns the squared length of the bounding box diagonal */
+ inline Scalar squaredDiagonal() const {
+ return sides().squaredNorm(); }
+
+ /** \returns the length of the bounding box diagonal */
+ inline FloatingPoint diagonal() const {
+ return ei_sqrt(FloatingPoint(squaredDiagonal())); }
+
+ /** \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
+ template<typename OtherVectorType>
+ inline bool contains(const OtherVectorType& p) const
{ 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 Os, int Od>
+ inline bool contains(const AlignedBox<Os,Od>& 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)
+ template<typename OtherVectorType>
+ inline AlignedBox& extend(const OtherVectorType& p)
{ 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 Os, int Od>
+ inline AlignedBox& extend(const AlignedBox<Os,Od>& 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 Os, int Od>
+ inline AlignedBox& clamp(const AlignedBox<Os,Od>& 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 Os, int Od>
+ inline AlignedBox intersection(const AlignedBox<Os,Od> &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 Os, int Od>
+ inline AlignedBox merged(const AlignedBox<Os,Od> &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)
+ template<typename OtherVectorType>
+ inline AlignedBox& translate(const OtherVectorType& t)
{ 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 OtherVectorType>
+ inline Scalar squaredExteriorDistance(const OtherVectorType& 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 Os,int Od>
+ inline Scalar squaredExteriorDistance(const AlignedBox<Os,Od>& 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 OtherVectorType>
+ inline FloatingPoint exteriorDistance(const OtherVectorType& 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 Os, int Od>
+ inline FloatingPoint exteriorDistance(const AlignedBox<Os,Od>& b) const
+ { return ei_sqrt(FloatingPoint(squaredExteriorDistance(b))); }
/** \returns \c *this with scalar type casted to \a NewScalarType
*
@@ -171,7 +300,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 +308,72 @@
VectorType m_min, m_max;
};
+
+template<typename Scalar>
+struct ei_computeBoxSideOp EIGEN_EMPTY_STRUCT
+{
+ typedef Scalar result_type;
+ inline
+ Scalar operator()(const Scalar& max, const Scalar& min) const
+ {
+ if( NumTraits<Scalar>::HasFloatingPoint ){ return max - min; }
+ else{ return max - min + 1; }
+ }
+};
+
+
template<typename Scalar,int AmbiantDim>
-inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const VectorType& p) const
+inline Scalar AlignedBox<Scalar,AmbiantDim>::side( size_t k ) const {
+ return ei_computeBoxSideOp<Scalar>()( m_max[k], m_min[k] ); }
+
+template<typename Scalar,int AmbiantDim>
+inline
+typename AlignedBox<Scalar,AmbiantDim>::VectorType
+AlignedBox<Scalar,AmbiantDim>::sides() const
+{
+ return m_max.binaryExpr(m_min, ei_computeBoxSideOp<Scalar>() );
+}
+
+template<typename Scalar,int AmbiantDim>
+template<typename OtherVectorType>
+inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const OtherVectorType& p) const
{
Scalar dist2 = 0.;
Scalar aux;
for (int k=0; k<dim(); ++k)
{
- if ((aux = (p[k]-m_min[k]))<0.)
- dist2 += aux*aux;
- else if ( (aux = (m_max[k]-p[k]))<0. )
- dist2 += aux*aux;
+ 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 AmbiantDim>
-inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const AlignedBox& b) const
+template<typename Os, int Od>
+inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const AlignedBox<Os,Od>& 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;
+ if( m_min[k] > b.m_max[k] )
+ {
+ aux = m_min[k] - b.m_max[k];
+ dist2 += aux*aux;
+ }
+ else if( b.m_min[k] > m_max[k] )
+ {
+ aux = b.m_min[k] - m_max[k];
+ dist2 += aux*aux;
+ }
}
return dist2;
}
# HG changeset patch
# User Manuel Yguel <manuel.yguel@xxxxxxxxx>
# Date 1261346995 -3600
# Node ID 1932912e81d040cae27a40f3f0627629d15df977
# Parent 8d31bd6406a8f9302d977172fb5a3d7d0172110e
Tests modified for the enhanced AlignedBox + (one bug fix in the tests for bounding boxes with one side that is less than 1).
diff -r 8d31bd6406a8 -r 1932912e81d0 test/geo_alignedbox.cpp
--- a/test/geo_alignedbox.cpp Sun Dec 20 23:08:19 2009 +0100
+++ b/test/geo_alignedbox.cpp Sun Dec 20 23:09:55 2009 +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,119 @@
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.sides() );
+ VERIFY_IS_APPROX(sides[1], box.side(1) );
+ VERIFY_IS_APPROX(sides[1], box.maxSide() );
+ VERIFY_IS_APPROX(sides[0], box.minSide() );
+
+ VERIFY_IS_APPROX( 14.0f, box.volume() );
+ VERIFY_IS_APPROX( 53.0f, box.squaredDiagonal() );
+ VERIFY_IS_APPROX( ei_sqrt( 53.0f ), box.diagonal() );
+
+ 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 + Vector3i(1,1,1);
+ VERIFY_IS_APPROX(sides, box.sides() );
+ VERIFY_IS_APPROX(sides[1], box.side(1) );
+ VERIFY_IS_APPROX(sides[1], box.maxSide() );
+ VERIFY_IS_APPROX(sides[0], box.minSide() );
+
+ VERIFY_IS_APPROX( 96, box.volume() );
+ VERIFY_IS_APPROX( 89, box.squaredDiagonal() );
+ VERIFY_IS_APPROX( ei_sqrt( 89.0 ), box.diagonal() );
+
+ VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeft ) );
+ VERIFY_IS_APPROX( M, box.corner( BoxType::TopRightCeil ) );
+ Vector3i bottomRight; bottomRight << M[0], m[1], m[2];
+ Vector3i topLeft; topLeft << m[0], M[1], m[2];
+ VERIFY_IS_APPROX( bottomRight, box.corner( BoxType::BottomRight ) );
+ VERIFY_IS_APPROX( topLeft, box.corner( BoxType::TopLeft ) );
+
+}
+
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() );
}