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