nist.cc 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2023 Google Inc. All rights reserved.
  3. // http://ceres-solver.org/
  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-Marquardt 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 <cstdlib>
  73. #include <fstream>
  74. #include <iostream>
  75. #include <iterator>
  76. #include <string>
  77. #include <vector>
  78. #include "Eigen/Core"
  79. #include "ceres/ceres.h"
  80. #include "ceres/tiny_solver.h"
  81. #include "ceres/tiny_solver_cost_function_adapter.h"
  82. #include "gflags/gflags.h"
  83. #include "glog/logging.h"
  84. DEFINE_bool(use_tiny_solver, false, "Use TinySolver instead of Ceres::Solver");
  85. DEFINE_string(nist_data_dir,
  86. "",
  87. "Directory containing the NIST non-linear regression examples");
  88. DEFINE_string(minimizer,
  89. "trust_region",
  90. "Minimizer type to use, choices are: line_search & trust_region");
  91. DEFINE_string(trust_region_strategy,
  92. "levenberg_marquardt",
  93. "Options are: levenberg_marquardt, dogleg");
  94. DEFINE_string(dogleg,
  95. "traditional_dogleg",
  96. "Options are: traditional_dogleg, subspace_dogleg");
  97. DEFINE_string(linear_solver,
  98. "dense_qr",
  99. "Options are: sparse_cholesky, dense_qr, dense_normal_cholesky "
  100. "and cgnr");
  101. DEFINE_string(dense_linear_algebra_library,
  102. "eigen",
  103. "Options are: eigen, lapack, and cuda.");
  104. DEFINE_string(preconditioner, "jacobi", "Options are: identity, jacobi");
  105. DEFINE_string(line_search,
  106. "wolfe",
  107. "Line search algorithm to use, choices are: armijo and wolfe.");
  108. DEFINE_string(line_search_direction,
  109. "lbfgs",
  110. "Line search direction algorithm to use, choices: lbfgs, bfgs");
  111. DEFINE_int32(max_line_search_iterations,
  112. 20,
  113. "Maximum number of iterations for each line search.");
  114. DEFINE_int32(max_line_search_restarts,
  115. 10,
  116. "Maximum number of restarts of line search direction algorithm.");
  117. DEFINE_string(line_search_interpolation,
  118. "cubic",
  119. "Degree of polynomial approximation in line search, choices are: "
  120. "bisection, quadratic & cubic.");
  121. DEFINE_int32(lbfgs_rank,
  122. 20,
  123. "Rank of L-BFGS inverse Hessian approximation in line search.");
  124. DEFINE_bool(approximate_eigenvalue_bfgs_scaling,
  125. false,
  126. "Use approximate eigenvalue scaling in (L)BFGS line search.");
  127. DEFINE_double(sufficient_decrease,
  128. 1.0e-4,
  129. "Line search Armijo sufficient (function) decrease factor.");
  130. DEFINE_double(sufficient_curvature_decrease,
  131. 0.9,
  132. "Line search Wolfe sufficient curvature decrease factor.");
  133. DEFINE_int32(num_iterations, 10000, "Number of iterations");
  134. DEFINE_bool(nonmonotonic_steps,
  135. false,
  136. "Trust region algorithm can use nonmonotic steps");
  137. DEFINE_double(initial_trust_region_radius, 1e4, "Initial trust region radius");
  138. DEFINE_bool(use_numeric_diff,
  139. false,
  140. "Use numeric differentiation instead of automatic "
  141. "differentiation.");
  142. DEFINE_string(numeric_diff_method,
  143. "ridders",
  144. "When using numeric differentiation, selects algorithm. Options "
  145. "are: central, forward, ridders.");
  146. DEFINE_double(ridders_step_size,
  147. 1e-9,
  148. "Initial step size for Ridders numeric differentiation.");
  149. DEFINE_int32(ridders_extrapolations,
  150. 3,
  151. "Maximal number of Ridders extrapolations.");
  152. namespace ceres::examples {
  153. namespace {
  154. using Eigen::Dynamic;
  155. using Eigen::RowMajor;
  156. using Vector = Eigen::Matrix<double, Dynamic, 1>;
  157. using Matrix = Eigen::Matrix<double, Dynamic, Dynamic, RowMajor>;
  158. void SplitStringUsingChar(const std::string& full,
  159. const char delim,
  160. std::vector<std::string>* result) {
  161. std::back_insert_iterator<std::vector<std::string>> it(*result);
  162. const char* p = full.data();
  163. const char* end = p + full.size();
  164. while (p != end) {
  165. if (*p == delim) {
  166. ++p;
  167. } else {
  168. const char* start = p;
  169. while (++p != end && *p != delim) {
  170. // Skip to the next occurrence of the delimiter.
  171. }
  172. *it++ = std::string(start, p - start);
  173. }
  174. }
  175. }
  176. bool GetAndSplitLine(std::ifstream& ifs, std::vector<std::string>* pieces) {
  177. pieces->clear();
  178. char buf[256];
  179. ifs.getline(buf, 256);
  180. SplitStringUsingChar(std::string(buf), ' ', pieces);
  181. return true;
  182. }
  183. void SkipLines(std::ifstream& ifs, int num_lines) {
  184. char buf[256];
  185. for (int i = 0; i < num_lines; ++i) {
  186. ifs.getline(buf, 256);
  187. }
  188. }
  189. class NISTProblem {
  190. public:
  191. explicit NISTProblem(const std::string& filename) {
  192. std::ifstream ifs(filename.c_str(), std::ifstream::in);
  193. CHECK(ifs) << "Unable to open : " << filename;
  194. std::vector<std::string> pieces;
  195. SkipLines(ifs, 24);
  196. GetAndSplitLine(ifs, &pieces);
  197. const int kNumResponses = std::atoi(pieces[1].c_str());
  198. GetAndSplitLine(ifs, &pieces);
  199. const int kNumPredictors = std::atoi(pieces[0].c_str());
  200. GetAndSplitLine(ifs, &pieces);
  201. const int kNumObservations = std::atoi(pieces[0].c_str());
  202. SkipLines(ifs, 4);
  203. GetAndSplitLine(ifs, &pieces);
  204. const int kNumParameters = std::atoi(pieces[0].c_str());
  205. SkipLines(ifs, 8);
  206. // Get the first line of initial and final parameter values to
  207. // determine the number of tries.
  208. GetAndSplitLine(ifs, &pieces);
  209. const int kNumTries = pieces.size() - 4;
  210. predictor_.resize(kNumObservations, kNumPredictors);
  211. response_.resize(kNumObservations, kNumResponses);
  212. initial_parameters_.resize(kNumTries, kNumParameters);
  213. final_parameters_.resize(1, kNumParameters);
  214. // Parse the line for parameter b1.
  215. int parameter_id = 0;
  216. for (int i = 0; i < kNumTries; ++i) {
  217. initial_parameters_(i, parameter_id) = std::atof(pieces[i + 2].c_str());
  218. }
  219. final_parameters_(0, parameter_id) =
  220. std::atof(pieces[2 + kNumTries].c_str());
  221. // Parse the remaining parameter lines.
  222. for (int parameter_id = 1; parameter_id < kNumParameters; ++parameter_id) {
  223. GetAndSplitLine(ifs, &pieces);
  224. // b2, b3, ....
  225. for (int i = 0; i < kNumTries; ++i) {
  226. initial_parameters_(i, parameter_id) = std::atof(pieces[i + 2].c_str());
  227. }
  228. final_parameters_(0, parameter_id) =
  229. std::atof(pieces[2 + kNumTries].c_str());
  230. }
  231. // Certified cost
  232. SkipLines(ifs, 1);
  233. GetAndSplitLine(ifs, &pieces);
  234. certified_cost_ = std::atof(pieces[4].c_str()) / 2.0;
  235. // Read the observations.
  236. SkipLines(ifs, 18 - kNumParameters);
  237. for (int i = 0; i < kNumObservations; ++i) {
  238. GetAndSplitLine(ifs, &pieces);
  239. // Response.
  240. for (int j = 0; j < kNumResponses; ++j) {
  241. response_(i, j) = std::atof(pieces[j].c_str());
  242. }
  243. // Predictor variables.
  244. for (int j = 0; j < kNumPredictors; ++j) {
  245. predictor_(i, j) = std::atof(pieces[j + kNumResponses].c_str());
  246. }
  247. }
  248. }
  249. Matrix initial_parameters(int start) const {
  250. return initial_parameters_.row(start);
  251. } // NOLINT
  252. Matrix final_parameters() const { return final_parameters_; }
  253. Matrix predictor() const { return predictor_; }
  254. Matrix response() const { return response_; }
  255. int predictor_size() const { return predictor_.cols(); }
  256. int num_observations() const { return predictor_.rows(); }
  257. int response_size() const { return response_.cols(); }
  258. int num_parameters() const { return initial_parameters_.cols(); }
  259. int num_starts() const { return initial_parameters_.rows(); }
  260. double certified_cost() const { return certified_cost_; }
  261. private:
  262. Matrix predictor_;
  263. Matrix response_;
  264. Matrix initial_parameters_;
  265. Matrix final_parameters_;
  266. double certified_cost_;
  267. };
  268. #define NIST_BEGIN(CostFunctionName) \
  269. struct CostFunctionName { \
  270. CostFunctionName(const double* const x, \
  271. const double* const y, \
  272. const int n) \
  273. : x_(x), y_(y), n_(n) {} \
  274. const double* x_; \
  275. const double* y_; \
  276. const int n_; \
  277. template <typename T> \
  278. bool operator()(const T* const b, T* residual) const { \
  279. for (int i = 0; i < n_; ++i) { \
  280. const T x(x_[i]); \
  281. residual[i] = y_[i] - (
  282. // clang-format off
  283. #define NIST_END ); } return true; }};
  284. // y = b1 * (b2+x)**(-1/b3) + e
  285. NIST_BEGIN(Bennet5)
  286. b[0] * pow(b[1] + x, -1.0 / b[2])
  287. NIST_END
  288. // y = b1*(1-exp[-b2*x]) + e
  289. NIST_BEGIN(BoxBOD)
  290. b[0] * (1.0 - exp(-b[1] * x))
  291. NIST_END
  292. // y = exp[-b1*x]/(b2+b3*x) + e
  293. NIST_BEGIN(Chwirut)
  294. exp(-b[0] * x) / (b[1] + b[2] * x)
  295. NIST_END
  296. // y = b1*x**b2 + e
  297. NIST_BEGIN(DanWood)
  298. b[0] * pow(x, b[1])
  299. NIST_END
  300. // y = b1*exp( -b2*x ) + b3*exp( -(x-b4)**2 / b5**2 )
  301. // + b6*exp( -(x-b7)**2 / b8**2 ) + e
  302. NIST_BEGIN(Gauss)
  303. b[0] * exp(-b[1] * x) +
  304. b[2] * exp(-pow((x - b[3])/b[4], 2)) +
  305. b[5] * exp(-pow((x - b[6])/b[7], 2))
  306. NIST_END
  307. // y = b1*exp(-b2*x) + b3*exp(-b4*x) + b5*exp(-b6*x) + e
  308. NIST_BEGIN(Lanczos)
  309. b[0] * exp(-b[1] * x) + b[2] * exp(-b[3] * x) + b[4] * exp(-b[5] * x)
  310. NIST_END
  311. // y = (b1+b2*x+b3*x**2+b4*x**3) /
  312. // (1+b5*x+b6*x**2+b7*x**3) + e
  313. NIST_BEGIN(Hahn1)
  314. (b[0] + b[1] * x + b[2] * x * x + b[3] * x * x * x) /
  315. (1.0 + b[4] * x + b[5] * x * x + b[6] * x * x * x)
  316. NIST_END
  317. // y = (b1 + b2*x + b3*x**2) /
  318. // (1 + b4*x + b5*x**2) + e
  319. NIST_BEGIN(Kirby2)
  320. (b[0] + b[1] * x + b[2] * x * x) /
  321. (1.0 + b[3] * x + b[4] * x * x)
  322. NIST_END
  323. // y = b1*(x**2+x*b2) / (x**2+x*b3+b4) + e
  324. NIST_BEGIN(MGH09)
  325. b[0] * (x * x + x * b[1]) / (x * x + x * b[2] + b[3])
  326. NIST_END
  327. // y = b1 * exp[b2/(x+b3)] + e
  328. NIST_BEGIN(MGH10)
  329. b[0] * exp(b[1] / (x + b[2]))
  330. NIST_END
  331. // y = b1 + b2*exp[-x*b4] + b3*exp[-x*b5]
  332. NIST_BEGIN(MGH17)
  333. b[0] + b[1] * exp(-x * b[3]) + b[2] * exp(-x * b[4])
  334. NIST_END
  335. // y = b1*(1-exp[-b2*x]) + e
  336. NIST_BEGIN(Misra1a)
  337. b[0] * (1.0 - exp(-b[1] * x))
  338. NIST_END
  339. // y = b1 * (1-(1+b2*x/2)**(-2)) + e
  340. NIST_BEGIN(Misra1b)
  341. b[0] * (1.0 - 1.0/ ((1.0 + b[1] * x / 2.0) * (1.0 + b[1] * x / 2.0))) // NOLINT
  342. NIST_END
  343. // y = b1 * (1-(1+2*b2*x)**(-.5)) + e
  344. NIST_BEGIN(Misra1c)
  345. b[0] * (1.0 - pow(1.0 + 2.0 * b[1] * x, -0.5))
  346. NIST_END
  347. // y = b1*b2*x*((1+b2*x)**(-1)) + e
  348. NIST_BEGIN(Misra1d)
  349. b[0] * b[1] * x / (1.0 + b[1] * x)
  350. NIST_END
  351. const double kPi = 3.141592653589793238462643383279;
  352. // pi = 3.141592653589793238462643383279E0
  353. // y = b1 - b2*x - arctan[b3/(x-b4)]/pi + e
  354. NIST_BEGIN(Roszman1)
  355. b[0] - b[1] * x - atan2(b[2], (x - b[3])) / kPi
  356. NIST_END
  357. // y = b1 / (1+exp[b2-b3*x]) + e
  358. NIST_BEGIN(Rat42)
  359. b[0] / (1.0 + exp(b[1] - b[2] * x))
  360. NIST_END
  361. // y = b1 / ((1+exp[b2-b3*x])**(1/b4)) + e
  362. NIST_BEGIN(Rat43)
  363. b[0] / pow(1.0 + exp(b[1] - b[2] * x), 1.0 / b[3])
  364. NIST_END
  365. // y = (b1 + b2*x + b3*x**2 + b4*x**3) /
  366. // (1 + b5*x + b6*x**2 + b7*x**3) + e
  367. NIST_BEGIN(Thurber)
  368. (b[0] + b[1] * x + b[2] * x * x + b[3] * x * x * x) /
  369. (1.0 + b[4] * x + b[5] * x * x + b[6] * x * x * x)
  370. NIST_END
  371. // y = b1 + b2*cos( 2*pi*x/12 ) + b3*sin( 2*pi*x/12 )
  372. // + b5*cos( 2*pi*x/b4 ) + b6*sin( 2*pi*x/b4 )
  373. // + b8*cos( 2*pi*x/b7 ) + b9*sin( 2*pi*x/b7 ) + e
  374. NIST_BEGIN(ENSO)
  375. b[0] + b[1] * cos(2.0 * kPi * x / 12.0) +
  376. b[2] * sin(2.0 * kPi * x / 12.0) +
  377. b[4] * cos(2.0 * kPi * x / b[3]) +
  378. b[5] * sin(2.0 * kPi * x / b[3]) +
  379. b[7] * cos(2.0 * kPi * x / b[6]) +
  380. b[8] * sin(2.0 * kPi * x / b[6])
  381. NIST_END
  382. // y = (b1/b2) * exp[-0.5*((x-b3)/b2)**2] + e
  383. NIST_BEGIN(Eckerle4)
  384. b[0] / b[1] * exp(-0.5 * pow((x - b[2])/b[1], 2))
  385. NIST_END
  386. struct Nelson {
  387. public:
  388. Nelson(const double* const x, const double* const y, const int n)
  389. : x_(x), y_(y), n_(n) {}
  390. template <typename T>
  391. bool operator()(const T* const b, T* residual) const {
  392. // log[y] = b1 - b2*x1 * exp[-b3*x2] + e
  393. for (int i = 0; i < n_; ++i) {
  394. residual[i] = log(y_[i]) - (b[0] - b[1] * x_[2 * i] * exp(-b[2] * x_[2 * i + 1]));
  395. }
  396. return true;
  397. }
  398. private:
  399. const double* x_;
  400. const double* y_;
  401. const int n_;
  402. };
  403. // clang-format on
  404. static void SetNumericDiffOptions(ceres::NumericDiffOptions* options) {
  405. options->max_num_ridders_extrapolations =
  406. CERES_GET_FLAG(FLAGS_ridders_extrapolations);
  407. options->ridders_relative_initial_step_size =
  408. CERES_GET_FLAG(FLAGS_ridders_step_size);
  409. }
  410. void SetMinimizerOptions(ceres::Solver::Options* options) {
  411. CHECK(ceres::StringToMinimizerType(CERES_GET_FLAG(FLAGS_minimizer),
  412. &options->minimizer_type));
  413. CHECK(ceres::StringToLinearSolverType(CERES_GET_FLAG(FLAGS_linear_solver),
  414. &options->linear_solver_type));
  415. CHECK(StringToDenseLinearAlgebraLibraryType(
  416. CERES_GET_FLAG(FLAGS_dense_linear_algebra_library),
  417. &options->dense_linear_algebra_library_type));
  418. CHECK(ceres::StringToPreconditionerType(CERES_GET_FLAG(FLAGS_preconditioner),
  419. &options->preconditioner_type));
  420. CHECK(ceres::StringToTrustRegionStrategyType(
  421. CERES_GET_FLAG(FLAGS_trust_region_strategy),
  422. &options->trust_region_strategy_type));
  423. CHECK(ceres::StringToDoglegType(CERES_GET_FLAG(FLAGS_dogleg),
  424. &options->dogleg_type));
  425. CHECK(ceres::StringToLineSearchDirectionType(
  426. CERES_GET_FLAG(FLAGS_line_search_direction),
  427. &options->line_search_direction_type));
  428. CHECK(ceres::StringToLineSearchType(CERES_GET_FLAG(FLAGS_line_search),
  429. &options->line_search_type));
  430. CHECK(ceres::StringToLineSearchInterpolationType(
  431. CERES_GET_FLAG(FLAGS_line_search_interpolation),
  432. &options->line_search_interpolation_type));
  433. options->max_num_iterations = CERES_GET_FLAG(FLAGS_num_iterations);
  434. options->use_nonmonotonic_steps = CERES_GET_FLAG(FLAGS_nonmonotonic_steps);
  435. options->initial_trust_region_radius =
  436. CERES_GET_FLAG(FLAGS_initial_trust_region_radius);
  437. options->max_lbfgs_rank = CERES_GET_FLAG(FLAGS_lbfgs_rank);
  438. options->line_search_sufficient_function_decrease =
  439. CERES_GET_FLAG(FLAGS_sufficient_decrease);
  440. options->line_search_sufficient_curvature_decrease =
  441. CERES_GET_FLAG(FLAGS_sufficient_curvature_decrease);
  442. options->max_num_line_search_step_size_iterations =
  443. CERES_GET_FLAG(FLAGS_max_line_search_iterations);
  444. options->max_num_line_search_direction_restarts =
  445. CERES_GET_FLAG(FLAGS_max_line_search_restarts);
  446. options->use_approximate_eigenvalue_bfgs_scaling =
  447. CERES_GET_FLAG(FLAGS_approximate_eigenvalue_bfgs_scaling);
  448. options->function_tolerance = std::numeric_limits<double>::epsilon();
  449. options->gradient_tolerance = std::numeric_limits<double>::epsilon();
  450. options->parameter_tolerance = std::numeric_limits<double>::epsilon();
  451. }
  452. std::string JoinPath(const std::string& dirname, const std::string& basename) {
  453. #ifdef _WIN32
  454. static const char separator = '\\';
  455. #else
  456. static const char separator = '/';
  457. #endif // _WIN32
  458. if ((!basename.empty() && basename[0] == separator) || dirname.empty()) {
  459. return basename;
  460. } else if (dirname[dirname.size() - 1] == separator) {
  461. return dirname + basename;
  462. } else {
  463. return dirname + std::string(&separator, 1) + basename;
  464. }
  465. }
  466. template <typename Model, int num_parameters>
  467. CostFunction* CreateCostFunction(const Matrix& predictor,
  468. const Matrix& response,
  469. const int num_observations) {
  470. auto* model = new Model(predictor.data(), response.data(), num_observations);
  471. ceres::CostFunction* cost_function = nullptr;
  472. if (CERES_GET_FLAG(FLAGS_use_numeric_diff)) {
  473. ceres::NumericDiffOptions options;
  474. SetNumericDiffOptions(&options);
  475. if (CERES_GET_FLAG(FLAGS_numeric_diff_method) == "central") {
  476. cost_function = new NumericDiffCostFunction<Model,
  477. ceres::CENTRAL,
  478. ceres::DYNAMIC,
  479. num_parameters>(
  480. model, ceres::TAKE_OWNERSHIP, num_observations, options);
  481. } else if (CERES_GET_FLAG(FLAGS_numeric_diff_method) == "forward") {
  482. cost_function = new NumericDiffCostFunction<Model,
  483. ceres::FORWARD,
  484. ceres::DYNAMIC,
  485. num_parameters>(
  486. model, ceres::TAKE_OWNERSHIP, num_observations, options);
  487. } else if (CERES_GET_FLAG(FLAGS_numeric_diff_method) == "ridders") {
  488. cost_function = new NumericDiffCostFunction<Model,
  489. ceres::RIDDERS,
  490. ceres::DYNAMIC,
  491. num_parameters>(
  492. model, ceres::TAKE_OWNERSHIP, num_observations, options);
  493. } else {
  494. LOG(ERROR) << "Invalid numeric diff method specified";
  495. return nullptr;
  496. }
  497. } else {
  498. cost_function =
  499. new ceres::AutoDiffCostFunction<Model, ceres::DYNAMIC, num_parameters>(
  500. model, num_observations);
  501. }
  502. return cost_function;
  503. }
  504. double ComputeLRE(const Matrix& expected, const Matrix& actual) {
  505. // Compute the LRE by comparing each component of the solution
  506. // with the ground truth, and taking the minimum.
  507. const double kMaxNumSignificantDigits = 11;
  508. double log_relative_error = kMaxNumSignificantDigits + 1;
  509. for (int i = 0; i < expected.cols(); ++i) {
  510. const double tmp_lre = -std::log10(std::fabs(expected(i) - actual(i)) /
  511. std::fabs(expected(i)));
  512. // The maximum LRE is capped at 11 - the precision at which the
  513. // ground truth is known.
  514. //
  515. // The minimum LRE is capped at 0 - no digits match between the
  516. // computed solution and the ground truth.
  517. log_relative_error =
  518. std::min(log_relative_error,
  519. std::max(0.0, std::min(kMaxNumSignificantDigits, tmp_lre)));
  520. }
  521. return log_relative_error;
  522. }
  523. template <typename Model, int num_parameters>
  524. int RegressionDriver(const std::string& filename) {
  525. NISTProblem nist_problem(
  526. JoinPath(CERES_GET_FLAG(FLAGS_nist_data_dir), filename));
  527. CHECK_EQ(num_parameters, nist_problem.num_parameters());
  528. Matrix predictor = nist_problem.predictor();
  529. Matrix response = nist_problem.response();
  530. Matrix final_parameters = nist_problem.final_parameters();
  531. printf("%s\n", filename.c_str());
  532. // Each NIST problem comes with multiple starting points, so we
  533. // construct the problem from scratch for each case and solve it.
  534. int num_success = 0;
  535. for (int start = 0; start < nist_problem.num_starts(); ++start) {
  536. Matrix initial_parameters = nist_problem.initial_parameters(start);
  537. ceres::CostFunction* cost_function =
  538. CreateCostFunction<Model, num_parameters>(
  539. predictor, response, nist_problem.num_observations());
  540. double initial_cost;
  541. double final_cost;
  542. if (!CERES_GET_FLAG(FLAGS_use_tiny_solver)) {
  543. ceres::Problem problem;
  544. problem.AddResidualBlock(
  545. cost_function, nullptr, initial_parameters.data());
  546. ceres::Solver::Summary summary;
  547. ceres::Solver::Options options;
  548. SetMinimizerOptions(&options);
  549. Solve(options, &problem, &summary);
  550. initial_cost = summary.initial_cost;
  551. final_cost = summary.final_cost;
  552. } else {
  553. ceres::TinySolverCostFunctionAdapter<Eigen::Dynamic, num_parameters> cfa(
  554. *cost_function);
  555. using Solver = ceres::TinySolver<
  556. ceres::TinySolverCostFunctionAdapter<Eigen::Dynamic, num_parameters>>;
  557. Solver solver;
  558. solver.options.max_num_iterations = CERES_GET_FLAG(FLAGS_num_iterations);
  559. solver.options.gradient_tolerance =
  560. std::numeric_limits<double>::epsilon();
  561. solver.options.parameter_tolerance =
  562. std::numeric_limits<double>::epsilon();
  563. solver.options.function_tolerance = 0.0;
  564. Eigen::Matrix<double, num_parameters, 1> x;
  565. x = initial_parameters.transpose();
  566. typename Solver::Summary summary = solver.Solve(cfa, &x);
  567. initial_parameters = x;
  568. initial_cost = summary.initial_cost;
  569. final_cost = summary.final_cost;
  570. delete cost_function;
  571. }
  572. const double log_relative_error =
  573. ComputeLRE(nist_problem.final_parameters(), initial_parameters);
  574. const int kMinNumMatchingDigits = 4;
  575. if (log_relative_error > kMinNumMatchingDigits) {
  576. ++num_success;
  577. }
  578. printf(
  579. "start: %d status: %s lre: %4.1f initial cost: %e final cost:%e "
  580. "certified cost: %e\n",
  581. start + 1,
  582. log_relative_error < kMinNumMatchingDigits ? "FAILURE" : "SUCCESS",
  583. log_relative_error,
  584. initial_cost,
  585. final_cost,
  586. nist_problem.certified_cost());
  587. }
  588. return num_success;
  589. }
  590. void SolveNISTProblems() {
  591. if (CERES_GET_FLAG(FLAGS_nist_data_dir).empty()) {
  592. LOG(FATAL) << "Must specify the directory containing the NIST problems";
  593. }
  594. std::cout << "Lower Difficulty\n";
  595. int easy_success = 0;
  596. easy_success += RegressionDriver<Misra1a, 2>("Misra1a.dat");
  597. easy_success += RegressionDriver<Chwirut, 3>("Chwirut1.dat");
  598. easy_success += RegressionDriver<Chwirut, 3>("Chwirut2.dat");
  599. easy_success += RegressionDriver<Lanczos, 6>("Lanczos3.dat");
  600. easy_success += RegressionDriver<Gauss, 8>("Gauss1.dat");
  601. easy_success += RegressionDriver<Gauss, 8>("Gauss2.dat");
  602. easy_success += RegressionDriver<DanWood, 2>("DanWood.dat");
  603. easy_success += RegressionDriver<Misra1b, 2>("Misra1b.dat");
  604. std::cout << "\nMedium Difficulty\n";
  605. int medium_success = 0;
  606. medium_success += RegressionDriver<Kirby2, 5>("Kirby2.dat");
  607. medium_success += RegressionDriver<Hahn1, 7>("Hahn1.dat");
  608. medium_success += RegressionDriver<Nelson, 3>("Nelson.dat");
  609. medium_success += RegressionDriver<MGH17, 5>("MGH17.dat");
  610. medium_success += RegressionDriver<Lanczos, 6>("Lanczos1.dat");
  611. medium_success += RegressionDriver<Lanczos, 6>("Lanczos2.dat");
  612. medium_success += RegressionDriver<Gauss, 8>("Gauss3.dat");
  613. medium_success += RegressionDriver<Misra1c, 2>("Misra1c.dat");
  614. medium_success += RegressionDriver<Misra1d, 2>("Misra1d.dat");
  615. medium_success += RegressionDriver<Roszman1, 4>("Roszman1.dat");
  616. medium_success += RegressionDriver<ENSO, 9>("ENSO.dat");
  617. std::cout << "\nHigher Difficulty\n";
  618. int hard_success = 0;
  619. hard_success += RegressionDriver<MGH09, 4>("MGH09.dat");
  620. hard_success += RegressionDriver<Thurber, 7>("Thurber.dat");
  621. hard_success += RegressionDriver<BoxBOD, 2>("BoxBOD.dat");
  622. hard_success += RegressionDriver<Rat42, 3>("Rat42.dat");
  623. hard_success += RegressionDriver<MGH10, 3>("MGH10.dat");
  624. hard_success += RegressionDriver<Eckerle4, 3>("Eckerle4.dat");
  625. hard_success += RegressionDriver<Rat43, 4>("Rat43.dat");
  626. hard_success += RegressionDriver<Bennet5, 3>("Bennett5.dat");
  627. std::cout << "\n";
  628. std::cout << "Easy : " << easy_success << "/16\n";
  629. std::cout << "Medium : " << medium_success << "/22\n";
  630. std::cout << "Hard : " << hard_success << "/16\n";
  631. std::cout << "Total : " << easy_success + medium_success + hard_success
  632. << "/54\n";
  633. }
  634. } // namespace
  635. } // namespace ceres::examples
  636. int main(int argc, char** argv) {
  637. GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
  638. google::InitGoogleLogging(argv[0]);
  639. ceres::examples::SolveNISTProblems();
  640. return 0;
  641. }