You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

178 lines
5.1 KiB

  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. #include "main.h"
  10. #include <Eigen/Geometry>
  11. #include <Eigen/LU>
  12. #include <Eigen/QR>
  13. #include<iostream>
  14. using namespace std;
  15. template<typename T> EIGEN_DONT_INLINE
  16. void kill_extra_precision(T& x) { eigen_assert(&x != 0); }
  17. template<typename BoxType> void alignedbox(const BoxType& _box)
  18. {
  19. /* this test covers the following files:
  20. AlignedBox.h
  21. */
  22. typedef typename BoxType::Index Index;
  23. typedef typename BoxType::Scalar Scalar;
  24. typedef typename NumTraits<Scalar>::Real RealScalar;
  25. typedef Matrix<Scalar, BoxType::AmbientDimAtCompileTime, 1> VectorType;
  26. const Index dim = _box.dim();
  27. VectorType p0 = VectorType::Random(dim);
  28. VectorType p1 = VectorType::Random(dim);
  29. while( p1 == p0 ){
  30. p1 = VectorType::Random(dim); }
  31. RealScalar s1 = internal::random<RealScalar>(0,1);
  32. BoxType b0(dim);
  33. BoxType b1(VectorType::Random(dim),VectorType::Random(dim));
  34. BoxType b2;
  35. kill_extra_precision(b1);
  36. kill_extra_precision(p0);
  37. kill_extra_precision(p1);
  38. b0.extend(p0);
  39. b0.extend(p1);
  40. VERIFY(b0.contains(p0*s1+(Scalar(1)-s1)*p1));
  41. (b2 = b0).extend(b1);
  42. VERIFY(b2.contains(b0));
  43. VERIFY(b2.contains(b1));
  44. VERIFY_IS_APPROX(b2.clamp(b0), b0);
  45. // alignment -- make sure there is no memory alignment assertion
  46. BoxType *bp0 = new BoxType(dim);
  47. BoxType *bp1 = new BoxType(dim);
  48. bp0->extend(*bp1);
  49. delete bp0;
  50. delete bp1;
  51. // sampling
  52. for( int i=0; i<10; ++i )
  53. {
  54. VectorType r = b0.sample();
  55. VERIFY(b0.contains(r));
  56. }
  57. }
  58. template<typename BoxType>
  59. void alignedboxCastTests(const BoxType& _box)
  60. {
  61. // casting
  62. typedef typename BoxType::Index Index;
  63. typedef typename BoxType::Scalar Scalar;
  64. typedef Matrix<Scalar, BoxType::AmbientDimAtCompileTime, 1> VectorType;
  65. const Index dim = _box.dim();
  66. VectorType p0 = VectorType::Random(dim);
  67. VectorType p1 = VectorType::Random(dim);
  68. BoxType b0(dim);
  69. b0.extend(p0);
  70. b0.extend(p1);
  71. const int Dim = BoxType::AmbientDimAtCompileTime;
  72. typedef typename GetDifferentType<Scalar>::type OtherScalar;
  73. AlignedBox<OtherScalar,Dim> hp1f = b0.template cast<OtherScalar>();
  74. VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),b0);
  75. AlignedBox<Scalar,Dim> hp1d = b0.template cast<Scalar>();
  76. VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),b0);
  77. }
  78. void specificTest1()
  79. {
  80. Vector2f m; m << -1.0f, -2.0f;
  81. Vector2f M; M << 1.0f, 5.0f;
  82. typedef AlignedBox2f BoxType;
  83. BoxType box( m, M );
  84. Vector2f sides = M-m;
  85. VERIFY_IS_APPROX(sides, box.sizes() );
  86. VERIFY_IS_APPROX(sides[1], box.sizes()[1] );
  87. VERIFY_IS_APPROX(sides[1], box.sizes().maxCoeff() );
  88. VERIFY_IS_APPROX(sides[0], box.sizes().minCoeff() );
  89. VERIFY_IS_APPROX( 14.0f, box.volume() );
  90. VERIFY_IS_APPROX( 53.0f, box.diagonal().squaredNorm() );
  91. VERIFY_IS_APPROX( std::sqrt( 53.0f ), box.diagonal().norm() );
  92. VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeft ) );
  93. VERIFY_IS_APPROX( M, box.corner( BoxType::TopRight ) );
  94. Vector2f bottomRight; bottomRight << M[0], m[1];
  95. Vector2f topLeft; topLeft << m[0], M[1];
  96. VERIFY_IS_APPROX( bottomRight, box.corner( BoxType::BottomRight ) );
  97. VERIFY_IS_APPROX( topLeft, box.corner( BoxType::TopLeft ) );
  98. }
  99. void specificTest2()
  100. {
  101. Vector3i m; m << -1, -2, 0;
  102. Vector3i M; M << 1, 5, 3;
  103. typedef AlignedBox3i BoxType;
  104. BoxType box( m, M );
  105. Vector3i sides = M-m;
  106. VERIFY_IS_APPROX(sides, box.sizes() );
  107. VERIFY_IS_APPROX(sides[1], box.sizes()[1] );
  108. VERIFY_IS_APPROX(sides[1], box.sizes().maxCoeff() );
  109. VERIFY_IS_APPROX(sides[0], box.sizes().minCoeff() );
  110. VERIFY_IS_APPROX( 42, box.volume() );
  111. VERIFY_IS_APPROX( 62, box.diagonal().squaredNorm() );
  112. VERIFY_IS_APPROX( m, box.corner( BoxType::BottomLeftFloor ) );
  113. VERIFY_IS_APPROX( M, box.corner( BoxType::TopRightCeil ) );
  114. Vector3i bottomRightFloor; bottomRightFloor << M[0], m[1], m[2];
  115. Vector3i topLeftFloor; topLeftFloor << m[0], M[1], m[2];
  116. VERIFY_IS_APPROX( bottomRightFloor, box.corner( BoxType::BottomRightFloor ) );
  117. VERIFY_IS_APPROX( topLeftFloor, box.corner( BoxType::TopLeftFloor ) );
  118. }
  119. void test_geo_alignedbox()
  120. {
  121. for(int i = 0; i < g_repeat; i++)
  122. {
  123. CALL_SUBTEST_1( alignedbox(AlignedBox2f()) );
  124. CALL_SUBTEST_2( alignedboxCastTests(AlignedBox2f()) );
  125. CALL_SUBTEST_3( alignedbox(AlignedBox3f()) );
  126. CALL_SUBTEST_4( alignedboxCastTests(AlignedBox3f()) );
  127. CALL_SUBTEST_5( alignedbox(AlignedBox4d()) );
  128. CALL_SUBTEST_6( alignedboxCastTests(AlignedBox4d()) );
  129. CALL_SUBTEST_7( alignedbox(AlignedBox1d()) );
  130. CALL_SUBTEST_8( alignedboxCastTests(AlignedBox1d()) );
  131. CALL_SUBTEST_9( alignedbox(AlignedBox1i()) );
  132. CALL_SUBTEST_10( alignedbox(AlignedBox2i()) );
  133. CALL_SUBTEST_11( alignedbox(AlignedBox3i()) );
  134. }
  135. CALL_SUBTEST_12( specificTest1() );
  136. CALL_SUBTEST_13( specificTest2() );
  137. }