blob: e7b154e511fd495db99565debed3a5d829b40e59 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 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
54#include <algorithm>
55#include <cmath>
56#include <cstdio>
57#include <cstdlib>
58#include <string>
59#include <vector>
60
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
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080067// clang-format makes the gflags definitions too verbose
68// clang-format off
69
Austin Schuh70cc9552019-01-21 19:46:48 -080070DEFINE_string(input, "", "Input File name");
71DEFINE_string(trust_region_strategy, "levenberg_marquardt",
72 "Options are: levenberg_marquardt, dogleg.");
73DEFINE_string(dogleg, "traditional_dogleg", "Options are: traditional_dogleg,"
74 "subspace_dogleg.");
75
76DEFINE_bool(inner_iterations, false, "Use inner iterations to non-linearly "
77 "refine each successful trust region step.");
78
79DEFINE_string(blocks_for_inner_iterations, "automatic", "Options are: "
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080080 "automatic, cameras, points, cameras,points, points,cameras");
Austin Schuh70cc9552019-01-21 19:46:48 -080081
82DEFINE_string(linear_solver, "sparse_schur", "Options are: "
83 "sparse_schur, dense_schur, iterative_schur, sparse_normal_cholesky, "
84 "dense_qr, dense_normal_cholesky and cgnr.");
85DEFINE_bool(explicit_schur_complement, false, "If using ITERATIVE_SCHUR "
86 "then explicitly compute the Schur complement.");
87DEFINE_string(preconditioner, "jacobi", "Options are: "
88 "identity, jacobi, schur_jacobi, cluster_jacobi, "
89 "cluster_tridiagonal.");
90DEFINE_string(visibility_clustering, "canonical_views",
91 "single_linkage, canonical_views");
92
93DEFINE_string(sparse_linear_algebra_library, "suite_sparse",
94 "Options are: suite_sparse and cx_sparse.");
95DEFINE_string(dense_linear_algebra_library, "eigen",
96 "Options are: eigen and lapack.");
97DEFINE_string(ordering, "automatic", "Options are: automatic, user.");
98
99DEFINE_bool(use_quaternions, false, "If true, uses quaternions to represent "
100 "rotations. If false, angle axis is used.");
101DEFINE_bool(use_local_parameterization, false, "For quaternions, use a local "
102 "parameterization.");
103DEFINE_bool(robustify, false, "Use a robust loss function.");
104
105DEFINE_double(eta, 1e-2, "Default value for eta. Eta determines the "
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800106 "accuracy of each linear solve of the truncated newton step. "
107 "Changing this parameter can affect solve performance.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800108
109DEFINE_int32(num_threads, 1, "Number of threads.");
110DEFINE_int32(num_iterations, 5, "Number of iterations.");
111DEFINE_double(max_solver_time, 1e32, "Maximum solve time in seconds.");
112DEFINE_bool(nonmonotonic_steps, false, "Trust region algorithm can use"
113 " nonmonotic steps.");
114
115DEFINE_double(rotation_sigma, 0.0, "Standard deviation of camera rotation "
116 "perturbation.");
117DEFINE_double(translation_sigma, 0.0, "Standard deviation of the camera "
118 "translation perturbation.");
119DEFINE_double(point_sigma, 0.0, "Standard deviation of the point "
120 "perturbation.");
121DEFINE_int32(random_seed, 38401, "Random seed used to set the state "
122 "of the pseudo random number generator used to generate "
123 "the pertubations.");
124DEFINE_bool(line_search, false, "Use a line search instead of trust region "
125 "algorithm.");
126DEFINE_bool(mixed_precision_solves, false, "Use mixed precision solves.");
127DEFINE_int32(max_num_refinement_iterations, 0, "Iterative refinement iterations");
128DEFINE_string(initial_ply, "", "Export the BAL file data as a PLY file.");
129DEFINE_string(final_ply, "", "Export the refined BAL file data as a PLY "
130 "file.");
131
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800132// clang-format on
133
Austin Schuh70cc9552019-01-21 19:46:48 -0800134namespace ceres {
135namespace examples {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800136namespace {
Austin Schuh70cc9552019-01-21 19:46:48 -0800137
138void SetLinearSolver(Solver::Options* options) {
139 CHECK(StringToLinearSolverType(FLAGS_linear_solver,
140 &options->linear_solver_type));
141 CHECK(StringToPreconditionerType(FLAGS_preconditioner,
142 &options->preconditioner_type));
143 CHECK(StringToVisibilityClusteringType(FLAGS_visibility_clustering,
144 &options->visibility_clustering_type));
145 CHECK(StringToSparseLinearAlgebraLibraryType(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800146 FLAGS_sparse_linear_algebra_library,
147 &options->sparse_linear_algebra_library_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800148 CHECK(StringToDenseLinearAlgebraLibraryType(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800149 FLAGS_dense_linear_algebra_library,
150 &options->dense_linear_algebra_library_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800151 options->use_explicit_schur_complement = FLAGS_explicit_schur_complement;
152 options->use_mixed_precision_solves = FLAGS_mixed_precision_solves;
153 options->max_num_refinement_iterations = FLAGS_max_num_refinement_iterations;
154}
155
156void SetOrdering(BALProblem* bal_problem, Solver::Options* options) {
157 const int num_points = bal_problem->num_points();
158 const int point_block_size = bal_problem->point_block_size();
159 double* points = bal_problem->mutable_points();
160
161 const int num_cameras = bal_problem->num_cameras();
162 const int camera_block_size = bal_problem->camera_block_size();
163 double* cameras = bal_problem->mutable_cameras();
164
165 if (options->use_inner_iterations) {
166 if (FLAGS_blocks_for_inner_iterations == "cameras") {
167 LOG(INFO) << "Camera blocks for inner iterations";
168 options->inner_iteration_ordering.reset(new ParameterBlockOrdering);
169 for (int i = 0; i < num_cameras; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800170 options->inner_iteration_ordering->AddElementToGroup(
171 cameras + camera_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800172 }
173 } else if (FLAGS_blocks_for_inner_iterations == "points") {
174 LOG(INFO) << "Point blocks for inner iterations";
175 options->inner_iteration_ordering.reset(new ParameterBlockOrdering);
176 for (int i = 0; i < num_points; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800177 options->inner_iteration_ordering->AddElementToGroup(
178 points + point_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800179 }
180 } else if (FLAGS_blocks_for_inner_iterations == "cameras,points") {
181 LOG(INFO) << "Camera followed by point blocks for inner iterations";
182 options->inner_iteration_ordering.reset(new ParameterBlockOrdering);
183 for (int i = 0; i < num_cameras; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800184 options->inner_iteration_ordering->AddElementToGroup(
185 cameras + camera_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800186 }
187 for (int i = 0; i < num_points; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800188 options->inner_iteration_ordering->AddElementToGroup(
189 points + point_block_size * i, 1);
Austin Schuh70cc9552019-01-21 19:46:48 -0800190 }
191 } else if (FLAGS_blocks_for_inner_iterations == "points,cameras") {
192 LOG(INFO) << "Point followed by camera blocks for inner iterations";
193 options->inner_iteration_ordering.reset(new ParameterBlockOrdering);
194 for (int i = 0; i < num_cameras; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800195 options->inner_iteration_ordering->AddElementToGroup(
196 cameras + camera_block_size * i, 1);
Austin Schuh70cc9552019-01-21 19:46:48 -0800197 }
198 for (int i = 0; i < num_points; ++i) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800199 options->inner_iteration_ordering->AddElementToGroup(
200 points + point_block_size * i, 0);
Austin Schuh70cc9552019-01-21 19:46:48 -0800201 }
202 } else if (FLAGS_blocks_for_inner_iterations == "automatic") {
203 LOG(INFO) << "Choosing automatic blocks for inner iterations";
204 } else {
205 LOG(FATAL) << "Unknown block type for inner iterations: "
206 << FLAGS_blocks_for_inner_iterations;
207 }
208 }
209
210 // Bundle adjustment problems have a sparsity structure that makes
211 // them amenable to more specialized and much more efficient
212 // solution strategies. The SPARSE_SCHUR, DENSE_SCHUR and
213 // ITERATIVE_SCHUR solvers make use of this specialized
214 // structure.
215 //
216 // This can either be done by specifying Options::ordering_type =
217 // ceres::SCHUR, in which case Ceres will automatically determine
218 // the right ParameterBlock ordering, or by manually specifying a
219 // suitable ordering vector and defining
220 // Options::num_eliminate_blocks.
221 if (FLAGS_ordering == "automatic") {
222 return;
223 }
224
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800225 ceres::ParameterBlockOrdering* ordering = new ceres::ParameterBlockOrdering;
Austin Schuh70cc9552019-01-21 19:46:48 -0800226
227 // The points come before the cameras.
228 for (int i = 0; i < num_points; ++i) {
229 ordering->AddElementToGroup(points + point_block_size * i, 0);
230 }
231
232 for (int i = 0; i < num_cameras; ++i) {
233 // When using axis-angle, there is a single parameter block for
234 // the entire camera.
235 ordering->AddElementToGroup(cameras + camera_block_size * i, 1);
236 }
237
238 options->linear_solver_ordering.reset(ordering);
239}
240
241void SetMinimizerOptions(Solver::Options* options) {
242 options->max_num_iterations = FLAGS_num_iterations;
243 options->minimizer_progress_to_stdout = true;
244 options->num_threads = FLAGS_num_threads;
245 options->eta = FLAGS_eta;
246 options->max_solver_time_in_seconds = FLAGS_max_solver_time;
247 options->use_nonmonotonic_steps = FLAGS_nonmonotonic_steps;
248 if (FLAGS_line_search) {
249 options->minimizer_type = ceres::LINE_SEARCH;
250 }
251
252 CHECK(StringToTrustRegionStrategyType(FLAGS_trust_region_strategy,
253 &options->trust_region_strategy_type));
254 CHECK(StringToDoglegType(FLAGS_dogleg, &options->dogleg_type));
255 options->use_inner_iterations = FLAGS_inner_iterations;
256}
257
258void SetSolverOptionsFromFlags(BALProblem* bal_problem,
259 Solver::Options* options) {
260 SetMinimizerOptions(options);
261 SetLinearSolver(options);
262 SetOrdering(bal_problem, options);
263}
264
265void BuildProblem(BALProblem* bal_problem, Problem* problem) {
266 const int point_block_size = bal_problem->point_block_size();
267 const int camera_block_size = bal_problem->camera_block_size();
268 double* points = bal_problem->mutable_points();
269 double* cameras = bal_problem->mutable_cameras();
270
271 // Observations is 2*num_observations long array observations =
272 // [u_1, u_2, ... , u_n], where each u_i is two dimensional, the x
273 // and y positions of the observation.
274 const double* observations = bal_problem->observations();
275 for (int i = 0; i < bal_problem->num_observations(); ++i) {
276 CostFunction* cost_function;
277 // Each Residual block takes a point and a camera as input and
278 // outputs a 2 dimensional residual.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800279 cost_function = (FLAGS_use_quaternions)
280 ? SnavelyReprojectionErrorWithQuaternions::Create(
281 observations[2 * i + 0], observations[2 * i + 1])
282 : SnavelyReprojectionError::Create(
283 observations[2 * i + 0], observations[2 * i + 1]);
Austin Schuh70cc9552019-01-21 19:46:48 -0800284
285 // If enabled use Huber's loss function.
286 LossFunction* loss_function = FLAGS_robustify ? new HuberLoss(1.0) : NULL;
287
288 // Each observation correponds to a pair of a camera and a point
289 // which are identified by camera_index()[i] and point_index()[i]
290 // respectively.
291 double* camera =
292 cameras + camera_block_size * bal_problem->camera_index()[i];
293 double* point = points + point_block_size * bal_problem->point_index()[i];
294 problem->AddResidualBlock(cost_function, loss_function, camera, point);
295 }
296
297 if (FLAGS_use_quaternions && FLAGS_use_local_parameterization) {
298 LocalParameterization* camera_parameterization =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800299 new ProductParameterization(new QuaternionParameterization(),
300 new IdentityParameterization(6));
Austin Schuh70cc9552019-01-21 19:46:48 -0800301 for (int i = 0; i < bal_problem->num_cameras(); ++i) {
302 problem->SetParameterization(cameras + camera_block_size * i,
303 camera_parameterization);
304 }
305 }
306}
307
308void SolveProblem(const char* filename) {
309 BALProblem bal_problem(filename, FLAGS_use_quaternions);
310
311 if (!FLAGS_initial_ply.empty()) {
312 bal_problem.WriteToPLYFile(FLAGS_initial_ply);
313 }
314
315 Problem problem;
316
317 srand(FLAGS_random_seed);
318 bal_problem.Normalize();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800319 bal_problem.Perturb(
320 FLAGS_rotation_sigma, FLAGS_translation_sigma, FLAGS_point_sigma);
Austin Schuh70cc9552019-01-21 19:46:48 -0800321
322 BuildProblem(&bal_problem, &problem);
323 Solver::Options options;
324 SetSolverOptionsFromFlags(&bal_problem, &options);
325 options.gradient_tolerance = 1e-16;
326 options.function_tolerance = 1e-16;
327 Solver::Summary summary;
328 Solve(options, &problem, &summary);
329 std::cout << summary.FullReport() << "\n";
330
331 if (!FLAGS_final_ply.empty()) {
332 bal_problem.WriteToPLYFile(FLAGS_final_ply);
333 }
334}
335
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800336} // namespace
Austin Schuh70cc9552019-01-21 19:46:48 -0800337} // namespace examples
338} // namespace ceres
339
340int main(int argc, char** argv) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800341 GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
Austin Schuh70cc9552019-01-21 19:46:48 -0800342 google::InitGoogleLogging(argv[0]);
343 if (FLAGS_input.empty()) {
344 LOG(ERROR) << "Usage: bundle_adjuster --input=bal_problem";
345 return 1;
346 }
347
348 CHECK(FLAGS_use_quaternions || !FLAGS_use_local_parameterization)
349 << "--use_local_parameterization can only be used with "
350 << "--use_quaternions.";
351 ceres::examples::SolveProblem(FLAGS_input.c_str());
352 return 0;
353}