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