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: keir (at) google.com (Keir Mierle) 30 // sameeragarwal (at) google.com (Sameer Agarwal) 31 // 32 // System level tests for Ceres. The current suite of two tests. The 33 // first test is a small test based on Powell's Function. It is a 34 // scalar problem with 4 variables. The second problem is a bundle 35 // adjustment problem with 16 cameras and two thousand cameras. The 36 // first problem is to test the sanity test the factorization based 37 // solvers. The second problem is used to test the various 38 // combinations of solvers, orderings, preconditioners and 39 // multithreading. 40 41 #include <cmath> 42 #include <cstdio> 43 #include <cstdlib> 44 #include <string> 45 46 #include "ceres/autodiff_cost_function.h" 47 #include "ceres/ordered_groups.h" 48 #include "ceres/problem.h" 49 #include "ceres/rotation.h" 50 #include "ceres/solver.h" 51 #include "ceres/stringprintf.h" 52 #include "ceres/test_util.h" 53 #include "ceres/types.h" 54 #include "gflags/gflags.h" 55 #include "glog/logging.h" 56 #include "gtest/gtest.h" 57 58 namespace ceres { 59 namespace internal { 60 61 const bool kAutomaticOrdering = true; 62 const bool kUserOrdering = false; 63 64 // Struct used for configuring the solver. 65 struct SolverConfig { 66 SolverConfig(LinearSolverType linear_solver_type, 67 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type, 68 bool use_automatic_ordering) 69 : linear_solver_type(linear_solver_type), 70 sparse_linear_algebra_library_type(sparse_linear_algebra_library_type), 71 use_automatic_ordering(use_automatic_ordering), 72 preconditioner_type(IDENTITY), 73 num_threads(1) { 74 } 75 76 SolverConfig(LinearSolverType linear_solver_type, 77 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type, 78 bool use_automatic_ordering, 79 PreconditionerType preconditioner_type) 80 : linear_solver_type(linear_solver_type), 81 sparse_linear_algebra_library_type(sparse_linear_algebra_library_type), 82 use_automatic_ordering(use_automatic_ordering), 83 preconditioner_type(preconditioner_type), 84 num_threads(1) { 85 } 86 87 string ToString() const { 88 return StringPrintf( 89 "(%s, %s, %s, %s, %d)", 90 LinearSolverTypeToString(linear_solver_type), 91 SparseLinearAlgebraLibraryTypeToString(sparse_linear_algebra_library_type), 92 use_automatic_ordering ? "AUTOMATIC" : "USER", 93 PreconditionerTypeToString(preconditioner_type), 94 num_threads); 95 } 96 97 LinearSolverType linear_solver_type; 98 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type; 99 bool use_automatic_ordering; 100 PreconditionerType preconditioner_type; 101 int num_threads; 102 }; 103 104 // Templated function that given a set of solver configurations, 105 // instantiates a new copy of SystemTestProblem for each configuration 106 // and solves it. The solutions are expected to have residuals with 107 // coordinate-wise maximum absolute difference less than or equal to 108 // max_abs_difference. 109 // 110 // The template parameter SystemTestProblem is expected to implement 111 // the following interface. 112 // 113 // class SystemTestProblem { 114 // public: 115 // SystemTestProblem(); 116 // Problem* mutable_problem(); 117 // Solver::Options* mutable_solver_options(); 118 // }; 119 template <typename SystemTestProblem> 120 void RunSolversAndCheckTheyMatch(const vector<SolverConfig>& configurations, 121 const double max_abs_difference) { 122 int num_configurations = configurations.size(); 123 vector<SystemTestProblem*> problems; 124 vector<vector<double> > final_residuals(num_configurations); 125 126 for (int i = 0; i < num_configurations; ++i) { 127 SystemTestProblem* system_test_problem = new SystemTestProblem(); 128 129 const SolverConfig& config = configurations[i]; 130 131 Solver::Options& options = *(system_test_problem->mutable_solver_options()); 132 options.linear_solver_type = config.linear_solver_type; 133 options.sparse_linear_algebra_library_type = 134 config.sparse_linear_algebra_library_type; 135 options.preconditioner_type = config.preconditioner_type; 136 options.num_threads = config.num_threads; 137 options.num_linear_solver_threads = config.num_threads; 138 139 if (config.use_automatic_ordering) { 140 delete options.linear_solver_ordering; 141 options.linear_solver_ordering = NULL; 142 } 143 144 LOG(INFO) << "Running solver configuration: " 145 << config.ToString(); 146 147 Solver::Summary summary; 148 Solve(options, 149 system_test_problem->mutable_problem(), 150 &summary); 151 152 system_test_problem 153 ->mutable_problem() 154 ->Evaluate(Problem::EvaluateOptions(), 155 NULL, 156 &final_residuals[i], 157 NULL, 158 NULL); 159 160 CHECK_NE(summary.termination_type, ceres::NUMERICAL_FAILURE) 161 << "Solver configuration " << i << " failed."; 162 problems.push_back(system_test_problem); 163 164 // Compare the resulting solutions to each other. Arbitrarily take 165 // SPARSE_NORMAL_CHOLESKY as the golden solve. We compare 166 // solutions by comparing their residual vectors. We do not 167 // compare parameter vectors because it is much more brittle and 168 // error prone to do so, since the same problem can have nearly 169 // the same residuals at two completely different positions in 170 // parameter space. 171 if (i > 0) { 172 const vector<double>& reference_residuals = final_residuals[0]; 173 const vector<double>& current_residuals = final_residuals[i]; 174 175 for (int j = 0; j < reference_residuals.size(); ++j) { 176 EXPECT_NEAR(current_residuals[j], 177 reference_residuals[j], 178 max_abs_difference) 179 << "Not close enough residual:" << j 180 << " reference " << reference_residuals[j] 181 << " current " << current_residuals[j]; 182 } 183 } 184 } 185 186 for (int i = 0; i < num_configurations; ++i) { 187 delete problems[i]; 188 } 189 } 190 191 // This class implements the SystemTestProblem interface and provides 192 // access to an implementation of Powell's singular function. 193 // 194 // F = 1/2 (f1^2 + f2^2 + f3^2 + f4^2) 195 // 196 // f1 = x1 + 10*x2; 197 // f2 = sqrt(5) * (x3 - x4) 198 // f3 = (x2 - 2*x3)^2 199 // f4 = sqrt(10) * (x1 - x4)^2 200 // 201 // The starting values are x1 = 3, x2 = -1, x3 = 0, x4 = 1. 202 // The minimum is 0 at (x1, x2, x3, x4) = 0. 203 // 204 // From: Testing Unconstrained Optimization Software by Jorge J. More, Burton S. 205 // Garbow and Kenneth E. Hillstrom in ACM Transactions on Mathematical Software, 206 // Vol 7(1), March 1981. 207 class PowellsFunction { 208 public: 209 PowellsFunction() { 210 x_[0] = 3.0; 211 x_[1] = -1.0; 212 x_[2] = 0.0; 213 x_[3] = 1.0; 214 215 problem_.AddResidualBlock( 216 new AutoDiffCostFunction<F1, 1, 1, 1>(new F1), NULL, &x_[0], &x_[1]); 217 problem_.AddResidualBlock( 218 new AutoDiffCostFunction<F2, 1, 1, 1>(new F2), NULL, &x_[2], &x_[3]); 219 problem_.AddResidualBlock( 220 new AutoDiffCostFunction<F3, 1, 1, 1>(new F3), NULL, &x_[1], &x_[2]); 221 problem_.AddResidualBlock( 222 new AutoDiffCostFunction<F4, 1, 1, 1>(new F4), NULL, &x_[0], &x_[3]); 223 224 options_.max_num_iterations = 10; 225 } 226 227 Problem* mutable_problem() { return &problem_; } 228 Solver::Options* mutable_solver_options() { return &options_; } 229 230 private: 231 // Templated functions used for automatically differentiated cost 232 // functions. 233 class F1 { 234 public: 235 template <typename T> bool operator()(const T* const x1, 236 const T* const x2, 237 T* residual) const { 238 // f1 = x1 + 10 * x2; 239 *residual = *x1 + T(10.0) * *x2; 240 return true; 241 } 242 }; 243 244 class F2 { 245 public: 246 template <typename T> bool operator()(const T* const x3, 247 const T* const x4, 248 T* residual) const { 249 // f2 = sqrt(5) (x3 - x4) 250 *residual = T(sqrt(5.0)) * (*x3 - *x4); 251 return true; 252 } 253 }; 254 255 class F3 { 256 public: 257 template <typename T> bool operator()(const T* const x2, 258 const T* const x4, 259 T* residual) const { 260 // f3 = (x2 - 2 x3)^2 261 residual[0] = (x2[0] - T(2.0) * x4[0]) * (x2[0] - T(2.0) * x4[0]); 262 return true; 263 } 264 }; 265 266 class F4 { 267 public: 268 template <typename T> bool operator()(const T* const x1, 269 const T* const x4, 270 T* residual) const { 271 // f4 = sqrt(10) (x1 - x4)^2 272 residual[0] = T(sqrt(10.0)) * (x1[0] - x4[0]) * (x1[0] - x4[0]); 273 return true; 274 } 275 }; 276 277 double x_[4]; 278 Problem problem_; 279 Solver::Options options_; 280 }; 281 282 TEST(SystemTest, PowellsFunction) { 283 vector<SolverConfig> configs; 284 #define CONFIGURE(linear_solver, sparse_linear_algebra_library_type, ordering) \ 285 configs.push_back(SolverConfig(linear_solver, \ 286 sparse_linear_algebra_library_type, \ 287 ordering)) 288 289 CONFIGURE(DENSE_QR, SUITE_SPARSE, kAutomaticOrdering); 290 CONFIGURE(DENSE_NORMAL_CHOLESKY, SUITE_SPARSE, kAutomaticOrdering); 291 CONFIGURE(DENSE_SCHUR, SUITE_SPARSE, kAutomaticOrdering); 292 293 #ifndef CERES_NO_SUITESPARSE 294 CONFIGURE(SPARSE_NORMAL_CHOLESKY, SUITE_SPARSE, kAutomaticOrdering); 295 #endif // CERES_NO_SUITESPARSE 296 297 #ifndef CERES_NO_CXSPARSE 298 CONFIGURE(SPARSE_NORMAL_CHOLESKY, CX_SPARSE, kAutomaticOrdering); 299 #endif // CERES_NO_CXSPARSE 300 301 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kAutomaticOrdering); 302 303 #undef CONFIGURE 304 305 const double kMaxAbsoluteDifference = 1e-8; 306 RunSolversAndCheckTheyMatch<PowellsFunction>(configs, kMaxAbsoluteDifference); 307 } 308 309 // This class implements the SystemTestProblem interface and provides 310 // access to a bundle adjustment problem. It is based on 311 // examples/bundle_adjustment_example.cc. Currently a small 16 camera 312 // problem is hard coded in the constructor. Going forward we may 313 // extend this to a larger number of problems. 314 class BundleAdjustmentProblem { 315 public: 316 BundleAdjustmentProblem() { 317 const string input_file = TestFileAbsolutePath("problem-16-22106-pre.txt"); 318 ReadData(input_file); 319 BuildProblem(); 320 } 321 322 ~BundleAdjustmentProblem() { 323 delete []point_index_; 324 delete []camera_index_; 325 delete []observations_; 326 delete []parameters_; 327 } 328 329 Problem* mutable_problem() { return &problem_; } 330 Solver::Options* mutable_solver_options() { return &options_; } 331 332 int num_cameras() const { return num_cameras_; } 333 int num_points() const { return num_points_; } 334 int num_observations() const { return num_observations_; } 335 const int* point_index() const { return point_index_; } 336 const int* camera_index() const { return camera_index_; } 337 const double* observations() const { return observations_; } 338 double* mutable_cameras() { return parameters_; } 339 double* mutable_points() { return parameters_ + 9 * num_cameras_; } 340 341 private: 342 void ReadData(const string& filename) { 343 FILE * fptr = fopen(filename.c_str(), "r"); 344 345 if (!fptr) { 346 LOG(FATAL) << "File Error: unable to open file " << filename; 347 }; 348 349 // This will die horribly on invalid files. Them's the breaks. 350 FscanfOrDie(fptr, "%d", &num_cameras_); 351 FscanfOrDie(fptr, "%d", &num_points_); 352 FscanfOrDie(fptr, "%d", &num_observations_); 353 354 VLOG(1) << "Header: " << num_cameras_ 355 << " " << num_points_ 356 << " " << num_observations_; 357 358 point_index_ = new int[num_observations_]; 359 camera_index_ = new int[num_observations_]; 360 observations_ = new double[2 * num_observations_]; 361 362 num_parameters_ = 9 * num_cameras_ + 3 * num_points_; 363 parameters_ = new double[num_parameters_]; 364 365 for (int i = 0; i < num_observations_; ++i) { 366 FscanfOrDie(fptr, "%d", camera_index_ + i); 367 FscanfOrDie(fptr, "%d", point_index_ + i); 368 for (int j = 0; j < 2; ++j) { 369 FscanfOrDie(fptr, "%lf", observations_ + 2*i + j); 370 } 371 } 372 373 for (int i = 0; i < num_parameters_; ++i) { 374 FscanfOrDie(fptr, "%lf", parameters_ + i); 375 } 376 } 377 378 void BuildProblem() { 379 double* points = mutable_points(); 380 double* cameras = mutable_cameras(); 381 382 for (int i = 0; i < num_observations(); ++i) { 383 // Each Residual block takes a point and a camera as input and 384 // outputs a 2 dimensional residual. 385 CostFunction* cost_function = 386 new AutoDiffCostFunction<BundlerResidual, 2, 9, 3>( 387 new BundlerResidual(observations_[2*i + 0], 388 observations_[2*i + 1])); 389 390 // Each observation correponds to a pair of a camera and a point 391 // which are identified by camera_index()[i] and 392 // point_index()[i] respectively. 393 double* camera = cameras + 9 * camera_index_[i]; 394 double* point = points + 3 * point_index()[i]; 395 problem_.AddResidualBlock(cost_function, NULL, camera, point); 396 } 397 398 options_.linear_solver_ordering = new ParameterBlockOrdering; 399 400 // The points come before the cameras. 401 for (int i = 0; i < num_points_; ++i) { 402 options_.linear_solver_ordering->AddElementToGroup(points + 3 * i, 0); 403 } 404 405 for (int i = 0; i < num_cameras_; ++i) { 406 options_.linear_solver_ordering->AddElementToGroup(cameras + 9 * i, 1); 407 } 408 409 options_.max_num_iterations = 25; 410 options_.function_tolerance = 1e-10; 411 options_.gradient_tolerance = 1e-10; 412 options_.parameter_tolerance = 1e-10; 413 } 414 415 template<typename T> 416 void FscanfOrDie(FILE *fptr, const char *format, T *value) { 417 int num_scanned = fscanf(fptr, format, value); 418 if (num_scanned != 1) { 419 LOG(FATAL) << "Invalid UW data file."; 420 } 421 } 422 423 // Templated pinhole camera model. The camera is parameterized 424 // using 9 parameters. 3 for rotation, 3 for translation, 1 for 425 // focal length and 2 for radial distortion. The principal point is 426 // not modeled (i.e. it is assumed be located at the image center). 427 struct BundlerResidual { 428 // (u, v): the position of the observation with respect to the image 429 // center point. 430 BundlerResidual(double u, double v): u(u), v(v) {} 431 432 template <typename T> 433 bool operator()(const T* const camera, 434 const T* const point, 435 T* residuals) const { 436 T p[3]; 437 AngleAxisRotatePoint(camera, point, p); 438 439 // Add the translation vector 440 p[0] += camera[3]; 441 p[1] += camera[4]; 442 p[2] += camera[5]; 443 444 const T& focal = camera[6]; 445 const T& l1 = camera[7]; 446 const T& l2 = camera[8]; 447 448 // Compute the center of distortion. The sign change comes from 449 // the camera model that Noah Snavely's Bundler assumes, whereby 450 // the camera coordinate system has a negative z axis. 451 T xp = - focal * p[0] / p[2]; 452 T yp = - focal * p[1] / p[2]; 453 454 // Apply second and fourth order radial distortion. 455 T r2 = xp*xp + yp*yp; 456 T distortion = T(1.0) + r2 * (l1 + l2 * r2); 457 458 residuals[0] = distortion * xp - T(u); 459 residuals[1] = distortion * yp - T(v); 460 461 return true; 462 } 463 464 double u; 465 double v; 466 }; 467 468 469 Problem problem_; 470 Solver::Options options_; 471 472 int num_cameras_; 473 int num_points_; 474 int num_observations_; 475 int num_parameters_; 476 477 int* point_index_; 478 int* camera_index_; 479 double* observations_; 480 // The parameter vector is laid out as follows 481 // [camera_1, ..., camera_n, point_1, ..., point_m] 482 double* parameters_; 483 }; 484 485 TEST(SystemTest, BundleAdjustmentProblem) { 486 vector<SolverConfig> configs; 487 488 #define CONFIGURE(linear_solver, sparse_linear_algebra_library_type, ordering, preconditioner) \ 489 configs.push_back(SolverConfig(linear_solver, \ 490 sparse_linear_algebra_library_type, \ 491 ordering, \ 492 preconditioner)) 493 494 #ifndef CERES_NO_SUITESPARSE 495 CONFIGURE(SPARSE_NORMAL_CHOLESKY, SUITE_SPARSE, kAutomaticOrdering, IDENTITY); 496 CONFIGURE(SPARSE_NORMAL_CHOLESKY, SUITE_SPARSE, kUserOrdering, IDENTITY); 497 498 CONFIGURE(SPARSE_SCHUR, SUITE_SPARSE, kAutomaticOrdering, IDENTITY); 499 CONFIGURE(SPARSE_SCHUR, SUITE_SPARSE, kUserOrdering, IDENTITY); 500 #endif // CERES_NO_SUITESPARSE 501 502 #ifndef CERES_NO_CXSPARSE 503 CONFIGURE(SPARSE_SCHUR, CX_SPARSE, kAutomaticOrdering, IDENTITY); 504 CONFIGURE(SPARSE_SCHUR, CX_SPARSE, kUserOrdering, IDENTITY); 505 #endif // CERES_NO_CXSPARSE 506 507 CONFIGURE(DENSE_SCHUR, SUITE_SPARSE, kAutomaticOrdering, IDENTITY); 508 CONFIGURE(DENSE_SCHUR, SUITE_SPARSE, kUserOrdering, IDENTITY); 509 510 CONFIGURE(CGNR, SUITE_SPARSE, kAutomaticOrdering, JACOBI); 511 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kUserOrdering, JACOBI); 512 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kUserOrdering, SCHUR_JACOBI); 513 514 #ifndef CERES_NO_SUITESPARSE 515 516 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kUserOrdering, CLUSTER_JACOBI); 517 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kUserOrdering, CLUSTER_TRIDIAGONAL); 518 #endif // CERES_NO_SUITESPARSE 519 520 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kAutomaticOrdering, JACOBI); 521 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kAutomaticOrdering, SCHUR_JACOBI); 522 523 #ifndef CERES_NO_SUITESPARSE 524 525 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kAutomaticOrdering, CLUSTER_JACOBI); 526 CONFIGURE(ITERATIVE_SCHUR, SUITE_SPARSE, kAutomaticOrdering, CLUSTER_TRIDIAGONAL); 527 #endif // CERES_NO_SUITESPARSE 528 529 #undef CONFIGURE 530 531 // Single threaded evaluators and linear solvers. 532 const double kMaxAbsoluteDifference = 1e-4; 533 RunSolversAndCheckTheyMatch<BundleAdjustmentProblem>(configs, 534 kMaxAbsoluteDifference); 535 536 #ifdef CERES_USE_OPENMP 537 // Multithreaded evaluators and linear solvers. 538 for (int i = 0; i < configs.size(); ++i) { 539 configs[i].num_threads = 2; 540 } 541 RunSolversAndCheckTheyMatch<BundleAdjustmentProblem>(configs, 542 kMaxAbsoluteDifference); 543 #endif // CERES_USE_OPENMP 544 } 545 546 } // namespace internal 547 } // namespace ceres 548