blob: b7c9eda8b52e0e136b722bff617ac2d2d9d81a3e [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// 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// Copyright (c) 2014 libmv authors.
30//
31// Permission is hereby granted, free of charge, to any person obtaining a copy
32// of this software and associated documentation files (the "Software"), to
33// deal in the Software without restriction, including without limitation the
34// rights to use, copy, modify, merge, publish, distribute, sublicense, and/or
35// sell copies of the Software, and to permit persons to whom the Software is
36// furnished to do so, subject to the following conditions:
37//
38// The above copyright notice and this permission notice shall be included in
39// all copies or substantial portions of the Software.
40//
41// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
42// IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
43// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
44// AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
45// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
46// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS
47// IN THE SOFTWARE.
48//
49// Author: sergey.vfx@gmail.com (Sergey Sharybin)
50//
51// This file demonstrates solving for a homography between two sets of points.
52// A homography describes a transformation between a sets of points on a plane,
53// perspectively projected into two images. The first step is to solve a
54// homogeneous system of equations via singular value decomposition, giving an
55// algebraic solution for the homography, then solving for a final solution by
56// minimizing the symmetric transfer error in image space with Ceres (called the
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080057// Gold Standard Solution in "Multiple View Geometry"). The routines are based
58// on the routines from the Libmv library.
Austin Schuh70cc9552019-01-21 19:46:48 -080059//
60// This example demonstrates custom exit criterion by having a callback check
61// for image-space error.
62
Austin Schuh3de38b02024-06-25 18:25:10 -070063#include <utility>
64
Austin Schuh70cc9552019-01-21 19:46:48 -080065#include "ceres/ceres.h"
66#include "glog/logging.h"
67
Austin Schuh3de38b02024-06-25 18:25:10 -070068using EigenDouble = Eigen::NumTraits<double>;
Austin Schuh70cc9552019-01-21 19:46:48 -080069
Austin Schuh3de38b02024-06-25 18:25:10 -070070using Mat = Eigen::MatrixXd;
71using Vec = Eigen::VectorXd;
72using Mat3 = Eigen::Matrix<double, 3, 3>;
73using Vec2 = Eigen::Matrix<double, 2, 1>;
74using MatX8 = Eigen::Matrix<double, Eigen::Dynamic, 8>;
75using Vec3 = Eigen::Vector3d;
Austin Schuh70cc9552019-01-21 19:46:48 -080076
77namespace {
78
79// This structure contains options that controls how the homography
80// estimation operates.
81//
82// Defaults should be suitable for a wide range of use cases, but
83// better performance and accuracy might require tweaking.
84struct EstimateHomographyOptions {
85 // Default settings for homography estimation which should be suitable
86 // for a wide range of use cases.
Austin Schuh3de38b02024-06-25 18:25:10 -070087 EstimateHomographyOptions() = default;
Austin Schuh70cc9552019-01-21 19:46:48 -080088
89 // Maximal number of iterations for the refinement step.
Austin Schuh3de38b02024-06-25 18:25:10 -070090 int max_num_iterations{50};
Austin Schuh70cc9552019-01-21 19:46:48 -080091
92 // Expected average of symmetric geometric distance between
93 // actual destination points and original ones transformed by
94 // estimated homography matrix.
95 //
96 // Refinement will finish as soon as average of symmetric
97 // geometric distance is less or equal to this value.
98 //
99 // This distance is measured in the same units as input points are.
Austin Schuh3de38b02024-06-25 18:25:10 -0700100 double expected_average_symmetric_distance{1e-16};
Austin Schuh70cc9552019-01-21 19:46:48 -0800101};
102
103// Calculate symmetric geometric cost terms:
104//
105// forward_error = D(H * x1, x2)
106// backward_error = D(H^-1 * x2, x1)
107//
108// Templated to be used with autodifferentiation.
109template <typename T>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800110void SymmetricGeometricDistanceTerms(const Eigen::Matrix<T, 3, 3>& H,
111 const Eigen::Matrix<T, 2, 1>& x1,
112 const Eigen::Matrix<T, 2, 1>& x2,
Austin Schuh70cc9552019-01-21 19:46:48 -0800113 T forward_error[2],
114 T backward_error[2]) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700115 using Vec3 = Eigen::Matrix<T, 3, 1>;
Austin Schuh70cc9552019-01-21 19:46:48 -0800116 Vec3 x(x1(0), x1(1), T(1.0));
117 Vec3 y(x2(0), x2(1), T(1.0));
118
119 Vec3 H_x = H * x;
120 Vec3 Hinv_y = H.inverse() * y;
121
122 H_x /= H_x(2);
123 Hinv_y /= Hinv_y(2);
124
125 forward_error[0] = H_x(0) - y(0);
126 forward_error[1] = H_x(1) - y(1);
127 backward_error[0] = Hinv_y(0) - x(0);
128 backward_error[1] = Hinv_y(1) - x(1);
129}
130
131// Calculate symmetric geometric cost:
132//
133// D(H * x1, x2)^2 + D(H^-1 * x2, x1)^2
134//
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800135double SymmetricGeometricDistance(const Mat3& H,
136 const Vec2& x1,
137 const Vec2& x2) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800138 Vec2 forward_error, backward_error;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800139 SymmetricGeometricDistanceTerms<double>(
140 H, x1, x2, forward_error.data(), backward_error.data());
141 return forward_error.squaredNorm() + backward_error.squaredNorm();
Austin Schuh70cc9552019-01-21 19:46:48 -0800142}
143
144// A parameterization of the 2D homography matrix that uses 8 parameters so
145// that the matrix is normalized (H(2,2) == 1).
146// The homography matrix H is built from a list of 8 parameters (a, b,...g, h)
147// as follows
148//
149// |a b c|
150// H = |d e f|
151// |g h 1|
152//
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800153template <typename T = double>
Austin Schuh70cc9552019-01-21 19:46:48 -0800154class Homography2DNormalizedParameterization {
155 public:
Austin Schuh3de38b02024-06-25 18:25:10 -0700156 using Parameters = Eigen::Matrix<T, 8, 1>; // a, b, ... g, h
157 using Parameterized = Eigen::Matrix<T, 3, 3>; // H
Austin Schuh70cc9552019-01-21 19:46:48 -0800158
159 // Convert from the 8 parameters to a H matrix.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800160 static void To(const Parameters& p, Parameterized* h) {
161 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800162 *h << p(0), p(1), p(2),
163 p(3), p(4), p(5),
164 p(6), p(7), 1.0;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800165 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800166 }
167
168 // Convert from a H matrix to the 8 parameters.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800169 static void From(const Parameterized& h, Parameters* p) {
170 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800171 *p << h(0, 0), h(0, 1), h(0, 2),
172 h(1, 0), h(1, 1), h(1, 2),
173 h(2, 0), h(2, 1);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800174 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800175 }
176};
177
178// 2D Homography transformation estimation in the case that points are in
179// euclidean coordinates.
180//
181// x = H y
182//
183// x and y vector must have the same direction, we could write
184//
185// crossproduct(|x|, * H * |y| ) = |0|
186//
187// | 0 -1 x2| |a b c| |y1| |0|
188// | 1 0 -x1| * |d e f| * |y2| = |0|
189// |-x2 x1 0| |g h 1| |1 | |0|
190//
191// That gives:
192//
193// (-d+x2*g)*y1 + (-e+x2*h)*y2 + -f+x2 |0|
194// (a-x1*g)*y1 + (b-x1*h)*y2 + c-x1 = |0|
195// (-x2*a+x1*d)*y1 + (-x2*b+x1*e)*y2 + -x2*c+x1*f |0|
196//
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800197bool Homography2DFromCorrespondencesLinearEuc(const Mat& x1,
198 const Mat& x2,
199 Mat3* H,
200 double expected_precision) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800201 assert(2 == x1.rows());
202 assert(4 <= x1.cols());
203 assert(x1.rows() == x2.rows());
204 assert(x1.cols() == x2.cols());
205
Austin Schuh3de38b02024-06-25 18:25:10 -0700206 const int64_t n = x1.cols();
Austin Schuh70cc9552019-01-21 19:46:48 -0800207 MatX8 L = Mat::Zero(n * 3, 8);
208 Mat b = Mat::Zero(n * 3, 1);
Austin Schuh3de38b02024-06-25 18:25:10 -0700209 for (int64_t i = 0; i < n; ++i) {
210 int64_t j = 3 * i;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800211 L(j, 0) = x1(0, i); // a
212 L(j, 1) = x1(1, i); // b
213 L(j, 2) = 1.0; // c
Austin Schuh70cc9552019-01-21 19:46:48 -0800214 L(j, 6) = -x2(0, i) * x1(0, i); // g
215 L(j, 7) = -x2(0, i) * x1(1, i); // h
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800216 b(j, 0) = x2(0, i); // i
Austin Schuh70cc9552019-01-21 19:46:48 -0800217
218 ++j;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800219 L(j, 3) = x1(0, i); // d
220 L(j, 4) = x1(1, i); // e
221 L(j, 5) = 1.0; // f
Austin Schuh70cc9552019-01-21 19:46:48 -0800222 L(j, 6) = -x2(1, i) * x1(0, i); // g
223 L(j, 7) = -x2(1, i) * x1(1, i); // h
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800224 b(j, 0) = x2(1, i); // i
Austin Schuh70cc9552019-01-21 19:46:48 -0800225
226 // This ensures better stability
227 // TODO(julien) make a lite version without this 3rd set
228 ++j;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800229 L(j, 0) = x2(1, i) * x1(0, i); // a
230 L(j, 1) = x2(1, i) * x1(1, i); // b
231 L(j, 2) = x2(1, i); // c
Austin Schuh70cc9552019-01-21 19:46:48 -0800232 L(j, 3) = -x2(0, i) * x1(0, i); // d
233 L(j, 4) = -x2(0, i) * x1(1, i); // e
234 L(j, 5) = -x2(0, i); // f
235 }
236 // Solve Lx=B
237 const Vec h = L.fullPivLu().solve(b);
238 Homography2DNormalizedParameterization<double>::To(h, H);
239 return (L * h).isApprox(b, expected_precision);
240}
241
242// Cost functor which computes symmetric geometric distance
243// used for homography matrix refinement.
244class HomographySymmetricGeometricCostFunctor {
245 public:
Austin Schuh3de38b02024-06-25 18:25:10 -0700246 HomographySymmetricGeometricCostFunctor(Vec2 x, Vec2 y)
247 : x_(std::move(x)), y_(std::move(y)) {}
Austin Schuh70cc9552019-01-21 19:46:48 -0800248
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800249 template <typename T>
Austin Schuh70cc9552019-01-21 19:46:48 -0800250 bool operator()(const T* homography_parameters, T* residuals) const {
Austin Schuh3de38b02024-06-25 18:25:10 -0700251 using Mat3 = Eigen::Matrix<T, 3, 3>;
252 using Vec2 = Eigen::Matrix<T, 2, 1>;
Austin Schuh70cc9552019-01-21 19:46:48 -0800253
254 Mat3 H(homography_parameters);
255 Vec2 x(T(x_(0)), T(x_(1)));
256 Vec2 y(T(y_(0)), T(y_(1)));
257
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800258 SymmetricGeometricDistanceTerms<T>(H, x, y, &residuals[0], &residuals[2]);
Austin Schuh70cc9552019-01-21 19:46:48 -0800259 return true;
260 }
261
262 const Vec2 x_;
263 const Vec2 y_;
264};
265
266// Termination checking callback. This is needed to finish the
267// optimization when an absolute error threshold is met, as opposed
268// to Ceres's function_tolerance, which provides for finishing when
269// successful steps reduce the cost function by a fractional amount.
270// In this case, the callback checks for the absolute average reprojection
271// error and terminates when it's below a threshold (for example all
272// points < 0.5px error).
273class TerminationCheckingCallback : public ceres::IterationCallback {
274 public:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800275 TerminationCheckingCallback(const Mat& x1,
276 const Mat& x2,
277 const EstimateHomographyOptions& options,
278 Mat3* H)
Austin Schuh70cc9552019-01-21 19:46:48 -0800279 : options_(options), x1_(x1), x2_(x2), H_(H) {}
280
Austin Schuh3de38b02024-06-25 18:25:10 -0700281 ceres::CallbackReturnType operator()(
282 const ceres::IterationSummary& summary) override {
Austin Schuh70cc9552019-01-21 19:46:48 -0800283 // If the step wasn't successful, there's nothing to do.
284 if (!summary.step_is_successful) {
285 return ceres::SOLVER_CONTINUE;
286 }
287
288 // Calculate average of symmetric geometric distance.
289 double average_distance = 0.0;
290 for (int i = 0; i < x1_.cols(); i++) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800291 average_distance +=
292 SymmetricGeometricDistance(*H_, x1_.col(i), x2_.col(i));
Austin Schuh70cc9552019-01-21 19:46:48 -0800293 }
294 average_distance /= x1_.cols();
295
296 if (average_distance <= options_.expected_average_symmetric_distance) {
297 return ceres::SOLVER_TERMINATE_SUCCESSFULLY;
298 }
299
300 return ceres::SOLVER_CONTINUE;
301 }
302
303 private:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800304 const EstimateHomographyOptions& options_;
305 const Mat& x1_;
306 const Mat& x2_;
307 Mat3* H_;
Austin Schuh70cc9552019-01-21 19:46:48 -0800308};
309
310bool EstimateHomography2DFromCorrespondences(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800311 const Mat& x1,
312 const Mat& x2,
313 const EstimateHomographyOptions& options,
314 Mat3* H) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800315 assert(2 == x1.rows());
316 assert(4 <= x1.cols());
317 assert(x1.rows() == x2.rows());
318 assert(x1.cols() == x2.cols());
319
320 // Step 1: Algebraic homography estimation.
321 // Assume algebraic estimation always succeeds.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800322 Homography2DFromCorrespondencesLinearEuc(
323 x1, x2, H, EigenDouble::dummy_precision());
Austin Schuh70cc9552019-01-21 19:46:48 -0800324
325 LOG(INFO) << "Estimated matrix after algebraic estimation:\n" << *H;
326
327 // Step 2: Refine matrix using Ceres minimizer.
328 ceres::Problem problem;
329 for (int i = 0; i < x1.cols(); i++) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800330 problem.AddResidualBlock(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800331 new ceres::AutoDiffCostFunction<HomographySymmetricGeometricCostFunctor,
332 4, // num_residuals
Austin Schuh3de38b02024-06-25 18:25:10 -0700333 9>(x1.col(i), x2.col(i)),
334 nullptr,
Austin Schuh70cc9552019-01-21 19:46:48 -0800335 H->data());
336 }
337
338 // Configure the solve.
339 ceres::Solver::Options solver_options;
340 solver_options.linear_solver_type = ceres::DENSE_QR;
341 solver_options.max_num_iterations = options.max_num_iterations;
342 solver_options.update_state_every_iteration = true;
343
344 // Terminate if the average symmetric distance is good enough.
345 TerminationCheckingCallback callback(x1, x2, options, H);
346 solver_options.callbacks.push_back(&callback);
347
348 // Run the solve.
349 ceres::Solver::Summary summary;
350 ceres::Solve(solver_options, &problem, &summary);
351
352 LOG(INFO) << "Summary:\n" << summary.FullReport();
353 LOG(INFO) << "Final refined matrix:\n" << *H;
354
355 return summary.IsSolutionUsable();
356}
357
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800358} // namespace
Austin Schuh70cc9552019-01-21 19:46:48 -0800359
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800360int main(int argc, char** argv) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800361 google::InitGoogleLogging(argv[0]);
362
363 Mat x1(2, 100);
364 for (int i = 0; i < x1.cols(); ++i) {
365 x1(0, i) = rand() % 1024;
366 x1(1, i) = rand() % 1024;
367 }
368
369 Mat3 homography_matrix;
370 // This matrix has been dumped from a Blender test file of plane tracking.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800371 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800372 homography_matrix << 1.243715, -0.461057, -111.964454,
373 0.0, 0.617589, -192.379252,
374 0.0, -0.000983, 1.0;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800375 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800376
377 Mat x2 = x1;
378 for (int i = 0; i < x2.cols(); ++i) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700379 Vec3 homogeneous_x1 = Vec3(x1(0, i), x1(1, i), 1.0);
380 Vec3 homogeneous_x2 = homography_matrix * homogeneous_x1;
381 x2(0, i) = homogeneous_x2(0) / homogeneous_x2(2);
382 x2(1, i) = homogeneous_x2(1) / homogeneous_x2(2);
Austin Schuh70cc9552019-01-21 19:46:48 -0800383
384 // Apply some noise so algebraic estimation is not good enough.
385 x2(0, i) += static_cast<double>(rand() % 1000) / 5000.0;
386 x2(1, i) += static_cast<double>(rand() % 1000) / 5000.0;
387 }
388
389 Mat3 estimated_matrix;
390
391 EstimateHomographyOptions options;
392 options.expected_average_symmetric_distance = 0.02;
393 EstimateHomography2DFromCorrespondences(x1, x2, options, &estimated_matrix);
394
395 // Normalize the matrix for easier comparison.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800396 estimated_matrix /= estimated_matrix(2, 2);
Austin Schuh70cc9552019-01-21 19:46:48 -0800397
398 std::cout << "Original matrix:\n" << homography_matrix << "\n";
399 std::cout << "Estimated matrix:\n" << estimated_matrix << "\n";
400
401 return EXIT_SUCCESS;
402}