You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

1861 lines
64 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 Thomas Capricelli <orzel@freehackers.org>
  5. #include <stdio.h>
  6. #include "main.h"
  7. #include <unsupported/Eigen/NonLinearOptimization>
  8. // This disables some useless Warnings on MSVC.
  9. // It is intended to be done for this test only.
  10. #include <Eigen/src/Core/util/DisableStupidWarnings.h>
  11. int fcn_chkder(const VectorXd &x, VectorXd &fvec, MatrixXd &fjac, int iflag)
  12. {
  13. /* subroutine fcn for chkder example. */
  14. int i;
  15. assert(15 == fvec.size());
  16. assert(3 == x.size());
  17. double tmp1, tmp2, tmp3, tmp4;
  18. static const double y[15]={1.4e-1, 1.8e-1, 2.2e-1, 2.5e-1, 2.9e-1, 3.2e-1, 3.5e-1,
  19. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  20. if (iflag == 0)
  21. return 0;
  22. if (iflag != 2)
  23. for (i=0; i<15; i++) {
  24. tmp1 = i+1;
  25. tmp2 = 16-i-1;
  26. tmp3 = tmp1;
  27. if (i >= 8) tmp3 = tmp2;
  28. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  29. }
  30. else {
  31. for (i = 0; i < 15; i++) {
  32. tmp1 = i+1;
  33. tmp2 = 16-i-1;
  34. /* error introduced into next statement for illustration. */
  35. /* corrected statement should read tmp3 = tmp1 . */
  36. tmp3 = tmp2;
  37. if (i >= 8) tmp3 = tmp2;
  38. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4=tmp4*tmp4;
  39. fjac(i,0) = -1.;
  40. fjac(i,1) = tmp1*tmp2/tmp4;
  41. fjac(i,2) = tmp1*tmp3/tmp4;
  42. }
  43. }
  44. return 0;
  45. }
  46. void testChkder()
  47. {
  48. const int m=15, n=3;
  49. VectorXd x(n), fvec(m), xp, fvecp(m), err;
  50. MatrixXd fjac(m,n);
  51. VectorXi ipvt;
  52. /* the following values should be suitable for */
  53. /* checking the jacobian matrix. */
  54. x << 9.2e-1, 1.3e-1, 5.4e-1;
  55. internal::chkder(x, fvec, fjac, xp, fvecp, 1, err);
  56. fcn_chkder(x, fvec, fjac, 1);
  57. fcn_chkder(x, fvec, fjac, 2);
  58. fcn_chkder(xp, fvecp, fjac, 1);
  59. internal::chkder(x, fvec, fjac, xp, fvecp, 2, err);
  60. fvecp -= fvec;
  61. // check those
  62. VectorXd fvec_ref(m), fvecp_ref(m), err_ref(m);
  63. fvec_ref <<
  64. -1.181606, -1.429655, -1.606344,
  65. -1.745269, -1.840654, -1.921586,
  66. -1.984141, -2.022537, -2.468977,
  67. -2.827562, -3.473582, -4.437612,
  68. -6.047662, -9.267761, -18.91806;
  69. fvecp_ref <<
  70. -7.724666e-09, -3.432406e-09, -2.034843e-10,
  71. 2.313685e-09, 4.331078e-09, 5.984096e-09,
  72. 7.363281e-09, 8.53147e-09, 1.488591e-08,
  73. 2.33585e-08, 3.522012e-08, 5.301255e-08,
  74. 8.26666e-08, 1.419747e-07, 3.19899e-07;
  75. err_ref <<
  76. 0.1141397, 0.09943516, 0.09674474,
  77. 0.09980447, 0.1073116, 0.1220445,
  78. 0.1526814, 1, 1,
  79. 1, 1, 1,
  80. 1, 1, 1;
  81. VERIFY_IS_APPROX(fvec, fvec_ref);
  82. VERIFY_IS_APPROX(fvecp, fvecp_ref);
  83. VERIFY_IS_APPROX(err, err_ref);
  84. }
  85. // Generic functor
  86. template<typename _Scalar, int NX=Dynamic, int NY=Dynamic>
  87. struct Functor
  88. {
  89. typedef _Scalar Scalar;
  90. enum {
  91. InputsAtCompileTime = NX,
  92. ValuesAtCompileTime = NY
  93. };
  94. typedef Matrix<Scalar,InputsAtCompileTime,1> InputType;
  95. typedef Matrix<Scalar,ValuesAtCompileTime,1> ValueType;
  96. typedef Matrix<Scalar,ValuesAtCompileTime,InputsAtCompileTime> JacobianType;
  97. const int m_inputs, m_values;
  98. Functor() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
  99. Functor(int inputs, int values) : m_inputs(inputs), m_values(values) {}
  100. int inputs() const { return m_inputs; }
  101. int values() const { return m_values; }
  102. // you should define that in the subclass :
  103. // void operator() (const InputType& x, ValueType* v, JacobianType* _j=0) const;
  104. };
  105. struct lmder_functor : Functor<double>
  106. {
  107. lmder_functor(void): Functor<double>(3,15) {}
  108. int operator()(const VectorXd &x, VectorXd &fvec) const
  109. {
  110. double tmp1, tmp2, tmp3;
  111. static const double y[15] = {1.4e-1, 1.8e-1, 2.2e-1, 2.5e-1, 2.9e-1, 3.2e-1, 3.5e-1,
  112. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  113. for (int i = 0; i < values(); i++)
  114. {
  115. tmp1 = i+1;
  116. tmp2 = 16 - i - 1;
  117. tmp3 = (i>=8)? tmp2 : tmp1;
  118. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  119. }
  120. return 0;
  121. }
  122. int df(const VectorXd &x, MatrixXd &fjac) const
  123. {
  124. double tmp1, tmp2, tmp3, tmp4;
  125. for (int i = 0; i < values(); i++)
  126. {
  127. tmp1 = i+1;
  128. tmp2 = 16 - i - 1;
  129. tmp3 = (i>=8)? tmp2 : tmp1;
  130. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4 = tmp4*tmp4;
  131. fjac(i,0) = -1;
  132. fjac(i,1) = tmp1*tmp2/tmp4;
  133. fjac(i,2) = tmp1*tmp3/tmp4;
  134. }
  135. return 0;
  136. }
  137. };
  138. void testLmder1()
  139. {
  140. int n=3, info;
  141. VectorXd x;
  142. /* the following starting values provide a rough fit. */
  143. x.setConstant(n, 1.);
  144. // do the computation
  145. lmder_functor functor;
  146. LevenbergMarquardt<lmder_functor> lm(functor);
  147. info = lm.lmder1(x);
  148. // check return value
  149. VERIFY_IS_EQUAL(info, 1);
  150. VERIFY_IS_EQUAL(lm.nfev, 6);
  151. VERIFY_IS_EQUAL(lm.njev, 5);
  152. // check norm
  153. VERIFY_IS_APPROX(lm.fvec.blueNorm(), 0.09063596);
  154. // check x
  155. VectorXd x_ref(n);
  156. x_ref << 0.08241058, 1.133037, 2.343695;
  157. VERIFY_IS_APPROX(x, x_ref);
  158. }
  159. void testLmder()
  160. {
  161. const int m=15, n=3;
  162. int info;
  163. double fnorm, covfac;
  164. VectorXd x;
  165. /* the following starting values provide a rough fit. */
  166. x.setConstant(n, 1.);
  167. // do the computation
  168. lmder_functor functor;
  169. LevenbergMarquardt<lmder_functor> lm(functor);
  170. info = lm.minimize(x);
  171. // check return values
  172. VERIFY_IS_EQUAL(info, 1);
  173. VERIFY_IS_EQUAL(lm.nfev, 6);
  174. VERIFY_IS_EQUAL(lm.njev, 5);
  175. // check norm
  176. fnorm = lm.fvec.blueNorm();
  177. VERIFY_IS_APPROX(fnorm, 0.09063596);
  178. // check x
  179. VectorXd x_ref(n);
  180. x_ref << 0.08241058, 1.133037, 2.343695;
  181. VERIFY_IS_APPROX(x, x_ref);
  182. // check covariance
  183. covfac = fnorm*fnorm/(m-n);
  184. internal::covar(lm.fjac, lm.permutation.indices()); // TODO : move this as a function of lm
  185. MatrixXd cov_ref(n,n);
  186. cov_ref <<
  187. 0.0001531202, 0.002869941, -0.002656662,
  188. 0.002869941, 0.09480935, -0.09098995,
  189. -0.002656662, -0.09098995, 0.08778727;
  190. // std::cout << fjac*covfac << std::endl;
  191. MatrixXd cov;
  192. cov = covfac*lm.fjac.topLeftCorner<n,n>();
  193. VERIFY_IS_APPROX( cov, cov_ref);
  194. // TODO: why isn't this allowed ? :
  195. // VERIFY_IS_APPROX( covfac*fjac.topLeftCorner<n,n>() , cov_ref);
  196. }
  197. struct hybrj_functor : Functor<double>
  198. {
  199. hybrj_functor(void) : Functor<double>(9,9) {}
  200. int operator()(const VectorXd &x, VectorXd &fvec)
  201. {
  202. double temp, temp1, temp2;
  203. const int n = x.size();
  204. assert(fvec.size()==n);
  205. for (int k = 0; k < n; k++)
  206. {
  207. temp = (3. - 2.*x[k])*x[k];
  208. temp1 = 0.;
  209. if (k) temp1 = x[k-1];
  210. temp2 = 0.;
  211. if (k != n-1) temp2 = x[k+1];
  212. fvec[k] = temp - temp1 - 2.*temp2 + 1.;
  213. }
  214. return 0;
  215. }
  216. int df(const VectorXd &x, MatrixXd &fjac)
  217. {
  218. const int n = x.size();
  219. assert(fjac.rows()==n);
  220. assert(fjac.cols()==n);
  221. for (int k = 0; k < n; k++)
  222. {
  223. for (int j = 0; j < n; j++)
  224. fjac(k,j) = 0.;
  225. fjac(k,k) = 3.- 4.*x[k];
  226. if (k) fjac(k,k-1) = -1.;
  227. if (k != n-1) fjac(k,k+1) = -2.;
  228. }
  229. return 0;
  230. }
  231. };
  232. void testHybrj1()
  233. {
  234. const int n=9;
  235. int info;
  236. VectorXd x(n);
  237. /* the following starting values provide a rough fit. */
  238. x.setConstant(n, -1.);
  239. // do the computation
  240. hybrj_functor functor;
  241. HybridNonLinearSolver<hybrj_functor> solver(functor);
  242. info = solver.hybrj1(x);
  243. // check return value
  244. VERIFY_IS_EQUAL(info, 1);
  245. VERIFY_IS_EQUAL(solver.nfev, 11);
  246. VERIFY_IS_EQUAL(solver.njev, 1);
  247. // check norm
  248. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  249. // check x
  250. VectorXd x_ref(n);
  251. x_ref <<
  252. -0.5706545, -0.6816283, -0.7017325,
  253. -0.7042129, -0.701369, -0.6918656,
  254. -0.665792, -0.5960342, -0.4164121;
  255. VERIFY_IS_APPROX(x, x_ref);
  256. }
  257. void testHybrj()
  258. {
  259. const int n=9;
  260. int info;
  261. VectorXd x(n);
  262. /* the following starting values provide a rough fit. */
  263. x.setConstant(n, -1.);
  264. // do the computation
  265. hybrj_functor functor;
  266. HybridNonLinearSolver<hybrj_functor> solver(functor);
  267. solver.diag.setConstant(n, 1.);
  268. solver.useExternalScaling = true;
  269. info = solver.solve(x);
  270. // check return value
  271. VERIFY_IS_EQUAL(info, 1);
  272. VERIFY_IS_EQUAL(solver.nfev, 11);
  273. VERIFY_IS_EQUAL(solver.njev, 1);
  274. // check norm
  275. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  276. // check x
  277. VectorXd x_ref(n);
  278. x_ref <<
  279. -0.5706545, -0.6816283, -0.7017325,
  280. -0.7042129, -0.701369, -0.6918656,
  281. -0.665792, -0.5960342, -0.4164121;
  282. VERIFY_IS_APPROX(x, x_ref);
  283. }
  284. struct hybrd_functor : Functor<double>
  285. {
  286. hybrd_functor(void) : Functor<double>(9,9) {}
  287. int operator()(const VectorXd &x, VectorXd &fvec) const
  288. {
  289. double temp, temp1, temp2;
  290. const int n = x.size();
  291. assert(fvec.size()==n);
  292. for (int k=0; k < n; k++)
  293. {
  294. temp = (3. - 2.*x[k])*x[k];
  295. temp1 = 0.;
  296. if (k) temp1 = x[k-1];
  297. temp2 = 0.;
  298. if (k != n-1) temp2 = x[k+1];
  299. fvec[k] = temp - temp1 - 2.*temp2 + 1.;
  300. }
  301. return 0;
  302. }
  303. };
  304. void testHybrd1()
  305. {
  306. int n=9, info;
  307. VectorXd x(n);
  308. /* the following starting values provide a rough solution. */
  309. x.setConstant(n, -1.);
  310. // do the computation
  311. hybrd_functor functor;
  312. HybridNonLinearSolver<hybrd_functor> solver(functor);
  313. info = solver.hybrd1(x);
  314. // check return value
  315. VERIFY_IS_EQUAL(info, 1);
  316. VERIFY_IS_EQUAL(solver.nfev, 20);
  317. // check norm
  318. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  319. // check x
  320. VectorXd x_ref(n);
  321. x_ref << -0.5706545, -0.6816283, -0.7017325, -0.7042129, -0.701369, -0.6918656, -0.665792, -0.5960342, -0.4164121;
  322. VERIFY_IS_APPROX(x, x_ref);
  323. }
  324. void testHybrd()
  325. {
  326. const int n=9;
  327. int info;
  328. VectorXd x;
  329. /* the following starting values provide a rough fit. */
  330. x.setConstant(n, -1.);
  331. // do the computation
  332. hybrd_functor functor;
  333. HybridNonLinearSolver<hybrd_functor> solver(functor);
  334. solver.parameters.nb_of_subdiagonals = 1;
  335. solver.parameters.nb_of_superdiagonals = 1;
  336. solver.diag.setConstant(n, 1.);
  337. solver.useExternalScaling = true;
  338. info = solver.solveNumericalDiff(x);
  339. // check return value
  340. VERIFY_IS_EQUAL(info, 1);
  341. VERIFY_IS_EQUAL(solver.nfev, 14);
  342. // check norm
  343. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  344. // check x
  345. VectorXd x_ref(n);
  346. x_ref <<
  347. -0.5706545, -0.6816283, -0.7017325,
  348. -0.7042129, -0.701369, -0.6918656,
  349. -0.665792, -0.5960342, -0.4164121;
  350. VERIFY_IS_APPROX(x, x_ref);
  351. }
  352. struct lmstr_functor : Functor<double>
  353. {
  354. lmstr_functor(void) : Functor<double>(3,15) {}
  355. int operator()(const VectorXd &x, VectorXd &fvec)
  356. {
  357. /* subroutine fcn for lmstr1 example. */
  358. double tmp1, tmp2, tmp3;
  359. static const double y[15]={1.4e-1, 1.8e-1, 2.2e-1, 2.5e-1, 2.9e-1, 3.2e-1, 3.5e-1,
  360. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  361. assert(15==fvec.size());
  362. assert(3==x.size());
  363. for (int i=0; i<15; i++)
  364. {
  365. tmp1 = i+1;
  366. tmp2 = 16 - i - 1;
  367. tmp3 = (i>=8)? tmp2 : tmp1;
  368. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  369. }
  370. return 0;
  371. }
  372. int df(const VectorXd &x, VectorXd &jac_row, VectorXd::Index rownb)
  373. {
  374. assert(x.size()==3);
  375. assert(jac_row.size()==x.size());
  376. double tmp1, tmp2, tmp3, tmp4;
  377. int i = rownb-2;
  378. tmp1 = i+1;
  379. tmp2 = 16 - i - 1;
  380. tmp3 = (i>=8)? tmp2 : tmp1;
  381. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4 = tmp4*tmp4;
  382. jac_row[0] = -1;
  383. jac_row[1] = tmp1*tmp2/tmp4;
  384. jac_row[2] = tmp1*tmp3/tmp4;
  385. return 0;
  386. }
  387. };
  388. void testLmstr1()
  389. {
  390. const int n=3;
  391. int info;
  392. VectorXd x(n);
  393. /* the following starting values provide a rough fit. */
  394. x.setConstant(n, 1.);
  395. // do the computation
  396. lmstr_functor functor;
  397. LevenbergMarquardt<lmstr_functor> lm(functor);
  398. info = lm.lmstr1(x);
  399. // check return value
  400. VERIFY_IS_EQUAL(info, 1);
  401. VERIFY_IS_EQUAL(lm.nfev, 6);
  402. VERIFY_IS_EQUAL(lm.njev, 5);
  403. // check norm
  404. VERIFY_IS_APPROX(lm.fvec.blueNorm(), 0.09063596);
  405. // check x
  406. VectorXd x_ref(n);
  407. x_ref << 0.08241058, 1.133037, 2.343695 ;
  408. VERIFY_IS_APPROX(x, x_ref);
  409. }
  410. void testLmstr()
  411. {
  412. const int n=3;
  413. int info;
  414. double fnorm;
  415. VectorXd x(n);
  416. /* the following starting values provide a rough fit. */
  417. x.setConstant(n, 1.);
  418. // do the computation
  419. lmstr_functor functor;
  420. LevenbergMarquardt<lmstr_functor> lm(functor);
  421. info = lm.minimizeOptimumStorage(x);
  422. // check return values
  423. VERIFY_IS_EQUAL(info, 1);
  424. VERIFY_IS_EQUAL(lm.nfev, 6);
  425. VERIFY_IS_EQUAL(lm.njev, 5);
  426. // check norm
  427. fnorm = lm.fvec.blueNorm();
  428. VERIFY_IS_APPROX(fnorm, 0.09063596);
  429. // check x
  430. VectorXd x_ref(n);
  431. x_ref << 0.08241058, 1.133037, 2.343695;
  432. VERIFY_IS_APPROX(x, x_ref);
  433. }
  434. struct lmdif_functor : Functor<double>
  435. {
  436. lmdif_functor(void) : Functor<double>(3,15) {}
  437. int operator()(const VectorXd &x, VectorXd &fvec) const
  438. {
  439. int i;
  440. double tmp1,tmp2,tmp3;
  441. static const double y[15]={1.4e-1,1.8e-1,2.2e-1,2.5e-1,2.9e-1,3.2e-1,3.5e-1,3.9e-1,
  442. 3.7e-1,5.8e-1,7.3e-1,9.6e-1,1.34e0,2.1e0,4.39e0};
  443. assert(x.size()==3);
  444. assert(fvec.size()==15);
  445. for (i=0; i<15; i++)
  446. {
  447. tmp1 = i+1;
  448. tmp2 = 15 - i;
  449. tmp3 = tmp1;
  450. if (i >= 8) tmp3 = tmp2;
  451. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  452. }
  453. return 0;
  454. }
  455. };
  456. void testLmdif1()
  457. {
  458. const int n=3;
  459. int info;
  460. VectorXd x(n), fvec(15);
  461. /* the following starting values provide a rough fit. */
  462. x.setConstant(n, 1.);
  463. // do the computation
  464. lmdif_functor functor;
  465. DenseIndex nfev;
  466. info = LevenbergMarquardt<lmdif_functor>::lmdif1(functor, x, &nfev);
  467. // check return value
  468. VERIFY_IS_EQUAL(info, 1);
  469. VERIFY_IS_EQUAL(nfev, 26);
  470. // check norm
  471. functor(x, fvec);
  472. VERIFY_IS_APPROX(fvec.blueNorm(), 0.09063596);
  473. // check x
  474. VectorXd x_ref(n);
  475. x_ref << 0.0824106, 1.1330366, 2.3436947;
  476. VERIFY_IS_APPROX(x, x_ref);
  477. }
  478. void testLmdif()
  479. {
  480. const int m=15, n=3;
  481. int info;
  482. double fnorm, covfac;
  483. VectorXd x(n);
  484. /* the following starting values provide a rough fit. */
  485. x.setConstant(n, 1.);
  486. // do the computation
  487. lmdif_functor functor;
  488. NumericalDiff<lmdif_functor> numDiff(functor);
  489. LevenbergMarquardt<NumericalDiff<lmdif_functor> > lm(numDiff);
  490. info = lm.minimize(x);
  491. // check return values
  492. VERIFY_IS_EQUAL(info, 1);
  493. VERIFY_IS_EQUAL(lm.nfev, 26);
  494. // check norm
  495. fnorm = lm.fvec.blueNorm();
  496. VERIFY_IS_APPROX(fnorm, 0.09063596);
  497. // check x
  498. VectorXd x_ref(n);
  499. x_ref << 0.08241058, 1.133037, 2.343695;
  500. VERIFY_IS_APPROX(x, x_ref);
  501. // check covariance
  502. covfac = fnorm*fnorm/(m-n);
  503. internal::covar(lm.fjac, lm.permutation.indices()); // TODO : move this as a function of lm
  504. MatrixXd cov_ref(n,n);
  505. cov_ref <<
  506. 0.0001531202, 0.002869942, -0.002656662,
  507. 0.002869942, 0.09480937, -0.09098997,
  508. -0.002656662, -0.09098997, 0.08778729;
  509. // std::cout << fjac*covfac << std::endl;
  510. MatrixXd cov;
  511. cov = covfac*lm.fjac.topLeftCorner<n,n>();
  512. VERIFY_IS_APPROX( cov, cov_ref);
  513. // TODO: why isn't this allowed ? :
  514. // VERIFY_IS_APPROX( covfac*fjac.topLeftCorner<n,n>() , cov_ref);
  515. }
  516. struct chwirut2_functor : Functor<double>
  517. {
  518. chwirut2_functor(void) : Functor<double>(3,54) {}
  519. static const double m_x[54];
  520. static const double m_y[54];
  521. int operator()(const VectorXd &b, VectorXd &fvec)
  522. {
  523. int i;
  524. assert(b.size()==3);
  525. assert(fvec.size()==54);
  526. for(i=0; i<54; i++) {
  527. double x = m_x[i];
  528. fvec[i] = exp(-b[0]*x)/(b[1]+b[2]*x) - m_y[i];
  529. }
  530. return 0;
  531. }
  532. int df(const VectorXd &b, MatrixXd &fjac)
  533. {
  534. assert(b.size()==3);
  535. assert(fjac.rows()==54);
  536. assert(fjac.cols()==3);
  537. for(int i=0; i<54; i++) {
  538. double x = m_x[i];
  539. double factor = 1./(b[1]+b[2]*x);
  540. double e = exp(-b[0]*x);
  541. fjac(i,0) = -x*e*factor;
  542. fjac(i,1) = -e*factor*factor;
  543. fjac(i,2) = -x*e*factor*factor;
  544. }
  545. return 0;
  546. }
  547. };
  548. const double chwirut2_functor::m_x[54] = { 0.500E0, 1.000E0, 1.750E0, 3.750E0, 5.750E0, 0.875E0, 2.250E0, 3.250E0, 5.250E0, 0.750E0, 1.750E0, 2.750E0, 4.750E0, 0.625E0, 1.250E0, 2.250E0, 4.250E0, .500E0, 3.000E0, .750E0, 3.000E0, 1.500E0, 6.000E0, 3.000E0, 6.000E0, 1.500E0, 3.000E0, .500E0, 2.000E0, 4.000E0, .750E0, 2.000E0, 5.000E0, .750E0, 2.250E0, 3.750E0, 5.750E0, 3.000E0, .750E0, 2.500E0, 4.000E0, .750E0, 2.500E0, 4.000E0, .750E0, 2.500E0, 4.000E0, .500E0, 6.000E0, 3.000E0, .500E0, 2.750E0, .500E0, 1.750E0};
  549. const double chwirut2_functor::m_y[54] = { 92.9000E0 ,57.1000E0 ,31.0500E0 ,11.5875E0 ,8.0250E0 ,63.6000E0 ,21.4000E0 ,14.2500E0 ,8.4750E0 ,63.8000E0 ,26.8000E0 ,16.4625E0 ,7.1250E0 ,67.3000E0 ,41.0000E0 ,21.1500E0 ,8.1750E0 ,81.5000E0 ,13.1200E0 ,59.9000E0 ,14.6200E0 ,32.9000E0 ,5.4400E0 ,12.5600E0 ,5.4400E0 ,32.0000E0 ,13.9500E0 ,75.8000E0 ,20.0000E0 ,10.4200E0 ,59.5000E0 ,21.6700E0 ,8.5500E0 ,62.0000E0 ,20.2000E0 ,7.7600E0 ,3.7500E0 ,11.8100E0 ,54.7000E0 ,23.7000E0 ,11.5500E0 ,61.3000E0 ,17.7000E0 ,8.7400E0 ,59.2000E0 ,16.3000E0 ,8.6200E0 ,81.0000E0 ,4.8700E0 ,14.6200E0 ,81.7000E0 ,17.1700E0 ,81.3000E0 ,28.9000E0 };
  550. // http://www.itl.nist.gov/div898/strd/nls/data/chwirut2.shtml
  551. void testNistChwirut2(void)
  552. {
  553. const int n=3;
  554. int info;
  555. VectorXd x(n);
  556. /*
  557. * First try
  558. */
  559. x<< 0.1, 0.01, 0.02;
  560. // do the computation
  561. chwirut2_functor functor;
  562. LevenbergMarquardt<chwirut2_functor> lm(functor);
  563. info = lm.minimize(x);
  564. // check return value
  565. VERIFY_IS_EQUAL(info, 1);
  566. VERIFY_IS_EQUAL(lm.nfev, 10);
  567. VERIFY_IS_EQUAL(lm.njev, 8);
  568. // check norm^2
  569. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.1304802941E+02);
  570. // check x
  571. VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
  572. VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
  573. VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
  574. /*
  575. * Second try
  576. */
  577. x<< 0.15, 0.008, 0.010;
  578. // do the computation
  579. lm.resetParameters();
  580. lm.parameters.ftol = 1.E6*NumTraits<double>::epsilon();
  581. lm.parameters.xtol = 1.E6*NumTraits<double>::epsilon();
  582. info = lm.minimize(x);
  583. // check return value
  584. VERIFY_IS_EQUAL(info, 1);
  585. VERIFY_IS_EQUAL(lm.nfev, 7);
  586. VERIFY_IS_EQUAL(lm.njev, 6);
  587. // check norm^2
  588. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.1304802941E+02);
  589. // check x
  590. VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
  591. VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
  592. VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
  593. }
  594. struct misra1a_functor : Functor<double>
  595. {
  596. misra1a_functor(void) : Functor<double>(2,14) {}
  597. static const double m_x[14];
  598. static const double m_y[14];
  599. int operator()(const VectorXd &b, VectorXd &fvec)
  600. {
  601. assert(b.size()==2);
  602. assert(fvec.size()==14);
  603. for(int i=0; i<14; i++) {
  604. fvec[i] = b[0]*(1.-exp(-b[1]*m_x[i])) - m_y[i] ;
  605. }
  606. return 0;
  607. }
  608. int df(const VectorXd &b, MatrixXd &fjac)
  609. {
  610. assert(b.size()==2);
  611. assert(fjac.rows()==14);
  612. assert(fjac.cols()==2);
  613. for(int i=0; i<14; i++) {
  614. fjac(i,0) = (1.-exp(-b[1]*m_x[i]));
  615. fjac(i,1) = (b[0]*m_x[i]*exp(-b[1]*m_x[i]));
  616. }
  617. return 0;
  618. }
  619. };
  620. const double misra1a_functor::m_x[14] = { 77.6E0, 114.9E0, 141.1E0, 190.8E0, 239.9E0, 289.0E0, 332.8E0, 378.4E0, 434.8E0, 477.3E0, 536.8E0, 593.1E0, 689.1E0, 760.0E0};
  621. const double misra1a_functor::m_y[14] = { 10.07E0, 14.73E0, 17.94E0, 23.93E0, 29.61E0, 35.18E0, 40.02E0, 44.82E0, 50.76E0, 55.05E0, 61.01E0, 66.40E0, 75.47E0, 81.78E0};
  622. // http://www.itl.nist.gov/div898/strd/nls/data/misra1a.shtml
  623. void testNistMisra1a(void)
  624. {
  625. const int n=2;
  626. int info;
  627. VectorXd x(n);
  628. /*
  629. * First try
  630. */
  631. x<< 500., 0.0001;
  632. // do the computation
  633. misra1a_functor functor;
  634. LevenbergMarquardt<misra1a_functor> lm(functor);
  635. info = lm.minimize(x);
  636. // check return value
  637. VERIFY_IS_EQUAL(info, 1);
  638. VERIFY_IS_EQUAL(lm.nfev, 19);
  639. VERIFY_IS_EQUAL(lm.njev, 15);
  640. // check norm^2
  641. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.2455138894E-01);
  642. // check x
  643. VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
  644. VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
  645. /*
  646. * Second try
  647. */
  648. x<< 250., 0.0005;
  649. // do the computation
  650. info = lm.minimize(x);
  651. // check return value
  652. VERIFY_IS_EQUAL(info, 1);
  653. VERIFY_IS_EQUAL(lm.nfev, 5);
  654. VERIFY_IS_EQUAL(lm.njev, 4);
  655. // check norm^2
  656. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.2455138894E-01);
  657. // check x
  658. VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
  659. VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
  660. }
  661. struct hahn1_functor : Functor<double>
  662. {
  663. hahn1_functor(void) : Functor<double>(7,236) {}
  664. static const double m_x[236];
  665. int operator()(const VectorXd &b, VectorXd &fvec)
  666. {
  667. static const double m_y[236] = { .591E0 , 1.547E0 , 2.902E0 , 2.894E0 , 4.703E0 , 6.307E0 , 7.03E0 , 7.898E0 , 9.470E0 , 9.484E0 , 10.072E0 , 10.163E0 , 11.615E0 , 12.005E0 , 12.478E0 , 12.982E0 , 12.970E0 , 13.926E0 , 14.452E0 , 14.404E0 , 15.190E0 , 15.550E0 , 15.528E0 , 15.499E0 , 16.131E0 , 16.438E0 , 16.387E0 , 16.549E0 , 16.872E0 , 16.830E0 , 16.926E0 , 16.907E0 , 16.966E0 , 17.060E0 , 17.122E0 , 17.311E0 , 17.355E0 , 17.668E0 , 17.767E0 , 17.803E0 , 17.765E0 , 17.768E0 , 17.736E0 , 17.858E0 , 17.877E0 , 17.912E0 , 18.046E0 , 18.085E0 , 18.291E0 , 18.357E0 , 18.426E0 , 18.584E0 , 18.610E0 , 18.870E0 , 18.795E0 , 19.111E0 , .367E0 , .796E0 , 0.892E0 , 1.903E0 , 2.150E0 , 3.697E0 , 5.870E0 , 6.421E0 , 7.422E0 , 9.944E0 , 11.023E0 , 11.87E0 , 12.786E0 , 14.067E0 , 13.974E0 , 14.462E0 , 14.464E0 , 15.381E0 , 15.483E0 , 15.59E0 , 16.075E0 , 16.347E0 , 16.181E0 , 16.915E0 , 17.003E0 , 16.978E0 , 17.756E0 , 17.808E0 , 17.868E0 , 18.481E0 , 18.486E0 , 19.090E0 , 16.062E0 , 16.337E0 , 16.345E0 , 16.388E0 , 17.159E0 , 17.116E0 , 17.164E0 , 17.123E0 , 17.979E0 , 17.974E0 , 18.007E0 , 17.993E0 , 18.523E0 , 18.669E0 , 18.617E0 , 19.371E0 , 19.330E0 , 0.080E0 , 0.248E0 , 1.089E0 , 1.418E0 , 2.278E0 , 3.624E0 , 4.574E0 , 5.556E0 , 7.267E0 , 7.695E0 , 9.136E0 , 9.959E0 , 9.957E0 , 11.600E0 , 13.138E0 , 13.564E0 , 13.871E0 , 13.994E0 , 14.947E0 , 15.473E0 , 15.379E0 , 15.455E0 , 15.908E0 , 16.114E0 , 17.071E0 , 17.135E0 , 17.282E0 , 17.368E0 , 17.483E0 , 17.764E0 , 18.185E0 , 18.271E0 , 18.236E0 , 18.237E0 , 18.523E0 , 18.627E0 , 18.665E0 , 19.086E0 , 0.214E0 , 0.943E0 , 1.429E0 , 2.241E0 , 2.951E0 , 3.782E0 , 4.757E0 , 5.602E0 , 7.169E0 , 8.920E0 , 10.055E0 , 12.035E0 , 12.861E0 , 13.436E0 , 14.167E0 , 14.755E0 , 15.168E0 , 15.651E0 , 15.746E0 , 16.216E0 , 16.445E0 , 16.965E0 , 17.121E0 , 17.206E0 , 17.250E0 , 17.339E0 , 17.793E0 , 18.123E0 , 18.49E0 , 18.566E0 , 18.645E0 , 18.706E0 , 18.924E0 , 19.1E0 , 0.375E0 , 0.471E0 , 1.504E0 , 2.204E0 , 2.813E0 , 4.765E0 , 9.835E0 , 10.040E0 , 11.946E0 , 12.596E0 , 13.303E0 , 13.922E0 , 14.440E0 , 14.951E0 , 15.627E0 , 15.639E0 , 15.814E0 , 16.315E0 , 16.334E0 , 16.430E0 , 16.423E0 , 17.024E0 , 17.009E0 , 17.165E0 , 17.134E0 , 17.349E0 , 17.576E0 , 17.848E0 , 18.090E0 , 18.276E0 , 18.404E0 , 18.519E0 , 19.133E0 , 19.074E0 , 19.239E0 , 19.280E0 , 19.101E0 , 19.398E0 , 19.252E0 , 19.89E0 , 20.007E0 , 19.929E0 , 19.268E0 , 19.324E0 , 20.049E0 , 20.107E0 , 20.062E0 , 20.065E0 , 19.286E0 , 19.972E0 , 20.088E0 , 20.743E0 , 20.83E0 , 20.935E0 , 21.035E0 , 20.93E0 , 21.074E0 , 21.085E0 , 20.935E0 };
  668. // int called=0; printf("call hahn1_functor with iflag=%d, called=%d\n", iflag, called); if (iflag==1) called++;
  669. assert(b.size()==7);
  670. assert(fvec.size()==236);
  671. for(int i=0; i<236; i++) {
  672. double x=m_x[i], xx=x*x, xxx=xx*x;
  673. fvec[i] = (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) / (1.+b[4]*x+b[5]*xx+b[6]*xxx) - m_y[i];
  674. }
  675. return 0;
  676. }
  677. int df(const VectorXd &b, MatrixXd &fjac)
  678. {
  679. assert(b.size()==7);
  680. assert(fjac.rows()==236);
  681. assert(fjac.cols()==7);
  682. for(int i=0; i<236; i++) {
  683. double x=m_x[i], xx=x*x, xxx=xx*x;
  684. double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
  685. fjac(i,0) = 1.*fact;
  686. fjac(i,1) = x*fact;
  687. fjac(i,2) = xx*fact;
  688. fjac(i,3) = xxx*fact;
  689. fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
  690. fjac(i,4) = x*fact;
  691. fjac(i,5) = xx*fact;
  692. fjac(i,6) = xxx*fact;
  693. }
  694. return 0;
  695. }
  696. };
  697. const double hahn1_functor::m_x[236] = { 24.41E0 , 34.82E0 , 44.09E0 , 45.07E0 , 54.98E0 , 65.51E0 , 70.53E0 , 75.70E0 , 89.57E0 , 91.14E0 , 96.40E0 , 97.19E0 , 114.26E0 , 120.25E0 , 127.08E0 , 133.55E0 , 133.61E0 , 158.67E0 , 172.74E0 , 171.31E0 , 202.14E0 , 220.55E0 , 221.05E0 , 221.39E0 , 250.99E0 , 268.99E0 , 271.80E0 , 271.97E0 , 321.31E0 , 321.69E0 , 330.14E0 , 333.03E0 , 333.47E0 , 340.77E0 , 345.65E0 , 373.11E0 , 373.79E0 , 411.82E0 , 419.51E0 , 421.59E0 , 422.02E0 , 422.47E0 , 422.61E0 , 441.75E0 , 447.41E0 , 448.7E0 , 472.89E0 , 476.69E0 , 522.47E0 , 522.62E0 , 524.43E0 , 546.75E0 , 549.53E0 , 575.29E0 , 576.00E0 , 625.55E0 , 20.15E0 , 28.78E0 , 29.57E0 , 37.41E0 , 39.12E0 , 50.24E0 , 61.38E0 , 66.25E0 , 73.42E0 , 95.52E0 , 107.32E0 , 122.04E0 , 134.03E0 , 163.19E0 , 163.48E0 , 175.70E0 , 179.86E0 , 211.27E0 , 217.78E0 , 219.14E0 , 262.52E0 , 268.01E0 , 268.62E0 , 336.25E0 , 337.23E0 , 339.33E0 , 427.38E0 , 428.58E0 , 432.68E0 , 528.99E0 , 531.08E0 , 628.34E0 , 253.24E0 , 273.13E0 , 273.66E0 , 282.10E0 , 346.62E0 , 347.19E0 , 348.78E0 , 351.18E0 , 450.10E0 , 450.35E0 , 451.92E0 , 455.56E0 , 552.22E0 , 553.56E0 , 555.74E0 , 652.59E0 , 656.20E0 , 14.13E0 , 20.41E0 , 31.30E0 , 33.84E0 , 39.70E0 , 48.83E0 , 54.50E0 , 60.41E0 , 72.77E0 , 75.25E0 , 86.84E0 , 94.88E0 , 96.40E0 , 117.37E0 , 139.08E0 , 147.73E0 , 158.63E0 , 161.84E0 , 192.11E0 , 206.76E0 , 209.07E0 , 213.32E0 , 226.44E0 , 237.12E0 , 330.90E0 , 358.72E0 , 370.77E0 , 372.72E0 , 396.24E0 , 416.59E0 , 484.02E0 , 495.47E0 , 514.78E0 , 515.65E0 , 519.47E0 , 544.47E0 , 560.11E0 , 620.77E0 , 18.97E0 , 28.93E0 , 33.91E0 , 40.03E0 , 44.66E0 , 49.87E0 , 55.16E0 , 60.90E0 , 72.08E0 , 85.15E0 , 97.06E0 , 119.63E0 , 133.27E0 , 143.84E0 , 161.91E0 , 180.67E0 , 198.44E0 , 226.86E0 , 229.65E0 , 258.27E0 , 273.77E0 , 339.15E0 , 350.13E0 , 362.75E0 , 371.03E0 , 393.32E0 , 448.53E0 , 473.78E0 , 511.12E0 , 524.70E0 , 548.75E0 , 551.64E0 , 574.02E0 , 623.86E0 , 21.46E0 , 24.33E0 , 33.43E0 , 39.22E0 , 44.18E0 , 55.02E0 , 94.33E0 , 96.44E0 , 118.82E0 , 128.48E0 , 141.94E0 , 156.92E0 , 171.65E0 , 190.00E0 , 223.26E0 , 223.88E0 , 231.50E0 , 265.05E0 , 269.44E0 , 271.78E0 , 273.46E0 , 334.61E0 , 339.79E0 , 349.52E0 , 358.18E0 , 377.98E0 , 394.77E0 , 429.66E0 , 468.22E0 , 487.27E0 , 519.54E0 , 523.03E0 , 612.99E0 , 638.59E0 , 641.36E0 , 622.05E0 , 631.50E0 , 663.97E0 , 646.9E0 , 748.29E0 , 749.21E0 , 750.14E0 , 647.04E0 , 646.89E0 , 746.9E0 , 748.43E0 , 747.35E0 , 749.27E0 , 647.61E0 , 747.78E0 , 750.51E0 , 851.37E0 , 845.97E0 , 847.54E0 , 849.93E0 , 851.61E0 , 849.75E0 , 850.98E0 , 848.23E0};
  698. // http://www.itl.nist.gov/div898/strd/nls/data/hahn1.shtml
  699. void testNistHahn1(void)
  700. {
  701. const int n=7;
  702. int info;
  703. VectorXd x(n);
  704. /*
  705. * First try
  706. */
  707. x<< 10., -1., .05, -.00001, -.05, .001, -.000001;
  708. // do the computation
  709. hahn1_functor functor;
  710. LevenbergMarquardt<hahn1_functor> lm(functor);
  711. info = lm.minimize(x);
  712. // check return value
  713. VERIFY_IS_EQUAL(info, 1);
  714. VERIFY_IS_EQUAL(lm.nfev, 11);
  715. VERIFY_IS_EQUAL(lm.njev, 10);
  716. // check norm^2
  717. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.5324382854E+00);
  718. // check x
  719. VERIFY_IS_APPROX(x[0], 1.0776351733E+00);
  720. VERIFY_IS_APPROX(x[1],-1.2269296921E-01);
  721. VERIFY_IS_APPROX(x[2], 4.0863750610E-03);
  722. VERIFY_IS_APPROX(x[3],-1.426264e-06); // shoulde be : -1.4262662514E-06
  723. VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
  724. VERIFY_IS_APPROX(x[5], 2.4053735503E-04);
  725. VERIFY_IS_APPROX(x[6],-1.2314450199E-07);
  726. /*
  727. * Second try
  728. */
  729. x<< .1, -.1, .005, -.000001, -.005, .0001, -.0000001;
  730. // do the computation
  731. info = lm.minimize(x);
  732. // check return value
  733. VERIFY_IS_EQUAL(info, 1);
  734. VERIFY_IS_EQUAL(lm.nfev, 11);
  735. VERIFY_IS_EQUAL(lm.njev, 10);
  736. // check norm^2
  737. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.5324382854E+00);
  738. // check x
  739. VERIFY_IS_APPROX(x[0], 1.077640); // should be : 1.0776351733E+00
  740. VERIFY_IS_APPROX(x[1], -0.1226933); // should be : -1.2269296921E-01
  741. VERIFY_IS_APPROX(x[2], 0.004086383); // should be : 4.0863750610E-03
  742. VERIFY_IS_APPROX(x[3], -1.426277e-06); // shoulde be : -1.4262662514E-06
  743. VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
  744. VERIFY_IS_APPROX(x[5], 0.00024053772); // should be : 2.4053735503E-04
  745. VERIFY_IS_APPROX(x[6], -1.231450e-07); // should be : -1.2314450199E-07
  746. }
  747. struct misra1d_functor : Functor<double>
  748. {
  749. misra1d_functor(void) : Functor<double>(2,14) {}
  750. static const double x[14];
  751. static const double y[14];
  752. int operator()(const VectorXd &b, VectorXd &fvec)
  753. {
  754. assert(b.size()==2);
  755. assert(fvec.size()==14);
  756. for(int i=0; i<14; i++) {
  757. fvec[i] = b[0]*b[1]*x[i]/(1.+b[1]*x[i]) - y[i];
  758. }
  759. return 0;
  760. }
  761. int df(const VectorXd &b, MatrixXd &fjac)
  762. {
  763. assert(b.size()==2);
  764. assert(fjac.rows()==14);
  765. assert(fjac.cols()==2);
  766. for(int i=0; i<14; i++) {
  767. double den = 1.+b[1]*x[i];
  768. fjac(i,0) = b[1]*x[i] / den;
  769. fjac(i,1) = b[0]*x[i]*(den-b[1]*x[i])/den/den;
  770. }
  771. return 0;
  772. }
  773. };
  774. const double misra1d_functor::x[14] = { 77.6E0, 114.9E0, 141.1E0, 190.8E0, 239.9E0, 289.0E0, 332.8E0, 378.4E0, 434.8E0, 477.3E0, 536.8E0, 593.1E0, 689.1E0, 760.0E0};
  775. const double misra1d_functor::y[14] = { 10.07E0, 14.73E0, 17.94E0, 23.93E0, 29.61E0, 35.18E0, 40.02E0, 44.82E0, 50.76E0, 55.05E0, 61.01E0, 66.40E0, 75.47E0, 81.78E0};
  776. // http://www.itl.nist.gov/div898/strd/nls/data/misra1d.shtml
  777. void testNistMisra1d(void)
  778. {
  779. const int n=2;
  780. int info;
  781. VectorXd x(n);
  782. /*
  783. * First try
  784. */
  785. x<< 500., 0.0001;
  786. // do the computation
  787. misra1d_functor functor;
  788. LevenbergMarquardt<misra1d_functor> lm(functor);
  789. info = lm.minimize(x);
  790. // check return value
  791. VERIFY_IS_EQUAL(info, 3);
  792. VERIFY_IS_EQUAL(lm.nfev, 9);
  793. VERIFY_IS_EQUAL(lm.njev, 7);
  794. // check norm^2
  795. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6419295283E-02);
  796. // check x
  797. VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
  798. VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
  799. /*
  800. * Second try
  801. */
  802. x<< 450., 0.0003;
  803. // do the computation
  804. info = lm.minimize(x);
  805. // check return value
  806. VERIFY_IS_EQUAL(info, 1);
  807. VERIFY_IS_EQUAL(lm.nfev, 4);
  808. VERIFY_IS_EQUAL(lm.njev, 3);
  809. // check norm^2
  810. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6419295283E-02);
  811. // check x
  812. VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
  813. VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
  814. }
  815. struct lanczos1_functor : Functor<double>
  816. {
  817. lanczos1_functor(void) : Functor<double>(6,24) {}
  818. static const double x[24];
  819. static const double y[24];
  820. int operator()(const VectorXd &b, VectorXd &fvec)
  821. {
  822. assert(b.size()==6);
  823. assert(fvec.size()==24);
  824. for(int i=0; i<24; i++)
  825. fvec[i] = b[0]*exp(-b[1]*x[i]) + b[2]*exp(-b[3]*x[i]) + b[4]*exp(-b[5]*x[i]) - y[i];
  826. return 0;
  827. }
  828. int df(const VectorXd &b, MatrixXd &fjac)
  829. {
  830. assert(b.size()==6);
  831. assert(fjac.rows()==24);
  832. assert(fjac.cols()==6);
  833. for(int i=0; i<24; i++) {
  834. fjac(i,0) = exp(-b[1]*x[i]);
  835. fjac(i,1) = -b[0]*x[i]*exp(-b[1]*x[i]);
  836. fjac(i,2) = exp(-b[3]*x[i]);
  837. fjac(i,3) = -b[2]*x[i]*exp(-b[3]*x[i]);
  838. fjac(i,4) = exp(-b[5]*x[i]);
  839. fjac(i,5) = -b[4]*x[i]*exp(-b[5]*x[i]);
  840. }
  841. return 0;
  842. }
  843. };
  844. const double lanczos1_functor::x[24] = { 0.000000000000E+00, 5.000000000000E-02, 1.000000000000E-01, 1.500000000000E-01, 2.000000000000E-01, 2.500000000000E-01, 3.000000000000E-01, 3.500000000000E-01, 4.000000000000E-01, 4.500000000000E-01, 5.000000000000E-01, 5.500000000000E-01, 6.000000000000E-01, 6.500000000000E-01, 7.000000000000E-01, 7.500000000000E-01, 8.000000000000E-01, 8.500000000000E-01, 9.000000000000E-01, 9.500000000000E-01, 1.000000000000E+00, 1.050000000000E+00, 1.100000000000E+00, 1.150000000000E+00 };
  845. const double lanczos1_functor::y[24] = { 2.513400000000E+00 ,2.044333373291E+00 ,1.668404436564E+00 ,1.366418021208E+00 ,1.123232487372E+00 ,9.268897180037E-01 ,7.679338563728E-01 ,6.388775523106E-01 ,5.337835317402E-01 ,4.479363617347E-01 ,3.775847884350E-01 ,3.197393199326E-01 ,2.720130773746E-01 ,2.324965529032E-01 ,1.996589546065E-01 ,1.722704126914E-01 ,1.493405660168E-01 ,1.300700206922E-01 ,1.138119324644E-01 ,1.000415587559E-01 ,8.833209084540E-02 ,7.833544019350E-02 ,6.976693743449E-02 ,6.239312536719E-02 };
  846. // http://www.itl.nist.gov/div898/strd/nls/data/lanczos1.shtml
  847. void testNistLanczos1(void)
  848. {
  849. const int n=6;
  850. int info;
  851. VectorXd x(n);
  852. /*
  853. * First try
  854. */
  855. x<< 1.2, 0.3, 5.6, 5.5, 6.5, 7.6;
  856. // do the computation
  857. lanczos1_functor functor;
  858. LevenbergMarquardt<lanczos1_functor> lm(functor);
  859. info = lm.minimize(x);
  860. // check return value
  861. VERIFY_IS_EQUAL(info, 2);
  862. VERIFY_IS_EQUAL(lm.nfev, 79);
  863. VERIFY_IS_EQUAL(lm.njev, 72);
  864. // check norm^2
  865. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.430899764097e-25); // should be 1.4307867721E-25, but nist results are on 128-bit floats
  866. // check x
  867. VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
  868. VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
  869. VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
  870. VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
  871. VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
  872. VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
  873. /*
  874. * Second try
  875. */
  876. x<< 0.5, 0.7, 3.6, 4.2, 4., 6.3;
  877. // do the computation
  878. info = lm.minimize(x);
  879. // check return value
  880. VERIFY_IS_EQUAL(info, 2);
  881. VERIFY_IS_EQUAL(lm.nfev, 9);
  882. VERIFY_IS_EQUAL(lm.njev, 8);
  883. // check norm^2
  884. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.428595533845e-25); // should be 1.4307867721E-25, but nist results are on 128-bit floats
  885. // check x
  886. VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
  887. VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
  888. VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
  889. VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
  890. VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
  891. VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
  892. }
  893. struct rat42_functor : Functor<double>
  894. {
  895. rat42_functor(void) : Functor<double>(3,9) {}
  896. static const double x[9];
  897. static const double y[9];
  898. int operator()(const VectorXd &b, VectorXd &fvec)
  899. {
  900. assert(b.size()==3);
  901. assert(fvec.size()==9);
  902. for(int i=0; i<9; i++) {
  903. fvec[i] = b[0] / (1.+exp(b[1]-b[2]*x[i])) - y[i];
  904. }
  905. return 0;
  906. }
  907. int df(const VectorXd &b, MatrixXd &fjac)
  908. {
  909. assert(b.size()==3);
  910. assert(fjac.rows()==9);
  911. assert(fjac.cols()==3);
  912. for(int i=0; i<9; i++) {
  913. double e = exp(b[1]-b[2]*x[i]);
  914. fjac(i,0) = 1./(1.+e);
  915. fjac(i,1) = -b[0]*e/(1.+e)/(1.+e);
  916. fjac(i,2) = +b[0]*e*x[i]/(1.+e)/(1.+e);
  917. }
  918. return 0;
  919. }
  920. };
  921. const double rat42_functor::x[9] = { 9.000E0, 14.000E0, 21.000E0, 28.000E0, 42.000E0, 57.000E0, 63.000E0, 70.000E0, 79.000E0 };
  922. const double rat42_functor::y[9] = { 8.930E0 ,10.800E0 ,18.590E0 ,22.330E0 ,39.350E0 ,56.110E0 ,61.730E0 ,64.620E0 ,67.080E0 };
  923. // http://www.itl.nist.gov/div898/strd/nls/data/ratkowsky2.shtml
  924. void testNistRat42(void)
  925. {
  926. const int n=3;
  927. int info;
  928. VectorXd x(n);
  929. /*
  930. * First try
  931. */
  932. x<< 100., 1., 0.1;
  933. // do the computation
  934. rat42_functor functor;
  935. LevenbergMarquardt<rat42_functor> lm(functor);
  936. info = lm.minimize(x);
  937. // check return value
  938. VERIFY_IS_EQUAL(info, 1);
  939. VERIFY_IS_EQUAL(lm.nfev, 10);
  940. VERIFY_IS_EQUAL(lm.njev, 8);
  941. // check norm^2
  942. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.0565229338E+00);
  943. // check x
  944. VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
  945. VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
  946. VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
  947. /*
  948. * Second try
  949. */
  950. x<< 75., 2.5, 0.07;
  951. // do the computation
  952. info = lm.minimize(x);
  953. // check return value
  954. VERIFY_IS_EQUAL(info, 1);
  955. VERIFY_IS_EQUAL(lm.nfev, 6);
  956. VERIFY_IS_EQUAL(lm.njev, 5);
  957. // check norm^2
  958. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.0565229338E+00);
  959. // check x
  960. VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
  961. VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
  962. VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
  963. }
  964. struct MGH10_functor : Functor<double>
  965. {
  966. MGH10_functor(void) : Functor<double>(3,16) {}
  967. static const double x[16];
  968. static const double y[16];
  969. int operator()(const VectorXd &b, VectorXd &fvec)
  970. {
  971. assert(b.size()==3);
  972. assert(fvec.size()==16);
  973. for(int i=0; i<16; i++)
  974. fvec[i] = b[0] * exp(b[1]/(x[i]+b[2])) - y[i];
  975. return 0;
  976. }
  977. int df(const VectorXd &b, MatrixXd &fjac)
  978. {
  979. assert(b.size()==3);
  980. assert(fjac.rows()==16);
  981. assert(fjac.cols()==3);
  982. for(int i=0; i<16; i++) {
  983. double factor = 1./(x[i]+b[2]);
  984. double e = exp(b[1]*factor);
  985. fjac(i,0) = e;
  986. fjac(i,1) = b[0]*factor*e;
  987. fjac(i,2) = -b[1]*b[0]*factor*factor*e;
  988. }
  989. return 0;
  990. }
  991. };
  992. const double MGH10_functor::x[16] = { 5.000000E+01, 5.500000E+01, 6.000000E+01, 6.500000E+01, 7.000000E+01, 7.500000E+01, 8.000000E+01, 8.500000E+01, 9.000000E+01, 9.500000E+01, 1.000000E+02, 1.050000E+02, 1.100000E+02, 1.150000E+02, 1.200000E+02, 1.250000E+02 };
  993. const double MGH10_functor::y[16] = { 3.478000E+04, 2.861000E+04, 2.365000E+04, 1.963000E+04, 1.637000E+04, 1.372000E+04, 1.154000E+04, 9.744000E+03, 8.261000E+03, 7.030000E+03, 6.005000E+03, 5.147000E+03, 4.427000E+03, 3.820000E+03, 3.307000E+03, 2.872000E+03 };
  994. // http://www.itl.nist.gov/div898/strd/nls/data/mgh10.shtml
  995. void testNistMGH10(void)
  996. {
  997. const int n=3;
  998. int info;
  999. VectorXd x(n);
  1000. /*
  1001. * First try
  1002. */
  1003. x<< 2., 400000., 25000.;
  1004. // do the computation
  1005. MGH10_functor functor;
  1006. LevenbergMarquardt<MGH10_functor> lm(functor);
  1007. info = lm.minimize(x);
  1008. // check return value
  1009. VERIFY_IS_EQUAL(info, 2);
  1010. VERIFY_IS_EQUAL(lm.nfev, 284 );
  1011. VERIFY_IS_EQUAL(lm.njev, 249 );
  1012. // check norm^2
  1013. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7945855171E+01);
  1014. // check x
  1015. VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
  1016. VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
  1017. VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
  1018. /*
  1019. * Second try
  1020. */
  1021. x<< 0.02, 4000., 250.;
  1022. // do the computation
  1023. info = lm.minimize(x);
  1024. // check return value
  1025. VERIFY_IS_EQUAL(info, 3);
  1026. VERIFY_IS_EQUAL(lm.nfev, 126);
  1027. VERIFY_IS_EQUAL(lm.njev, 116);
  1028. // check norm^2
  1029. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7945855171E+01);
  1030. // check x
  1031. VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
  1032. VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
  1033. VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
  1034. }
  1035. struct BoxBOD_functor : Functor<double>
  1036. {
  1037. BoxBOD_functor(void) : Functor<double>(2,6) {}
  1038. static const double x[6];
  1039. int operator()(const VectorXd &b, VectorXd &fvec)
  1040. {
  1041. static const double y[6] = { 109., 149., 149., 191., 213., 224. };
  1042. assert(b.size()==2);
  1043. assert(fvec.size()==6);
  1044. for(int i=0; i<6; i++)
  1045. fvec[i] = b[0]*(1.-exp(-b[1]*x[i])) - y[i];
  1046. return 0;
  1047. }
  1048. int df(const VectorXd &b, MatrixXd &fjac)
  1049. {
  1050. assert(b.size()==2);
  1051. assert(fjac.rows()==6);
  1052. assert(fjac.cols()==2);
  1053. for(int i=0; i<6; i++) {
  1054. double e = exp(-b[1]*x[i]);
  1055. fjac(i,0) = 1.-e;
  1056. fjac(i,1) = b[0]*x[i]*e;
  1057. }
  1058. return 0;
  1059. }
  1060. };
  1061. const double BoxBOD_functor::x[6] = { 1., 2., 3., 5., 7., 10. };
  1062. // http://www.itl.nist.gov/div898/strd/nls/data/boxbod.shtml
  1063. void testNistBoxBOD(void)
  1064. {
  1065. const int n=2;
  1066. int info;
  1067. VectorXd x(n);
  1068. /*
  1069. * First try
  1070. */
  1071. x<< 1., 1.;
  1072. // do the computation
  1073. BoxBOD_functor functor;
  1074. LevenbergMarquardt<BoxBOD_functor> lm(functor);
  1075. lm.parameters.ftol = 1.E6*NumTraits<double>::epsilon();
  1076. lm.parameters.xtol = 1.E6*NumTraits<double>::epsilon();
  1077. lm.parameters.factor = 10.;
  1078. info = lm.minimize(x);
  1079. // check return value
  1080. VERIFY_IS_EQUAL(info, 1);
  1081. VERIFY_IS_EQUAL(lm.nfev, 31);
  1082. VERIFY_IS_EQUAL(lm.njev, 25);
  1083. // check norm^2
  1084. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.1680088766E+03);
  1085. // check x
  1086. VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
  1087. VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
  1088. /*
  1089. * Second try
  1090. */
  1091. x<< 100., 0.75;
  1092. // do the computation
  1093. lm.resetParameters();
  1094. lm.parameters.ftol = NumTraits<double>::epsilon();
  1095. lm.parameters.xtol = NumTraits<double>::epsilon();
  1096. info = lm.minimize(x);
  1097. // check return value
  1098. VERIFY_IS_EQUAL(info, 1);
  1099. VERIFY_IS_EQUAL(lm.nfev, 15 );
  1100. VERIFY_IS_EQUAL(lm.njev, 14 );
  1101. // check norm^2
  1102. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.1680088766E+03);
  1103. // check x
  1104. VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
  1105. VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
  1106. }
  1107. struct MGH17_functor : Functor<double>
  1108. {
  1109. MGH17_functor(void) : Functor<double>(5,33) {}
  1110. static const double x[33];
  1111. static const double y[33];
  1112. int operator()(const VectorXd &b, VectorXd &fvec)
  1113. {
  1114. assert(b.size()==5);
  1115. assert(fvec.size()==33);
  1116. for(int i=0; i<33; i++)
  1117. fvec[i] = b[0] + b[1]*exp(-b[3]*x[i]) + b[2]*exp(-b[4]*x[i]) - y[i];
  1118. return 0;
  1119. }
  1120. int df(const VectorXd &b, MatrixXd &fjac)
  1121. {
  1122. assert(b.size()==5);
  1123. assert(fjac.rows()==33);
  1124. assert(fjac.cols()==5);
  1125. for(int i=0; i<33; i++) {
  1126. fjac(i,0) = 1.;
  1127. fjac(i,1) = exp(-b[3]*x[i]);
  1128. fjac(i,2) = exp(-b[4]*x[i]);
  1129. fjac(i,3) = -x[i]*b[1]*exp(-b[3]*x[i]);
  1130. fjac(i,4) = -x[i]*b[2]*exp(-b[4]*x[i]);
  1131. }
  1132. return 0;
  1133. }
  1134. };
  1135. const double MGH17_functor::x[33] = { 0.000000E+00, 1.000000E+01, 2.000000E+01, 3.000000E+01, 4.000000E+01, 5.000000E+01, 6.000000E+01, 7.000000E+01, 8.000000E+01, 9.000000E+01, 1.000000E+02, 1.100000E+02, 1.200000E+02, 1.300000E+02, 1.400000E+02, 1.500000E+02, 1.600000E+02, 1.700000E+02, 1.800000E+02, 1.900000E+02, 2.000000E+02, 2.100000E+02, 2.200000E+02, 2.300000E+02, 2.400000E+02, 2.500000E+02, 2.600000E+02, 2.700000E+02, 2.800000E+02, 2.900000E+02, 3.000000E+02, 3.100000E+02, 3.200000E+02 };
  1136. const double MGH17_functor::y[33] = { 8.440000E-01, 9.080000E-01, 9.320000E-01, 9.360000E-01, 9.250000E-01, 9.080000E-01, 8.810000E-01, 8.500000E-01, 8.180000E-01, 7.840000E-01, 7.510000E-01, 7.180000E-01, 6.850000E-01, 6.580000E-01, 6.280000E-01, 6.030000E-01, 5.800000E-01, 5.580000E-01, 5.380000E-01, 5.220000E-01, 5.060000E-01, 4.900000E-01, 4.780000E-01, 4.670000E-01, 4.570000E-01, 4.480000E-01, 4.380000E-01, 4.310000E-01, 4.240000E-01, 4.200000E-01, 4.140000E-01, 4.110000E-01, 4.060000E-01 };
  1137. // http://www.itl.nist.gov/div898/strd/nls/data/mgh17.shtml
  1138. void testNistMGH17(void)
  1139. {
  1140. const int n=5;
  1141. int info;
  1142. VectorXd x(n);
  1143. /*
  1144. * First try
  1145. */
  1146. x<< 50., 150., -100., 1., 2.;
  1147. // do the computation
  1148. MGH17_functor functor;
  1149. LevenbergMarquardt<MGH17_functor> lm(functor);
  1150. lm.parameters.ftol = NumTraits<double>::epsilon();
  1151. lm.parameters.xtol = NumTraits<double>::epsilon();
  1152. lm.parameters.maxfev = 1000;
  1153. info = lm.minimize(x);
  1154. // check return value
  1155. VERIFY_IS_EQUAL(info, 2);
  1156. VERIFY_IS_EQUAL(lm.nfev, 602 );
  1157. VERIFY_IS_EQUAL(lm.njev, 545 );
  1158. // check norm^2
  1159. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.4648946975E-05);
  1160. // check x
  1161. VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
  1162. VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
  1163. VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
  1164. VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
  1165. VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
  1166. /*
  1167. * Second try
  1168. */
  1169. x<< 0.5 ,1.5 ,-1 ,0.01 ,0.02;
  1170. // do the computation
  1171. lm.resetParameters();
  1172. info = lm.minimize(x);
  1173. // check return value
  1174. VERIFY_IS_EQUAL(info, 1);
  1175. VERIFY_IS_EQUAL(lm.nfev, 18);
  1176. VERIFY_IS_EQUAL(lm.njev, 15);
  1177. // check norm^2
  1178. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.4648946975E-05);
  1179. // check x
  1180. VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
  1181. VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
  1182. VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
  1183. VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
  1184. VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
  1185. }
  1186. struct MGH09_functor : Functor<double>
  1187. {
  1188. MGH09_functor(void) : Functor<double>(4,11) {}
  1189. static const double _x[11];
  1190. static const double y[11];
  1191. int operator()(const VectorXd &b, VectorXd &fvec)
  1192. {
  1193. assert(b.size()==4);
  1194. assert(fvec.size()==11);
  1195. for(int i=0; i<11; i++) {
  1196. double x = _x[i], xx=x*x;
  1197. fvec[i] = b[0]*(xx+x*b[1])/(xx+x*b[2]+b[3]) - y[i];
  1198. }
  1199. return 0;
  1200. }
  1201. int df(const VectorXd &b, MatrixXd &fjac)
  1202. {
  1203. assert(b.size()==4);
  1204. assert(fjac.rows()==11);
  1205. assert(fjac.cols()==4);
  1206. for(int i=0; i<11; i++) {
  1207. double x = _x[i], xx=x*x;
  1208. double factor = 1./(xx+x*b[2]+b[3]);
  1209. fjac(i,0) = (xx+x*b[1]) * factor;
  1210. fjac(i,1) = b[0]*x* factor;
  1211. fjac(i,2) = - b[0]*(xx+x*b[1]) * x * factor * factor;
  1212. fjac(i,3) = - b[0]*(xx+x*b[1]) * factor * factor;
  1213. }
  1214. return 0;
  1215. }
  1216. };
  1217. const double MGH09_functor::_x[11] = { 4., 2., 1., 5.E-1 , 2.5E-01, 1.670000E-01, 1.250000E-01, 1.E-01, 8.330000E-02, 7.140000E-02, 6.250000E-02 };
  1218. const double MGH09_functor::y[11] = { 1.957000E-01, 1.947000E-01, 1.735000E-01, 1.600000E-01, 8.440000E-02, 6.270000E-02, 4.560000E-02, 3.420000E-02, 3.230000E-02, 2.350000E-02, 2.460000E-02 };
  1219. // http://www.itl.nist.gov/div898/strd/nls/data/mgh09.shtml
  1220. void testNistMGH09(void)
  1221. {
  1222. const int n=4;
  1223. int info;
  1224. VectorXd x(n);
  1225. /*
  1226. * First try
  1227. */
  1228. x<< 25., 39, 41.5, 39.;
  1229. // do the computation
  1230. MGH09_functor functor;
  1231. LevenbergMarquardt<MGH09_functor> lm(functor);
  1232. lm.parameters.maxfev = 1000;
  1233. info = lm.minimize(x);
  1234. // check return value
  1235. VERIFY_IS_EQUAL(info, 1);
  1236. VERIFY_IS_EQUAL(lm.nfev, 490 );
  1237. VERIFY_IS_EQUAL(lm.njev, 376 );
  1238. // check norm^2
  1239. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 3.0750560385E-04);
  1240. // check x
  1241. VERIFY_IS_APPROX(x[0], 0.1928077089); // should be 1.9280693458E-01
  1242. VERIFY_IS_APPROX(x[1], 0.19126423573); // should be 1.9128232873E-01
  1243. VERIFY_IS_APPROX(x[2], 0.12305309914); // should be 1.2305650693E-01
  1244. VERIFY_IS_APPROX(x[3], 0.13605395375); // should be 1.3606233068E-01
  1245. /*
  1246. * Second try
  1247. */
  1248. x<< 0.25, 0.39, 0.415, 0.39;
  1249. // do the computation
  1250. lm.resetParameters();
  1251. info = lm.minimize(x);
  1252. // check return value
  1253. VERIFY_IS_EQUAL(info, 1);
  1254. VERIFY_IS_EQUAL(lm.nfev, 18);
  1255. VERIFY_IS_EQUAL(lm.njev, 16);
  1256. // check norm^2
  1257. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 3.0750560385E-04);
  1258. // check x
  1259. VERIFY_IS_APPROX(x[0], 0.19280781); // should be 1.9280693458E-01
  1260. VERIFY_IS_APPROX(x[1], 0.19126265); // should be 1.9128232873E-01
  1261. VERIFY_IS_APPROX(x[2], 0.12305280); // should be 1.2305650693E-01
  1262. VERIFY_IS_APPROX(x[3], 0.13605322); // should be 1.3606233068E-01
  1263. }
  1264. struct Bennett5_functor : Functor<double>
  1265. {
  1266. Bennett5_functor(void) : Functor<double>(3,154) {}
  1267. static const double x[154];
  1268. static const double y[154];
  1269. int operator()(const VectorXd &b, VectorXd &fvec)
  1270. {
  1271. assert(b.size()==3);
  1272. assert(fvec.size()==154);
  1273. for(int i=0; i<154; i++)
  1274. fvec[i] = b[0]* pow(b[1]+x[i],-1./b[2]) - y[i];
  1275. return 0;
  1276. }
  1277. int df(const VectorXd &b, MatrixXd &fjac)
  1278. {
  1279. assert(b.size()==3);
  1280. assert(fjac.rows()==154);
  1281. assert(fjac.cols()==3);
  1282. for(int i=0; i<154; i++) {
  1283. double e = pow(b[1]+x[i],-1./b[2]);
  1284. fjac(i,0) = e;
  1285. fjac(i,1) = - b[0]*e/b[2]/(b[1]+x[i]);
  1286. fjac(i,2) = b[0]*e*log(b[1]+x[i])/b[2]/b[2];
  1287. }
  1288. return 0;
  1289. }
  1290. };
  1291. const double Bennett5_functor::x[154] = { 7.447168E0, 8.102586E0, 8.452547E0, 8.711278E0, 8.916774E0, 9.087155E0, 9.232590E0, 9.359535E0, 9.472166E0, 9.573384E0, 9.665293E0, 9.749461E0, 9.827092E0, 9.899128E0, 9.966321E0, 10.029280E0, 10.088510E0, 10.144430E0, 10.197380E0, 10.247670E0, 10.295560E0, 10.341250E0, 10.384950E0, 10.426820E0, 10.467000E0, 10.505640E0, 10.542830E0, 10.578690E0, 10.613310E0, 10.646780E0, 10.679150E0, 10.710520E0, 10.740920E0, 10.770440E0, 10.799100E0, 10.826970E0, 10.854080E0, 10.880470E0, 10.906190E0, 10.931260E0, 10.955720E0, 10.979590E0, 11.002910E0, 11.025700E0, 11.047980E0, 11.069770E0, 11.091100E0, 11.111980E0, 11.132440E0, 11.152480E0, 11.172130E0, 11.191410E0, 11.210310E0, 11.228870E0, 11.247090E0, 11.264980E0, 11.282560E0, 11.299840E0, 11.316820E0, 11.333520E0, 11.349940E0, 11.366100E0, 11.382000E0, 11.397660E0, 11.413070E0, 11.428240E0, 11.443200E0, 11.457930E0, 11.472440E0, 11.486750E0, 11.500860E0, 11.514770E0, 11.528490E0, 11.542020E0, 11.555380E0, 11.568550E0, 11.581560E0, 11.594420E0, 11.607121E0, 11.619640E0, 11.632000E0, 11.644210E0, 11.656280E0, 11.668200E0, 11.679980E0, 11.691620E0, 11.703130E0, 11.714510E0, 11.725760E0, 11.736880E0, 11.747890E0, 11.758780E0, 11.769550E0, 11.780200E0, 11.790730E0, 11.801160E0, 11.811480E0, 11.821700E0, 11.831810E0, 11.841820E0, 11.851730E0, 11.861550E0, 11.871270E0, 11.880890E0, 11.890420E0, 11.899870E0, 11.909220E0, 11.918490E0, 11.927680E0, 11.936780E0, 11.945790E0, 11.954730E0, 11.963590E0, 11.972370E0, 11.981070E0, 11.989700E0, 11.998260E0, 12.006740E0, 12.015150E0, 12.023490E0, 12.031760E0, 12.039970E0, 12.048100E0, 12.056170E0, 12.064180E0, 12.072120E0, 12.080010E0, 12.087820E0, 12.095580E0, 12.103280E0, 12.110920E0, 12.118500E0, 12.126030E0, 12.133500E0, 12.140910E0, 12.148270E0, 12.155570E0, 12.162830E0, 12.170030E0, 12.177170E0, 12.184270E0, 12.191320E0, 12.198320E0, 12.205270E0, 12.212170E0, 12.219030E0, 12.225840E0, 12.232600E0, 12.239320E0, 12.245990E0, 12.252620E0, 12.259200E0, 12.265750E0, 12.272240E0 };
  1292. const double Bennett5_functor::y[154] = { -34.834702E0 ,-34.393200E0 ,-34.152901E0 ,-33.979099E0 ,-33.845901E0 ,-33.732899E0 ,-33.640301E0 ,-33.559200E0 ,-33.486801E0 ,-33.423100E0 ,-33.365101E0 ,-33.313000E0 ,-33.260899E0 ,-33.217400E0 ,-33.176899E0 ,-33.139198E0 ,-33.101601E0 ,-33.066799E0 ,-33.035000E0 ,-33.003101E0 ,-32.971298E0 ,-32.942299E0 ,-32.916302E0 ,-32.890202E0 ,-32.864101E0 ,-32.841000E0 ,-32.817799E0 ,-32.797501E0 ,-32.774300E0 ,-32.757000E0 ,-32.733799E0 ,-32.716400E0 ,-32.699100E0 ,-32.678799E0 ,-32.661400E0 ,-32.644001E0 ,-32.626701E0 ,-32.612202E0 ,-32.597698E0 ,-32.583199E0 ,-32.568699E0 ,-32.554298E0 ,-32.539799E0 ,-32.525299E0 ,-32.510799E0 ,-32.499199E0 ,-32.487598E0 ,-32.473202E0 ,-32.461601E0 ,-32.435501E0 ,-32.435501E0 ,-32.426800E0 ,-32.412300E0 ,-32.400799E0 ,-32.392101E0 ,-32.380501E0 ,-32.366001E0 ,-32.357300E0 ,-32.348598E0 ,-32.339901E0 ,-32.328400E0 ,-32.319698E0 ,-32.311001E0 ,-32.299400E0 ,-32.290699E0 ,-32.282001E0 ,-32.273300E0 ,-32.264599E0 ,-32.256001E0 ,-32.247299E0 ,-32.238602E0 ,-32.229900E0 ,-32.224098E0 ,-32.215401E0 ,-32.203800E0 ,-32.198002E0 ,-32.189400E0 ,-32.183601E0 ,-32.174900E0 ,-32.169102E0 ,-32.163300E0 ,-32.154598E0 ,-32.145901E0 ,-32.140099E0 ,-32.131401E0 ,-32.125599E0 ,-32.119801E0 ,-32.111198E0 ,-32.105400E0 ,-32.096699E0 ,-32.090900E0 ,-32.088001E0 ,-32.079300E0 ,-32.073502E0 ,-32.067699E0 ,-32.061901E0 ,-32.056099E0 ,-32.050301E0 ,-32.044498E0 ,-32.038799E0 ,-32.033001E0 ,-32.027199E0 ,-32.024300E0 ,-32.018501E0 ,-32.012699E0 ,-32.004002E0 ,-32.001099E0 ,-31.995300E0 ,-31.989500E0 ,-31.983700E0 ,-31.977900E0 ,-31.972099E0 ,-31.969299E0 ,-31.963501E0 ,-31.957701E0 ,-31.951900E0 ,-31.946100E0 ,-31.940300E0 ,-31.937401E0 ,-31.931601E0 ,-31.925800E0 ,-31.922899E0 ,-31.917101E0 ,-31.911301E0 ,-31.908400E0 ,-31.902599E0 ,-31.896900E0 ,-31.893999E0 ,-31.888201E0 ,-31.885300E0 ,-31.882401E0 ,-31.876600E0 ,-31.873699E0 ,-31.867901E0 ,-31.862101E0 ,-31.859200E0 ,-31.856300E0 ,-31.850500E0 ,-31.844700E0 ,-31.841801E0 ,-31.838900E0 ,-31.833099E0 ,-31.830200E0 ,-31.827299E0 ,-31.821600E0 ,-31.818701E0 ,-31.812901E0 ,-31.809999E0 ,-31.807100E0 ,-31.801300E0 ,-31.798401E0 ,-31.795500E0 ,-31.789700E0 ,-31.786800E0 };
  1293. // http://www.itl.nist.gov/div898/strd/nls/data/bennett5.shtml
  1294. void testNistBennett5(void)
  1295. {
  1296. const int n=3;
  1297. int info;
  1298. VectorXd x(n);
  1299. /*
  1300. * First try
  1301. */
  1302. x<< -2000., 50., 0.8;
  1303. // do the computation
  1304. Bennett5_functor functor;
  1305. LevenbergMarquardt<Bennett5_functor> lm(functor);
  1306. lm.parameters.maxfev = 1000;
  1307. info = lm.minimize(x);
  1308. // check return value
  1309. VERIFY_IS_EQUAL(info, 1);
  1310. VERIFY_IS_EQUAL(lm.nfev, 758);
  1311. VERIFY_IS_EQUAL(lm.njev, 744);
  1312. // check norm^2
  1313. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.2404744073E-04);
  1314. // check x
  1315. VERIFY_IS_APPROX(x[0], -2.5235058043E+03);
  1316. VERIFY_IS_APPROX(x[1], 4.6736564644E+01);
  1317. VERIFY_IS_APPROX(x[2], 9.3218483193E-01);
  1318. /*
  1319. * Second try
  1320. */
  1321. x<< -1500., 45., 0.85;
  1322. // do the computation
  1323. lm.resetParameters();
  1324. info = lm.minimize(x);
  1325. // check return value
  1326. VERIFY_IS_EQUAL(info, 1);
  1327. VERIFY_IS_EQUAL(lm.nfev, 203);
  1328. VERIFY_IS_EQUAL(lm.njev, 192);
  1329. // check norm^2
  1330. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.2404744073E-04);
  1331. // check x
  1332. VERIFY_IS_APPROX(x[0], -2523.3007865); // should be -2.5235058043E+03
  1333. VERIFY_IS_APPROX(x[1], 46.735705771); // should be 4.6736564644E+01);
  1334. VERIFY_IS_APPROX(x[2], 0.93219881891); // should be 9.3218483193E-01);
  1335. }
  1336. struct thurber_functor : Functor<double>
  1337. {
  1338. thurber_functor(void) : Functor<double>(7,37) {}
  1339. static const double _x[37];
  1340. static const double _y[37];
  1341. int operator()(const VectorXd &b, VectorXd &fvec)
  1342. {
  1343. // int called=0; printf("call hahn1_functor with iflag=%d, called=%d\n", iflag, called); if (iflag==1) called++;
  1344. assert(b.size()==7);
  1345. assert(fvec.size()==37);
  1346. for(int i=0; i<37; i++) {
  1347. double x=_x[i], xx=x*x, xxx=xx*x;
  1348. fvec[i] = (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) / (1.+b[4]*x+b[5]*xx+b[6]*xxx) - _y[i];
  1349. }
  1350. return 0;
  1351. }
  1352. int df(const VectorXd &b, MatrixXd &fjac)
  1353. {
  1354. assert(b.size()==7);
  1355. assert(fjac.rows()==37);
  1356. assert(fjac.cols()==7);
  1357. for(int i=0; i<37; i++) {
  1358. double x=_x[i], xx=x*x, xxx=xx*x;
  1359. double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
  1360. fjac(i,0) = 1.*fact;
  1361. fjac(i,1) = x*fact;
  1362. fjac(i,2) = xx*fact;
  1363. fjac(i,3) = xxx*fact;
  1364. fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
  1365. fjac(i,4) = x*fact;
  1366. fjac(i,5) = xx*fact;
  1367. fjac(i,6) = xxx*fact;
  1368. }
  1369. return 0;
  1370. }
  1371. };
  1372. const double thurber_functor::_x[37] = { -3.067E0, -2.981E0, -2.921E0, -2.912E0, -2.840E0, -2.797E0, -2.702E0, -2.699E0, -2.633E0, -2.481E0, -2.363E0, -2.322E0, -1.501E0, -1.460E0, -1.274E0, -1.212E0, -1.100E0, -1.046E0, -0.915E0, -0.714E0, -0.566E0, -0.545E0, -0.400E0, -0.309E0, -0.109E0, -0.103E0, 0.010E0, 0.119E0, 0.377E0, 0.790E0, 0.963E0, 1.006E0, 1.115E0, 1.572E0, 1.841E0, 2.047E0, 2.200E0 };
  1373. const double thurber_functor::_y[37] = { 80.574E0, 84.248E0, 87.264E0, 87.195E0, 89.076E0, 89.608E0, 89.868E0, 90.101E0, 92.405E0, 95.854E0, 100.696E0, 101.060E0, 401.672E0, 390.724E0, 567.534E0, 635.316E0, 733.054E0, 759.087E0, 894.206E0, 990.785E0, 1090.109E0, 1080.914E0, 1122.643E0, 1178.351E0, 1260.531E0, 1273.514E0, 1288.339E0, 1327.543E0, 1353.863E0, 1414.509E0, 1425.208E0, 1421.384E0, 1442.962E0, 1464.350E0, 1468.705E0, 1447.894E0, 1457.628E0};
  1374. // http://www.itl.nist.gov/div898/strd/nls/data/thurber.shtml
  1375. void testNistThurber(void)
  1376. {
  1377. const int n=7;
  1378. int info;
  1379. VectorXd x(n);
  1380. /*
  1381. * First try
  1382. */
  1383. x<< 1000 ,1000 ,400 ,40 ,0.7,0.3,0.0 ;
  1384. // do the computation
  1385. thurber_functor functor;
  1386. LevenbergMarquardt<thurber_functor> lm(functor);
  1387. lm.parameters.ftol = 1.E4*NumTraits<double>::epsilon();
  1388. lm.parameters.xtol = 1.E4*NumTraits<double>::epsilon();
  1389. info = lm.minimize(x);
  1390. // check return value
  1391. VERIFY_IS_EQUAL(info, 1);
  1392. VERIFY_IS_EQUAL(lm.nfev, 39);
  1393. VERIFY_IS_EQUAL(lm.njev, 36);
  1394. // check norm^2
  1395. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6427082397E+03);
  1396. // check x
  1397. VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
  1398. VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
  1399. VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
  1400. VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
  1401. VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
  1402. VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
  1403. VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
  1404. /*
  1405. * Second try
  1406. */
  1407. x<< 1300 ,1500 ,500 ,75 ,1 ,0.4 ,0.05 ;
  1408. // do the computation
  1409. lm.resetParameters();
  1410. lm.parameters.ftol = 1.E4*NumTraits<double>::epsilon();
  1411. lm.parameters.xtol = 1.E4*NumTraits<double>::epsilon();
  1412. info = lm.minimize(x);
  1413. // check return value
  1414. VERIFY_IS_EQUAL(info, 1);
  1415. VERIFY_IS_EQUAL(lm.nfev, 29);
  1416. VERIFY_IS_EQUAL(lm.njev, 28);
  1417. // check norm^2
  1418. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6427082397E+03);
  1419. // check x
  1420. VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
  1421. VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
  1422. VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
  1423. VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
  1424. VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
  1425. VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
  1426. VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
  1427. }
  1428. struct rat43_functor : Functor<double>
  1429. {
  1430. rat43_functor(void) : Functor<double>(4,15) {}
  1431. static const double x[15];
  1432. static const double y[15];
  1433. int operator()(const VectorXd &b, VectorXd &fvec)
  1434. {
  1435. assert(b.size()==4);
  1436. assert(fvec.size()==15);
  1437. for(int i=0; i<15; i++)
  1438. fvec[i] = b[0] * pow(1.+exp(b[1]-b[2]*x[i]),-1./b[3]) - y[i];
  1439. return 0;
  1440. }
  1441. int df(const VectorXd &b, MatrixXd &fjac)
  1442. {
  1443. assert(b.size()==4);
  1444. assert(fjac.rows()==15);
  1445. assert(fjac.cols()==4);
  1446. for(int i=0; i<15; i++) {
  1447. double e = exp(b[1]-b[2]*x[i]);
  1448. double power = -1./b[3];
  1449. fjac(i,0) = pow(1.+e, power);
  1450. fjac(i,1) = power*b[0]*e*pow(1.+e, power-1.);
  1451. fjac(i,2) = -power*b[0]*e*x[i]*pow(1.+e, power-1.);
  1452. fjac(i,3) = b[0]*power*power*log(1.+e)*pow(1.+e, power);
  1453. }
  1454. return 0;
  1455. }
  1456. };
  1457. const double rat43_functor::x[15] = { 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14., 15. };
  1458. const double rat43_functor::y[15] = { 16.08, 33.83, 65.80, 97.20, 191.55, 326.20, 386.87, 520.53, 590.03, 651.92, 724.93, 699.56, 689.96, 637.56, 717.41 };
  1459. // http://www.itl.nist.gov/div898/strd/nls/data/ratkowsky3.shtml
  1460. void testNistRat43(void)
  1461. {
  1462. const int n=4;
  1463. int info;
  1464. VectorXd x(n);
  1465. /*
  1466. * First try
  1467. */
  1468. x<< 100., 10., 1., 1.;
  1469. // do the computation
  1470. rat43_functor functor;
  1471. LevenbergMarquardt<rat43_functor> lm(functor);
  1472. lm.parameters.ftol = 1.E6*NumTraits<double>::epsilon();
  1473. lm.parameters.xtol = 1.E6*NumTraits<double>::epsilon();
  1474. info = lm.minimize(x);
  1475. // check return value
  1476. VERIFY_IS_EQUAL(info, 1);
  1477. VERIFY_IS_EQUAL(lm.nfev, 27);
  1478. VERIFY_IS_EQUAL(lm.njev, 20);
  1479. // check norm^2
  1480. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7864049080E+03);
  1481. // check x
  1482. VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
  1483. VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
  1484. VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
  1485. VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
  1486. /*
  1487. * Second try
  1488. */
  1489. x<< 700., 5., 0.75, 1.3;
  1490. // do the computation
  1491. lm.resetParameters();
  1492. lm.parameters.ftol = 1.E5*NumTraits<double>::epsilon();
  1493. lm.parameters.xtol = 1.E5*NumTraits<double>::epsilon();
  1494. info = lm.minimize(x);
  1495. // check return value
  1496. VERIFY_IS_EQUAL(info, 1);
  1497. VERIFY_IS_EQUAL(lm.nfev, 9);
  1498. VERIFY_IS_EQUAL(lm.njev, 8);
  1499. // check norm^2
  1500. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7864049080E+03);
  1501. // check x
  1502. VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
  1503. VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
  1504. VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
  1505. VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
  1506. }
  1507. struct eckerle4_functor : Functor<double>
  1508. {
  1509. eckerle4_functor(void) : Functor<double>(3,35) {}
  1510. static const double x[35];
  1511. static const double y[35];
  1512. int operator()(const VectorXd &b, VectorXd &fvec)
  1513. {
  1514. assert(b.size()==3);
  1515. assert(fvec.size()==35);
  1516. for(int i=0; i<35; i++)
  1517. fvec[i] = b[0]/b[1] * exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/(b[1]*b[1])) - y[i];
  1518. return 0;
  1519. }
  1520. int df(const VectorXd &b, MatrixXd &fjac)
  1521. {
  1522. assert(b.size()==3);
  1523. assert(fjac.rows()==35);
  1524. assert(fjac.cols()==3);
  1525. for(int i=0; i<35; i++) {
  1526. double b12 = b[1]*b[1];
  1527. double e = exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/b12);
  1528. fjac(i,0) = e / b[1];
  1529. fjac(i,1) = ((x[i]-b[2])*(x[i]-b[2])/b12-1.) * b[0]*e/b12;
  1530. fjac(i,2) = (x[i]-b[2])*e*b[0]/b[1]/b12;
  1531. }
  1532. return 0;
  1533. }
  1534. };
  1535. const double eckerle4_functor::x[35] = { 400.0, 405.0, 410.0, 415.0, 420.0, 425.0, 430.0, 435.0, 436.5, 438.0, 439.5, 441.0, 442.5, 444.0, 445.5, 447.0, 448.5, 450.0, 451.5, 453.0, 454.5, 456.0, 457.5, 459.0, 460.5, 462.0, 463.5, 465.0, 470.0, 475.0, 480.0, 485.0, 490.0, 495.0, 500.0};
  1536. const double eckerle4_functor::y[35] = { 0.0001575, 0.0001699, 0.0002350, 0.0003102, 0.0004917, 0.0008710, 0.0017418, 0.0046400, 0.0065895, 0.0097302, 0.0149002, 0.0237310, 0.0401683, 0.0712559, 0.1264458, 0.2073413, 0.2902366, 0.3445623, 0.3698049, 0.3668534, 0.3106727, 0.2078154, 0.1164354, 0.0616764, 0.0337200, 0.0194023, 0.0117831, 0.0074357, 0.0022732, 0.0008800, 0.0004579, 0.0002345, 0.0001586, 0.0001143, 0.0000710 };
  1537. // http://www.itl.nist.gov/div898/strd/nls/data/eckerle4.shtml
  1538. void testNistEckerle4(void)
  1539. {
  1540. const int n=3;
  1541. int info;
  1542. VectorXd x(n);
  1543. /*
  1544. * First try
  1545. */
  1546. x<< 1., 10., 500.;
  1547. // do the computation
  1548. eckerle4_functor functor;
  1549. LevenbergMarquardt<eckerle4_functor> lm(functor);
  1550. info = lm.minimize(x);
  1551. // check return value
  1552. VERIFY_IS_EQUAL(info, 1);
  1553. VERIFY_IS_EQUAL(lm.nfev, 18);
  1554. VERIFY_IS_EQUAL(lm.njev, 15);
  1555. // check norm^2
  1556. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.4635887487E-03);
  1557. // check x
  1558. VERIFY_IS_APPROX(x[0], 1.5543827178);
  1559. VERIFY_IS_APPROX(x[1], 4.0888321754);
  1560. VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
  1561. /*
  1562. * Second try
  1563. */
  1564. x<< 1.5, 5., 450.;
  1565. // do the computation
  1566. info = lm.minimize(x);
  1567. // check return value
  1568. VERIFY_IS_EQUAL(info, 1);
  1569. VERIFY_IS_EQUAL(lm.nfev, 7);
  1570. VERIFY_IS_EQUAL(lm.njev, 6);
  1571. // check norm^2
  1572. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.4635887487E-03);
  1573. // check x
  1574. VERIFY_IS_APPROX(x[0], 1.5543827178);
  1575. VERIFY_IS_APPROX(x[1], 4.0888321754);
  1576. VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
  1577. }
  1578. void test_NonLinearOptimization()
  1579. {
  1580. // Tests using the examples provided by (c)minpack
  1581. CALL_SUBTEST/*_1*/(testChkder());
  1582. CALL_SUBTEST/*_1*/(testLmder1());
  1583. CALL_SUBTEST/*_1*/(testLmder());
  1584. CALL_SUBTEST/*_2*/(testHybrj1());
  1585. CALL_SUBTEST/*_2*/(testHybrj());
  1586. CALL_SUBTEST/*_2*/(testHybrd1());
  1587. CALL_SUBTEST/*_2*/(testHybrd());
  1588. CALL_SUBTEST/*_3*/(testLmstr1());
  1589. CALL_SUBTEST/*_3*/(testLmstr());
  1590. CALL_SUBTEST/*_3*/(testLmdif1());
  1591. CALL_SUBTEST/*_3*/(testLmdif());
  1592. // NIST tests, level of difficulty = "Lower"
  1593. CALL_SUBTEST/*_4*/(testNistMisra1a());
  1594. CALL_SUBTEST/*_4*/(testNistChwirut2());
  1595. // NIST tests, level of difficulty = "Average"
  1596. CALL_SUBTEST/*_5*/(testNistHahn1());
  1597. CALL_SUBTEST/*_6*/(testNistMisra1d());
  1598. CALL_SUBTEST/*_7*/(testNistMGH17());
  1599. CALL_SUBTEST/*_8*/(testNistLanczos1());
  1600. // NIST tests, level of difficulty = "Higher"
  1601. CALL_SUBTEST/*_9*/(testNistRat42());
  1602. CALL_SUBTEST/*_10*/(testNistMGH10());
  1603. CALL_SUBTEST/*_11*/(testNistBoxBOD());
  1604. CALL_SUBTEST/*_12*/(testNistMGH09());
  1605. CALL_SUBTEST/*_13*/(testNistBennett5());
  1606. CALL_SUBTEST/*_14*/(testNistThurber());
  1607. CALL_SUBTEST/*_15*/(testNistRat43());
  1608. CALL_SUBTEST/*_16*/(testNistEckerle4());
  1609. }
  1610. /*
  1611. * Can be useful for debugging...
  1612. printf("info, nfev : %d, %d\n", info, lm.nfev);
  1613. printf("info, nfev, njev : %d, %d, %d\n", info, solver.nfev, solver.njev);
  1614. printf("info, nfev : %d, %d\n", info, solver.nfev);
  1615. printf("x[0] : %.32g\n", x[0]);
  1616. printf("x[1] : %.32g\n", x[1]);
  1617. printf("x[2] : %.32g\n", x[2]);
  1618. printf("x[3] : %.32g\n", x[3]);
  1619. printf("fvec.blueNorm() : %.32g\n", solver.fvec.blueNorm());
  1620. printf("fvec.blueNorm() : %.32g\n", lm.fvec.blueNorm());
  1621. printf("info, nfev, njev : %d, %d, %d\n", info, lm.nfev, lm.njev);
  1622. printf("fvec.squaredNorm() : %.13g\n", lm.fvec.squaredNorm());
  1623. std::cout << x << std::endl;
  1624. std::cout.precision(9);
  1625. std::cout << x[0] << std::endl;
  1626. std::cout << x[1] << std::endl;
  1627. std::cout << x[2] << std::endl;
  1628. std::cout << x[3] << std::endl;
  1629. */