blob: 87e00676d34348b063ae42993df2dae50eba5b3c [file] [log] [blame]
Austin Schuh1d1e6ea2020-12-23 21:56:30 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2020 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// Author: darius.rueckert@fau.de (Darius Rueckert)
30
31#include <memory>
32#include <random>
33#include <utility>
34
35#include "benchmark/benchmark.h"
36#include "ceres/autodiff_benchmarks/brdf_cost_function.h"
37#include "ceres/autodiff_benchmarks/constant_cost_function.h"
38#include "ceres/autodiff_benchmarks/linear_cost_functions.h"
39#include "ceres/autodiff_benchmarks/photometric_error.h"
40#include "ceres/autodiff_benchmarks/relative_pose_error.h"
41#include "ceres/autodiff_benchmarks/snavely_reprojection_error.h"
42#include "ceres/ceres.h"
43
44namespace ceres {
45
46enum Dynamic { kNotDynamic, kDynamic };
47
48// Transforms a static functor into a dynamic one.
49template <typename CostFunctionType, int kNumParameterBlocks>
50class ToDynamic {
51 public:
52 template <typename... _Args>
53 explicit ToDynamic(_Args&&... __args)
54 : cost_function_(std::forward<_Args>(__args)...) {}
55
56 template <typename T>
57 bool operator()(const T* const* parameters, T* residuals) const {
58 return Apply(
59 parameters, residuals, std::make_index_sequence<kNumParameterBlocks>());
60 }
61
62 private:
63 template <typename T, size_t... Indices>
64 bool Apply(const T* const* parameters,
65 T* residuals,
66 std::index_sequence<Indices...>) const {
67 return cost_function_(parameters[Indices]..., residuals);
68 }
69
70 CostFunctionType cost_function_;
71};
72
73template <int kParameterBlockSize>
74static void BM_ConstantAnalytic(benchmark::State& state) {
75 constexpr int num_residuals = 1;
76 std::array<double, kParameterBlockSize> parameters_values;
77 std::iota(parameters_values.begin(), parameters_values.end(), 0);
78 double* parameters[] = {parameters_values.data()};
79
80 std::array<double, num_residuals> residuals;
81
82 std::array<double, num_residuals * kParameterBlockSize> jacobian_values;
83 double* jacobians[] = {jacobian_values.data()};
84
85 std::unique_ptr<ceres::CostFunction> cost_function(
86 new AnalyticConstantCostFunction<kParameterBlockSize>());
87
88 for (auto _ : state) {
89 cost_function->Evaluate(parameters, residuals.data(), jacobians);
90 }
91}
92
93// Helpers for CostFunctionFactory.
94template <typename DynamicCostFunctionType>
95void AddParameterBlocks(DynamicCostFunctionType*) {}
96
97template <int HeadN, int... TailNs, typename DynamicCostFunctionType>
98void AddParameterBlocks(DynamicCostFunctionType* dynamic_function) {
99 dynamic_function->AddParameterBlock(HeadN);
100 AddParameterBlocks<TailNs...>(dynamic_function);
101}
102
103// Creates an autodiff cost function wrapping `CostFunctor`, with
104// `kNumResiduals` residuals and parameter blocks with sized `Ns..`.
105// Depending on `kIsDynamic`, either a static or dynamic cost function is
106// created.
107// `args` are forwarded to the `CostFunctor` constructor.
108template <Dynamic kIsDynamic>
109struct CostFunctionFactory {};
110
111template <>
112struct CostFunctionFactory<kNotDynamic> {
113 template <typename CostFunctor,
114 int kNumResiduals,
115 int... Ns,
116 typename... Args>
117 static std::unique_ptr<ceres::CostFunction> Create(Args&&... args) {
118 return std::make_unique<
119 ceres::AutoDiffCostFunction<CostFunctor, kNumResiduals, Ns...>>(
120 new CostFunctor(std::forward<Args>(args)...));
121 }
122};
123
124template <>
125struct CostFunctionFactory<kDynamic> {
126 template <typename CostFunctor,
127 int kNumResiduals,
128 int... Ns,
129 typename... Args>
130 static std::unique_ptr<ceres::CostFunction> Create(Args&&... args) {
131 constexpr const int kNumParameterBlocks = sizeof...(Ns);
132 auto dynamic_function = std::make_unique<ceres::DynamicAutoDiffCostFunction<
133 ToDynamic<CostFunctor, kNumParameterBlocks>>>(
134 new ToDynamic<CostFunctor, kNumParameterBlocks>(
135 std::forward<Args>(args)...));
136 dynamic_function->SetNumResiduals(kNumResiduals);
137 AddParameterBlocks<Ns...>(dynamic_function.get());
138 return dynamic_function;
139 }
140};
141
142template <int kParameterBlockSize, Dynamic kIsDynamic>
143static void BM_ConstantAutodiff(benchmark::State& state) {
144 constexpr int num_residuals = 1;
145 std::array<double, kParameterBlockSize> parameters_values;
146 std::iota(parameters_values.begin(), parameters_values.end(), 0);
147 double* parameters[] = {parameters_values.data()};
148
149 std::array<double, num_residuals> residuals;
150
151 std::array<double, num_residuals * kParameterBlockSize> jacobian_values;
152 double* jacobians[] = {jacobian_values.data()};
153
154 std::unique_ptr<ceres::CostFunction> cost_function =
155 CostFunctionFactory<kIsDynamic>::
156 template Create<ConstantCostFunction<kParameterBlockSize>, 1, 1>();
157
158 for (auto _ : state) {
159 cost_function->Evaluate(parameters, residuals.data(), jacobians);
160 }
161}
162
163BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 1);
164BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 1, kNotDynamic);
165BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 1, kDynamic);
166BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 10);
167BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 10, kNotDynamic);
168BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 10, kDynamic);
169BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 20);
170BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 20, kNotDynamic);
171BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 20, kDynamic);
172BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 30);
173BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 30, kNotDynamic);
174BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 30, kDynamic);
175BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 40);
176BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 40, kNotDynamic);
177BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 40, kDynamic);
178BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 50);
179BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 50, kNotDynamic);
180BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 50, kDynamic);
181BENCHMARK_TEMPLATE(BM_ConstantAnalytic, 60);
182BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 60, kNotDynamic);
183BENCHMARK_TEMPLATE(BM_ConstantAutodiff, 60, kDynamic);
184
185template <Dynamic kIsDynamic>
186static void BM_Linear1AutoDiff(benchmark::State& state) {
187 double parameter_block1[] = {1.};
188 double* parameters[] = {parameter_block1};
189
190 double jacobian1[1];
191 double residuals[1];
192 double* jacobians[] = {jacobian1};
193
194 std::unique_ptr<ceres::CostFunction> cost_function = CostFunctionFactory<
195 kIsDynamic>::template Create<Linear1CostFunction, 1, 1>();
196
197 for (auto _ : state) {
198 cost_function->Evaluate(
199 parameters, residuals, state.range(0) ? jacobians : nullptr);
200 }
201}
202BENCHMARK_TEMPLATE(BM_Linear1AutoDiff, kNotDynamic)->Arg(0)->Arg(1);
203BENCHMARK_TEMPLATE(BM_Linear1AutoDiff, kDynamic)->Arg(0)->Arg(1);
204
205template <Dynamic kIsDynamic>
206static void BM_Linear10AutoDiff(benchmark::State& state) {
207 double parameter_block1[] = {1., 2., 3., 4., 5., 6., 7., 8., 9., 10.};
208 double* parameters[] = {parameter_block1};
209
210 double jacobian1[10 * 10];
211 double residuals[10];
212 double* jacobians[] = {jacobian1};
213
214 std::unique_ptr<ceres::CostFunction> cost_function = CostFunctionFactory<
215 kIsDynamic>::template Create<Linear10CostFunction, 10, 10>();
216
217 for (auto _ : state) {
218 cost_function->Evaluate(
219 parameters, residuals, state.range(0) ? jacobians : nullptr);
220 }
221}
222BENCHMARK_TEMPLATE(BM_Linear10AutoDiff, kNotDynamic)->Arg(0)->Arg(1);
223BENCHMARK_TEMPLATE(BM_Linear10AutoDiff, kDynamic)->Arg(0)->Arg(1);
224
225// From the NIST problem collection.
226struct Rat43CostFunctor {
227 Rat43CostFunctor(const double x, const double y) : x_(x), y_(y) {}
228
229 template <typename T>
230 inline bool operator()(const T* parameters, T* residuals) const {
231 const T& b1 = parameters[0];
232 const T& b2 = parameters[1];
233 const T& b3 = parameters[2];
234 const T& b4 = parameters[3];
235 residuals[0] = b1 * pow(1.0 + exp(b2 - b3 * x_), -1.0 / b4) - y_;
236 return true;
237 }
238
239 static constexpr int kNumParameterBlocks = 1;
240
241 private:
242 const double x_;
243 const double y_;
244};
245
246template <Dynamic kIsDynamic>
247static void BM_Rat43AutoDiff(benchmark::State& state) {
248 double parameter_block1[] = {1., 2., 3., 4.};
249 double* parameters[] = {parameter_block1};
250
251 double jacobian1[] = {0.0, 0.0, 0.0, 0.0};
252 double residuals;
253 double* jacobians[] = {jacobian1};
254 const double x = 0.2;
255 const double y = 0.3;
256 std::unique_ptr<ceres::CostFunction> cost_function =
257 CostFunctionFactory<kIsDynamic>::template Create<Rat43CostFunctor, 1, 4>(
258 x, y);
259
260 for (auto _ : state) {
261 cost_function->Evaluate(
262 parameters, &residuals, state.range(0) ? jacobians : nullptr);
263 }
264}
265BENCHMARK_TEMPLATE(BM_Rat43AutoDiff, kNotDynamic)->Arg(0)->Arg(1);
266BENCHMARK_TEMPLATE(BM_Rat43AutoDiff, kDynamic)->Arg(0)->Arg(1);
267
268template <Dynamic kIsDynamic>
269static void BM_SnavelyReprojectionAutoDiff(benchmark::State& state) {
270 double parameter_block1[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
271 double parameter_block2[] = {1., 2., 3.};
272 double* parameters[] = {parameter_block1, parameter_block2};
273
274 double jacobian1[2 * 9];
275 double jacobian2[2 * 3];
276 double residuals[2];
277 double* jacobians[] = {jacobian1, jacobian2};
278
279 const double x = 0.2;
280 const double y = 0.3;
281 std::unique_ptr<ceres::CostFunction> cost_function = CostFunctionFactory<
282 kIsDynamic>::template Create<SnavelyReprojectionError, 2, 9, 3>(x, y);
283
284 for (auto _ : state) {
285 cost_function->Evaluate(
286 parameters, residuals, state.range(0) ? jacobians : nullptr);
287 }
288}
289
290BENCHMARK_TEMPLATE(BM_SnavelyReprojectionAutoDiff, kNotDynamic)->Arg(0)->Arg(1);
291BENCHMARK_TEMPLATE(BM_SnavelyReprojectionAutoDiff, kDynamic)->Arg(0)->Arg(1);
292
293template <Dynamic kIsDynamic>
294static void BM_PhotometricAutoDiff(benchmark::State& state) {
295 constexpr int PATCH_SIZE = 8;
296
297 using FunctorType = PhotometricError<PATCH_SIZE>;
298 using ImageType = Eigen::Matrix<uint8_t, 128, 128, Eigen::RowMajor>;
299
300 // Prepare parameter / residual / jacobian blocks.
301 double parameter_block1[] = {1., 2., 3., 4., 5., 6., 7.};
302 double parameter_block2[] = {1.1, 2.1, 3.1, 4.1, 5.1, 6.1, 7.1};
303 double parameter_block3[] = {1.};
304 double* parameters[] = {parameter_block1, parameter_block2, parameter_block3};
305
306 Eigen::Map<Eigen::Quaterniond>(parameter_block1).normalize();
307 Eigen::Map<Eigen::Quaterniond>(parameter_block2).normalize();
308
309 double jacobian1[FunctorType::PATCH_SIZE * FunctorType::POSE_SIZE];
310 double jacobian2[FunctorType::PATCH_SIZE * FunctorType::POSE_SIZE];
311 double jacobian3[FunctorType::PATCH_SIZE * FunctorType::POINT_SIZE];
312 double residuals[FunctorType::PATCH_SIZE];
313 double* jacobians[] = {jacobian1, jacobian2, jacobian3};
314
315 // Prepare data (fixed seed for repeatability).
316 std::mt19937::result_type seed = 42;
317 std::mt19937 gen(seed);
318 std::uniform_real_distribution<double> uniform01(0.0, 1.0);
319 std::uniform_int_distribution<unsigned int> uniform0255(0, 255);
320
321 FunctorType::Patch<double> intensities_host =
322 FunctorType::Patch<double>::NullaryExpr(
323 [&]() { return uniform0255(gen); });
324
325 // Set bearing vector's z component to 1, i.e. pointing away from the camera,
326 // to ensure they are (likely) in the domain of the projection function (given
327 // a small rotation between host and target frame).
328 FunctorType::PatchVectors<double> bearings_host =
329 FunctorType::PatchVectors<double>::NullaryExpr(
330 [&]() { return uniform01(gen); });
331 bearings_host.row(2).array() = 1;
332 bearings_host.colwise().normalize();
333
334 ImageType image = ImageType::NullaryExpr(
335 [&]() { return static_cast<uint8_t>(uniform0255(gen)); });
336 FunctorType::Grid grid(image.data(), 0, image.rows(), 0, image.cols());
337 FunctorType::Interpolator image_target(grid);
338
339 FunctorType::Intrinsics intrinsics;
340 intrinsics << 128, 128, 1, -1, 0.5, 0.5;
341
342 std::unique_ptr<ceres::CostFunction> cost_function =
343 CostFunctionFactory<kIsDynamic>::template Create<FunctorType,
344 FunctorType::PATCH_SIZE,
345 FunctorType::POSE_SIZE,
346 FunctorType::POSE_SIZE,
347 FunctorType::POINT_SIZE>(
348 intensities_host, bearings_host, image_target, intrinsics);
349
350 for (auto _ : state) {
351 cost_function->Evaluate(
352 parameters, residuals, state.range(0) ? jacobians : nullptr);
353 }
354}
355
356BENCHMARK_TEMPLATE(BM_PhotometricAutoDiff, kNotDynamic)->Arg(0)->Arg(1);
357BENCHMARK_TEMPLATE(BM_PhotometricAutoDiff, kDynamic)->Arg(0)->Arg(1);
358
359template <Dynamic kIsDynamic>
360static void BM_RelativePoseAutoDiff(benchmark::State& state) {
361 using FunctorType = RelativePoseError;
362
363 double parameter_block1[] = {1., 2., 3., 4., 5., 6., 7.};
364 double parameter_block2[] = {1.1, 2.1, 3.1, 4.1, 5.1, 6.1, 7.1};
365 double* parameters[] = {parameter_block1, parameter_block2};
366
367 Eigen::Map<Eigen::Quaterniond>(parameter_block1).normalize();
368 Eigen::Map<Eigen::Quaterniond>(parameter_block2).normalize();
369
370 double jacobian1[6 * 7];
371 double jacobian2[6 * 7];
372 double residuals[6];
373 double* jacobians[] = {jacobian1, jacobian2};
374
375 Eigen::Quaterniond q_i_j = Eigen::Quaterniond(1, 2, 3, 4).normalized();
376 Eigen::Vector3d t_i_j(1, 2, 3);
377
378 std::unique_ptr<ceres::CostFunction> cost_function =
379 CostFunctionFactory<kIsDynamic>::template Create<FunctorType, 6, 7, 7>(
380 q_i_j, t_i_j);
381
382 for (auto _ : state) {
383 cost_function->Evaluate(
384 parameters, residuals, state.range(0) ? jacobians : nullptr);
385 }
386}
387
388BENCHMARK_TEMPLATE(BM_RelativePoseAutoDiff, kNotDynamic)->Arg(0)->Arg(1);
389BENCHMARK_TEMPLATE(BM_RelativePoseAutoDiff, kDynamic)->Arg(0)->Arg(1);
390
391template <Dynamic kIsDynamic>
392static void BM_BrdfAutoDiff(benchmark::State& state) {
393 using FunctorType = Brdf;
394
395 double material[] = {1., 2., 3., 4., 5., 6., 7., 8., 9., 10.};
396 auto c = Eigen::Vector3d(0.1, 0.2, 0.3);
397 auto n = Eigen::Vector3d(-0.1, 0.5, 0.2).normalized();
398 auto v = Eigen::Vector3d(0.5, -0.2, 0.9).normalized();
399 auto l = Eigen::Vector3d(-0.3, 0.4, -0.3).normalized();
400 auto x = Eigen::Vector3d(0.5, 0.7, -0.1).normalized();
401 auto y = Eigen::Vector3d(0.2, -0.2, -0.2).normalized();
402
403 double* parameters[7] = {
404 material, c.data(), n.data(), v.data(), l.data(), x.data(), y.data()};
405
406 double jacobian[(10 + 6 * 3) * 3];
407 double residuals[3];
408 // clang-format off
409 double* jacobians[7] = {
410 jacobian + 0, jacobian + 10 * 3, jacobian + 13 * 3,
411 jacobian + 16 * 3, jacobian + 19 * 3, jacobian + 22 * 3,
412 jacobian + 25 * 3,
413 };
414 // clang-format on
415
416 std::unique_ptr<ceres::CostFunction> cost_function = CostFunctionFactory<
417 kIsDynamic>::template Create<FunctorType, 3, 10, 3, 3, 3, 3, 3, 3>();
418
419 for (auto _ : state) {
420 cost_function->Evaluate(
421 parameters, residuals, state.range(0) ? jacobians : nullptr);
422 }
423}
424
425BENCHMARK_TEMPLATE(BM_BrdfAutoDiff, kNotDynamic)->Arg(0)->Arg(1);
426BENCHMARK_TEMPLATE(BM_BrdfAutoDiff, kDynamic)->Arg(0)->Arg(1);
427
428} // namespace ceres
429
430BENCHMARK_MAIN();