autodiff_benchmarks.cc 16 KB

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