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.

126 lines
4.8 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra. Eigen itself is part of the KDE project.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
  5. // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #include "main.h"
  11. #include <Eigen/Geometry>
  12. #include <Eigen/LU>
  13. #include <Eigen/QR>
  14. template<typename HyperplaneType> void hyperplane(const HyperplaneType& _plane)
  15. {
  16. /* this test covers the following files:
  17. Hyperplane.h
  18. */
  19. const int dim = _plane.dim();
  20. typedef typename HyperplaneType::Scalar Scalar;
  21. typedef typename NumTraits<Scalar>::Real RealScalar;
  22. typedef Matrix<Scalar, HyperplaneType::AmbientDimAtCompileTime, 1> VectorType;
  23. typedef Matrix<Scalar, HyperplaneType::AmbientDimAtCompileTime,
  24. HyperplaneType::AmbientDimAtCompileTime> MatrixType;
  25. VectorType p0 = VectorType::Random(dim);
  26. VectorType p1 = VectorType::Random(dim);
  27. VectorType n0 = VectorType::Random(dim).normalized();
  28. VectorType n1 = VectorType::Random(dim).normalized();
  29. HyperplaneType pl0(n0, p0);
  30. HyperplaneType pl1(n1, p1);
  31. HyperplaneType pl2 = pl1;
  32. Scalar s0 = ei_random<Scalar>();
  33. Scalar s1 = ei_random<Scalar>();
  34. VERIFY_IS_APPROX( n1.eigen2_dot(n1), Scalar(1) );
  35. VERIFY_IS_MUCH_SMALLER_THAN( pl0.absDistance(p0), Scalar(1) );
  36. VERIFY_IS_APPROX( pl1.signedDistance(p1 + n1 * s0), s0 );
  37. VERIFY_IS_MUCH_SMALLER_THAN( pl1.signedDistance(pl1.projection(p0)), Scalar(1) );
  38. VERIFY_IS_MUCH_SMALLER_THAN( pl1.absDistance(p1 + pl1.normal().unitOrthogonal() * s1), Scalar(1) );
  39. // transform
  40. if (!NumTraits<Scalar>::IsComplex)
  41. {
  42. MatrixType rot = MatrixType::Random(dim,dim).qr().matrixQ();
  43. Scaling<Scalar,HyperplaneType::AmbientDimAtCompileTime> scaling(VectorType::Random());
  44. Translation<Scalar,HyperplaneType::AmbientDimAtCompileTime> translation(VectorType::Random());
  45. pl2 = pl1;
  46. VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot).absDistance(rot * p1), Scalar(1) );
  47. pl2 = pl1;
  48. VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot,Isometry).absDistance(rot * p1), Scalar(1) );
  49. pl2 = pl1;
  50. VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot*scaling).absDistance((rot*scaling) * p1), Scalar(1) );
  51. pl2 = pl1;
  52. VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot*scaling*translation)
  53. .absDistance((rot*scaling*translation) * p1), Scalar(1) );
  54. pl2 = pl1;
  55. VERIFY_IS_MUCH_SMALLER_THAN( pl2.transform(rot*translation,Isometry)
  56. .absDistance((rot*translation) * p1), Scalar(1) );
  57. }
  58. // casting
  59. const int Dim = HyperplaneType::AmbientDimAtCompileTime;
  60. typedef typename GetDifferentType<Scalar>::type OtherScalar;
  61. Hyperplane<OtherScalar,Dim> hp1f = pl1.template cast<OtherScalar>();
  62. VERIFY_IS_APPROX(hp1f.template cast<Scalar>(),pl1);
  63. Hyperplane<Scalar,Dim> hp1d = pl1.template cast<Scalar>();
  64. VERIFY_IS_APPROX(hp1d.template cast<Scalar>(),pl1);
  65. }
  66. template<typename Scalar> void lines()
  67. {
  68. typedef Hyperplane<Scalar, 2> HLine;
  69. typedef ParametrizedLine<Scalar, 2> PLine;
  70. typedef Matrix<Scalar,2,1> Vector;
  71. typedef Matrix<Scalar,3,1> CoeffsType;
  72. for(int i = 0; i < 10; i++)
  73. {
  74. Vector center = Vector::Random();
  75. Vector u = Vector::Random();
  76. Vector v = Vector::Random();
  77. Scalar a = ei_random<Scalar>();
  78. while (ei_abs(a-1) < 1e-4) a = ei_random<Scalar>();
  79. while (u.norm() < 1e-4) u = Vector::Random();
  80. while (v.norm() < 1e-4) v = Vector::Random();
  81. HLine line_u = HLine::Through(center + u, center + a*u);
  82. HLine line_v = HLine::Through(center + v, center + a*v);
  83. // the line equations should be normalized so that a^2+b^2=1
  84. VERIFY_IS_APPROX(line_u.normal().norm(), Scalar(1));
  85. VERIFY_IS_APPROX(line_v.normal().norm(), Scalar(1));
  86. Vector result = line_u.intersection(line_v);
  87. // the lines should intersect at the point we called "center"
  88. VERIFY_IS_APPROX(result, center);
  89. // check conversions between two types of lines
  90. PLine pl(line_u); // gcc 3.3 will commit suicide if we don't name this variable
  91. CoeffsType converted_coeffs(HLine(pl).coeffs());
  92. converted_coeffs *= line_u.coeffs()(0)/converted_coeffs(0);
  93. VERIFY(line_u.coeffs().isApprox(converted_coeffs));
  94. }
  95. }
  96. void test_eigen2_hyperplane()
  97. {
  98. for(int i = 0; i < g_repeat; i++) {
  99. CALL_SUBTEST_1( hyperplane(Hyperplane<float,2>()) );
  100. CALL_SUBTEST_2( hyperplane(Hyperplane<float,3>()) );
  101. CALL_SUBTEST_3( hyperplane(Hyperplane<double,4>()) );
  102. CALL_SUBTEST_4( hyperplane(Hyperplane<std::complex<double>,5>()) );
  103. CALL_SUBTEST_5( lines<float>() );
  104. CALL_SUBTEST_6( lines<double>() );
  105. }
  106. }