blob: 68458a0dfdd272d24b20bb06624ff70b14a54bc4 [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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>
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -070057#include <cstdlib>
Keir Mierle8ebb0732012-04-30 23:09:08 -070058#include <string>
59#include <vector>
60
Keir Mierle8ebb0732012-04-30 23:09:08 -070061#include "bal_problem.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070062#include "ceres/ceres.h"
Sameer Agarwal5476df52012-08-09 21:46:19 -070063#include "ceres/random.h"
Sameer Agarwal4b040432012-08-13 14:38:41 -070064#include "gflags/gflags.h"
65#include "glog/logging.h"
66#include "snavely_reprojection_error.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070067
68DEFINE_string(input, "", "Input File name");
Sameer Agarwal5476df52012-08-09 21:46:19 -070069DEFINE_bool(use_quaternions, false, "If true, uses quaternions to represent "
70 "rotations. If false, angle axis is used");
71DEFINE_bool(use_local_parameterization, false, "For quaternions, use a local "
72 "parameterization.");
73DEFINE_bool(robustify, false, "Use a robust loss function");
74
75DEFINE_string(trust_region_strategy, "lm", "Options are: lm, dogleg");
76DEFINE_double(eta, 1e-2, "Default value for eta. Eta determines the "
77 "accuracy of each linear solve of the truncated newton step. "
78 "Changing this parameter can affect solve performance ");
Keir Mierle8ebb0732012-04-30 23:09:08 -070079DEFINE_string(solver_type, "sparse_schur", "Options are: "
Sameer Agarwalb9f15a52012-08-18 13:06:19 -070080 "sparse_schur, dense_schur, iterative_schur, sparse_cholesky, "
81 "dense_qr, dense_cholesky and conjugate_gradients");
Keir Mierle8ebb0732012-04-30 23:09:08 -070082DEFINE_string(preconditioner_type, "jacobi", "Options are: "
83 "identity, jacobi, schur_jacobi, cluster_jacobi, "
84 "cluster_tridiagonal");
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -070085DEFINE_string(sparse_linear_algebra_library, "suitesparse",
86 "Options are: suitesparse and cxsparse");
Sameer Agarwal5476df52012-08-09 21:46:19 -070087
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -070088DEFINE_string(ordering_type, "schur", "Options are: schur, user, natural");
Markus Moll51cf7cb2012-08-20 20:10:20 +020089DEFINE_string(dogleg_type, "traditional", "Options are: traditional, subspace");
Sameer Agarwal5476df52012-08-09 21:46:19 -070090DEFINE_bool(use_block_amd, true, "Use a block oriented fill reducing "
91 "ordering.");
92
93DEFINE_int32(num_threads, 1, "Number of threads");
94DEFINE_int32(num_iterations, 5, "Number of iterations");
Sameer Agarwalfa015192012-06-11 14:21:42 -070095DEFINE_double(max_solver_time, 1e32, "Maximum solve time in seconds.");
Sameer Agarwala8f87d72012-08-08 10:38:31 -070096DEFINE_bool(nonmonotonic_steps, false, "Trust region algorithm can use"
97 " nonmonotic steps");
Keir Mierle8ebb0732012-04-30 23:09:08 -070098
Sameer Agarwal5476df52012-08-09 21:46:19 -070099DEFINE_double(rotation_sigma, 0.0, "Standard deviation of camera rotation "
100 "perturbation.");
101DEFINE_double(translation_sigma, 0.0, "Standard deviation of the camera "
102 "translation perturbation.");
103DEFINE_double(point_sigma, 0.0, "Standard deviation of the point "
104 "perturbation");
105DEFINE_int32(random_seed, 38401, "Random seed used to set the state "
106 "of the pseudo random number generator used to generate "
107 "the pertubations.");
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700108DEFINE_string(solver_log, "", "File to record the solver execution to.");
Sameer Agarwal5476df52012-08-09 21:46:19 -0700109
Keir Mierle8ebb0732012-04-30 23:09:08 -0700110namespace ceres {
111namespace examples {
112
113void SetLinearSolver(Solver::Options* options) {
114 if (FLAGS_solver_type == "sparse_schur") {
115 options->linear_solver_type = ceres::SPARSE_SCHUR;
116 } else if (FLAGS_solver_type == "dense_schur") {
117 options->linear_solver_type = ceres::DENSE_SCHUR;
118 } else if (FLAGS_solver_type == "iterative_schur") {
119 options->linear_solver_type = ceres::ITERATIVE_SCHUR;
Sameer Agarwalb9f15a52012-08-18 13:06:19 -0700120 } else if (FLAGS_solver_type == "sparse_cholesky") {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700121 options->linear_solver_type = ceres::SPARSE_NORMAL_CHOLESKY;
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700122 } else if (FLAGS_solver_type == "cgnr") {
123 options->linear_solver_type = ceres::CGNR;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700124 } else if (FLAGS_solver_type == "dense_qr") {
125 // DENSE_QR is included here for completeness, but actually using
Keir Mierle211812a2012-05-07 04:33:50 -0700126 // this option is a bad idea due to the amount of memory needed
Keir Mierle8ebb0732012-04-30 23:09:08 -0700127 // to store even the smallest of the bundle adjustment jacobian
128 // arrays
129 options->linear_solver_type = ceres::DENSE_QR;
Sameer Agarwalb9f15a52012-08-18 13:06:19 -0700130 } else if (FLAGS_solver_type == "dense_cholesky") {
131 // DENSE_NORMAL_CHOLESKY is included here for completeness, but
132 // actually using this option is a bad idea due to the amount of
133 // memory needed to store even the smallest of the bundle
134 // adjustment jacobian arrays
135 options->linear_solver_type = ceres::DENSE_NORMAL_CHOLESKY;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700136 } else {
137 LOG(FATAL) << "Unknown ceres solver type: "
138 << FLAGS_solver_type;
139 }
140
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700141 if (options->linear_solver_type == ceres::CGNR) {
Keir Mierle211812a2012-05-07 04:33:50 -0700142 options->linear_solver_min_num_iterations = 5;
143 if (FLAGS_preconditioner_type == "identity") {
144 options->preconditioner_type = ceres::IDENTITY;
145 } else if (FLAGS_preconditioner_type == "jacobi") {
146 options->preconditioner_type = ceres::JACOBI;
147 } else {
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700148 LOG(FATAL) << "For CGNR, only identity and jacobian "
Keir Mierle211812a2012-05-07 04:33:50 -0700149 << "preconditioners are supported. Got: "
150 << FLAGS_preconditioner_type;
151 }
152 }
153
154 if (options->linear_solver_type == ceres::ITERATIVE_SCHUR) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700155 options->linear_solver_min_num_iterations = 5;
156 if (FLAGS_preconditioner_type == "identity") {
157 options->preconditioner_type = ceres::IDENTITY;
158 } else if (FLAGS_preconditioner_type == "jacobi") {
159 options->preconditioner_type = ceres::JACOBI;
160 } else if (FLAGS_preconditioner_type == "schur_jacobi") {
161 options->preconditioner_type = ceres::SCHUR_JACOBI;
162 } else if (FLAGS_preconditioner_type == "cluster_jacobi") {
163 options->preconditioner_type = ceres::CLUSTER_JACOBI;
164 } else if (FLAGS_preconditioner_type == "cluster_tridiagonal") {
165 options->preconditioner_type = ceres::CLUSTER_TRIDIAGONAL;
166 } else {
167 LOG(FATAL) << "Unknown ceres preconditioner type: "
Keir Mierle211812a2012-05-07 04:33:50 -0700168 << FLAGS_preconditioner_type;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700169 }
170 }
171
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700172 if (FLAGS_sparse_linear_algebra_library == "suitesparse") {
173 options->sparse_linear_algebra_library = SUITE_SPARSE;
174 } else if (FLAGS_sparse_linear_algebra_library == "cxsparse") {
175 options->sparse_linear_algebra_library = CX_SPARSE;
176 } else {
177 LOG(FATAL) << "Unknown sparse linear algebra library type.";
178 }
179
Keir Mierle8ebb0732012-04-30 23:09:08 -0700180 options->num_linear_solver_threads = FLAGS_num_threads;
181}
182
183void SetOrdering(BALProblem* bal_problem, Solver::Options* options) {
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700184 options->use_block_amd = FLAGS_use_block_amd;
185
186 // Only non-Schur solvers support the natural ordering for this
187 // problem.
188 if (FLAGS_ordering_type == "natural") {
189 if (options->linear_solver_type == SPARSE_SCHUR ||
190 options->linear_solver_type == DENSE_SCHUR ||
191 options->linear_solver_type == ITERATIVE_SCHUR) {
192 LOG(FATAL) << "Natural ordering with Schur type solver does not work.";
193 }
194 return;
195 }
196
Keir Mierle8ebb0732012-04-30 23:09:08 -0700197 // Bundle adjustment problems have a sparsity structure that makes
198 // them amenable to more specialized and much more efficient
199 // solution strategies. The SPARSE_SCHUR, DENSE_SCHUR and
200 // ITERATIVE_SCHUR solvers make use of this specialized
201 // structure. Using them however requires that the ParameterBlocks
202 // are in a particular order (points before cameras) and
203 // Solver::Options::num_eliminate_blocks is set to the number of
204 // points.
205 //
206 // This can either be done by specifying Options::ordering_type =
207 // ceres::SCHUR, in which case Ceres will automatically determine
208 // the right ParameterBlock ordering, or by manually specifying a
209 // suitable ordering vector and defining
210 // Options::num_eliminate_blocks.
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700211 if (FLAGS_ordering_type == "schur") {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700212 options->ordering_type = ceres::SCHUR;
213 return;
214 }
215
216 options->ordering_type = ceres::USER;
217 const int num_points = bal_problem->num_points();
218 const int point_block_size = bal_problem->point_block_size();
219 double* points = bal_problem->mutable_points();
220 const int num_cameras = bal_problem->num_cameras();
221 const int camera_block_size = bal_problem->camera_block_size();
222 double* cameras = bal_problem->mutable_cameras();
223
224 // The points come before the cameras.
225 for (int i = 0; i < num_points; ++i) {
226 options->ordering.push_back(points + point_block_size * i);
227 }
228
229 for (int i = 0; i < num_cameras; ++i) {
230 // When using axis-angle, there is a single parameter block for
231 // the entire camera.
232 options->ordering.push_back(cameras + camera_block_size * i);
233
234 // If quaternions are used, there are two blocks, so add the
235 // second block to the ordering.
236 if (FLAGS_use_quaternions) {
237 options->ordering.push_back(cameras + camera_block_size * i + 4);
238 }
239 }
240
241 options->num_eliminate_blocks = num_points;
242}
243
244void SetMinimizerOptions(Solver::Options* options) {
245 options->max_num_iterations = FLAGS_num_iterations;
246 options->minimizer_progress_to_stdout = true;
247 options->num_threads = FLAGS_num_threads;
Keir Mierlef7898fb2012-05-05 20:55:08 -0700248 options->eta = FLAGS_eta;
Sameer Agarwalfa015192012-06-11 14:21:42 -0700249 options->max_solver_time_in_seconds = FLAGS_max_solver_time;
Sameer Agarwala8f87d72012-08-08 10:38:31 -0700250 options->use_nonmonotonic_steps = FLAGS_nonmonotonic_steps;
Sameer Agarwalfa015192012-06-11 14:21:42 -0700251 if (FLAGS_trust_region_strategy == "lm") {
252 options->trust_region_strategy_type = LEVENBERG_MARQUARDT;
253 } else if (FLAGS_trust_region_strategy == "dogleg") {
254 options->trust_region_strategy_type = DOGLEG;
255 } else {
256 LOG(FATAL) << "Unknown trust region strategy: "
257 << FLAGS_trust_region_strategy;
258 }
Markus Moll51cf7cb2012-08-20 20:10:20 +0200259 if (FLAGS_dogleg_type == "traditional") {
260 options->dogleg_type = TRADITIONAL_DOGLEG;
261 } else if (FLAGS_dogleg_type == "subspace") {
262 options->dogleg_type = SUBSPACE_DOGLEG;
263 } else {
264 LOG(FATAL) << "Unknown dogleg type: "
265 << FLAGS_dogleg_type;
266 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700267}
268
269void SetSolverOptionsFromFlags(BALProblem* bal_problem,
270 Solver::Options* options) {
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700271 SetMinimizerOptions(options);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700272 SetLinearSolver(options);
273 SetOrdering(bal_problem, options);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700274}
275
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700276// Uniform random numbers between 0 and 1.
277double UniformRandom() {
278 return static_cast<double>(random()) / static_cast<double>(RAND_MAX);
279}
280
281// Normal random numbers using the Box-Mueller algorithm. Its a bit
282// wasteful, as it generates two but only returns one.
283double RandNormal() {
284 double x1, x2, w, y1, y2;
285 do {
286 x1 = 2.0 * UniformRandom() - 1.0;
287 x2 = 2.0 * UniformRandom() - 1.0;
288 w = x1 * x1 + x2 * x2;
289 } while ( w >= 1.0 );
290
291 w = sqrt((-2.0 * log(w)) / w);
292 y1 = x1 * w;
293 y2 = x2 * w;
294 return y1;
295}
296
Keir Mierle8ebb0732012-04-30 23:09:08 -0700297void BuildProblem(BALProblem* bal_problem, Problem* problem) {
298 const int point_block_size = bal_problem->point_block_size();
299 const int camera_block_size = bal_problem->camera_block_size();
300 double* points = bal_problem->mutable_points();
301 double* cameras = bal_problem->mutable_cameras();
302
303 // Observations is 2*num_observations long array observations =
304 // [u_1, u_2, ... , u_n], where each u_i is two dimensional, the x
305 // and y positions of the observation.
306 const double* observations = bal_problem->observations();
307
308 for (int i = 0; i < bal_problem->num_observations(); ++i) {
309 CostFunction* cost_function;
310 // Each Residual block takes a point and a camera as input and
311 // outputs a 2 dimensional residual.
312 if (FLAGS_use_quaternions) {
313 cost_function = new AutoDiffCostFunction<
Sameer Agarwal5476df52012-08-09 21:46:19 -0700314 SnavelyReprojectionErrorWithQuaternions, 2, 4, 6, 3>(
315 new SnavelyReprojectionErrorWithQuaternions(
Keir Mierle8ebb0732012-04-30 23:09:08 -0700316 observations[2 * i + 0],
317 observations[2 * i + 1]));
318 } else {
319 cost_function =
320 new AutoDiffCostFunction<SnavelyReprojectionError, 2, 9, 3>(
321 new SnavelyReprojectionError(observations[2 * i + 0],
322 observations[2 * i + 1]));
323 }
324
325 // If enabled use Huber's loss function.
326 LossFunction* loss_function = FLAGS_robustify ? new HuberLoss(1.0) : NULL;
327
328 // Each observation correponds to a pair of a camera and a point
329 // which are identified by camera_index()[i] and point_index()[i]
330 // respectively.
331 double* camera =
332 cameras + camera_block_size * bal_problem->camera_index()[i];
333 double* point = points + point_block_size * bal_problem->point_index()[i];
334
335 if (FLAGS_use_quaternions) {
336 // When using quaternions, we split the camera into two
337 // parameter blocks. One of size 4 for the quaternion and the
338 // other of size 6 containing the translation, focal length and
339 // the radial distortion parameters.
340 problem->AddResidualBlock(cost_function,
341 loss_function,
342 camera,
343 camera + 4,
344 point);
345 } else {
346 problem->AddResidualBlock(cost_function, loss_function, camera, point);
347 }
348 }
349
350 if (FLAGS_use_quaternions && FLAGS_use_local_parameterization) {
351 LocalParameterization* quaternion_parameterization =
352 new QuaternionParameterization;
353 for (int i = 0; i < bal_problem->num_cameras(); ++i) {
354 problem->SetParameterization(cameras + camera_block_size * i,
355 quaternion_parameterization);
356 }
357 }
358}
359
360void SolveProblem(const char* filename) {
361 BALProblem bal_problem(filename, FLAGS_use_quaternions);
362 Problem problem;
Sameer Agarwal5476df52012-08-09 21:46:19 -0700363
364 SetRandomState(FLAGS_random_seed);
Sameer Agarwal36a33092012-08-10 19:52:09 -0700365 bal_problem.Normalize();
Sameer Agarwal5476df52012-08-09 21:46:19 -0700366 bal_problem.Perturb(FLAGS_rotation_sigma,
367 FLAGS_translation_sigma,
368 FLAGS_point_sigma);
369
Keir Mierle8ebb0732012-04-30 23:09:08 -0700370 BuildProblem(&bal_problem, &problem);
371 Solver::Options options;
372 SetSolverOptionsFromFlags(&bal_problem, &options);
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700373 options.solver_log = FLAGS_solver_log;
Sameer Agarwal36a33092012-08-10 19:52:09 -0700374 options.gradient_tolerance *= 1e-3;
375
Keir Mierle8ebb0732012-04-30 23:09:08 -0700376 Solver::Summary summary;
377 Solve(options, &problem, &summary);
378 std::cout << summary.FullReport() << "\n";
379}
380
381} // namespace examples
382} // namespace ceres
383
384int main(int argc, char** argv) {
385 google::ParseCommandLineFlags(&argc, &argv, true);
386 google::InitGoogleLogging(argv[0]);
387 if (FLAGS_input.empty()) {
388 LOG(ERROR) << "Usage: bundle_adjustment_example --input=bal_problem";
389 return 1;
390 }
391
392 CHECK(FLAGS_use_quaternions || !FLAGS_use_local_parameterization)
393 << "--use_local_parameterization can only be used with "
394 << "--use_quaternions.";
395 ceres::examples::SolveProblem(FLAGS_input.c_str());
396 return 0;
397}