bicubic_interpolation.cc 5.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155
  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. // Bicubic interpolation with automatic differentiation
  30. //
  31. // We will use estimation of 2d shift as a sample problem for bicubic
  32. // interpolation.
  33. //
  34. // Let us define f(x, y) = x * x - y * x + y * y
  35. // And optimize cost function sum_i [f(x_i + s_x, y_i + s_y) - v_i]^2
  36. //
  37. // Bicubic interpolation of f(x, y) will be exact, thus we can expect close to
  38. // perfect convergence
  39. #include <utility>
  40. #include "ceres/ceres.h"
  41. #include "ceres/cubic_interpolation.h"
  42. #include "glog/logging.h"
  43. using Grid = ceres::Grid2D<double>;
  44. using Interpolator = ceres::BiCubicInterpolator<Grid>;
  45. // Cost-function using autodiff interface of BiCubicInterpolator
  46. struct AutoDiffBiCubicCost {
  47. EIGEN_MAKE_ALIGNED_OPERATOR_NEW;
  48. template <typename T>
  49. bool operator()(const T* s, T* residual) const {
  50. using Vector2T = Eigen::Matrix<T, 2, 1>;
  51. Eigen::Map<const Vector2T> shift(s);
  52. const Vector2T point = point_ + shift;
  53. T v;
  54. interpolator_.Evaluate(point.y(), point.x(), &v);
  55. *residual = v - value_;
  56. return true;
  57. }
  58. AutoDiffBiCubicCost(const Interpolator& interpolator,
  59. Eigen::Vector2d point,
  60. double value)
  61. : point_(std::move(point)), value_(value), interpolator_(interpolator) {}
  62. static ceres::CostFunction* Create(const Interpolator& interpolator,
  63. const Eigen::Vector2d& point,
  64. double value) {
  65. return new ceres::AutoDiffCostFunction<AutoDiffBiCubicCost, 1, 2>(
  66. new AutoDiffBiCubicCost(interpolator, point, value));
  67. }
  68. const Eigen::Vector2d point_;
  69. const double value_;
  70. const Interpolator& interpolator_;
  71. };
  72. // Function for input data generation
  73. static double f(const double& x, const double& y) {
  74. return x * x - y * x + y * y;
  75. }
  76. int main(int argc, char** argv) {
  77. google::InitGoogleLogging(argv[0]);
  78. // Problem sizes
  79. const int kGridRowsHalf = 9;
  80. const int kGridColsHalf = 11;
  81. const int kGridRows = 2 * kGridRowsHalf + 1;
  82. const int kGridCols = 2 * kGridColsHalf + 1;
  83. const int kPoints = 4;
  84. const Eigen::Vector2d shift(1.234, 2.345);
  85. const std::array<Eigen::Vector2d, kPoints> points = {
  86. Eigen::Vector2d{-2., -3.},
  87. Eigen::Vector2d{-2., 3.},
  88. Eigen::Vector2d{2., 3.},
  89. Eigen::Vector2d{2., -3.}};
  90. // Data is a row-major array of kGridRows x kGridCols values of function
  91. // f(x, y) on the grid, with x in {-kGridColsHalf, ..., +kGridColsHalf},
  92. // and y in {-kGridRowsHalf, ..., +kGridRowsHalf}
  93. double data[kGridRows * kGridCols];
  94. for (int i = 0; i < kGridRows; ++i) {
  95. for (int j = 0; j < kGridCols; ++j) {
  96. // Using row-major order
  97. int index = i * kGridCols + j;
  98. double y = i - kGridRowsHalf;
  99. double x = j - kGridColsHalf;
  100. data[index] = f(x, y);
  101. }
  102. }
  103. const Grid grid(data,
  104. -kGridRowsHalf,
  105. kGridRowsHalf + 1,
  106. -kGridColsHalf,
  107. kGridColsHalf + 1);
  108. const Interpolator interpolator(grid);
  109. Eigen::Vector2d shift_estimate(3.1415, 1.337);
  110. ceres::Problem problem;
  111. problem.AddParameterBlock(shift_estimate.data(), 2);
  112. for (const auto& p : points) {
  113. const Eigen::Vector2d shifted = p + shift;
  114. const double v = f(shifted.x(), shifted.y());
  115. problem.AddResidualBlock(AutoDiffBiCubicCost::Create(interpolator, p, v),
  116. nullptr,
  117. shift_estimate.data());
  118. }
  119. ceres::Solver::Options options;
  120. options.minimizer_progress_to_stdout = true;
  121. ceres::Solver::Summary summary;
  122. ceres::Solve(options, &problem, &summary);
  123. std::cout << summary.BriefReport() << '\n';
  124. std::cout << "Bicubic interpolation with automatic derivatives:\n";
  125. std::cout << "Estimated shift: " << shift_estimate.transpose()
  126. << ", ground-truth: " << shift.transpose()
  127. << " (error: " << (shift_estimate - shift).transpose() << ")"
  128. << std::endl;
  129. CHECK_LT((shift_estimate - shift).norm(), 1e-9);
  130. return 0;
  131. }