blob: c54247503a795fbe1980bd5ed9b6b9a31b08139f [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 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//
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080051// [1] http://www.mathworks.com/matlabcentral/fileexchange/5557-circle-fit/content/circfit.m // NOLINT
Austin Schuh70cc9552019-01-21 19:46:48 -080052
53#include <cstdio>
54#include <vector>
55
56#include "ceres/ceres.h"
57#include "gflags/gflags.h"
58#include "glog/logging.h"
59
60using ceres::AutoDiffCostFunction;
61using ceres::CauchyLoss;
62using ceres::CostFunction;
63using ceres::LossFunction;
64using ceres::Problem;
65using ceres::Solve;
66using ceres::Solver;
67
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080068DEFINE_double(robust_threshold,
69 0.0,
70 "Robust loss parameter. Set to 0 for normal squared error (no "
71 "robustification).");
Austin Schuh70cc9552019-01-21 19:46:48 -080072
73// The cost for a single sample. The returned residual is related to the
74// distance of the point from the circle (passed in as x, y, m parameters).
75//
76// Note that the radius is parameterized as r = m^2 to constrain the radius to
77// positive values.
78class DistanceFromCircleCost {
79 public:
80 DistanceFromCircleCost(double xx, double yy) : xx_(xx), yy_(yy) {}
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080081 template <typename T>
82 bool operator()(const T* const x,
83 const T* const y,
84 const T* const m, // r = m^2
85 T* residual) const {
Austin Schuh70cc9552019-01-21 19:46:48 -080086 // Since the radius is parameterized as m^2, unpack m to get r.
87 T r = *m * *m;
88
89 // Get the position of the sample in the circle's coordinate system.
90 T xp = xx_ - *x;
91 T yp = yy_ - *y;
92
93 // It is tempting to use the following cost:
94 //
95 // residual[0] = r - sqrt(xp*xp + yp*yp);
96 //
97 // which is the distance of the sample from the circle. This works
98 // reasonably well, but the sqrt() adds strong nonlinearities to the cost
99 // function. Instead, a different cost is used, which while not strictly a
100 // distance in the metric sense (it has units distance^2) it produces more
101 // robust fits when there are outliers. This is because the cost surface is
102 // more convex.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800103 residual[0] = r * r - xp * xp - yp * yp;
Austin Schuh70cc9552019-01-21 19:46:48 -0800104 return true;
105 }
106
107 private:
108 // The measured x,y coordinate that should be on the circle.
109 double xx_, yy_;
110};
111
112int main(int argc, char** argv) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800113 GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
Austin Schuh70cc9552019-01-21 19:46:48 -0800114 google::InitGoogleLogging(argv[0]);
115
116 double x, y, r;
117 if (scanf("%lg %lg %lg", &x, &y, &r) != 3) {
118 fprintf(stderr, "Couldn't read first line.\n");
119 return 1;
120 }
121 fprintf(stderr, "Got x, y, r %lg, %lg, %lg\n", x, y, r);
122
123 // Save initial values for comparison.
124 double initial_x = x;
125 double initial_y = y;
126 double initial_r = r;
127
128 // Parameterize r as m^2 so that it can't be negative.
129 double m = sqrt(r);
130
131 Problem problem;
132
133 // Configure the loss function.
134 LossFunction* loss = NULL;
135 if (FLAGS_robust_threshold) {
136 loss = new CauchyLoss(FLAGS_robust_threshold);
137 }
138
139 // Add the residuals.
140 double xx, yy;
141 int num_points = 0;
142 while (scanf("%lf %lf\n", &xx, &yy) == 2) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800143 CostFunction* cost =
Austin Schuh70cc9552019-01-21 19:46:48 -0800144 new AutoDiffCostFunction<DistanceFromCircleCost, 1, 1, 1, 1>(
145 new DistanceFromCircleCost(xx, yy));
146 problem.AddResidualBlock(cost, loss, &x, &y, &m);
147 num_points++;
148 }
149
150 std::cout << "Got " << num_points << " points.\n";
151
152 // Build and solve the problem.
153 Solver::Options options;
154 options.max_num_iterations = 500;
155 options.linear_solver_type = ceres::DENSE_QR;
156 Solver::Summary summary;
157 Solve(options, &problem, &summary);
158
159 // Recover r from m.
160 r = m * m;
161
162 std::cout << summary.BriefReport() << "\n";
163 std::cout << "x : " << initial_x << " -> " << x << "\n";
164 std::cout << "y : " << initial_y << " -> " << y << "\n";
165 std::cout << "r : " << initial_r << " -> " << r << "\n";
166 return 0;
167}