GoogleGit

blob: 7c177e930464060bc9adb5f42de7e6bcebe26e80 [file] [log] [blame]
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 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: sameeragarwal@google.com (Sameer Agarwal)
  30. //
  31. // The National Institute of Standards and Technology has released a
  32. // set of problems to test non-linear least squares solvers.
  33. //
  34. // More information about the background on these problems and
  35. // suggested evaluation methodology can be found at:
  36. //
  37. // http://www.itl.nist.gov/div898/strd/nls/nls_info.shtml
  38. //
  39. // The problem data themselves can be found at
  40. //
  41. // http://www.itl.nist.gov/div898/strd/nls/nls_main.shtml
  42. //
  43. // The problems are divided into three levels of difficulty, Easy,
  44. // Medium and Hard. For each problem there are two starting guesses,
  45. // the first one far away from the global minimum and the second
  46. // closer to it.
  47. //
  48. // A problem is considered successfully solved, if every components of
  49. // the solution matches the globally optimal solution in at least 4
  50. // digits or more.
  51. //
  52. // This dataset was used for an evaluation of Non-linear least squares
  53. // solvers:
  54. //
  55. // P. F. Mondragon & B. Borchers, A Comparison of Nonlinear Regression
  56. // Codes, Journal of Modern Applied Statistical Methods, 4(1):343-351,
  57. // 2005.
  58. //
  59. // The results from Mondragon & Borchers can be summarized as
  60. // Excel Gnuplot GaussFit HBN MinPack
  61. // Average LRE 2.3 4.3 4.0 6.8 4.4
  62. // Winner 1 5 12 29 12
  63. //
  64. // Where the row Winner counts, the number of problems for which the
  65. // solver had the highest LRE.
  66. // In this file, we implement the same evaluation methodology using
  67. // Ceres. Currently using Levenberg-Marquard with DENSE_QR, we get
  68. //
  69. // Excel Gnuplot GaussFit HBN MinPack Ceres
  70. // Average LRE 2.3 4.3 4.0 6.8 4.4 9.4
  71. // Winner 0 0 5 11 2 41
  72. #include <iostream>
  73. #include <iterator>
  74. #include <fstream>
  75. #include "ceres/ceres.h"
  76. #include "gflags/gflags.h"
  77. #include "glog/logging.h"
  78. #include "Eigen/Core"
  79. DEFINE_string(nist_data_dir, "", "Directory containing the NIST non-linear"
  80. "regression examples");
  81. DEFINE_string(minimizer, "trust_region",
  82. "Minimizer type to use, choices are: line_search & trust_region");
  83. DEFINE_string(trust_region_strategy, "levenberg_marquardt",
  84. "Options are: levenberg_marquardt, dogleg");
  85. DEFINE_string(dogleg, "traditional_dogleg",
  86. "Options are: traditional_dogleg, subspace_dogleg");
  87. DEFINE_string(linear_solver, "dense_qr", "Options are: "
  88. "sparse_cholesky, dense_qr, dense_normal_cholesky and"
  89. "cgnr");
  90. DEFINE_string(preconditioner, "jacobi", "Options are: "
  91. "identity, jacobi");
  92. DEFINE_string(line_search, "wolfe",
  93. "Line search algorithm to use, choices are: armijo and wolfe.");
  94. DEFINE_string(line_search_direction, "lbfgs",
  95. "Line search direction algorithm to use, choices: lbfgs, bfgs");
  96. DEFINE_int32(max_line_search_iterations, 20,
  97. "Maximum number of iterations for each line search.");
  98. DEFINE_int32(max_line_search_restarts, 10,
  99. "Maximum number of restarts of line search direction algorithm.");
  100. DEFINE_string(line_search_interpolation, "cubic",
  101. "Degree of polynomial aproximation in line search, "
  102. "choices are: bisection, quadratic & cubic.");
  103. DEFINE_int32(lbfgs_rank, 20,
  104. "Rank of L-BFGS inverse Hessian approximation in line search.");
  105. DEFINE_bool(approximate_eigenvalue_bfgs_scaling, false,
  106. "Use approximate eigenvalue scaling in (L)BFGS line search.");
  107. DEFINE_double(sufficient_decrease, 1.0e-4,
  108. "Line search Armijo sufficient (function) decrease factor.");
  109. DEFINE_double(sufficient_curvature_decrease, 0.9,
  110. "Line search Wolfe sufficient curvature decrease factor.");
  111. DEFINE_int32(num_iterations, 10000, "Number of iterations");
  112. DEFINE_bool(nonmonotonic_steps, false, "Trust region algorithm can use"
  113. " nonmonotic steps");
  114. DEFINE_double(initial_trust_region_radius, 1e4, "Initial trust region radius");
  115. namespace ceres {
  116. namespace examples {
  117. using Eigen::Dynamic;
  118. using Eigen::RowMajor;
  119. typedef Eigen::Matrix<double, Dynamic, 1> Vector;
  120. typedef Eigen::Matrix<double, Dynamic, Dynamic, RowMajor> Matrix;
  121. using std::atof;
  122. using std::atoi;
  123. using std::cout;
  124. using std::ifstream;
  125. using std::string;
  126. using std::vector;
  127. void SplitStringUsingChar(const string& full,
  128. const char delim,
  129. vector<string>* result) {
  130. std::back_insert_iterator< vector<string> > it(*result);
  131. const char* p = full.data();
  132. const char* end = p + full.size();
  133. while (p != end) {
  134. if (*p == delim) {
  135. ++p;
  136. } else {
  137. const char* start = p;
  138. while (++p != end && *p != delim) {
  139. // Skip to the next occurence of the delimiter.
  140. }
  141. *it++ = string(start, p - start);
  142. }
  143. }
  144. }
  145. bool GetAndSplitLine(ifstream& ifs, vector<string>* pieces) {
  146. pieces->clear();
  147. char buf[256];
  148. ifs.getline(buf, 256);
  149. SplitStringUsingChar(string(buf), ' ', pieces);
  150. return true;
  151. }
  152. void SkipLines(ifstream& ifs, int num_lines) {
  153. char buf[256];
  154. for (int i = 0; i < num_lines; ++i) {
  155. ifs.getline(buf, 256);
  156. }
  157. }
  158. class NISTProblem {
  159. public:
  160. explicit NISTProblem(const string& filename) {
  161. ifstream ifs(filename.c_str(), ifstream::in);
  162. vector<string> pieces;
  163. SkipLines(ifs, 24);
  164. GetAndSplitLine(ifs, &pieces);
  165. const int kNumResponses = atoi(pieces[1].c_str());
  166. GetAndSplitLine(ifs, &pieces);
  167. const int kNumPredictors = atoi(pieces[0].c_str());
  168. GetAndSplitLine(ifs, &pieces);
  169. const int kNumObservations = atoi(pieces[0].c_str());
  170. SkipLines(ifs, 4);
  171. GetAndSplitLine(ifs, &pieces);
  172. const int kNumParameters = atoi(pieces[0].c_str());
  173. SkipLines(ifs, 8);
  174. // Get the first line of initial and final parameter values to
  175. // determine the number of tries.
  176. GetAndSplitLine(ifs, &pieces);
  177. const int kNumTries = pieces.size() - 4;
  178. predictor_.resize(kNumObservations, kNumPredictors);
  179. response_.resize(kNumObservations, kNumResponses);
  180. initial_parameters_.resize(kNumTries, kNumParameters);
  181. final_parameters_.resize(1, kNumParameters);
  182. // Parse the line for parameter b1.
  183. int parameter_id = 0;
  184. for (int i = 0; i < kNumTries; ++i) {
  185. initial_parameters_(i, parameter_id) = atof(pieces[i + 2].c_str());
  186. }
  187. final_parameters_(0, parameter_id) = atof(pieces[2 + kNumTries].c_str());
  188. // Parse the remaining parameter lines.
  189. for (int parameter_id = 1; parameter_id < kNumParameters; ++parameter_id) {
  190. GetAndSplitLine(ifs, &pieces);
  191. // b2, b3, ....
  192. for (int i = 0; i < kNumTries; ++i) {
  193. initial_parameters_(i, parameter_id) = atof(pieces[i + 2].c_str());
  194. }
  195. final_parameters_(0, parameter_id) = atof(pieces[2 + kNumTries].c_str());
  196. }
  197. // Certfied cost
  198. SkipLines(ifs, 1);
  199. GetAndSplitLine(ifs, &pieces);
  200. certified_cost_ = atof(pieces[4].c_str()) / 2.0;
  201. // Read the observations.
  202. SkipLines(ifs, 18 - kNumParameters);
  203. for (int i = 0; i < kNumObservations; ++i) {
  204. GetAndSplitLine(ifs, &pieces);
  205. // Response.
  206. for (int j = 0; j < kNumResponses; ++j) {
  207. response_(i, j) = atof(pieces[j].c_str());
  208. }
  209. // Predictor variables.
  210. for (int j = 0; j < kNumPredictors; ++j) {
  211. predictor_(i, j) = atof(pieces[j + kNumResponses].c_str());
  212. }
  213. }
  214. }
  215. Matrix initial_parameters(int start) const { return initial_parameters_.row(start); }
  216. Matrix final_parameters() const { return final_parameters_; }
  217. Matrix predictor() const { return predictor_; }
  218. Matrix response() const { return response_; }
  219. int predictor_size() const { return predictor_.cols(); }
  220. int num_observations() const { return predictor_.rows(); }
  221. int response_size() const { return response_.cols(); }
  222. int num_parameters() const { return initial_parameters_.cols(); }
  223. int num_starts() const { return initial_parameters_.rows(); }
  224. double certified_cost() const { return certified_cost_; }
  225. private:
  226. Matrix predictor_;
  227. Matrix response_;
  228. Matrix initial_parameters_;
  229. Matrix final_parameters_;
  230. double certified_cost_;
  231. };
  232. #define NIST_BEGIN(CostFunctionName) \
  233. struct CostFunctionName { \
  234. CostFunctionName(const double* const x, \
  235. const double* const y) \
  236. : x_(*x), y_(*y) {} \
  237. double x_; \
  238. double y_; \
  239. template <typename T> \
  240. bool operator()(const T* const b, T* residual) const { \
  241. const T y(y_); \
  242. const T x(x_); \
  243. residual[0] = y - (
  244. #define NIST_END ); return true; }};
  245. // y = b1 * (b2+x)**(-1/b3) + e
  246. NIST_BEGIN(Bennet5)
  247. b[0] * pow(b[1] + x, T(-1.0) / b[2])
  248. NIST_END
  249. // y = b1*(1-exp[-b2*x]) + e
  250. NIST_BEGIN(BoxBOD)
  251. b[0] * (T(1.0) - exp(-b[1] * x))
  252. NIST_END
  253. // y = exp[-b1*x]/(b2+b3*x) + e
  254. NIST_BEGIN(Chwirut)
  255. exp(-b[0] * x) / (b[1] + b[2] * x)
  256. NIST_END
  257. // y = b1*x**b2 + e
  258. NIST_BEGIN(DanWood)
  259. b[0] * pow(x, b[1])
  260. NIST_END
  261. // y = b1*exp( -b2*x ) + b3*exp( -(x-b4)**2 / b5**2 )
  262. // + b6*exp( -(x-b7)**2 / b8**2 ) + e
  263. NIST_BEGIN(Gauss)
  264. b[0] * exp(-b[1] * x) +
  265. b[2] * exp(-pow((x - b[3])/b[4], 2)) +
  266. b[5] * exp(-pow((x - b[6])/b[7],2))
  267. NIST_END
  268. // y = b1*exp(-b2*x) + b3*exp(-b4*x) + b5*exp(-b6*x) + e
  269. NIST_BEGIN(Lanczos)
  270. b[0] * exp(-b[1] * x) + b[2] * exp(-b[3] * x) + b[4] * exp(-b[5] * x)
  271. NIST_END
  272. // y = (b1+b2*x+b3*x**2+b4*x**3) /
  273. // (1+b5*x+b6*x**2+b7*x**3) + e
  274. NIST_BEGIN(Hahn1)
  275. (b[0] + b[1] * x + b[2] * x * x + b[3] * x * x * x) /
  276. (T(1.0) + b[4] * x + b[5] * x * x + b[6] * x * x * x)
  277. NIST_END
  278. // y = (b1 + b2*x + b3*x**2) /
  279. // (1 + b4*x + b5*x**2) + e
  280. NIST_BEGIN(Kirby2)
  281. (b[0] + b[1] * x + b[2] * x * x) /
  282. (T(1.0) + b[3] * x + b[4] * x * x)
  283. NIST_END
  284. // y = b1*(x**2+x*b2) / (x**2+x*b3+b4) + e
  285. NIST_BEGIN(MGH09)
  286. b[0] * (x * x + x * b[1]) / (x * x + x * b[2] + b[3])
  287. NIST_END
  288. // y = b1 * exp[b2/(x+b3)] + e
  289. NIST_BEGIN(MGH10)
  290. b[0] * exp(b[1] / (x + b[2]))
  291. NIST_END
  292. // y = b1 + b2*exp[-x*b4] + b3*exp[-x*b5]
  293. NIST_BEGIN(MGH17)
  294. b[0] + b[1] * exp(-x * b[3]) + b[2] * exp(-x * b[4])
  295. NIST_END
  296. // y = b1*(1-exp[-b2*x]) + e
  297. NIST_BEGIN(Misra1a)
  298. b[0] * (T(1.0) - exp(-b[1] * x))
  299. NIST_END
  300. // y = b1 * (1-(1+b2*x/2)**(-2)) + e
  301. NIST_BEGIN(Misra1b)
  302. b[0] * (T(1.0) - T(1.0)/ ((T(1.0) + b[1] * x / 2.0) * (T(1.0) + b[1] * x / 2.0)))
  303. NIST_END
  304. // y = b1 * (1-(1+2*b2*x)**(-.5)) + e
  305. NIST_BEGIN(Misra1c)
  306. b[0] * (T(1.0) - pow(T(1.0) + T(2.0) * b[1] * x, -0.5))
  307. NIST_END
  308. // y = b1*b2*x*((1+b2*x)**(-1)) + e
  309. NIST_BEGIN(Misra1d)
  310. b[0] * b[1] * x / (T(1.0) + b[1] * x)
  311. NIST_END
  312. const double kPi = 3.141592653589793238462643383279;
  313. // pi = 3.141592653589793238462643383279E0
  314. // y = b1 - b2*x - arctan[b3/(x-b4)]/pi + e
  315. NIST_BEGIN(Roszman1)
  316. b[0] - b[1] * x - atan2(b[2], (x - b[3]))/T(kPi)
  317. NIST_END
  318. // y = b1 / (1+exp[b2-b3*x]) + e
  319. NIST_BEGIN(Rat42)
  320. b[0] / (T(1.0) + exp(b[1] - b[2] * x))
  321. NIST_END
  322. // y = b1 / ((1+exp[b2-b3*x])**(1/b4)) + e
  323. NIST_BEGIN(Rat43)
  324. b[0] / pow(T(1.0) + exp(b[1] - b[2] * x), T(1.0) / b[3])
  325. NIST_END
  326. // y = (b1 + b2*x + b3*x**2 + b4*x**3) /
  327. // (1 + b5*x + b6*x**2 + b7*x**3) + e
  328. NIST_BEGIN(Thurber)
  329. (b[0] + b[1] * x + b[2] * x * x + b[3] * x * x * x) /
  330. (T(1.0) + b[4] * x + b[5] * x * x + b[6] * x * x * x)
  331. NIST_END
  332. // y = b1 + b2*cos( 2*pi*x/12 ) + b3*sin( 2*pi*x/12 )
  333. // + b5*cos( 2*pi*x/b4 ) + b6*sin( 2*pi*x/b4 )
  334. // + b8*cos( 2*pi*x/b7 ) + b9*sin( 2*pi*x/b7 ) + e
  335. NIST_BEGIN(ENSO)
  336. b[0] + b[1] * cos(T(2.0 * kPi) * x / T(12.0)) +
  337. b[2] * sin(T(2.0 * kPi) * x / T(12.0)) +
  338. b[4] * cos(T(2.0 * kPi) * x / b[3]) +
  339. b[5] * sin(T(2.0 * kPi) * x / b[3]) +
  340. b[7] * cos(T(2.0 * kPi) * x / b[6]) +
  341. b[8] * sin(T(2.0 * kPi) * x / b[6])
  342. NIST_END
  343. // y = (b1/b2) * exp[-0.5*((x-b3)/b2)**2] + e
  344. NIST_BEGIN(Eckerle4)
  345. b[0] / b[1] * exp(T(-0.5) * pow((x - b[2])/b[1], 2))
  346. NIST_END
  347. struct Nelson {
  348. public:
  349. Nelson(const double* const x, const double* const y)
  350. : x1_(x[0]), x2_(x[1]), y_(y[0]) {}
  351. template <typename T>
  352. bool operator()(const T* const b, T* residual) const {
  353. // log[y] = b1 - b2*x1 * exp[-b3*x2] + e
  354. residual[0] = T(log(y_)) - (b[0] - b[1] * T(x1_) * exp(-b[2] * T(x2_)));
  355. return true;
  356. }
  357. private:
  358. double x1_;
  359. double x2_;
  360. double y_;
  361. };
  362. template <typename Model, int num_residuals, int num_parameters>
  363. int RegressionDriver(const string& filename,
  364. const ceres::Solver::Options& options) {
  365. NISTProblem nist_problem(FLAGS_nist_data_dir + filename);
  366. CHECK_EQ(num_residuals, nist_problem.response_size());
  367. CHECK_EQ(num_parameters, nist_problem.num_parameters());
  368. Matrix predictor = nist_problem.predictor();
  369. Matrix response = nist_problem.response();
  370. Matrix final_parameters = nist_problem.final_parameters();
  371. printf("%s\n", filename.c_str());
  372. // Each NIST problem comes with multiple starting points, so we
  373. // construct the problem from scratch for each case and solve it.
  374. int num_success = 0;
  375. for (int start = 0; start < nist_problem.num_starts(); ++start) {
  376. Matrix initial_parameters = nist_problem.initial_parameters(start);
  377. ceres::Problem problem;
  378. for (int i = 0; i < nist_problem.num_observations(); ++i) {
  379. problem.AddResidualBlock(
  380. new ceres::AutoDiffCostFunction<Model, num_residuals, num_parameters>(
  381. new Model(predictor.data() + nist_problem.predictor_size() * i,
  382. response.data() + nist_problem.response_size() * i)),
  383. NULL,
  384. initial_parameters.data());
  385. }
  386. ceres::Solver::Summary summary;
  387. Solve(options, &problem, &summary);
  388. // Compute the LRE by comparing each component of the solution
  389. // with the ground truth, and taking the minimum.
  390. Matrix final_parameters = nist_problem.final_parameters();
  391. const double kMaxNumSignificantDigits = 11;
  392. double log_relative_error = kMaxNumSignificantDigits + 1;
  393. for (int i = 0; i < num_parameters; ++i) {
  394. const double tmp_lre =
  395. -std::log10(std::fabs(final_parameters(i) - initial_parameters(i)) /
  396. std::fabs(final_parameters(i)));
  397. // The maximum LRE is capped at 11 - the precision at which the
  398. // ground truth is known.
  399. //
  400. // The minimum LRE is capped at 0 - no digits match between the
  401. // computed solution and the ground truth.
  402. log_relative_error =
  403. std::min(log_relative_error,
  404. std::max(0.0, std::min(kMaxNumSignificantDigits, tmp_lre)));
  405. }
  406. const int kMinNumMatchingDigits = 4;
  407. if (log_relative_error >= kMinNumMatchingDigits) {
  408. ++num_success;
  409. }
  410. printf("start: %d status: %s lre: %4.1f initial cost: %e final cost:%e "
  411. "certified cost: %e total iterations: %d\n",
  412. start + 1,
  413. log_relative_error < kMinNumMatchingDigits ? "FAILURE" : "SUCCESS",
  414. log_relative_error,
  415. summary.initial_cost,
  416. summary.final_cost,
  417. nist_problem.certified_cost(),
  418. (summary.num_successful_steps + summary.num_unsuccessful_steps));
  419. }
  420. return num_success;
  421. }
  422. void SetMinimizerOptions(ceres::Solver::Options* options) {
  423. CHECK(ceres::StringToMinimizerType(FLAGS_minimizer,
  424. &options->minimizer_type));
  425. CHECK(ceres::StringToLinearSolverType(FLAGS_linear_solver,
  426. &options->linear_solver_type));
  427. CHECK(ceres::StringToPreconditionerType(FLAGS_preconditioner,
  428. &options->preconditioner_type));
  429. CHECK(ceres::StringToTrustRegionStrategyType(
  430. FLAGS_trust_region_strategy,
  431. &options->trust_region_strategy_type));
  432. CHECK(ceres::StringToDoglegType(FLAGS_dogleg, &options->dogleg_type));
  433. CHECK(ceres::StringToLineSearchDirectionType(
  434. FLAGS_line_search_direction,
  435. &options->line_search_direction_type));
  436. CHECK(ceres::StringToLineSearchType(FLAGS_line_search,
  437. &options->line_search_type));
  438. CHECK(ceres::StringToLineSearchInterpolationType(
  439. FLAGS_line_search_interpolation,
  440. &options->line_search_interpolation_type));
  441. options->max_num_iterations = FLAGS_num_iterations;
  442. options->use_nonmonotonic_steps = FLAGS_nonmonotonic_steps;
  443. options->initial_trust_region_radius = FLAGS_initial_trust_region_radius;
  444. options->max_lbfgs_rank = FLAGS_lbfgs_rank;
  445. options->line_search_sufficient_function_decrease = FLAGS_sufficient_decrease;
  446. options->line_search_sufficient_curvature_decrease =
  447. FLAGS_sufficient_curvature_decrease;
  448. options->max_num_line_search_step_size_iterations =
  449. FLAGS_max_line_search_iterations;
  450. options->max_num_line_search_direction_restarts =
  451. FLAGS_max_line_search_restarts;
  452. options->use_approximate_eigenvalue_bfgs_scaling =
  453. FLAGS_approximate_eigenvalue_bfgs_scaling;
  454. options->function_tolerance = 1e-18;
  455. options->gradient_tolerance = 1e-18;
  456. options->parameter_tolerance = 1e-18;
  457. }
  458. void SolveNISTProblems() {
  459. if (FLAGS_nist_data_dir.empty()) {
  460. LOG(FATAL) << "Must specify the directory containing the NIST problems";
  461. }
  462. ceres::Solver::Options options;
  463. SetMinimizerOptions(&options);
  464. cout << "Lower Difficulty\n";
  465. int easy_success = 0;
  466. easy_success += RegressionDriver<Misra1a, 1, 2>("Misra1a.dat", options);
  467. easy_success += RegressionDriver<Chwirut, 1, 3>("Chwirut1.dat", options);
  468. easy_success += RegressionDriver<Chwirut, 1, 3>("Chwirut2.dat", options);
  469. easy_success += RegressionDriver<Lanczos, 1, 6>("Lanczos3.dat", options);
  470. easy_success += RegressionDriver<Gauss, 1, 8>("Gauss1.dat", options);
  471. easy_success += RegressionDriver<Gauss, 1, 8>("Gauss2.dat", options);
  472. easy_success += RegressionDriver<DanWood, 1, 2>("DanWood.dat", options);
  473. easy_success += RegressionDriver<Misra1b, 1, 2>("Misra1b.dat", options);
  474. cout << "\nMedium Difficulty\n";
  475. int medium_success = 0;
  476. medium_success += RegressionDriver<Kirby2, 1, 5>("Kirby2.dat", options);
  477. medium_success += RegressionDriver<Hahn1, 1, 7>("Hahn1.dat", options);
  478. medium_success += RegressionDriver<Nelson, 1, 3>("Nelson.dat", options);
  479. medium_success += RegressionDriver<MGH17, 1, 5>("MGH17.dat", options);
  480. medium_success += RegressionDriver<Lanczos, 1, 6>("Lanczos1.dat", options);
  481. medium_success += RegressionDriver<Lanczos, 1, 6>("Lanczos2.dat", options);
  482. medium_success += RegressionDriver<Gauss, 1, 8>("Gauss3.dat", options);
  483. medium_success += RegressionDriver<Misra1c, 1, 2>("Misra1c.dat", options);
  484. medium_success += RegressionDriver<Misra1d, 1, 2>("Misra1d.dat", options);
  485. medium_success += RegressionDriver<Roszman1, 1, 4>("Roszman1.dat", options);
  486. medium_success += RegressionDriver<ENSO, 1, 9>("ENSO.dat", options);
  487. cout << "\nHigher Difficulty\n";
  488. int hard_success = 0;
  489. hard_success += RegressionDriver<MGH09, 1, 4>("MGH09.dat", options);
  490. hard_success += RegressionDriver<Thurber, 1, 7>("Thurber.dat", options);
  491. hard_success += RegressionDriver<BoxBOD, 1, 2>("BoxBOD.dat", options);
  492. hard_success += RegressionDriver<Rat42, 1, 3>("Rat42.dat", options);
  493. hard_success += RegressionDriver<MGH10, 1, 3>("MGH10.dat", options);
  494. hard_success += RegressionDriver<Eckerle4, 1, 3>("Eckerle4.dat", options);
  495. hard_success += RegressionDriver<Rat43, 1, 4>("Rat43.dat", options);
  496. hard_success += RegressionDriver<Bennet5, 1, 3>("Bennett5.dat", options);
  497. cout << "\n";
  498. cout << "Easy : " << easy_success << "/16\n";
  499. cout << "Medium : " << medium_success << "/22\n";
  500. cout << "Hard : " << hard_success << "/16\n";
  501. cout << "Total : " << easy_success + medium_success + hard_success << "/54\n";
  502. }
  503. } // namespace examples
  504. } // namespace ceres
  505. int main(int argc, char** argv) {
  506. CERES_GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
  507. google::InitGoogleLogging(argv[0]);
  508. ceres::examples::SolveNISTProblems();
  509. return 0;
  510. };