bicubic_interpolation_analytic.cc 5.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2021 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 analytic 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 analytic interface of BiCubicInterpolator
  46. struct AnalyticBiCubicCost : public ceres::CostFunction {
  47. EIGEN_MAKE_ALIGNED_OPERATOR_NEW;
  48. bool Evaluate(double const* const* parameters,
  49. double* residuals,
  50. double** jacobians) const override {
  51. Eigen::Map<const Eigen::Vector2d> shift(parameters[0]);
  52. const Eigen::Vector2d point = point_ + shift;
  53. double* f = residuals;
  54. double* dfdr = nullptr;
  55. double* dfdc = nullptr;
  56. if (jacobians && jacobians[0]) {
  57. dfdc = jacobians[0];
  58. dfdr = dfdc + 1;
  59. }
  60. interpolator_.Evaluate(point.y(), point.x(), f, dfdr, dfdc);
  61. if (residuals) {
  62. *f -= value_;
  63. }
  64. return true;
  65. }
  66. AnalyticBiCubicCost(const Interpolator& interpolator,
  67. Eigen::Vector2d point,
  68. double value)
  69. : point_(std::move(point)), value_(value), interpolator_(interpolator) {
  70. set_num_residuals(1);
  71. *mutable_parameter_block_sizes() = {2};
  72. }
  73. static ceres::CostFunction* Create(const Interpolator& interpolator,
  74. const Eigen::Vector2d& point,
  75. double value) {
  76. return new AnalyticBiCubicCost(interpolator, point, value);
  77. }
  78. const Eigen::Vector2d point_;
  79. const double value_;
  80. const Interpolator& interpolator_;
  81. };
  82. // Function for input data generation
  83. static double f(const double& x, const double& y) {
  84. return x * x - y * x + y * y;
  85. }
  86. int main(int argc, char** argv) {
  87. google::InitGoogleLogging(argv[0]);
  88. // Problem sizes
  89. const int kGridRowsHalf = 9;
  90. const int kGridColsHalf = 11;
  91. const int kGridRows = 2 * kGridRowsHalf + 1;
  92. const int kGridCols = 2 * kGridColsHalf + 1;
  93. const int kPoints = 4;
  94. const Eigen::Vector2d shift(1.234, 2.345);
  95. const std::array<Eigen::Vector2d, kPoints> points = {
  96. Eigen::Vector2d{-2., -3.},
  97. Eigen::Vector2d{-2., 3.},
  98. Eigen::Vector2d{2., 3.},
  99. Eigen::Vector2d{2., -3.}};
  100. // Data is a row-major array of kGridRows x kGridCols values of function
  101. // f(x, y) on the grid, with x in {-kGridColsHalf, ..., +kGridColsHalf},
  102. // and y in {-kGridRowsHalf, ..., +kGridRowsHalf}
  103. double data[kGridRows * kGridCols];
  104. for (int i = 0; i < kGridRows; ++i) {
  105. for (int j = 0; j < kGridCols; ++j) {
  106. // Using row-major order
  107. int index = i * kGridCols + j;
  108. double y = i - kGridRowsHalf;
  109. double x = j - kGridColsHalf;
  110. data[index] = f(x, y);
  111. }
  112. }
  113. const Grid grid(data,
  114. -kGridRowsHalf,
  115. kGridRowsHalf + 1,
  116. -kGridColsHalf,
  117. kGridColsHalf + 1);
  118. const Interpolator interpolator(grid);
  119. Eigen::Vector2d shift_estimate(3.1415, 1.337);
  120. ceres::Problem problem;
  121. problem.AddParameterBlock(shift_estimate.data(), 2);
  122. for (const auto& p : points) {
  123. const Eigen::Vector2d shifted = p + shift;
  124. const double v = f(shifted.x(), shifted.y());
  125. problem.AddResidualBlock(AnalyticBiCubicCost::Create(interpolator, p, v),
  126. nullptr,
  127. shift_estimate.data());
  128. }
  129. ceres::Solver::Options options;
  130. options.minimizer_progress_to_stdout = true;
  131. ceres::Solver::Summary summary;
  132. ceres::Solve(options, &problem, &summary);
  133. std::cout << summary.BriefReport() << '\n';
  134. std::cout << "Bicubic interpolation with analytic derivatives:\n";
  135. std::cout << "Estimated shift: " << shift_estimate.transpose()
  136. << ", ground-truth: " << shift.transpose()
  137. << " (error: " << (shift_estimate - shift).transpose() << ")"
  138. << std::endl;
  139. CHECK_LT((shift_estimate - shift).norm(), 1e-9);
  140. return 0;
  141. }