Austin Schuh | 3de38b0 | 2024-06-25 18:25:10 -0700 | [diff] [blame] | 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 | |
| 40 | #include <utility> |
| 41 | |
| 42 | #include "ceres/ceres.h" |
| 43 | #include "ceres/cubic_interpolation.h" |
| 44 | #include "glog/logging.h" |
| 45 | |
| 46 | using Grid = ceres::Grid2D<double>; |
| 47 | using Interpolator = ceres::BiCubicInterpolator<Grid>; |
| 48 | |
| 49 | // Cost-function using autodiff interface of BiCubicInterpolator |
| 50 | struct AutoDiffBiCubicCost { |
| 51 | EIGEN_MAKE_ALIGNED_OPERATOR_NEW; |
| 52 | |
| 53 | template <typename T> |
| 54 | bool operator()(const T* s, T* residual) const { |
| 55 | using Vector2T = Eigen::Matrix<T, 2, 1>; |
| 56 | Eigen::Map<const Vector2T> shift(s); |
| 57 | |
| 58 | const Vector2T point = point_ + shift; |
| 59 | |
| 60 | T v; |
| 61 | interpolator_.Evaluate(point.y(), point.x(), &v); |
| 62 | |
| 63 | *residual = v - value_; |
| 64 | return true; |
| 65 | } |
| 66 | |
| 67 | AutoDiffBiCubicCost(const Interpolator& interpolator, |
| 68 | Eigen::Vector2d point, |
| 69 | double value) |
| 70 | : point_(std::move(point)), value_(value), interpolator_(interpolator) {} |
| 71 | |
| 72 | static ceres::CostFunction* Create(const Interpolator& interpolator, |
| 73 | const Eigen::Vector2d& point, |
| 74 | double value) { |
| 75 | return new ceres::AutoDiffCostFunction<AutoDiffBiCubicCost, 1, 2>( |
| 76 | interpolator, point, value); |
| 77 | } |
| 78 | |
| 79 | const Eigen::Vector2d point_; |
| 80 | const double value_; |
| 81 | const Interpolator& interpolator_; |
| 82 | }; |
| 83 | |
| 84 | // Function for input data generation |
| 85 | static double f(const double& x, const double& y) { |
| 86 | return x * x - y * x + y * y; |
| 87 | } |
| 88 | |
| 89 | int main(int argc, char** argv) { |
| 90 | google::InitGoogleLogging(argv[0]); |
| 91 | // Problem sizes |
| 92 | const int kGridRowsHalf = 9; |
| 93 | const int kGridColsHalf = 11; |
| 94 | const int kGridRows = 2 * kGridRowsHalf + 1; |
| 95 | const int kGridCols = 2 * kGridColsHalf + 1; |
| 96 | const int kPoints = 4; |
| 97 | |
| 98 | const Eigen::Vector2d shift(1.234, 2.345); |
| 99 | const std::array<Eigen::Vector2d, kPoints> points = { |
| 100 | Eigen::Vector2d{-2., -3.}, |
| 101 | Eigen::Vector2d{-2., 3.}, |
| 102 | Eigen::Vector2d{2., 3.}, |
| 103 | Eigen::Vector2d{2., -3.}}; |
| 104 | |
| 105 | // Data is a row-major array of kGridRows x kGridCols values of function |
| 106 | // f(x, y) on the grid, with x in {-kGridColsHalf, ..., +kGridColsHalf}, |
| 107 | // and y in {-kGridRowsHalf, ..., +kGridRowsHalf} |
| 108 | double data[kGridRows * kGridCols]; |
| 109 | for (int i = 0; i < kGridRows; ++i) { |
| 110 | for (int j = 0; j < kGridCols; ++j) { |
| 111 | // Using row-major order |
| 112 | int index = i * kGridCols + j; |
| 113 | double y = i - kGridRowsHalf; |
| 114 | double x = j - kGridColsHalf; |
| 115 | |
| 116 | data[index] = f(x, y); |
| 117 | } |
| 118 | } |
| 119 | const Grid grid(data, |
| 120 | -kGridRowsHalf, |
| 121 | kGridRowsHalf + 1, |
| 122 | -kGridColsHalf, |
| 123 | kGridColsHalf + 1); |
| 124 | const Interpolator interpolator(grid); |
| 125 | |
| 126 | Eigen::Vector2d shift_estimate(3.1415, 1.337); |
| 127 | |
| 128 | ceres::Problem problem; |
| 129 | problem.AddParameterBlock(shift_estimate.data(), 2); |
| 130 | |
| 131 | for (const auto& p : points) { |
| 132 | const Eigen::Vector2d shifted = p + shift; |
| 133 | |
| 134 | const double v = f(shifted.x(), shifted.y()); |
| 135 | problem.AddResidualBlock(AutoDiffBiCubicCost::Create(interpolator, p, v), |
| 136 | nullptr, |
| 137 | shift_estimate.data()); |
| 138 | } |
| 139 | |
| 140 | ceres::Solver::Options options; |
| 141 | options.minimizer_progress_to_stdout = true; |
| 142 | |
| 143 | ceres::Solver::Summary summary; |
| 144 | ceres::Solve(options, &problem, &summary); |
| 145 | std::cout << summary.BriefReport() << '\n'; |
| 146 | |
| 147 | std::cout << "Bicubic interpolation with automatic derivatives:\n"; |
| 148 | std::cout << "Estimated shift: " << shift_estimate.transpose() |
| 149 | << ", ground-truth: " << shift.transpose() |
| 150 | << " (error: " << (shift_estimate - shift).transpose() << ")" |
| 151 | << std::endl; |
| 152 | |
| 153 | CHECK_LT((shift_estimate - shift).norm(), 1e-9); |
| 154 | return 0; |
| 155 | } |