blob: c02951ecdfb619f39492ebd608996184c04741c1 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// 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// Author: sameeragarwal@google.com (Sameer Agarwal)
30
31#include "ceres/cubic_interpolation.h"
32
33#include <memory>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080034
Austin Schuh70cc9552019-01-21 19:46:48 -080035#include "ceres/jet.h"
36#include "glog/logging.h"
37#include "gtest/gtest.h"
38
Austin Schuh3de38b02024-06-25 18:25:10 -070039namespace ceres::internal {
Austin Schuh70cc9552019-01-21 19:46:48 -080040
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080041static constexpr double kTolerance = 1e-12;
Austin Schuh70cc9552019-01-21 19:46:48 -080042
43TEST(Grid1D, OneDataDimension) {
44 int x[] = {1, 2, 3};
45 Grid1D<int, 1> grid(x, 0, 3);
46 for (int i = 0; i < 3; ++i) {
47 double value;
48 grid.GetValue(i, &value);
49 EXPECT_EQ(value, static_cast<double>(i + 1));
50 }
51}
52
53TEST(Grid1D, OneDataDimensionOutOfBounds) {
54 int x[] = {1, 2, 3};
55 Grid1D<int, 1> grid(x, 0, 3);
56 double value;
57 grid.GetValue(-1, &value);
58 EXPECT_EQ(value, x[0]);
59 grid.GetValue(-2, &value);
60 EXPECT_EQ(value, x[0]);
61 grid.GetValue(3, &value);
62 EXPECT_EQ(value, x[2]);
63 grid.GetValue(4, &value);
64 EXPECT_EQ(value, x[2]);
65}
66
67TEST(Grid1D, TwoDataDimensionIntegerDataInterleaved) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080068 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -080069 int x[] = {1, 5,
70 2, 6,
71 3, 7};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080072 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -080073
74 Grid1D<int, 2, true> grid(x, 0, 3);
75 for (int i = 0; i < 3; ++i) {
76 double value[2];
77 grid.GetValue(i, value);
78 EXPECT_EQ(value[0], static_cast<double>(i + 1));
79 EXPECT_EQ(value[1], static_cast<double>(i + 5));
80 }
81}
82
Austin Schuh70cc9552019-01-21 19:46:48 -080083TEST(Grid1D, TwoDataDimensionIntegerDataStacked) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080084 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -080085 int x[] = {1, 2, 3,
86 5, 6, 7};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080087 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -080088
89 Grid1D<int, 2, false> grid(x, 0, 3);
90 for (int i = 0; i < 3; ++i) {
91 double value[2];
92 grid.GetValue(i, value);
93 EXPECT_EQ(value[0], static_cast<double>(i + 1));
94 EXPECT_EQ(value[1], static_cast<double>(i + 5));
95 }
96}
97
98TEST(Grid2D, OneDataDimensionRowMajor) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080099 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800100 int x[] = {1, 2, 3,
101 2, 3, 4};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800102 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800103 Grid2D<int, 1, true, true> grid(x, 0, 2, 0, 3);
104 for (int r = 0; r < 2; ++r) {
105 for (int c = 0; c < 3; ++c) {
106 double value;
107 grid.GetValue(r, c, &value);
108 EXPECT_EQ(value, static_cast<double>(r + c + 1));
109 }
110 }
111}
112
113TEST(Grid2D, OneDataDimensionRowMajorOutOfBounds) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800114 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800115 int x[] = {1, 2, 3,
116 2, 3, 4};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800117 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800118 Grid2D<int, 1, true, true> grid(x, 0, 2, 0, 3);
119 double value;
120 grid.GetValue(-1, -1, &value);
121 EXPECT_EQ(value, x[0]);
122 grid.GetValue(-1, 0, &value);
123 EXPECT_EQ(value, x[0]);
124 grid.GetValue(-1, 1, &value);
125 EXPECT_EQ(value, x[1]);
126 grid.GetValue(-1, 2, &value);
127 EXPECT_EQ(value, x[2]);
128 grid.GetValue(-1, 3, &value);
129 EXPECT_EQ(value, x[2]);
130 grid.GetValue(0, 3, &value);
131 EXPECT_EQ(value, x[2]);
132 grid.GetValue(1, 3, &value);
133 EXPECT_EQ(value, x[5]);
134 grid.GetValue(2, 3, &value);
135 EXPECT_EQ(value, x[5]);
136 grid.GetValue(2, 2, &value);
137 EXPECT_EQ(value, x[5]);
138 grid.GetValue(2, 1, &value);
139 EXPECT_EQ(value, x[4]);
140 grid.GetValue(2, 0, &value);
141 EXPECT_EQ(value, x[3]);
142 grid.GetValue(2, -1, &value);
143 EXPECT_EQ(value, x[3]);
144 grid.GetValue(1, -1, &value);
145 EXPECT_EQ(value, x[3]);
146 grid.GetValue(0, -1, &value);
147 EXPECT_EQ(value, x[0]);
148}
149
150TEST(Grid2D, TwoDataDimensionRowMajorInterleaved) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800151 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800152 int x[] = {1, 4, 2, 8, 3, 12,
153 2, 8, 3, 12, 4, 16};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800154 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800155 Grid2D<int, 2, true, true> grid(x, 0, 2, 0, 3);
156 for (int r = 0; r < 2; ++r) {
157 for (int c = 0; c < 3; ++c) {
158 double value[2];
159 grid.GetValue(r, c, value);
160 EXPECT_EQ(value[0], static_cast<double>(r + c + 1));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800161 EXPECT_EQ(value[1], static_cast<double>(4 * (r + c + 1)));
Austin Schuh70cc9552019-01-21 19:46:48 -0800162 }
163 }
164}
165
166TEST(Grid2D, TwoDataDimensionRowMajorStacked) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800167 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800168 int x[] = {1, 2, 3,
169 2, 3, 4,
170 4, 8, 12,
171 8, 12, 16};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800172 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800173 Grid2D<int, 2, true, false> grid(x, 0, 2, 0, 3);
174 for (int r = 0; r < 2; ++r) {
175 for (int c = 0; c < 3; ++c) {
176 double value[2];
177 grid.GetValue(r, c, value);
178 EXPECT_EQ(value[0], static_cast<double>(r + c + 1));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800179 EXPECT_EQ(value[1], static_cast<double>(4 * (r + c + 1)));
Austin Schuh70cc9552019-01-21 19:46:48 -0800180 }
181 }
182}
183
184TEST(Grid2D, TwoDataDimensionColMajorInterleaved) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800185 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800186 int x[] = { 1, 4, 2, 8,
187 2, 8, 3, 12,
188 3, 12, 4, 16};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800189 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800190 Grid2D<int, 2, false, true> grid(x, 0, 2, 0, 3);
191 for (int r = 0; r < 2; ++r) {
192 for (int c = 0; c < 3; ++c) {
193 double value[2];
194 grid.GetValue(r, c, value);
195 EXPECT_EQ(value[0], static_cast<double>(r + c + 1));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800196 EXPECT_EQ(value[1], static_cast<double>(4 * (r + c + 1)));
Austin Schuh70cc9552019-01-21 19:46:48 -0800197 }
198 }
199}
200
201TEST(Grid2D, TwoDataDimensionColMajorStacked) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800202 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800203 int x[] = {1, 2,
204 2, 3,
205 3, 4,
206 4, 8,
207 8, 12,
208 12, 16};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800209 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800210 Grid2D<int, 2, false, false> grid(x, 0, 2, 0, 3);
211 for (int r = 0; r < 2; ++r) {
212 for (int c = 0; c < 3; ++c) {
213 double value[2];
214 grid.GetValue(r, c, value);
215 EXPECT_EQ(value[0], static_cast<double>(r + c + 1));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800216 EXPECT_EQ(value[1], static_cast<double>(4 * (r + c + 1)));
Austin Schuh70cc9552019-01-21 19:46:48 -0800217 }
218 }
219}
220
221class CubicInterpolatorTest : public ::testing::Test {
222 public:
223 template <int kDataDimension>
224 void RunPolynomialInterpolationTest(const double a,
225 const double b,
226 const double c,
227 const double d) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700228 values_ = std::make_unique<double[]>(kDataDimension * kNumSamples);
Austin Schuh70cc9552019-01-21 19:46:48 -0800229
230 for (int x = 0; x < kNumSamples; ++x) {
231 for (int dim = 0; dim < kDataDimension; ++dim) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800232 values_[x * kDataDimension + dim] =
233 (dim * dim + 1) * (a * x * x * x + b * x * x + c * x + d);
Austin Schuh70cc9552019-01-21 19:46:48 -0800234 }
235 }
236
237 Grid1D<double, kDataDimension> grid(values_.get(), 0, kNumSamples);
238 CubicInterpolator<Grid1D<double, kDataDimension>> interpolator(grid);
239
240 // Check values in the all the cells but the first and the last
241 // ones. In these cells, the interpolated function values should
242 // match exactly the values of the function being interpolated.
243 //
244 // On the boundary, we extrapolate the values of the function on
245 // the basis of its first derivative, so we do not expect the
246 // function values and its derivatives not to match.
247 for (int j = 0; j < kNumTestSamples; ++j) {
248 const double x = 1.0 + 7.0 / (kNumTestSamples - 1) * j;
249 double expected_f[kDataDimension], expected_dfdx[kDataDimension];
250 double f[kDataDimension], dfdx[kDataDimension];
251
252 for (int dim = 0; dim < kDataDimension; ++dim) {
253 expected_f[dim] =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800254 (dim * dim + 1) * (a * x * x * x + b * x * x + c * x + d);
255 expected_dfdx[dim] =
256 (dim * dim + 1) * (3.0 * a * x * x + 2.0 * b * x + c);
Austin Schuh70cc9552019-01-21 19:46:48 -0800257 }
258
259 interpolator.Evaluate(x, f, dfdx);
260 for (int dim = 0; dim < kDataDimension; ++dim) {
261 EXPECT_NEAR(f[dim], expected_f[dim], kTolerance)
262 << "x: " << x << " dim: " << dim
263 << " actual f(x): " << expected_f[dim]
264 << " estimated f(x): " << f[dim];
265 EXPECT_NEAR(dfdx[dim], expected_dfdx[dim], kTolerance)
266 << "x: " << x << " dim: " << dim
267 << " actual df(x)/dx: " << expected_dfdx[dim]
268 << " estimated df(x)/dx: " << dfdx[dim];
269 }
270 }
271 }
272
273 private:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800274 static constexpr int kNumSamples = 10;
275 static constexpr int kNumTestSamples = 100;
Austin Schuh70cc9552019-01-21 19:46:48 -0800276 std::unique_ptr<double[]> values_;
277};
278
279TEST_F(CubicInterpolatorTest, ConstantFunction) {
280 RunPolynomialInterpolationTest<1>(0.0, 0.0, 0.0, 0.5);
281 RunPolynomialInterpolationTest<2>(0.0, 0.0, 0.0, 0.5);
282 RunPolynomialInterpolationTest<3>(0.0, 0.0, 0.0, 0.5);
283}
284
285TEST_F(CubicInterpolatorTest, LinearFunction) {
286 RunPolynomialInterpolationTest<1>(0.0, 0.0, 1.0, 0.5);
287 RunPolynomialInterpolationTest<2>(0.0, 0.0, 1.0, 0.5);
288 RunPolynomialInterpolationTest<3>(0.0, 0.0, 1.0, 0.5);
289}
290
291TEST_F(CubicInterpolatorTest, QuadraticFunction) {
292 RunPolynomialInterpolationTest<1>(0.0, 0.4, 1.0, 0.5);
293 RunPolynomialInterpolationTest<2>(0.0, 0.4, 1.0, 0.5);
294 RunPolynomialInterpolationTest<3>(0.0, 0.4, 1.0, 0.5);
295}
296
Austin Schuh70cc9552019-01-21 19:46:48 -0800297TEST(CubicInterpolator, JetEvaluation) {
298 const double values[] = {1.0, 2.0, 2.0, 5.0, 3.0, 9.0, 2.0, 7.0};
299
300 Grid1D<double, 2, true> grid(values, 0, 4);
301 CubicInterpolator<Grid1D<double, 2, true>> interpolator(grid);
302
303 double f[2], dfdx[2];
304 const double x = 2.5;
305 interpolator.Evaluate(x, f, dfdx);
306
307 // Create a Jet with the same scalar part as x, so that the output
308 // Jet will be evaluated at x.
309 Jet<double, 4> x_jet;
310 x_jet.a = x;
311 x_jet.v(0) = 1.0;
312 x_jet.v(1) = 1.1;
313 x_jet.v(2) = 1.2;
314 x_jet.v(3) = 1.3;
315
316 Jet<double, 4> f_jets[2];
317 interpolator.Evaluate(x_jet, f_jets);
318
319 // Check that the scalar part of the Jet is f(x).
320 EXPECT_EQ(f_jets[0].a, f[0]);
321 EXPECT_EQ(f_jets[1].a, f[1]);
322
323 // Check that the derivative part of the Jet is dfdx * x_jet.v
324 // by the chain rule.
325 EXPECT_NEAR((f_jets[0].v - dfdx[0] * x_jet.v).norm(), 0.0, kTolerance);
326 EXPECT_NEAR((f_jets[1].v - dfdx[1] * x_jet.v).norm(), 0.0, kTolerance);
327}
328
329class BiCubicInterpolatorTest : public ::testing::Test {
330 public:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800331 // This class needs to have an Eigen aligned operator new as it contains
332 // fixed-size Eigen types.
333 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
334
Austin Schuh70cc9552019-01-21 19:46:48 -0800335 template <int kDataDimension>
336 void RunPolynomialInterpolationTest(const Eigen::Matrix3d& coeff) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700337 values_ = std::make_unique<double[]>(kNumRows * kNumCols * kDataDimension);
Austin Schuh70cc9552019-01-21 19:46:48 -0800338 coeff_ = coeff;
339 double* v = values_.get();
340 for (int r = 0; r < kNumRows; ++r) {
341 for (int c = 0; c < kNumCols; ++c) {
342 for (int dim = 0; dim < kDataDimension; ++dim) {
343 *v++ = (dim * dim + 1) * EvaluateF(r, c);
344 }
345 }
346 }
347
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800348 Grid2D<double, kDataDimension> grid(
349 values_.get(), 0, kNumRows, 0, kNumCols);
Austin Schuh70cc9552019-01-21 19:46:48 -0800350 BiCubicInterpolator<Grid2D<double, kDataDimension>> interpolator(grid);
351
352 for (int j = 0; j < kNumRowSamples; ++j) {
353 const double r = 1.0 + 7.0 / (kNumRowSamples - 1) * j;
354 for (int k = 0; k < kNumColSamples; ++k) {
355 const double c = 1.0 + 7.0 / (kNumColSamples - 1) * k;
356 double f[kDataDimension], dfdr[kDataDimension], dfdc[kDataDimension];
357 interpolator.Evaluate(r, c, f, dfdr, dfdc);
358 for (int dim = 0; dim < kDataDimension; ++dim) {
359 EXPECT_NEAR(f[dim], (dim * dim + 1) * EvaluateF(r, c), kTolerance);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800360 EXPECT_NEAR(
361 dfdr[dim], (dim * dim + 1) * EvaluatedFdr(r, c), kTolerance);
362 EXPECT_NEAR(
363 dfdc[dim], (dim * dim + 1) * EvaluatedFdc(r, c), kTolerance);
Austin Schuh70cc9552019-01-21 19:46:48 -0800364 }
365 }
366 }
367 }
368
369 private:
370 double EvaluateF(double r, double c) {
371 Eigen::Vector3d x;
372 x(0) = r;
373 x(1) = c;
374 x(2) = 1;
375 return x.transpose() * coeff_ * x;
376 }
377
378 double EvaluatedFdr(double r, double c) {
379 Eigen::Vector3d x;
380 x(0) = r;
381 x(1) = c;
382 x(2) = 1;
383 return (coeff_.row(0) + coeff_.col(0).transpose()) * x;
384 }
385
386 double EvaluatedFdc(double r, double c) {
387 Eigen::Vector3d x;
388 x(0) = r;
389 x(1) = c;
390 x(2) = 1;
391 return (coeff_.row(1) + coeff_.col(1).transpose()) * x;
392 }
393
Austin Schuh70cc9552019-01-21 19:46:48 -0800394 Eigen::Matrix3d coeff_;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800395 static constexpr int kNumRows = 10;
396 static constexpr int kNumCols = 10;
397 static constexpr int kNumRowSamples = 100;
398 static constexpr int kNumColSamples = 100;
Austin Schuh70cc9552019-01-21 19:46:48 -0800399 std::unique_ptr<double[]> values_;
400};
401
402TEST_F(BiCubicInterpolatorTest, ZeroFunction) {
403 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
404 RunPolynomialInterpolationTest<1>(coeff);
405 RunPolynomialInterpolationTest<2>(coeff);
406 RunPolynomialInterpolationTest<3>(coeff);
407}
408
409TEST_F(BiCubicInterpolatorTest, Degree00Function) {
410 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
411 coeff(2, 2) = 1.0;
412 RunPolynomialInterpolationTest<1>(coeff);
413 RunPolynomialInterpolationTest<2>(coeff);
414 RunPolynomialInterpolationTest<3>(coeff);
415}
416
417TEST_F(BiCubicInterpolatorTest, Degree01Function) {
418 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
419 coeff(2, 2) = 1.0;
420 coeff(0, 2) = 0.1;
421 coeff(2, 0) = 0.1;
422 RunPolynomialInterpolationTest<1>(coeff);
423 RunPolynomialInterpolationTest<2>(coeff);
424 RunPolynomialInterpolationTest<3>(coeff);
425}
426
427TEST_F(BiCubicInterpolatorTest, Degree10Function) {
428 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
429 coeff(2, 2) = 1.0;
430 coeff(0, 1) = 0.1;
431 coeff(1, 0) = 0.1;
432 RunPolynomialInterpolationTest<1>(coeff);
433 RunPolynomialInterpolationTest<2>(coeff);
434 RunPolynomialInterpolationTest<3>(coeff);
435}
436
437TEST_F(BiCubicInterpolatorTest, Degree11Function) {
438 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
439 coeff(2, 2) = 1.0;
440 coeff(0, 1) = 0.1;
441 coeff(1, 0) = 0.1;
442 coeff(0, 2) = 0.2;
443 coeff(2, 0) = 0.2;
444 RunPolynomialInterpolationTest<1>(coeff);
445 RunPolynomialInterpolationTest<2>(coeff);
446 RunPolynomialInterpolationTest<3>(coeff);
447}
448
449TEST_F(BiCubicInterpolatorTest, Degree12Function) {
450 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
451 coeff(2, 2) = 1.0;
452 coeff(0, 1) = 0.1;
453 coeff(1, 0) = 0.1;
454 coeff(0, 2) = 0.2;
455 coeff(2, 0) = 0.2;
456 coeff(1, 1) = 0.3;
457 RunPolynomialInterpolationTest<1>(coeff);
458 RunPolynomialInterpolationTest<2>(coeff);
459 RunPolynomialInterpolationTest<3>(coeff);
460}
461
462TEST_F(BiCubicInterpolatorTest, Degree21Function) {
463 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
464 coeff(2, 2) = 1.0;
465 coeff(0, 1) = 0.1;
466 coeff(1, 0) = 0.1;
467 coeff(0, 2) = 0.2;
468 coeff(2, 0) = 0.2;
469 coeff(0, 0) = 0.3;
470 RunPolynomialInterpolationTest<1>(coeff);
471 RunPolynomialInterpolationTest<2>(coeff);
472 RunPolynomialInterpolationTest<3>(coeff);
473}
474
475TEST_F(BiCubicInterpolatorTest, Degree22Function) {
476 Eigen::Matrix3d coeff = Eigen::Matrix3d::Zero();
477 coeff(2, 2) = 1.0;
478 coeff(0, 1) = 0.1;
479 coeff(1, 0) = 0.1;
480 coeff(0, 2) = 0.2;
481 coeff(2, 0) = 0.2;
482 coeff(0, 0) = 0.3;
483 coeff(0, 1) = -0.4;
484 coeff(1, 0) = -0.4;
485 RunPolynomialInterpolationTest<1>(coeff);
486 RunPolynomialInterpolationTest<2>(coeff);
487 RunPolynomialInterpolationTest<3>(coeff);
488}
489
490TEST(BiCubicInterpolator, JetEvaluation) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800491 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800492 const double values[] = {1.0, 5.0, 2.0, 10.0, 2.0, 6.0, 3.0, 5.0,
493 1.0, 2.0, 2.0, 2.0, 2.0, 2.0, 3.0, 1.0};
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800494 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800495
496 Grid2D<double, 2> grid(values, 0, 2, 0, 4);
497 BiCubicInterpolator<Grid2D<double, 2>> interpolator(grid);
498
499 double f[2], dfdr[2], dfdc[2];
500 const double r = 0.5;
501 const double c = 2.5;
502 interpolator.Evaluate(r, c, f, dfdr, dfdc);
503
504 // Create a Jet with the same scalar part as x, so that the output
505 // Jet will be evaluated at x.
506 Jet<double, 4> r_jet;
507 r_jet.a = r;
508 r_jet.v(0) = 1.0;
509 r_jet.v(1) = 1.1;
510 r_jet.v(2) = 1.2;
511 r_jet.v(3) = 1.3;
512
513 Jet<double, 4> c_jet;
514 c_jet.a = c;
515 c_jet.v(0) = 2.0;
516 c_jet.v(1) = 3.1;
517 c_jet.v(2) = 4.2;
518 c_jet.v(3) = 5.3;
519
520 Jet<double, 4> f_jets[2];
521 interpolator.Evaluate(r_jet, c_jet, f_jets);
522 EXPECT_EQ(f_jets[0].a, f[0]);
523 EXPECT_EQ(f_jets[1].a, f[1]);
524 EXPECT_NEAR((f_jets[0].v - dfdr[0] * r_jet.v - dfdc[0] * c_jet.v).norm(),
525 0.0,
526 kTolerance);
527 EXPECT_NEAR((f_jets[1].v - dfdr[1] * r_jet.v - dfdc[1] * c_jet.v).norm(),
528 0.0,
529 kTolerance);
530}
531
Austin Schuh3de38b02024-06-25 18:25:10 -0700532} // namespace ceres::internal