• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #include "main.h"
11 #include <Eigen/Geometry>
12 #include <Eigen/LU>
13 #include <Eigen/QR>
14 
15 #include<iostream>
16 using namespace std;
17 
18 template<typename T> EIGEN_DONT_INLINE
kill_extra_precision(T & x)19 void kill_extra_precision(T& x) { eigen_assert((void*)(&x) != (void*)0); }
20 
21 
alignedbox(const BoxType & _box)22 template<typename BoxType> void alignedbox(const BoxType& _box)
23 {
24   /* this test covers the following files:
25      AlignedBox.h
26   */
27   typedef typename BoxType::Index Index;
28   typedef typename BoxType::Scalar Scalar;
29   typedef typename NumTraits<Scalar>::Real RealScalar;
30   typedef Matrix<Scalar, BoxType::AmbientDimAtCompileTime, 1> VectorType;
31 
32   const Index dim = _box.dim();
33 
34   VectorType p0 = VectorType::Random(dim);
35   VectorType p1 = VectorType::Random(dim);
36   while( p1 == p0 ){
37       p1 =  VectorType::Random(dim); }
38   RealScalar s1 = internal::random<RealScalar>(0,1);
39 
40   BoxType b0(dim);
41   BoxType b1(VectorType::Random(dim),VectorType::Random(dim));
42   BoxType b2;
43 
44   kill_extra_precision(b1);
45   kill_extra_precision(p0);
46   kill_extra_precision(p1);
47 
48   b0.extend(p0);
49   b0.extend(p1);
50   VERIFY(b0.contains(p0*s1+(Scalar(1)-s1)*p1));
51   VERIFY(b0.contains(b0.center()));
52   VERIFY_IS_APPROX(b0.center(),(p0+p1)/Scalar(2));
53 
54   (b2 = b0).extend(b1);
55   VERIFY(b2.contains(b0));
56   VERIFY(b2.contains(b1));
57   VERIFY_IS_APPROX(b2.clamp(b0), b0);
58 
59   // intersection
60   BoxType box1(VectorType::Random(dim));
61   box1.extend(VectorType::Random(dim));
62   BoxType box2(VectorType::Random(dim));
63   box2.extend(VectorType::Random(dim));
64 
65   VERIFY(box1.intersects(box2) == !box1.intersection(box2).isEmpty());
66 
67   // alignment -- make sure there is no memory alignment assertion
68   BoxType *bp0 = new BoxType(dim);
69   BoxType *bp1 = new BoxType(dim);
70   bp0->extend(*bp1);
71   delete bp0;
72   delete bp1;
73 
74   // sampling
75   for( int i=0; i<10; ++i )
76   {
77       VectorType r = b0.sample();
78       VERIFY(b0.contains(r));
79   }
80 
81 }
82 
83 
84 
85 template<typename BoxType>
alignedboxCastTests(const BoxType & _box)86 void alignedboxCastTests(const BoxType& _box)
87 {
88   // casting
89   typedef typename BoxType::Index Index;
90   typedef typename BoxType::Scalar Scalar;
91   typedef Matrix<Scalar, BoxType::AmbientDimAtCompileTime, 1> VectorType;
92 
93   const Index dim = _box.dim();
94 
95   VectorType p0 = VectorType::Random(dim);
96   VectorType p1 = VectorType::Random(dim);
97 
98   BoxType b0(dim);
99 
100   b0.extend(p0);
101   b0.extend(p1);
102 
103   const int Dim = BoxType::AmbientDimAtCompileTime;
104   typedef typename GetDifferentType<Scalar>::type OtherScalar;
105   AlignedBox<OtherScalar,Dim> hp1f = b0.template cast<OtherScalar>();
106   VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),b0);
107   AlignedBox<Scalar,Dim> hp1d = b0.template cast<Scalar>();
108   VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),b0);
109 }
110 
111 
specificTest1()112 void specificTest1()
113 {
114     Vector2f m; m << -1.0f, -2.0f;
115     Vector2f M; M <<  1.0f,  5.0f;
116 
117     typedef AlignedBox2f  BoxType;
118     BoxType box( m, M );
119 
120     Vector2f sides = M-m;
121     VERIFY_IS_APPROX(sides, box.sizes() );
122     VERIFY_IS_APPROX(sides[1], box.sizes()[1] );
123     VERIFY_IS_APPROX(sides[1], box.sizes().maxCoeff() );
124     VERIFY_IS_APPROX(sides[0], box.sizes().minCoeff() );
125 
126     VERIFY_IS_APPROX( 14.0f, box.volume() );
127     VERIFY_IS_APPROX( 53.0f, box.diagonal().squaredNorm() );
128     VERIFY_IS_APPROX( std::sqrt( 53.0f ), box.diagonal().norm() );
129 
130     VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeft ) );
131     VERIFY_IS_APPROX( M, box.corner( BoxType::TopRight ) );
132     Vector2f bottomRight; bottomRight << M[0], m[1];
133     Vector2f topLeft; topLeft << m[0], M[1];
134     VERIFY_IS_APPROX( bottomRight, box.corner( BoxType::BottomRight ) );
135     VERIFY_IS_APPROX( topLeft, box.corner( BoxType::TopLeft ) );
136 }
137 
138 
specificTest2()139 void specificTest2()
140 {
141     Vector3i m; m << -1, -2, 0;
142     Vector3i M; M <<  1,  5, 3;
143 
144     typedef AlignedBox3i  BoxType;
145     BoxType box( m, M );
146 
147     Vector3i sides = M-m;
148     VERIFY_IS_APPROX(sides, box.sizes() );
149     VERIFY_IS_APPROX(sides[1], box.sizes()[1] );
150     VERIFY_IS_APPROX(sides[1], box.sizes().maxCoeff() );
151     VERIFY_IS_APPROX(sides[0], box.sizes().minCoeff() );
152 
153     VERIFY_IS_APPROX( 42, box.volume() );
154     VERIFY_IS_APPROX( 62, box.diagonal().squaredNorm() );
155 
156     VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeftFloor ) );
157     VERIFY_IS_APPROX( M, box.corner( BoxType::TopRightCeil ) );
158     Vector3i bottomRightFloor; bottomRightFloor << M[0], m[1], m[2];
159     Vector3i topLeftFloor; topLeftFloor << m[0], M[1], m[2];
160     VERIFY_IS_APPROX( bottomRightFloor, box.corner( BoxType::BottomRightFloor ) );
161     VERIFY_IS_APPROX( topLeftFloor, box.corner( BoxType::TopLeftFloor ) );
162 }
163 
164 
test_geo_alignedbox()165 void test_geo_alignedbox()
166 {
167   for(int i = 0; i < g_repeat; i++)
168   {
169     CALL_SUBTEST_1( alignedbox(AlignedBox2f()) );
170     CALL_SUBTEST_2( alignedboxCastTests(AlignedBox2f()) );
171 
172     CALL_SUBTEST_3( alignedbox(AlignedBox3f()) );
173     CALL_SUBTEST_4( alignedboxCastTests(AlignedBox3f()) );
174 
175     CALL_SUBTEST_5( alignedbox(AlignedBox4d()) );
176     CALL_SUBTEST_6( alignedboxCastTests(AlignedBox4d()) );
177 
178     CALL_SUBTEST_7( alignedbox(AlignedBox1d()) );
179     CALL_SUBTEST_8( alignedboxCastTests(AlignedBox1d()) );
180 
181     CALL_SUBTEST_9( alignedbox(AlignedBox1i()) );
182     CALL_SUBTEST_10( alignedbox(AlignedBox2i()) );
183     CALL_SUBTEST_11( alignedbox(AlignedBox3i()) );
184 
185     CALL_SUBTEST_14( alignedbox(AlignedBox<double,Dynamic>(4)) );
186   }
187   CALL_SUBTEST_12( specificTest1() );
188   CALL_SUBTEST_13( specificTest2() );
189 }
190