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