blob: 582ae2e752d009f4ea744fe4552a6574d4280457 [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// 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
54#include <algorithm>
55#include <cmath>
56#include <cstdio>
57#include <cstdlib>
Austin Schuh3de38b02024-06-25 18:25:10 -070058#include <memory>
Austin Schuh70cc9552019-01-21 19:46:48 -080059#include <string>
Austin Schuh3de38b02024-06-25 18:25:10 -070060#include <thread>
Austin Schuh70cc9552019-01-21 19:46:48 -080061#include <vector>
62
63#include "bal_problem.h"
64#include "ceres/ceres.h"
65#include "gflags/gflags.h"
66#include "glog/logging.h"
67#include "snavely_reprojection_error.h"
68
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080069// clang-format makes the gflags definitions too verbose
70// clang-format off
71
Austin Schuh70cc9552019-01-21 19:46:48 -080072DEFINE_string(input, "", "Input File name");
73DEFINE_string(trust_region_strategy, "levenberg_marquardt",
74 "Options are: levenberg_marquardt, dogleg.");
75DEFINE_string(dogleg, "traditional_dogleg", "Options are: traditional_dogleg,"
76 "subspace_dogleg.");
77
78DEFINE_bool(inner_iterations, false, "Use inner iterations to non-linearly "
79 "refine each successful trust region step.");
80
81DEFINE_string(blocks_for_inner_iterations, "automatic", "Options are: "
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080082 "automatic, cameras, points, cameras,points, points,cameras");
Austin Schuh70cc9552019-01-21 19:46:48 -080083
84DEFINE_string(linear_solver, "sparse_schur", "Options are: "
Austin Schuh3de38b02024-06-25 18:25:10 -070085 "sparse_schur, dense_schur, iterative_schur, "
86 "sparse_normal_cholesky, dense_qr, dense_normal_cholesky, "
87 "and cgnr.");
Austin Schuh70cc9552019-01-21 19:46:48 -080088DEFINE_bool(explicit_schur_complement, false, "If using ITERATIVE_SCHUR "
89 "then explicitly compute the Schur complement.");
90DEFINE_string(preconditioner, "jacobi", "Options are: "
Austin Schuh3de38b02024-06-25 18:25:10 -070091 "identity, jacobi, schur_jacobi, schur_power_series_expansion, cluster_jacobi, "
Austin Schuh70cc9552019-01-21 19:46:48 -080092 "cluster_tridiagonal.");
93DEFINE_string(visibility_clustering, "canonical_views",
94 "single_linkage, canonical_views");
Austin Schuh3de38b02024-06-25 18:25:10 -070095DEFINE_bool(use_spse_initialization, false,
96 "Use power series expansion to initialize the solution in ITERATIVE_SCHUR linear solver.");
Austin Schuh70cc9552019-01-21 19:46:48 -080097
98DEFINE_string(sparse_linear_algebra_library, "suite_sparse",
Austin Schuh3de38b02024-06-25 18:25:10 -070099 "Options are: suite_sparse, accelerate_sparse, eigen_sparse, and "
100 "cuda_sparse.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800101DEFINE_string(dense_linear_algebra_library, "eigen",
Austin Schuh3de38b02024-06-25 18:25:10 -0700102 "Options are: eigen, lapack, and cuda");
103DEFINE_string(ordering_type, "amd", "Options are: amd, nesdis");
104DEFINE_string(linear_solver_ordering, "user",
105 "Options are: automatic and user");
Austin Schuh70cc9552019-01-21 19:46:48 -0800106
107DEFINE_bool(use_quaternions, false, "If true, uses quaternions to represent "
108 "rotations. If false, angle axis is used.");
Austin Schuh3de38b02024-06-25 18:25:10 -0700109DEFINE_bool(use_manifolds, false, "For quaternions, use a manifold.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800110DEFINE_bool(robustify, false, "Use a robust loss function.");
111
112DEFINE_double(eta, 1e-2, "Default value for eta. Eta determines the "
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800113 "accuracy of each linear solve of the truncated newton step. "
114 "Changing this parameter can affect solve performance.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800115
Austin Schuh3de38b02024-06-25 18:25:10 -0700116DEFINE_int32(num_threads, -1, "Number of threads. -1 = std::thread::hardware_concurrency.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800117DEFINE_int32(num_iterations, 5, "Number of iterations.");
Austin Schuh3de38b02024-06-25 18:25:10 -0700118DEFINE_int32(max_linear_solver_iterations, 500, "Maximum number of iterations"
119 " for solution of linear system.");
120DEFINE_double(spse_tolerance, 0.1,
121 "Tolerance to reach during the iterations of power series expansion initialization or preconditioning.");
122DEFINE_int32(max_num_spse_iterations, 5,
123 "Maximum number of iterations for power series expansion initialization or preconditioning.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800124DEFINE_double(max_solver_time, 1e32, "Maximum solve time in seconds.");
125DEFINE_bool(nonmonotonic_steps, false, "Trust region algorithm can use"
126 " nonmonotic steps.");
127
128DEFINE_double(rotation_sigma, 0.0, "Standard deviation of camera rotation "
129 "perturbation.");
130DEFINE_double(translation_sigma, 0.0, "Standard deviation of the camera "
131 "translation perturbation.");
132DEFINE_double(point_sigma, 0.0, "Standard deviation of the point "
133 "perturbation.");
134DEFINE_int32(random_seed, 38401, "Random seed used to set the state "
135 "of the pseudo random number generator used to generate "
Austin Schuh3de38b02024-06-25 18:25:10 -0700136 "the perturbations.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800137DEFINE_bool(line_search, false, "Use a line search instead of trust region "
138 "algorithm.");
139DEFINE_bool(mixed_precision_solves, false, "Use mixed precision solves.");
140DEFINE_int32(max_num_refinement_iterations, 0, "Iterative refinement iterations");
141DEFINE_string(initial_ply, "", "Export the BAL file data as a PLY file.");
142DEFINE_string(final_ply, "", "Export the refined BAL file data as a PLY "
143 "file.");
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800144// clang-format on
145
Austin Schuh3de38b02024-06-25 18:25:10 -0700146namespace ceres::examples {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800147namespace {
Austin Schuh70cc9552019-01-21 19:46:48 -0800148
149void SetLinearSolver(Solver::Options* options) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700150 CHECK(StringToLinearSolverType(CERES_GET_FLAG(FLAGS_linear_solver),
Austin Schuh70cc9552019-01-21 19:46:48 -0800151 &options->linear_solver_type));
Austin Schuh3de38b02024-06-25 18:25:10 -0700152 CHECK(StringToPreconditionerType(CERES_GET_FLAG(FLAGS_preconditioner),
Austin Schuh70cc9552019-01-21 19:46:48 -0800153 &options->preconditioner_type));
Austin Schuh3de38b02024-06-25 18:25:10 -0700154 CHECK(StringToVisibilityClusteringType(
155 CERES_GET_FLAG(FLAGS_visibility_clustering),
156 &options->visibility_clustering_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800157 CHECK(StringToSparseLinearAlgebraLibraryType(
Austin Schuh3de38b02024-06-25 18:25:10 -0700158 CERES_GET_FLAG(FLAGS_sparse_linear_algebra_library),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800159 &options->sparse_linear_algebra_library_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800160 CHECK(StringToDenseLinearAlgebraLibraryType(
Austin Schuh3de38b02024-06-25 18:25:10 -0700161 CERES_GET_FLAG(FLAGS_dense_linear_algebra_library),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800162 &options->dense_linear_algebra_library_type));
Austin Schuh3de38b02024-06-25 18:25:10 -0700163 CHECK(
164 StringToLinearSolverOrderingType(CERES_GET_FLAG(FLAGS_ordering_type),
165 &options->linear_solver_ordering_type));
166 options->use_explicit_schur_complement =
167 CERES_GET_FLAG(FLAGS_explicit_schur_complement);
168 options->use_mixed_precision_solves =
169 CERES_GET_FLAG(FLAGS_mixed_precision_solves);
170 options->max_num_refinement_iterations =
171 CERES_GET_FLAG(FLAGS_max_num_refinement_iterations);
172 options->max_linear_solver_iterations =
173 CERES_GET_FLAG(FLAGS_max_linear_solver_iterations);
174 options->use_spse_initialization =
175 CERES_GET_FLAG(FLAGS_use_spse_initialization);
176 options->spse_tolerance = CERES_GET_FLAG(FLAGS_spse_tolerance);
177 options->max_num_spse_iterations =
178 CERES_GET_FLAG(FLAGS_max_num_spse_iterations);
Austin Schuh70cc9552019-01-21 19:46:48 -0800179}
180
181void SetOrdering(BALProblem* bal_problem, Solver::Options* options) {
182 const int num_points = bal_problem->num_points();
183 const int point_block_size = bal_problem->point_block_size();
184 double* points = bal_problem->mutable_points();
185
186 const int num_cameras = bal_problem->num_cameras();
187 const int camera_block_size = bal_problem->camera_block_size();
188 double* cameras = bal_problem->mutable_cameras();
189
190 if (options->use_inner_iterations) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700191 if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) == "cameras") {
Austin Schuh70cc9552019-01-21 19:46:48 -0800192 LOG(INFO) << "Camera blocks for inner iterations";
Austin Schuh3de38b02024-06-25 18:25:10 -0700193 options->inner_iteration_ordering =
194 std::make_shared<ParameterBlockOrdering>();
Austin Schuh70cc9552019-01-21 19:46:48 -0800195 for (int i = 0; i < num_cameras; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800196 options->inner_iteration_ordering->AddElementToGroup(
197 cameras + camera_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800198 }
Austin Schuh3de38b02024-06-25 18:25:10 -0700199 } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) == "points") {
Austin Schuh70cc9552019-01-21 19:46:48 -0800200 LOG(INFO) << "Point blocks for inner iterations";
Austin Schuh3de38b02024-06-25 18:25:10 -0700201 options->inner_iteration_ordering =
202 std::make_shared<ParameterBlockOrdering>();
Austin Schuh70cc9552019-01-21 19:46:48 -0800203 for (int i = 0; i < num_points; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800204 options->inner_iteration_ordering->AddElementToGroup(
205 points + point_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800206 }
Austin Schuh3de38b02024-06-25 18:25:10 -0700207 } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) ==
208 "cameras,points") {
Austin Schuh70cc9552019-01-21 19:46:48 -0800209 LOG(INFO) << "Camera followed by point blocks for inner iterations";
Austin Schuh3de38b02024-06-25 18:25:10 -0700210 options->inner_iteration_ordering =
211 std::make_shared<ParameterBlockOrdering>();
Austin Schuh70cc9552019-01-21 19:46:48 -0800212 for (int i = 0; i < num_cameras; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800213 options->inner_iteration_ordering->AddElementToGroup(
214 cameras + camera_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800215 }
216 for (int i = 0; i < num_points; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800217 options->inner_iteration_ordering->AddElementToGroup(
218 points + point_block_size * i, 1);
Austin Schuh70cc9552019-01-21 19:46:48 -0800219 }
Austin Schuh3de38b02024-06-25 18:25:10 -0700220 } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) ==
221 "points,cameras") {
Austin Schuh70cc9552019-01-21 19:46:48 -0800222 LOG(INFO) << "Point followed by camera blocks for inner iterations";
Austin Schuh3de38b02024-06-25 18:25:10 -0700223 options->inner_iteration_ordering =
224 std::make_shared<ParameterBlockOrdering>();
Austin Schuh70cc9552019-01-21 19:46:48 -0800225 for (int i = 0; i < num_cameras; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800226 options->inner_iteration_ordering->AddElementToGroup(
227 cameras + camera_block_size * i, 1);
Austin Schuh70cc9552019-01-21 19:46:48 -0800228 }
229 for (int i = 0; i < num_points; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800230 options->inner_iteration_ordering->AddElementToGroup(
231 points + point_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800232 }
Austin Schuh3de38b02024-06-25 18:25:10 -0700233 } else if (CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations) ==
234 "automatic") {
Austin Schuh70cc9552019-01-21 19:46:48 -0800235 LOG(INFO) << "Choosing automatic blocks for inner iterations";
236 } else {
237 LOG(FATAL) << "Unknown block type for inner iterations: "
Austin Schuh3de38b02024-06-25 18:25:10 -0700238 << CERES_GET_FLAG(FLAGS_blocks_for_inner_iterations);
Austin Schuh70cc9552019-01-21 19:46:48 -0800239 }
240 }
241
242 // Bundle adjustment problems have a sparsity structure that makes
243 // them amenable to more specialized and much more efficient
244 // solution strategies. The SPARSE_SCHUR, DENSE_SCHUR and
245 // ITERATIVE_SCHUR solvers make use of this specialized
246 // structure.
247 //
Austin Schuh3de38b02024-06-25 18:25:10 -0700248 // This can either be done by specifying a
249 // Options::linear_solver_ordering or having Ceres figure it out
250 // automatically using a greedy maximum independent set algorithm.
251 if (CERES_GET_FLAG(FLAGS_linear_solver_ordering) == "user") {
252 auto* ordering = new ceres::ParameterBlockOrdering;
253
254 // The points come before the cameras.
255 for (int i = 0; i < num_points; ++i) {
256 ordering->AddElementToGroup(points + point_block_size * i, 0);
257 }
258
259 for (int i = 0; i < num_cameras; ++i) {
260 // When using axis-angle, there is a single parameter block for
261 // the entire camera.
262 ordering->AddElementToGroup(cameras + camera_block_size * i, 1);
263 }
264
265 options->linear_solver_ordering.reset(ordering);
Austin Schuh70cc9552019-01-21 19:46:48 -0800266 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800267}
268
269void SetMinimizerOptions(Solver::Options* options) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700270 options->max_num_iterations = CERES_GET_FLAG(FLAGS_num_iterations);
Austin Schuh70cc9552019-01-21 19:46:48 -0800271 options->minimizer_progress_to_stdout = true;
Austin Schuh3de38b02024-06-25 18:25:10 -0700272 if (CERES_GET_FLAG(FLAGS_num_threads) == -1) {
273 const int num_available_threads =
274 static_cast<int>(std::thread::hardware_concurrency());
275 if (num_available_threads > 0) {
276 options->num_threads = num_available_threads;
277 }
278 } else {
279 options->num_threads = CERES_GET_FLAG(FLAGS_num_threads);
280 }
281 CHECK_GE(options->num_threads, 1);
282
283 options->eta = CERES_GET_FLAG(FLAGS_eta);
284 options->max_solver_time_in_seconds = CERES_GET_FLAG(FLAGS_max_solver_time);
285 options->use_nonmonotonic_steps = CERES_GET_FLAG(FLAGS_nonmonotonic_steps);
286 if (CERES_GET_FLAG(FLAGS_line_search)) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800287 options->minimizer_type = ceres::LINE_SEARCH;
288 }
289
Austin Schuh3de38b02024-06-25 18:25:10 -0700290 CHECK(StringToTrustRegionStrategyType(
291 CERES_GET_FLAG(FLAGS_trust_region_strategy),
292 &options->trust_region_strategy_type));
293 CHECK(
294 StringToDoglegType(CERES_GET_FLAG(FLAGS_dogleg), &options->dogleg_type));
295 options->use_inner_iterations = CERES_GET_FLAG(FLAGS_inner_iterations);
Austin Schuh70cc9552019-01-21 19:46:48 -0800296}
297
298void SetSolverOptionsFromFlags(BALProblem* bal_problem,
299 Solver::Options* options) {
300 SetMinimizerOptions(options);
301 SetLinearSolver(options);
302 SetOrdering(bal_problem, options);
303}
304
305void BuildProblem(BALProblem* bal_problem, Problem* problem) {
306 const int point_block_size = bal_problem->point_block_size();
307 const int camera_block_size = bal_problem->camera_block_size();
308 double* points = bal_problem->mutable_points();
309 double* cameras = bal_problem->mutable_cameras();
310
311 // Observations is 2*num_observations long array observations =
312 // [u_1, u_2, ... , u_n], where each u_i is two dimensional, the x
313 // and y positions of the observation.
314 const double* observations = bal_problem->observations();
315 for (int i = 0; i < bal_problem->num_observations(); ++i) {
316 CostFunction* cost_function;
317 // Each Residual block takes a point and a camera as input and
318 // outputs a 2 dimensional residual.
Austin Schuh3de38b02024-06-25 18:25:10 -0700319 cost_function = (CERES_GET_FLAG(FLAGS_use_quaternions))
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800320 ? SnavelyReprojectionErrorWithQuaternions::Create(
321 observations[2 * i + 0], observations[2 * i + 1])
322 : SnavelyReprojectionError::Create(
323 observations[2 * i + 0], observations[2 * i + 1]);
Austin Schuh70cc9552019-01-21 19:46:48 -0800324
325 // If enabled use Huber's loss function.
Austin Schuh3de38b02024-06-25 18:25:10 -0700326 LossFunction* loss_function =
327 CERES_GET_FLAG(FLAGS_robustify) ? new HuberLoss(1.0) : nullptr;
Austin Schuh70cc9552019-01-21 19:46:48 -0800328
Austin Schuh3de38b02024-06-25 18:25:10 -0700329 // Each observation corresponds to a pair of a camera and a point
Austin Schuh70cc9552019-01-21 19:46:48 -0800330 // which are identified by camera_index()[i] and point_index()[i]
331 // respectively.
332 double* camera =
333 cameras + camera_block_size * bal_problem->camera_index()[i];
334 double* point = points + point_block_size * bal_problem->point_index()[i];
335 problem->AddResidualBlock(cost_function, loss_function, camera, point);
336 }
337
Austin Schuh3de38b02024-06-25 18:25:10 -0700338 if (CERES_GET_FLAG(FLAGS_use_quaternions) &&
339 CERES_GET_FLAG(FLAGS_use_manifolds)) {
340 Manifold* camera_manifold =
341 new ProductManifold<QuaternionManifold, EuclideanManifold<6>>{};
Austin Schuh70cc9552019-01-21 19:46:48 -0800342 for (int i = 0; i < bal_problem->num_cameras(); ++i) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700343 problem->SetManifold(cameras + camera_block_size * i, camera_manifold);
Austin Schuh70cc9552019-01-21 19:46:48 -0800344 }
345 }
346}
347
348void SolveProblem(const char* filename) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700349 BALProblem bal_problem(filename, CERES_GET_FLAG(FLAGS_use_quaternions));
Austin Schuh70cc9552019-01-21 19:46:48 -0800350
Austin Schuh3de38b02024-06-25 18:25:10 -0700351 if (!CERES_GET_FLAG(FLAGS_initial_ply).empty()) {
352 bal_problem.WriteToPLYFile(CERES_GET_FLAG(FLAGS_initial_ply));
Austin Schuh70cc9552019-01-21 19:46:48 -0800353 }
354
355 Problem problem;
356
Austin Schuh3de38b02024-06-25 18:25:10 -0700357 srand(CERES_GET_FLAG(FLAGS_random_seed));
Austin Schuh70cc9552019-01-21 19:46:48 -0800358 bal_problem.Normalize();
Austin Schuh3de38b02024-06-25 18:25:10 -0700359 bal_problem.Perturb(CERES_GET_FLAG(FLAGS_rotation_sigma),
360 CERES_GET_FLAG(FLAGS_translation_sigma),
361 CERES_GET_FLAG(FLAGS_point_sigma));
Austin Schuh70cc9552019-01-21 19:46:48 -0800362
363 BuildProblem(&bal_problem, &problem);
364 Solver::Options options;
365 SetSolverOptionsFromFlags(&bal_problem, &options);
366 options.gradient_tolerance = 1e-16;
367 options.function_tolerance = 1e-16;
Austin Schuh3de38b02024-06-25 18:25:10 -0700368 options.parameter_tolerance = 1e-16;
Austin Schuh70cc9552019-01-21 19:46:48 -0800369 Solver::Summary summary;
370 Solve(options, &problem, &summary);
371 std::cout << summary.FullReport() << "\n";
372
Austin Schuh3de38b02024-06-25 18:25:10 -0700373 if (!CERES_GET_FLAG(FLAGS_final_ply).empty()) {
374 bal_problem.WriteToPLYFile(CERES_GET_FLAG(FLAGS_final_ply));
Austin Schuh70cc9552019-01-21 19:46:48 -0800375 }
376}
377
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800378} // namespace
Austin Schuh3de38b02024-06-25 18:25:10 -0700379} // namespace ceres::examples
Austin Schuh70cc9552019-01-21 19:46:48 -0800380
381int main(int argc, char** argv) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800382 GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
Austin Schuh70cc9552019-01-21 19:46:48 -0800383 google::InitGoogleLogging(argv[0]);
Austin Schuh3de38b02024-06-25 18:25:10 -0700384 if (CERES_GET_FLAG(FLAGS_input).empty()) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800385 LOG(ERROR) << "Usage: bundle_adjuster --input=bal_problem";
386 return 1;
387 }
388
Austin Schuh3de38b02024-06-25 18:25:10 -0700389 CHECK(CERES_GET_FLAG(FLAGS_use_quaternions) ||
390 !CERES_GET_FLAG(FLAGS_use_manifolds))
391 << "--use_manifolds can only be used with --use_quaternions.";
392 ceres::examples::SolveProblem(CERES_GET_FLAG(FLAGS_input).c_str());
Austin Schuh70cc9552019-01-21 19:46:48 -0800393 return 0;
394}