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.

136 lines
5.3 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 Benoit Jacob <jacob.benoit.1@gmail.com>
  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/LeastSquares>
  11. template<typename VectorType,
  12. typename HyperplaneType>
  13. void makeNoisyCohyperplanarPoints(int numPoints,
  14. VectorType **points,
  15. HyperplaneType *hyperplane,
  16. typename VectorType::Scalar noiseAmplitude)
  17. {
  18. typedef typename VectorType::Scalar Scalar;
  19. const int size = points[0]->size();
  20. // pick a random hyperplane, store the coefficients of its equation
  21. hyperplane->coeffs().resize(size + 1);
  22. for(int j = 0; j < size + 1; j++)
  23. {
  24. do {
  25. hyperplane->coeffs().coeffRef(j) = ei_random<Scalar>();
  26. } while(ei_abs(hyperplane->coeffs().coeff(j)) < 0.5);
  27. }
  28. // now pick numPoints random points on this hyperplane
  29. for(int i = 0; i < numPoints; i++)
  30. {
  31. VectorType& cur_point = *(points[i]);
  32. do
  33. {
  34. cur_point = VectorType::Random(size)/*.normalized()*/;
  35. // project cur_point onto the hyperplane
  36. Scalar x = - (hyperplane->coeffs().start(size).cwise()*cur_point).sum();
  37. cur_point *= hyperplane->coeffs().coeff(size) / x;
  38. } while( cur_point.norm() < 0.5
  39. || cur_point.norm() > 2.0 );
  40. }
  41. // add some noise to these points
  42. for(int i = 0; i < numPoints; i++ )
  43. *(points[i]) += noiseAmplitude * VectorType::Random(size);
  44. }
  45. template<typename VectorType>
  46. void check_linearRegression(int numPoints,
  47. VectorType **points,
  48. const VectorType& original,
  49. typename VectorType::Scalar tolerance)
  50. {
  51. int size = points[0]->size();
  52. assert(size==2);
  53. VectorType result(size);
  54. linearRegression(numPoints, points, &result, 1);
  55. typename VectorType::Scalar error = (result - original).norm() / original.norm();
  56. VERIFY(ei_abs(error) < ei_abs(tolerance));
  57. }
  58. template<typename VectorType,
  59. typename HyperplaneType>
  60. void check_fitHyperplane(int numPoints,
  61. VectorType **points,
  62. const HyperplaneType& original,
  63. typename VectorType::Scalar tolerance)
  64. {
  65. int size = points[0]->size();
  66. HyperplaneType result(size);
  67. fitHyperplane(numPoints, points, &result);
  68. result.coeffs() *= original.coeffs().coeff(size)/result.coeffs().coeff(size);
  69. typename VectorType::Scalar error = (result.coeffs() - original.coeffs()).norm() / original.coeffs().norm();
  70. std::cout << ei_abs(error) << " xxx " << ei_abs(tolerance) << std::endl;
  71. VERIFY(ei_abs(error) < ei_abs(tolerance));
  72. }
  73. void test_eigen2_regression()
  74. {
  75. for(int i = 0; i < g_repeat; i++)
  76. {
  77. #ifdef EIGEN_TEST_PART_1
  78. {
  79. Vector2f points2f [1000];
  80. Vector2f *points2f_ptrs [1000];
  81. for(int i = 0; i < 1000; i++) points2f_ptrs[i] = &(points2f[i]);
  82. Vector2f coeffs2f;
  83. Hyperplane<float,2> coeffs3f;
  84. makeNoisyCohyperplanarPoints(1000, points2f_ptrs, &coeffs3f, 0.01f);
  85. coeffs2f[0] = -coeffs3f.coeffs()[0]/coeffs3f.coeffs()[1];
  86. coeffs2f[1] = -coeffs3f.coeffs()[2]/coeffs3f.coeffs()[1];
  87. CALL_SUBTEST(check_linearRegression(10, points2f_ptrs, coeffs2f, 0.05f));
  88. CALL_SUBTEST(check_linearRegression(100, points2f_ptrs, coeffs2f, 0.01f));
  89. CALL_SUBTEST(check_linearRegression(1000, points2f_ptrs, coeffs2f, 0.002f));
  90. }
  91. #endif
  92. #ifdef EIGEN_TEST_PART_2
  93. {
  94. Vector2f points2f [1000];
  95. Vector2f *points2f_ptrs [1000];
  96. for(int i = 0; i < 1000; i++) points2f_ptrs[i] = &(points2f[i]);
  97. Hyperplane<float,2> coeffs3f;
  98. makeNoisyCohyperplanarPoints(1000, points2f_ptrs, &coeffs3f, 0.01f);
  99. CALL_SUBTEST(check_fitHyperplane(10, points2f_ptrs, coeffs3f, 0.05f));
  100. CALL_SUBTEST(check_fitHyperplane(100, points2f_ptrs, coeffs3f, 0.01f));
  101. CALL_SUBTEST(check_fitHyperplane(1000, points2f_ptrs, coeffs3f, 0.002f));
  102. }
  103. #endif
  104. #ifdef EIGEN_TEST_PART_3
  105. {
  106. Vector4d points4d [1000];
  107. Vector4d *points4d_ptrs [1000];
  108. for(int i = 0; i < 1000; i++) points4d_ptrs[i] = &(points4d[i]);
  109. Hyperplane<double,4> coeffs5d;
  110. makeNoisyCohyperplanarPoints(1000, points4d_ptrs, &coeffs5d, 0.01);
  111. CALL_SUBTEST(check_fitHyperplane(10, points4d_ptrs, coeffs5d, 0.05));
  112. CALL_SUBTEST(check_fitHyperplane(100, points4d_ptrs, coeffs5d, 0.01));
  113. CALL_SUBTEST(check_fitHyperplane(1000, points4d_ptrs, coeffs5d, 0.002));
  114. }
  115. #endif
  116. #ifdef EIGEN_TEST_PART_4
  117. {
  118. VectorXcd *points11cd_ptrs[1000];
  119. for(int i = 0; i < 1000; i++) points11cd_ptrs[i] = new VectorXcd(11);
  120. Hyperplane<std::complex<double>,Dynamic> *coeffs12cd = new Hyperplane<std::complex<double>,Dynamic>(11);
  121. makeNoisyCohyperplanarPoints(1000, points11cd_ptrs, coeffs12cd, 0.01);
  122. CALL_SUBTEST(check_fitHyperplane(100, points11cd_ptrs, *coeffs12cd, 0.025));
  123. CALL_SUBTEST(check_fitHyperplane(1000, points11cd_ptrs, *coeffs12cd, 0.006));
  124. delete coeffs12cd;
  125. for(int i = 0; i < 1000; i++) delete points11cd_ptrs[i];
  126. }
  127. #endif
  128. }
  129. }