bundle_adjuster.cc 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394
  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. // Author: sameeragarwal@google.com (Sameer Agarwal)
  30. //
  31. // An example of solving a dynamically sized problem with various
  32. // solvers and loss functions.
  33. //
  34. // For a simpler bare bones example of doing bundle adjustment with
  35. // Ceres, please see simple_bundle_adjuster.cc.
  36. //
  37. // NOTE: This example will not compile without gflags and SuiteSparse.
  38. //
  39. // The problem being solved here is known as a Bundle Adjustment
  40. // problem in computer vision. Given a set of 3d points X_1, ..., X_n,
  41. // a set of cameras P_1, ..., P_m. If the point X_i is visible in
  42. // image j, then there is a 2D observation u_ij that is the expected
  43. // projection of X_i using P_j. The aim of this optimization is to
  44. // find values of X_i and P_j such that the reprojection error
  45. //
  46. // E(X,P) = sum_ij |u_ij - P_j X_i|^2
  47. //
  48. // is minimized.
  49. //
  50. // The problem used here comes from a collection of bundle adjustment
  51. // problems published at University of Washington.
  52. // http://grail.cs.washington.edu/projects/bal
  53. #include <algorithm>
  54. #include <cmath>
  55. #include <cstdio>
  56. #include <cstdlib>
  57. #include <memory>
  58. #include <string>
  59. #include <thread>
  60. #include <vector>
  61. #include "bal_problem.h"
  62. #include "ceres/ceres.h"
  63. #include "gflags/gflags.h"
  64. #include "glog/logging.h"
  65. #include "snavely_reprojection_error.h"
  66. // clang-format makes the gflags definitions too verbose
  67. // clang-format off
  68. DEFINE_string(input, "", "Input File name");
  69. DEFINE_string(trust_region_strategy, "levenberg_marquardt",
  70. "Options are: levenberg_marquardt, dogleg.");
  71. DEFINE_string(dogleg, "traditional_dogleg", "Options are: traditional_dogleg,"
  72. "subspace_dogleg.");
  73. DEFINE_bool(inner_iterations, false, "Use inner iterations to non-linearly "
  74. "refine each successful trust region step.");
  75. DEFINE_string(blocks_for_inner_iterations, "automatic", "Options are: "
  76. "automatic, cameras, points, cameras,points, points,cameras");
  77. DEFINE_string(linear_solver, "sparse_schur", "Options are: "
  78. "sparse_schur, dense_schur, iterative_schur, "
  79. "sparse_normal_cholesky, dense_qr, dense_normal_cholesky, "
  80. "and cgnr.");
  81. DEFINE_bool(explicit_schur_complement, false, "If using ITERATIVE_SCHUR "
  82. "then explicitly compute the Schur complement.");
  83. DEFINE_string(preconditioner, "jacobi", "Options are: "
  84. "identity, jacobi, schur_jacobi, schur_power_series_expansion, cluster_jacobi, "
  85. "cluster_tridiagonal.");
  86. DEFINE_string(visibility_clustering, "canonical_views",
  87. "single_linkage, canonical_views");
  88. DEFINE_bool(use_spse_initialization, false,
  89. "Use power series expansion to initialize the solution in ITERATIVE_SCHUR linear solver.");
  90. DEFINE_string(sparse_linear_algebra_library, "suite_sparse",
  91. "Options are: suite_sparse, accelerate_sparse, eigen_sparse, and "
  92. "cuda_sparse.");
  93. DEFINE_string(dense_linear_algebra_library, "eigen",
  94. "Options are: eigen, lapack, and cuda");
  95. DEFINE_string(ordering_type, "amd", "Options are: amd, nesdis");
  96. DEFINE_string(linear_solver_ordering, "user",
  97. "Options are: automatic and user");
  98. DEFINE_bool(use_quaternions, false, "If true, uses quaternions to represent "
  99. "rotations. If false, angle axis is used.");
  100. DEFINE_bool(use_manifolds, false, "For quaternions, use a manifold.");
  101. DEFINE_bool(robustify, false, "Use a robust loss function.");
  102. DEFINE_double(eta, 1e-2, "Default value for eta. Eta determines the "
  103. "accuracy of each linear solve of the truncated newton step. "
  104. "Changing this parameter can affect solve performance.");
  105. DEFINE_int32(num_threads, -1, "Number of threads. -1 = std::thread::hardware_concurrency.");
  106. DEFINE_int32(num_iterations, 5, "Number of iterations.");
  107. DEFINE_int32(max_linear_solver_iterations, 500, "Maximum number of iterations"
  108. " for solution of linear system.");
  109. DEFINE_double(spse_tolerance, 0.1,
  110. "Tolerance to reach during the iterations of power series expansion initialization or preconditioning.");
  111. DEFINE_int32(max_num_spse_iterations, 5,
  112. "Maximum number of iterations for power series expansion initialization or preconditioning.");
  113. DEFINE_double(max_solver_time, 1e32, "Maximum solve time in seconds.");
  114. DEFINE_bool(nonmonotonic_steps, false, "Trust region algorithm can use"
  115. " nonmonotic steps.");
  116. DEFINE_double(rotation_sigma, 0.0, "Standard deviation of camera rotation "
  117. "perturbation.");
  118. DEFINE_double(translation_sigma, 0.0, "Standard deviation of the camera "
  119. "translation perturbation.");
  120. DEFINE_double(point_sigma, 0.0, "Standard deviation of the point "
  121. "perturbation.");
  122. DEFINE_int32(random_seed, 38401, "Random seed used to set the state "
  123. "of the pseudo random number generator used to generate "
  124. "the perturbations.");
  125. DEFINE_bool(line_search, false, "Use a line search instead of trust region "
  126. "algorithm.");
  127. DEFINE_bool(mixed_precision_solves, false, "Use mixed precision solves.");
  128. DEFINE_int32(max_num_refinement_iterations, 0, "Iterative refinement iterations");
  129. DEFINE_string(initial_ply, "", "Export the BAL file data as a PLY file.");
  130. DEFINE_string(final_ply, "", "Export the refined BAL file data as a PLY "
  131. "file.");
  132. // clang-format on
  133. namespace ceres::examples {
  134. namespace {
  135. void SetLinearSolver(Solver::Options* options) {
  136. CHECK(StringToLinearSolverType(CERES_GET_FLAG(FLAGS_linear_solver),
  137. &options->linear_solver_type));
  138. CHECK(StringToPreconditionerType(CERES_GET_FLAG(FLAGS_preconditioner),
  139. &options->preconditioner_type));
  140. CHECK(StringToVisibilityClusteringType(
  141. CERES_GET_FLAG(FLAGS_visibility_clustering),
  142. &options->visibility_clustering_type));
  143. CHECK(StringToSparseLinearAlgebraLibraryType(
  144. CERES_GET_FLAG(FLAGS_sparse_linear_algebra_library),
  145. &options->sparse_linear_algebra_library_type));
  146. CHECK(StringToDenseLinearAlgebraLibraryType(
  147. CERES_GET_FLAG(FLAGS_dense_linear_algebra_library),
  148. &options->dense_linear_algebra_library_type));
  149. CHECK(
  150. StringToLinearSolverOrderingType(CERES_GET_FLAG(FLAGS_ordering_type),
  151. &options->linear_solver_ordering_type));
  152. options->use_explicit_schur_complement =
  153. CERES_GET_FLAG(FLAGS_explicit_schur_complement);
  154. options->use_mixed_precision_solves =
  155. CERES_GET_FLAG(FLAGS_mixed_precision_solves);
  156. options->max_num_refinement_iterations =
  157. CERES_GET_FLAG(FLAGS_max_num_refinement_iterations);
  158. options->max_linear_solver_iterations =
  159. CERES_GET_FLAG(FLAGS_max_linear_solver_iterations);
  160. options->use_spse_initialization =
  161. CERES_GET_FLAG(FLAGS_use_spse_initialization);
  162. options->spse_tolerance = CERES_GET_FLAG(FLAGS_spse_tolerance);
  163. options->max_num_spse_iterations =
  164. CERES_GET_FLAG(FLAGS_max_num_spse_iterations);
  165. }
  166. void SetOrdering(BALProblem* bal_problem, Solver::Options* options) {
  167. const int num_points = bal_problem->num_points();
  168. const int point_block_size = bal_problem->point_block_size();
  169. double* points = bal_problem->mutable_points();
  170. const int num_cameras = bal_problem->num_cameras();
  171. const int camera_block_size = bal_problem->camera_block_size();
  172. double* cameras = bal_problem->mutable_cameras();
  173. if (options->use_inner_iterations) {
  174. if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) == "cameras") {
  175. LOG(INFO) << "Camera blocks for inner iterations";
  176. options->inner_iteration_ordering =
  177. std::make_shared<ParameterBlockOrdering>();
  178. for (int i = 0; i < num_cameras; ++i) {
  179. options->inner_iteration_ordering->AddElementToGroup(
  180. cameras + camera_block_size * i, 0);
  181. }
  182. } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) == "points") {
  183. LOG(INFO) << "Point blocks for inner iterations";
  184. options->inner_iteration_ordering =
  185. std::make_shared<ParameterBlockOrdering>();
  186. for (int i = 0; i < num_points; ++i) {
  187. options->inner_iteration_ordering->AddElementToGroup(
  188. points + point_block_size * i, 0);
  189. }
  190. } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) ==
  191. "cameras,points") {
  192. LOG(INFO) << "Camera followed by point blocks for inner iterations";
  193. options->inner_iteration_ordering =
  194. std::make_shared<ParameterBlockOrdering>();
  195. for (int i = 0; i < num_cameras; ++i) {
  196. options->inner_iteration_ordering->AddElementToGroup(
  197. cameras + camera_block_size * i, 0);
  198. }
  199. for (int i = 0; i < num_points; ++i) {
  200. options->inner_iteration_ordering->AddElementToGroup(
  201. points + point_block_size * i, 1);
  202. }
  203. } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) ==
  204. "points,cameras") {
  205. LOG(INFO) << "Point followed by camera blocks for inner iterations";
  206. options->inner_iteration_ordering =
  207. std::make_shared<ParameterBlockOrdering>();
  208. for (int i = 0; i < num_cameras; ++i) {
  209. options->inner_iteration_ordering->AddElementToGroup(
  210. cameras + camera_block_size * i, 1);
  211. }
  212. for (int i = 0; i < num_points; ++i) {
  213. options->inner_iteration_ordering->AddElementToGroup(
  214. points + point_block_size * i, 0);
  215. }
  216. } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) ==
  217. "automatic") {
  218. LOG(INFO) << "Choosing automatic blocks for inner iterations";
  219. } else {
  220. LOG(FATAL) << "Unknown block type for inner iterations: "
  221. << CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations);
  222. }
  223. }
  224. // Bundle adjustment problems have a sparsity structure that makes
  225. // them amenable to more specialized and much more efficient
  226. // solution strategies. The SPARSE_SCHUR, DENSE_SCHUR and
  227. // ITERATIVE_SCHUR solvers make use of this specialized
  228. // structure.
  229. //
  230. // This can either be done by specifying a
  231. // Options::linear_solver_ordering or having Ceres figure it out
  232. // automatically using a greedy maximum independent set algorithm.
  233. if (CERES_GET_FLAG(FLAGS_linear_solver_ordering) == "user") {
  234. auto* ordering = new ceres::ParameterBlockOrdering;
  235. // The points come before the cameras.
  236. for (int i = 0; i < num_points; ++i) {
  237. ordering->AddElementToGroup(points + point_block_size * i, 0);
  238. }
  239. for (int i = 0; i < num_cameras; ++i) {
  240. // When using axis-angle, there is a single parameter block for
  241. // the entire camera.
  242. ordering->AddElementToGroup(cameras + camera_block_size * i, 1);
  243. }
  244. options->linear_solver_ordering.reset(ordering);
  245. }
  246. }
  247. void SetMinimizerOptions(Solver::Options* options) {
  248. options->max_num_iterations = CERES_GET_FLAG(FLAGS_num_iterations);
  249. options->minimizer_progress_to_stdout = true;
  250. if (CERES_GET_FLAG(FLAGS_num_threads) == -1) {
  251. const int num_available_threads =
  252. static_cast<int>(std::thread::hardware_concurrency());
  253. if (num_available_threads > 0) {
  254. options->num_threads = num_available_threads;
  255. }
  256. } else {
  257. options->num_threads = CERES_GET_FLAG(FLAGS_num_threads);
  258. }
  259. CHECK_GE(options->num_threads, 1);
  260. options->eta = CERES_GET_FLAG(FLAGS_eta);
  261. options->max_solver_time_in_seconds = CERES_GET_FLAG(FLAGS_max_solver_time);
  262. options->use_nonmonotonic_steps = CERES_GET_FLAG(FLAGS_nonmonotonic_steps);
  263. if (CERES_GET_FLAG(FLAGS_line_search)) {
  264. options->minimizer_type = ceres::LINE_SEARCH;
  265. }
  266. CHECK(StringToTrustRegionStrategyType(
  267. CERES_GET_FLAG(FLAGS_trust_region_strategy),
  268. &options->trust_region_strategy_type));
  269. CHECK(
  270. StringToDoglegType(CERES_GET_FLAG(FLAGS_dogleg), &options->dogleg_type));
  271. options->use_inner_iterations = CERES_GET_FLAG(FLAGS_inner_iterations);
  272. }
  273. void SetSolverOptionsFromFlags(BALProblem* bal_problem,
  274. Solver::Options* options) {
  275. SetMinimizerOptions(options);
  276. SetLinearSolver(options);
  277. SetOrdering(bal_problem, options);
  278. }
  279. void BuildProblem(BALProblem* bal_problem, Problem* problem) {
  280. const int point_block_size = bal_problem->point_block_size();
  281. const int camera_block_size = bal_problem->camera_block_size();
  282. double* points = bal_problem->mutable_points();
  283. double* cameras = bal_problem->mutable_cameras();
  284. // Observations is 2*num_observations long array observations =
  285. // [u_1, u_2, ... , u_n], where each u_i is two dimensional, the x
  286. // and y positions of the observation.
  287. const double* observations = bal_problem->observations();
  288. for (int i = 0; i < bal_problem->num_observations(); ++i) {
  289. CostFunction* cost_function;
  290. // Each Residual block takes a point and a camera as input and
  291. // outputs a 2 dimensional residual.
  292. cost_function = (CERES_GET_FLAG(FLAGS_use_quaternions))
  293. ? SnavelyReprojectionErrorWithQuaternions::Create(
  294. observations[2 * i + 0], observations[2 * i + 1])
  295. : SnavelyReprojectionError::Create(
  296. observations[2 * i + 0], observations[2 * i + 1]);
  297. // If enabled use Huber's loss function.
  298. LossFunction* loss_function =
  299. CERES_GET_FLAG(FLAGS_robustify) ? new HuberLoss(1.0) : nullptr;
  300. // Each observation corresponds to a pair of a camera and a point
  301. // which are identified by camera_index()[i] and point_index()[i]
  302. // respectively.
  303. double* camera =
  304. cameras + camera_block_size * bal_problem->camera_index()[i];
  305. double* point = points + point_block_size * bal_problem->point_index()[i];
  306. problem->AddResidualBlock(cost_function, loss_function, camera, point);
  307. }
  308. if (CERES_GET_FLAG(FLAGS_use_quaternions) &&
  309. CERES_GET_FLAG(FLAGS_use_manifolds)) {
  310. Manifold* camera_manifold =
  311. new ProductManifold<QuaternionManifold, EuclideanManifold<6>>{};
  312. for (int i = 0; i < bal_problem->num_cameras(); ++i) {
  313. problem->SetManifold(cameras + camera_block_size * i, camera_manifold);
  314. }
  315. }
  316. }
  317. void SolveProblem(const char* filename) {
  318. BALProblem bal_problem(filename, CERES_GET_FLAG(FLAGS_use_quaternions));
  319. if (!CERES_GET_FLAG(FLAGS_initial_ply).empty()) {
  320. bal_problem.WriteToPLYFile(CERES_GET_FLAG(FLAGS_initial_ply));
  321. }
  322. Problem problem;
  323. srand(CERES_GET_FLAG(FLAGS_random_seed));
  324. bal_problem.Normalize();
  325. bal_problem.Perturb(CERES_GET_FLAG(FLAGS_rotation_sigma),
  326. CERES_GET_FLAG(FLAGS_translation_sigma),
  327. CERES_GET_FLAG(FLAGS_point_sigma));
  328. BuildProblem(&bal_problem, &problem);
  329. Solver::Options options;
  330. SetSolverOptionsFromFlags(&bal_problem, &options);
  331. options.gradient_tolerance = 1e-16;
  332. options.function_tolerance = 1e-16;
  333. options.parameter_tolerance = 1e-16;
  334. Solver::Summary summary;
  335. Solve(options, &problem, &summary);
  336. std::cout << summary.FullReport() << "\n";
  337. if (!CERES_GET_FLAG(FLAGS_final_ply).empty()) {
  338. bal_problem.WriteToPLYFile(CERES_GET_FLAG(FLAGS_final_ply));
  339. }
  340. }
  341. } // namespace
  342. } // namespace ceres::examples
  343. int main(int argc, char** argv) {
  344. GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
  345. google::InitGoogleLogging(argv[0]);
  346. if (CERES_GET_FLAG(FLAGS_input).empty()) {
  347. LOG(ERROR) << "Usage: bundle_adjuster --input=bal_problem";
  348. return 1;
  349. }
  350. CHECK(CERES_GET_FLAG(FLAGS_use_quaternions) ||
  351. !CERES_GET_FLAG(FLAGS_use_manifolds))
  352. << "--use_manifolds can only be used with --use_quaternions.";
  353. ceres::examples::SolveProblem(CERES_GET_FLAG(FLAGS_input).c_str());
  354. return 0;
  355. }