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