libmv_homography.cc 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406
  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. // 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
  57. // Gold Standard Solution in "Multiple View Geometry"). The routines are based
  58. // on the routines from the Libmv library.
  59. //
  60. // This example demonstrates custom exit criterion by having a callback check
  61. // for image-space error.
  62. #include <utility>
  63. #include "ceres/ceres.h"
  64. #include "glog/logging.h"
  65. using EigenDouble = Eigen::NumTraits<double>;
  66. using Mat = Eigen::MatrixXd;
  67. using Vec = Eigen::VectorXd;
  68. using Mat3 = Eigen::Matrix<double, 3, 3>;
  69. using Vec2 = Eigen::Matrix<double, 2, 1>;
  70. using MatX8 = Eigen::Matrix<double, Eigen::Dynamic, 8>;
  71. using Vec3 = Eigen::Vector3d;
  72. namespace {
  73. // This structure contains options that controls how the homography
  74. // estimation operates.
  75. //
  76. // Defaults should be suitable for a wide range of use cases, but
  77. // better performance and accuracy might require tweaking.
  78. struct EstimateHomographyOptions {
  79. // Default settings for homography estimation which should be suitable
  80. // for a wide range of use cases.
  81. EstimateHomographyOptions() = default;
  82. // Maximal number of iterations for the refinement step.
  83. int max_num_iterations{50};
  84. // Expected average of symmetric geometric distance between
  85. // actual destination points and original ones transformed by
  86. // estimated homography matrix.
  87. //
  88. // Refinement will finish as soon as average of symmetric
  89. // geometric distance is less or equal to this value.
  90. //
  91. // This distance is measured in the same units as input points are.
  92. double expected_average_symmetric_distance{1e-16};
  93. };
  94. // Calculate symmetric geometric cost terms:
  95. //
  96. // forward_error = D(H * x1, x2)
  97. // backward_error = D(H^-1 * x2, x1)
  98. //
  99. // Templated to be used with autodifferentiation.
  100. template <typename T>
  101. void SymmetricGeometricDistanceTerms(const Eigen::Matrix<T, 3, 3>& H,
  102. const Eigen::Matrix<T, 2, 1>& x1,
  103. const Eigen::Matrix<T, 2, 1>& x2,
  104. T forward_error[2],
  105. T backward_error[2]) {
  106. using Vec3 = Eigen::Matrix<T, 3, 1>;
  107. Vec3 x(x1(0), x1(1), T(1.0));
  108. Vec3 y(x2(0), x2(1), T(1.0));
  109. Vec3 H_x = H * x;
  110. Vec3 Hinv_y = H.inverse() * y;
  111. H_x /= H_x(2);
  112. Hinv_y /= Hinv_y(2);
  113. forward_error[0] = H_x(0) - y(0);
  114. forward_error[1] = H_x(1) - y(1);
  115. backward_error[0] = Hinv_y(0) - x(0);
  116. backward_error[1] = Hinv_y(1) - x(1);
  117. }
  118. // Calculate symmetric geometric cost:
  119. //
  120. // D(H * x1, x2)^2 + D(H^-1 * x2, x1)^2
  121. //
  122. double SymmetricGeometricDistance(const Mat3& H,
  123. const Vec2& x1,
  124. const Vec2& x2) {
  125. Vec2 forward_error, backward_error;
  126. SymmetricGeometricDistanceTerms<double>(
  127. H, x1, x2, forward_error.data(), backward_error.data());
  128. return forward_error.squaredNorm() + backward_error.squaredNorm();
  129. }
  130. // A parameterization of the 2D homography matrix that uses 8 parameters so
  131. // that the matrix is normalized (H(2,2) == 1).
  132. // The homography matrix H is built from a list of 8 parameters (a, b,...g, h)
  133. // as follows
  134. //
  135. // |a b c|
  136. // H = |d e f|
  137. // |g h 1|
  138. //
  139. template <typename T = double>
  140. class Homography2DNormalizedParameterization {
  141. public:
  142. using Parameters = Eigen::Matrix<T, 8, 1>; // a, b, ... g, h
  143. using Parameterized = Eigen::Matrix<T, 3, 3>; // H
  144. // Convert from the 8 parameters to a H matrix.
  145. static void To(const Parameters& p, Parameterized* h) {
  146. // clang-format off
  147. *h << p(0), p(1), p(2),
  148. p(3), p(4), p(5),
  149. p(6), p(7), 1.0;
  150. // clang-format on
  151. }
  152. // Convert from a H matrix to the 8 parameters.
  153. static void From(const Parameterized& h, Parameters* p) {
  154. // clang-format off
  155. *p << h(0, 0), h(0, 1), h(0, 2),
  156. h(1, 0), h(1, 1), h(1, 2),
  157. h(2, 0), h(2, 1);
  158. // clang-format on
  159. }
  160. };
  161. // 2D Homography transformation estimation in the case that points are in
  162. // euclidean coordinates.
  163. //
  164. // x = H y
  165. //
  166. // x and y vector must have the same direction, we could write
  167. //
  168. // crossproduct(|x|, * H * |y| ) = |0|
  169. //
  170. // | 0 -1 x2| |a b c| |y1| |0|
  171. // | 1 0 -x1| * |d e f| * |y2| = |0|
  172. // |-x2 x1 0| |g h 1| |1 | |0|
  173. //
  174. // That gives:
  175. //
  176. // (-d+x2*g)*y1 + (-e+x2*h)*y2 + -f+x2 |0|
  177. // (a-x1*g)*y1 + (b-x1*h)*y2 + c-x1 = |0|
  178. // (-x2*a+x1*d)*y1 + (-x2*b+x1*e)*y2 + -x2*c+x1*f |0|
  179. //
  180. bool Homography2DFromCorrespondencesLinearEuc(const Mat& x1,
  181. const Mat& x2,
  182. Mat3* H,
  183. double expected_precision) {
  184. assert(2 == x1.rows());
  185. assert(4 <= x1.cols());
  186. assert(x1.rows() == x2.rows());
  187. assert(x1.cols() == x2.cols());
  188. const int64_t n = x1.cols();
  189. MatX8 L = Mat::Zero(n * 3, 8);
  190. Mat b = Mat::Zero(n * 3, 1);
  191. for (int64_t i = 0; i < n; ++i) {
  192. int64_t j = 3 * i;
  193. L(j, 0) = x1(0, i); // a
  194. L(j, 1) = x1(1, i); // b
  195. L(j, 2) = 1.0; // c
  196. L(j, 6) = -x2(0, i) * x1(0, i); // g
  197. L(j, 7) = -x2(0, i) * x1(1, i); // h
  198. b(j, 0) = x2(0, i); // i
  199. ++j;
  200. L(j, 3) = x1(0, i); // d
  201. L(j, 4) = x1(1, i); // e
  202. L(j, 5) = 1.0; // f
  203. L(j, 6) = -x2(1, i) * x1(0, i); // g
  204. L(j, 7) = -x2(1, i) * x1(1, i); // h
  205. b(j, 0) = x2(1, i); // i
  206. // This ensures better stability
  207. // TODO(julien) make a lite version without this 3rd set
  208. ++j;
  209. L(j, 0) = x2(1, i) * x1(0, i); // a
  210. L(j, 1) = x2(1, i) * x1(1, i); // b
  211. L(j, 2) = x2(1, i); // c
  212. L(j, 3) = -x2(0, i) * x1(0, i); // d
  213. L(j, 4) = -x2(0, i) * x1(1, i); // e
  214. L(j, 5) = -x2(0, i); // f
  215. }
  216. // Solve Lx=B
  217. const Vec h = L.fullPivLu().solve(b);
  218. Homography2DNormalizedParameterization<double>::To(h, H);
  219. return (L * h).isApprox(b, expected_precision);
  220. }
  221. // Cost functor which computes symmetric geometric distance
  222. // used for homography matrix refinement.
  223. class HomographySymmetricGeometricCostFunctor {
  224. public:
  225. HomographySymmetricGeometricCostFunctor(Vec2 x, Vec2 y)
  226. : x_(std::move(x)), y_(std::move(y)) {}
  227. template <typename T>
  228. bool operator()(const T* homography_parameters, T* residuals) const {
  229. using Mat3 = Eigen::Matrix<T, 3, 3>;
  230. using Vec2 = Eigen::Matrix<T, 2, 1>;
  231. Mat3 H(homography_parameters);
  232. Vec2 x(T(x_(0)), T(x_(1)));
  233. Vec2 y(T(y_(0)), T(y_(1)));
  234. SymmetricGeometricDistanceTerms<T>(H, x, y, &residuals[0], &residuals[2]);
  235. return true;
  236. }
  237. const Vec2 x_;
  238. const Vec2 y_;
  239. };
  240. // Termination checking callback. This is needed to finish the
  241. // optimization when an absolute error threshold is met, as opposed
  242. // to Ceres's function_tolerance, which provides for finishing when
  243. // successful steps reduce the cost function by a fractional amount.
  244. // In this case, the callback checks for the absolute average reprojection
  245. // error and terminates when it's below a threshold (for example all
  246. // points < 0.5px error).
  247. class TerminationCheckingCallback : public ceres::IterationCallback {
  248. public:
  249. TerminationCheckingCallback(const Mat& x1,
  250. const Mat& x2,
  251. const EstimateHomographyOptions& options,
  252. Mat3* H)
  253. : options_(options), x1_(x1), x2_(x2), H_(H) {}
  254. ceres::CallbackReturnType operator()(
  255. const ceres::IterationSummary& summary) override {
  256. // If the step wasn't successful, there's nothing to do.
  257. if (!summary.step_is_successful) {
  258. return ceres::SOLVER_CONTINUE;
  259. }
  260. // Calculate average of symmetric geometric distance.
  261. double average_distance = 0.0;
  262. for (int i = 0; i < x1_.cols(); i++) {
  263. average_distance +=
  264. SymmetricGeometricDistance(*H_, x1_.col(i), x2_.col(i));
  265. }
  266. average_distance /= x1_.cols();
  267. if (average_distance <= options_.expected_average_symmetric_distance) {
  268. return ceres::SOLVER_TERMINATE_SUCCESSFULLY;
  269. }
  270. return ceres::SOLVER_CONTINUE;
  271. }
  272. private:
  273. const EstimateHomographyOptions& options_;
  274. const Mat& x1_;
  275. const Mat& x2_;
  276. Mat3* H_;
  277. };
  278. bool EstimateHomography2DFromCorrespondences(
  279. const Mat& x1,
  280. const Mat& x2,
  281. const EstimateHomographyOptions& options,
  282. Mat3* H) {
  283. assert(2 == x1.rows());
  284. assert(4 <= x1.cols());
  285. assert(x1.rows() == x2.rows());
  286. assert(x1.cols() == x2.cols());
  287. // Step 1: Algebraic homography estimation.
  288. // Assume algebraic estimation always succeeds.
  289. Homography2DFromCorrespondencesLinearEuc(
  290. x1, x2, H, EigenDouble::dummy_precision());
  291. LOG(INFO) << "Estimated matrix after algebraic estimation:\n" << *H;
  292. // Step 2: Refine matrix using Ceres minimizer.
  293. ceres::Problem problem;
  294. for (int i = 0; i < x1.cols(); i++) {
  295. auto* homography_symmetric_geometric_cost_function =
  296. new HomographySymmetricGeometricCostFunctor(x1.col(i), x2.col(i));
  297. problem.AddResidualBlock(
  298. new ceres::AutoDiffCostFunction<HomographySymmetricGeometricCostFunctor,
  299. 4, // num_residuals
  300. 9>(
  301. homography_symmetric_geometric_cost_function),
  302. nullptr,
  303. H->data());
  304. }
  305. // Configure the solve.
  306. ceres::Solver::Options solver_options;
  307. solver_options.linear_solver_type = ceres::DENSE_QR;
  308. solver_options.max_num_iterations = options.max_num_iterations;
  309. solver_options.update_state_every_iteration = true;
  310. // Terminate if the average symmetric distance is good enough.
  311. TerminationCheckingCallback callback(x1, x2, options, H);
  312. solver_options.callbacks.push_back(&callback);
  313. // Run the solve.
  314. ceres::Solver::Summary summary;
  315. ceres::Solve(solver_options, &problem, &summary);
  316. LOG(INFO) << "Summary:\n" << summary.FullReport();
  317. LOG(INFO) << "Final refined matrix:\n" << *H;
  318. return summary.IsSolutionUsable();
  319. }
  320. } // namespace
  321. int main(int argc, char** argv) {
  322. google::InitGoogleLogging(argv[0]);
  323. Mat x1(2, 100);
  324. for (int i = 0; i < x1.cols(); ++i) {
  325. x1(0, i) = rand() % 1024;
  326. x1(1, i) = rand() % 1024;
  327. }
  328. Mat3 homography_matrix;
  329. // This matrix has been dumped from a Blender test file of plane tracking.
  330. // clang-format off
  331. homography_matrix << 1.243715, -0.461057, -111.964454,
  332. 0.0, 0.617589, -192.379252,
  333. 0.0, -0.000983, 1.0;
  334. // clang-format on
  335. Mat x2 = x1;
  336. for (int i = 0; i < x2.cols(); ++i) {
  337. Vec3 homogeneous_x1 = Vec3(x1(0, i), x1(1, i), 1.0);
  338. Vec3 homogeneous_x2 = homography_matrix * homogeneous_x1;
  339. x2(0, i) = homogeneous_x2(0) / homogeneous_x2(2);
  340. x2(1, i) = homogeneous_x2(1) / homogeneous_x2(2);
  341. // Apply some noise so algebraic estimation is not good enough.
  342. x2(0, i) += static_cast<double>(rand() % 1000) / 5000.0;
  343. x2(1, i) += static_cast<double>(rand() % 1000) / 5000.0;
  344. }
  345. Mat3 estimated_matrix;
  346. EstimateHomographyOptions options;
  347. options.expected_average_symmetric_distance = 0.02;
  348. EstimateHomography2DFromCorrespondences(x1, x2, options, &estimated_matrix);
  349. // Normalize the matrix for easier comparison.
  350. estimated_matrix /= estimated_matrix(2, 2);
  351. std::cout << "Original matrix:\n" << homography_matrix << "\n";
  352. std::cout << "Estimated matrix:\n" << estimated_matrix << "\n";
  353. return EXIT_SUCCESS;
  354. }