circle_fit.cc 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2023 Google Inc. All rights reserved.
  3. // http://ceres-solver.org/
  4. //
  5. // Redistribution and use in source and binary forms, with or without
  6. // modification, are permitted provided that the following conditions are met:
  7. //
  8. // * Redistributions of source code must retain the above copyright notice,
  9. // this list of conditions and the following disclaimer.
  10. // * Redistributions in binary form must reproduce the above copyright notice,
  11. // this list of conditions and the following disclaimer in the documentation
  12. // and/or other materials provided with the distribution.
  13. // * Neither the name of Google Inc. nor the names of its contributors may be
  14. // used to endorse or promote products derived from this software without
  15. // specific prior written permission.
  16. //
  17. // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  18. // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  19. // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  20. // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
  21. // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
  22. // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
  23. // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
  24. // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
  25. // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
  26. // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
  27. // POSSIBILITY OF SUCH DAMAGE.
  28. //
  29. // Author: keir@google.com (Keir Mierle)
  30. //
  31. // This fits circles to a collection of points, where the error is related to
  32. // the distance of a point from the circle. This uses auto-differentiation to
  33. // take the derivatives.
  34. //
  35. // The input format is simple text. Feed on standard in:
  36. //
  37. // x_initial y_initial r_initial
  38. // x1 y1
  39. // x2 y2
  40. // y3 y3
  41. // ...
  42. //
  43. // And the result after solving will be printed to stdout:
  44. //
  45. // x y r
  46. //
  47. // There are closed form solutions [1] to this problem which you may want to
  48. // consider instead of using this one. If you already have a decent guess, Ceres
  49. // can squeeze down the last bit of error.
  50. //
  51. // [1] http://www.mathworks.com/matlabcentral/fileexchange/5557-circle-fit/content/circfit.m // NOLINT
  52. #include <cstdio>
  53. #include <vector>
  54. #include "ceres/ceres.h"
  55. #include "gflags/gflags.h"
  56. #include "glog/logging.h"
  57. DEFINE_double(robust_threshold,
  58. 0.0,
  59. "Robust loss parameter. Set to 0 for normal squared error (no "
  60. "robustification).");
  61. // The cost for a single sample. The returned residual is related to the
  62. // distance of the point from the circle (passed in as x, y, m parameters).
  63. //
  64. // Note that the radius is parameterized as r = m^2 to constrain the radius to
  65. // positive values.
  66. class DistanceFromCircleCost {
  67. public:
  68. DistanceFromCircleCost(double xx, double yy) : xx_(xx), yy_(yy) {}
  69. template <typename T>
  70. bool operator()(const T* const x,
  71. const T* const y,
  72. const T* const m, // r = m^2
  73. T* residual) const {
  74. // Since the radius is parameterized as m^2, unpack m to get r.
  75. T r = *m * *m;
  76. // Get the position of the sample in the circle's coordinate system.
  77. T xp = xx_ - *x;
  78. T yp = yy_ - *y;
  79. // It is tempting to use the following cost:
  80. //
  81. // residual[0] = r - sqrt(xp*xp + yp*yp);
  82. //
  83. // which is the distance of the sample from the circle. This works
  84. // reasonably well, but the sqrt() adds strong nonlinearities to the cost
  85. // function. Instead, a different cost is used, which while not strictly a
  86. // distance in the metric sense (it has units distance^2) it produces more
  87. // robust fits when there are outliers. This is because the cost surface is
  88. // more convex.
  89. residual[0] = r * r - xp * xp - yp * yp;
  90. return true;
  91. }
  92. private:
  93. // The measured x,y coordinate that should be on the circle.
  94. double xx_, yy_;
  95. };
  96. int main(int argc, char** argv) {
  97. GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
  98. google::InitGoogleLogging(argv[0]);
  99. double x, y, r;
  100. if (scanf("%lg %lg %lg", &x, &y, &r) != 3) {
  101. fprintf(stderr, "Couldn't read first line.\n");
  102. return 1;
  103. }
  104. fprintf(stderr, "Got x, y, r %lg, %lg, %lg\n", x, y, r);
  105. // Save initial values for comparison.
  106. double initial_x = x;
  107. double initial_y = y;
  108. double initial_r = r;
  109. // Parameterize r as m^2 so that it can't be negative.
  110. double m = sqrt(r);
  111. ceres::Problem problem;
  112. // Configure the loss function.
  113. ceres::LossFunction* loss = nullptr;
  114. if (CERES_GET_FLAG(FLAGS_robust_threshold)) {
  115. loss = new ceres::CauchyLoss(CERES_GET_FLAG(FLAGS_robust_threshold));
  116. }
  117. // Add the residuals.
  118. double xx, yy;
  119. int num_points = 0;
  120. while (scanf("%lf %lf\n", &xx, &yy) == 2) {
  121. ceres::CostFunction* cost =
  122. new ceres::AutoDiffCostFunction<DistanceFromCircleCost, 1, 1, 1, 1>(
  123. new DistanceFromCircleCost(xx, yy));
  124. problem.AddResidualBlock(cost, loss, &x, &y, &m);
  125. num_points++;
  126. }
  127. std::cout << "Got " << num_points << " points.\n";
  128. // Build and solve the problem.
  129. ceres::Solver::Options options;
  130. options.max_num_iterations = 500;
  131. options.linear_solver_type = ceres::DENSE_QR;
  132. ceres::Solver::Summary summary;
  133. ceres::Solve(options, &problem, &summary);
  134. // Recover r from m.
  135. r = m * m;
  136. std::cout << summary.BriefReport() << "\n";
  137. std::cout << "x : " << initial_x << " -> " << x << "\n";
  138. std::cout << "y : " << initial_y << " -> " << y << "\n";
  139. std::cout << "r : " << initial_r << " -> " << r << "\n";
  140. return 0;
  141. }