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