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