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