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.

1872 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/StormEigen/NonLinearOptimization>
  8. // This disables some useless Warnings on MSVC.
  9. // It is intended to be done for this test only.
  10. #include <StormEigen/src/Core/util/DisableStupidWarnings.h>
  11. using std::sqrt;
  12. int fcn_chkder(const VectorXd &x, VectorXd &fvec, MatrixXd &fjac, int iflag)
  13. {
  14. /* subroutine fcn for chkder example. */
  15. int i;
  16. assert(15 == fvec.size());
  17. assert(3 == x.size());
  18. double tmp1, tmp2, tmp3, tmp4;
  19. 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,
  20. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  21. if (iflag == 0)
  22. return 0;
  23. if (iflag != 2)
  24. for (i=0; i<15; i++) {
  25. tmp1 = i+1;
  26. tmp2 = 16-i-1;
  27. tmp3 = tmp1;
  28. if (i >= 8) tmp3 = tmp2;
  29. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  30. }
  31. else {
  32. for (i = 0; i < 15; i++) {
  33. tmp1 = i+1;
  34. tmp2 = 16-i-1;
  35. /* error introduced into next statement for illustration. */
  36. /* corrected statement should read tmp3 = tmp1 . */
  37. tmp3 = tmp2;
  38. if (i >= 8) tmp3 = tmp2;
  39. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4=tmp4*tmp4;
  40. fjac(i,0) = -1.;
  41. fjac(i,1) = tmp1*tmp2/tmp4;
  42. fjac(i,2) = tmp1*tmp3/tmp4;
  43. }
  44. }
  45. return 0;
  46. }
  47. void testChkder()
  48. {
  49. const int m=15, n=3;
  50. VectorXd x(n), fvec(m), xp, fvecp(m), err;
  51. MatrixXd fjac(m,n);
  52. VectorXi ipvt;
  53. /* the following values should be suitable for */
  54. /* checking the jacobian matrix. */
  55. x << 9.2e-1, 1.3e-1, 5.4e-1;
  56. internal::chkder(x, fvec, fjac, xp, fvecp, 1, err);
  57. fcn_chkder(x, fvec, fjac, 1);
  58. fcn_chkder(x, fvec, fjac, 2);
  59. fcn_chkder(xp, fvecp, fjac, 1);
  60. internal::chkder(x, fvec, fjac, xp, fvecp, 2, err);
  61. fvecp -= fvec;
  62. // check those
  63. VectorXd fvec_ref(m), fvecp_ref(m), err_ref(m);
  64. fvec_ref <<
  65. -1.181606, -1.429655, -1.606344,
  66. -1.745269, -1.840654, -1.921586,
  67. -1.984141, -2.022537, -2.468977,
  68. -2.827562, -3.473582, -4.437612,
  69. -6.047662, -9.267761, -18.91806;
  70. fvecp_ref <<
  71. -7.724666e-09, -3.432406e-09, -2.034843e-10,
  72. 2.313685e-09, 4.331078e-09, 5.984096e-09,
  73. 7.363281e-09, 8.53147e-09, 1.488591e-08,
  74. 2.33585e-08, 3.522012e-08, 5.301255e-08,
  75. 8.26666e-08, 1.419747e-07, 3.19899e-07;
  76. err_ref <<
  77. 0.1141397, 0.09943516, 0.09674474,
  78. 0.09980447, 0.1073116, 0.1220445,
  79. 0.1526814, 1, 1,
  80. 1, 1, 1,
  81. 1, 1, 1;
  82. VERIFY_IS_APPROX(fvec, fvec_ref);
  83. VERIFY_IS_APPROX(fvecp, fvecp_ref);
  84. VERIFY_IS_APPROX(err, err_ref);
  85. }
  86. // Generic functor
  87. template<typename _Scalar, int NX=Dynamic, int NY=Dynamic>
  88. struct Functor
  89. {
  90. typedef _Scalar Scalar;
  91. enum {
  92. InputsAtCompileTime = NX,
  93. ValuesAtCompileTime = NY
  94. };
  95. typedef Matrix<Scalar,InputsAtCompileTime,1> InputType;
  96. typedef Matrix<Scalar,ValuesAtCompileTime,1> ValueType;
  97. typedef Matrix<Scalar,ValuesAtCompileTime,InputsAtCompileTime> JacobianType;
  98. const int m_inputs, m_values;
  99. Functor() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
  100. Functor(int inputs, int values) : m_inputs(inputs), m_values(values) {}
  101. int inputs() const { return m_inputs; }
  102. int values() const { return m_values; }
  103. // you should define that in the subclass :
  104. // void operator() (const InputType& x, ValueType* v, JacobianType* _j=0) const;
  105. };
  106. struct lmder_functor : Functor<double>
  107. {
  108. lmder_functor(void): Functor<double>(3,15) {}
  109. int operator()(const VectorXd &x, VectorXd &fvec) const
  110. {
  111. double tmp1, tmp2, tmp3;
  112. 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,
  113. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  114. for (int i = 0; i < values(); i++)
  115. {
  116. tmp1 = i+1;
  117. tmp2 = 16 - i - 1;
  118. tmp3 = (i>=8)? tmp2 : tmp1;
  119. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  120. }
  121. return 0;
  122. }
  123. int df(const VectorXd &x, MatrixXd &fjac) const
  124. {
  125. double tmp1, tmp2, tmp3, tmp4;
  126. for (int i = 0; i < values(); i++)
  127. {
  128. tmp1 = i+1;
  129. tmp2 = 16 - i - 1;
  130. tmp3 = (i>=8)? tmp2 : tmp1;
  131. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4 = tmp4*tmp4;
  132. fjac(i,0) = -1;
  133. fjac(i,1) = tmp1*tmp2/tmp4;
  134. fjac(i,2) = tmp1*tmp3/tmp4;
  135. }
  136. return 0;
  137. }
  138. };
  139. void testLmder1()
  140. {
  141. int n=3, info;
  142. VectorXd x;
  143. /* the following starting values provide a rough fit. */
  144. x.setConstant(n, 1.);
  145. // do the computation
  146. lmder_functor functor;
  147. LevenbergMarquardt<lmder_functor> lm(functor);
  148. info = lm.lmder1(x);
  149. // check return value
  150. VERIFY_IS_EQUAL(info, 1);
  151. VERIFY_IS_EQUAL(lm.nfev, 6);
  152. VERIFY_IS_EQUAL(lm.njev, 5);
  153. // check norm
  154. VERIFY_IS_APPROX(lm.fvec.blueNorm(), 0.09063596);
  155. // check x
  156. VectorXd x_ref(n);
  157. x_ref << 0.08241058, 1.133037, 2.343695;
  158. VERIFY_IS_APPROX(x, x_ref);
  159. }
  160. void testLmder()
  161. {
  162. const int m=15, n=3;
  163. int info;
  164. double fnorm, covfac;
  165. VectorXd x;
  166. /* the following starting values provide a rough fit. */
  167. x.setConstant(n, 1.);
  168. // do the computation
  169. lmder_functor functor;
  170. LevenbergMarquardt<lmder_functor> lm(functor);
  171. info = lm.minimize(x);
  172. // check return values
  173. VERIFY_IS_EQUAL(info, 1);
  174. VERIFY_IS_EQUAL(lm.nfev, 6);
  175. VERIFY_IS_EQUAL(lm.njev, 5);
  176. // check norm
  177. fnorm = lm.fvec.blueNorm();
  178. VERIFY_IS_APPROX(fnorm, 0.09063596);
  179. // check x
  180. VectorXd x_ref(n);
  181. x_ref << 0.08241058, 1.133037, 2.343695;
  182. VERIFY_IS_APPROX(x, x_ref);
  183. // check covariance
  184. covfac = fnorm*fnorm/(m-n);
  185. internal::covar(lm.fjac, lm.permutation.indices()); // TODO : move this as a function of lm
  186. MatrixXd cov_ref(n,n);
  187. cov_ref <<
  188. 0.0001531202, 0.002869941, -0.002656662,
  189. 0.002869941, 0.09480935, -0.09098995,
  190. -0.002656662, -0.09098995, 0.08778727;
  191. // std::cout << fjac*covfac << std::endl;
  192. MatrixXd cov;
  193. cov = covfac*lm.fjac.topLeftCorner<n,n>();
  194. VERIFY_IS_APPROX( cov, cov_ref);
  195. // TODO: why isn't this allowed ? :
  196. // VERIFY_IS_APPROX( covfac*fjac.topLeftCorner<n,n>() , cov_ref);
  197. }
  198. struct hybrj_functor : Functor<double>
  199. {
  200. hybrj_functor(void) : Functor<double>(9,9) {}
  201. int operator()(const VectorXd &x, VectorXd &fvec)
  202. {
  203. double temp, temp1, temp2;
  204. const VectorXd::Index n = x.size();
  205. assert(fvec.size()==n);
  206. for (VectorXd::Index k = 0; k < n; k++)
  207. {
  208. temp = (3. - 2.*x[k])*x[k];
  209. temp1 = 0.;
  210. if (k) temp1 = x[k-1];
  211. temp2 = 0.;
  212. if (k != n-1) temp2 = x[k+1];
  213. fvec[k] = temp - temp1 - 2.*temp2 + 1.;
  214. }
  215. return 0;
  216. }
  217. int df(const VectorXd &x, MatrixXd &fjac)
  218. {
  219. const VectorXd::Index n = x.size();
  220. assert(fjac.rows()==n);
  221. assert(fjac.cols()==n);
  222. for (VectorXd::Index k = 0; k < n; k++)
  223. {
  224. for (VectorXd::Index j = 0; j < n; j++)
  225. fjac(k,j) = 0.;
  226. fjac(k,k) = 3.- 4.*x[k];
  227. if (k) fjac(k,k-1) = -1.;
  228. if (k != n-1) fjac(k,k+1) = -2.;
  229. }
  230. return 0;
  231. }
  232. };
  233. void testHybrj1()
  234. {
  235. const int n=9;
  236. int info;
  237. VectorXd x(n);
  238. /* the following starting values provide a rough fit. */
  239. x.setConstant(n, -1.);
  240. // do the computation
  241. hybrj_functor functor;
  242. HybridNonLinearSolver<hybrj_functor> solver(functor);
  243. info = solver.hybrj1(x);
  244. // check return value
  245. VERIFY_IS_EQUAL(info, 1);
  246. VERIFY_IS_EQUAL(solver.nfev, 11);
  247. VERIFY_IS_EQUAL(solver.njev, 1);
  248. // check norm
  249. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  250. // check x
  251. VectorXd x_ref(n);
  252. x_ref <<
  253. -0.5706545, -0.6816283, -0.7017325,
  254. -0.7042129, -0.701369, -0.6918656,
  255. -0.665792, -0.5960342, -0.4164121;
  256. VERIFY_IS_APPROX(x, x_ref);
  257. }
  258. void testHybrj()
  259. {
  260. const int n=9;
  261. int info;
  262. VectorXd x(n);
  263. /* the following starting values provide a rough fit. */
  264. x.setConstant(n, -1.);
  265. // do the computation
  266. hybrj_functor functor;
  267. HybridNonLinearSolver<hybrj_functor> solver(functor);
  268. solver.diag.setConstant(n, 1.);
  269. solver.useExternalScaling = true;
  270. info = solver.solve(x);
  271. // check return value
  272. VERIFY_IS_EQUAL(info, 1);
  273. VERIFY_IS_EQUAL(solver.nfev, 11);
  274. VERIFY_IS_EQUAL(solver.njev, 1);
  275. // check norm
  276. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  277. // check x
  278. VectorXd x_ref(n);
  279. x_ref <<
  280. -0.5706545, -0.6816283, -0.7017325,
  281. -0.7042129, -0.701369, -0.6918656,
  282. -0.665792, -0.5960342, -0.4164121;
  283. VERIFY_IS_APPROX(x, x_ref);
  284. }
  285. struct hybrd_functor : Functor<double>
  286. {
  287. hybrd_functor(void) : Functor<double>(9,9) {}
  288. int operator()(const VectorXd &x, VectorXd &fvec) const
  289. {
  290. double temp, temp1, temp2;
  291. const VectorXd::Index n = x.size();
  292. assert(fvec.size()==n);
  293. for (VectorXd::Index k=0; k < n; k++)
  294. {
  295. temp = (3. - 2.*x[k])*x[k];
  296. temp1 = 0.;
  297. if (k) temp1 = x[k-1];
  298. temp2 = 0.;
  299. if (k != n-1) temp2 = x[k+1];
  300. fvec[k] = temp - temp1 - 2.*temp2 + 1.;
  301. }
  302. return 0;
  303. }
  304. };
  305. void testHybrd1()
  306. {
  307. int n=9, info;
  308. VectorXd x(n);
  309. /* the following starting values provide a rough solution. */
  310. x.setConstant(n, -1.);
  311. // do the computation
  312. hybrd_functor functor;
  313. HybridNonLinearSolver<hybrd_functor> solver(functor);
  314. info = solver.hybrd1(x);
  315. // check return value
  316. VERIFY_IS_EQUAL(info, 1);
  317. VERIFY_IS_EQUAL(solver.nfev, 20);
  318. // check norm
  319. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  320. // check x
  321. VectorXd x_ref(n);
  322. x_ref << -0.5706545, -0.6816283, -0.7017325, -0.7042129, -0.701369, -0.6918656, -0.665792, -0.5960342, -0.4164121;
  323. VERIFY_IS_APPROX(x, x_ref);
  324. }
  325. void testHybrd()
  326. {
  327. const int n=9;
  328. int info;
  329. VectorXd x;
  330. /* the following starting values provide a rough fit. */
  331. x.setConstant(n, -1.);
  332. // do the computation
  333. hybrd_functor functor;
  334. HybridNonLinearSolver<hybrd_functor> solver(functor);
  335. solver.parameters.nb_of_subdiagonals = 1;
  336. solver.parameters.nb_of_superdiagonals = 1;
  337. solver.diag.setConstant(n, 1.);
  338. solver.useExternalScaling = true;
  339. info = solver.solveNumericalDiff(x);
  340. // check return value
  341. VERIFY_IS_EQUAL(info, 1);
  342. VERIFY_IS_EQUAL(solver.nfev, 14);
  343. // check norm
  344. VERIFY_IS_APPROX(solver.fvec.blueNorm(), 1.192636e-08);
  345. // check x
  346. VectorXd x_ref(n);
  347. x_ref <<
  348. -0.5706545, -0.6816283, -0.7017325,
  349. -0.7042129, -0.701369, -0.6918656,
  350. -0.665792, -0.5960342, -0.4164121;
  351. VERIFY_IS_APPROX(x, x_ref);
  352. }
  353. struct lmstr_functor : Functor<double>
  354. {
  355. lmstr_functor(void) : Functor<double>(3,15) {}
  356. int operator()(const VectorXd &x, VectorXd &fvec)
  357. {
  358. /* subroutine fcn for lmstr1 example. */
  359. double tmp1, tmp2, tmp3;
  360. 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,
  361. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  362. assert(15==fvec.size());
  363. assert(3==x.size());
  364. for (int i=0; i<15; i++)
  365. {
  366. tmp1 = i+1;
  367. tmp2 = 16 - i - 1;
  368. tmp3 = (i>=8)? tmp2 : tmp1;
  369. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  370. }
  371. return 0;
  372. }
  373. int df(const VectorXd &x, VectorXd &jac_row, VectorXd::Index rownb)
  374. {
  375. assert(x.size()==3);
  376. assert(jac_row.size()==x.size());
  377. double tmp1, tmp2, tmp3, tmp4;
  378. VectorXd::Index i = rownb-2;
  379. tmp1 = i+1;
  380. tmp2 = 16 - i - 1;
  381. tmp3 = (i>=8)? tmp2 : tmp1;
  382. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4 = tmp4*tmp4;
  383. jac_row[0] = -1;
  384. jac_row[1] = tmp1*tmp2/tmp4;
  385. jac_row[2] = tmp1*tmp3/tmp4;
  386. return 0;
  387. }
  388. };
  389. void testLmstr1()
  390. {
  391. const int n=3;
  392. int info;
  393. VectorXd x(n);
  394. /* the following starting values provide a rough fit. */
  395. x.setConstant(n, 1.);
  396. // do the computation
  397. lmstr_functor functor;
  398. LevenbergMarquardt<lmstr_functor> lm(functor);
  399. info = lm.lmstr1(x);
  400. // check return value
  401. VERIFY_IS_EQUAL(info, 1);
  402. VERIFY_IS_EQUAL(lm.nfev, 6);
  403. VERIFY_IS_EQUAL(lm.njev, 5);
  404. // check norm
  405. VERIFY_IS_APPROX(lm.fvec.blueNorm(), 0.09063596);
  406. // check x
  407. VectorXd x_ref(n);
  408. x_ref << 0.08241058, 1.133037, 2.343695 ;
  409. VERIFY_IS_APPROX(x, x_ref);
  410. }
  411. void testLmstr()
  412. {
  413. const int n=3;
  414. int info;
  415. double fnorm;
  416. VectorXd x(n);
  417. /* the following starting values provide a rough fit. */
  418. x.setConstant(n, 1.);
  419. // do the computation
  420. lmstr_functor functor;
  421. LevenbergMarquardt<lmstr_functor> lm(functor);
  422. info = lm.minimizeOptimumStorage(x);
  423. // check return values
  424. VERIFY_IS_EQUAL(info, 1);
  425. VERIFY_IS_EQUAL(lm.nfev, 6);
  426. VERIFY_IS_EQUAL(lm.njev, 5);
  427. // check norm
  428. fnorm = lm.fvec.blueNorm();
  429. VERIFY_IS_APPROX(fnorm, 0.09063596);
  430. // check x
  431. VectorXd x_ref(n);
  432. x_ref << 0.08241058, 1.133037, 2.343695;
  433. VERIFY_IS_APPROX(x, x_ref);
  434. }
  435. struct lmdif_functor : Functor<double>
  436. {
  437. lmdif_functor(void) : Functor<double>(3,15) {}
  438. int operator()(const VectorXd &x, VectorXd &fvec) const
  439. {
  440. int i;
  441. double tmp1,tmp2,tmp3;
  442. 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,
  443. 3.7e-1,5.8e-1,7.3e-1,9.6e-1,1.34e0,2.1e0,4.39e0};
  444. assert(x.size()==3);
  445. assert(fvec.size()==15);
  446. for (i=0; i<15; i++)
  447. {
  448. tmp1 = i+1;
  449. tmp2 = 15 - i;
  450. tmp3 = tmp1;
  451. if (i >= 8) tmp3 = tmp2;
  452. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  453. }
  454. return 0;
  455. }
  456. };
  457. void testLmdif1()
  458. {
  459. const int n=3;
  460. int info;
  461. VectorXd x(n), fvec(15);
  462. /* the following starting values provide a rough fit. */
  463. x.setConstant(n, 1.);
  464. // do the computation
  465. lmdif_functor functor;
  466. DenseIndex nfev;
  467. info = LevenbergMarquardt<lmdif_functor>::lmdif1(functor, x, &nfev);
  468. // check return value
  469. VERIFY_IS_EQUAL(info, 1);
  470. VERIFY_IS_EQUAL(nfev, 26);
  471. // check norm
  472. functor(x, fvec);
  473. VERIFY_IS_APPROX(fvec.blueNorm(), 0.09063596);
  474. // check x
  475. VectorXd x_ref(n);
  476. x_ref << 0.0824106, 1.1330366, 2.3436947;
  477. VERIFY_IS_APPROX(x, x_ref);
  478. }
  479. void testLmdif()
  480. {
  481. const int m=15, n=3;
  482. int info;
  483. double fnorm, covfac;
  484. VectorXd x(n);
  485. /* the following starting values provide a rough fit. */
  486. x.setConstant(n, 1.);
  487. // do the computation
  488. lmdif_functor functor;
  489. NumericalDiff<lmdif_functor> numDiff(functor);
  490. LevenbergMarquardt<NumericalDiff<lmdif_functor> > lm(numDiff);
  491. info = lm.minimize(x);
  492. // check return values
  493. VERIFY_IS_EQUAL(info, 1);
  494. VERIFY_IS_EQUAL(lm.nfev, 26);
  495. // check norm
  496. fnorm = lm.fvec.blueNorm();
  497. VERIFY_IS_APPROX(fnorm, 0.09063596);
  498. // check x
  499. VectorXd x_ref(n);
  500. x_ref << 0.08241058, 1.133037, 2.343695;
  501. VERIFY_IS_APPROX(x, x_ref);
  502. // check covariance
  503. covfac = fnorm*fnorm/(m-n);
  504. internal::covar(lm.fjac, lm.permutation.indices()); // TODO : move this as a function of lm
  505. MatrixXd cov_ref(n,n);
  506. cov_ref <<
  507. 0.0001531202, 0.002869942, -0.002656662,
  508. 0.002869942, 0.09480937, -0.09098997,
  509. -0.002656662, -0.09098997, 0.08778729;
  510. // std::cout << fjac*covfac << std::endl;
  511. MatrixXd cov;
  512. cov = covfac*lm.fjac.topLeftCorner<n,n>();
  513. VERIFY_IS_APPROX( cov, cov_ref);
  514. // TODO: why isn't this allowed ? :
  515. // VERIFY_IS_APPROX( covfac*fjac.topLeftCorner<n,n>() , cov_ref);
  516. }
  517. struct chwirut2_functor : Functor<double>
  518. {
  519. chwirut2_functor(void) : Functor<double>(3,54) {}
  520. static const double m_x[54];
  521. static const double m_y[54];
  522. int operator()(const VectorXd &b, VectorXd &fvec)
  523. {
  524. int i;
  525. assert(b.size()==3);
  526. assert(fvec.size()==54);
  527. for(i=0; i<54; i++) {
  528. double x = m_x[i];
  529. fvec[i] = exp(-b[0]*x)/(b[1]+b[2]*x) - m_y[i];
  530. }
  531. return 0;
  532. }
  533. int df(const VectorXd &b, MatrixXd &fjac)
  534. {
  535. assert(b.size()==3);
  536. assert(fjac.rows()==54);
  537. assert(fjac.cols()==3);
  538. for(int i=0; i<54; i++) {
  539. double x = m_x[i];
  540. double factor = 1./(b[1]+b[2]*x);
  541. double e = exp(-b[0]*x);
  542. fjac(i,0) = -x*e*factor;
  543. fjac(i,1) = -e*factor*factor;
  544. fjac(i,2) = -x*e*factor*factor;
  545. }
  546. return 0;
  547. }
  548. };
  549. 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};
  550. 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 };
  551. // http://www.itl.nist.gov/div898/strd/nls/data/chwirut2.shtml
  552. void testNistChwirut2(void)
  553. {
  554. const int n=3;
  555. int info;
  556. VectorXd x(n);
  557. /*
  558. * First try
  559. */
  560. x<< 0.1, 0.01, 0.02;
  561. // do the computation
  562. chwirut2_functor functor;
  563. LevenbergMarquardt<chwirut2_functor> lm(functor);
  564. info = lm.minimize(x);
  565. // check return value
  566. VERIFY_IS_EQUAL(info, 1);
  567. VERIFY_IS_EQUAL(lm.nfev, 10);
  568. VERIFY_IS_EQUAL(lm.njev, 8);
  569. // check norm^2
  570. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.1304802941E+02);
  571. // check x
  572. VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
  573. VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
  574. VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
  575. /*
  576. * Second try
  577. */
  578. x<< 0.15, 0.008, 0.010;
  579. // do the computation
  580. lm.resetParameters();
  581. lm.parameters.ftol = 1.E6*NumTraits<double>::epsilon();
  582. lm.parameters.xtol = 1.E6*NumTraits<double>::epsilon();
  583. info = lm.minimize(x);
  584. // check return value
  585. VERIFY_IS_EQUAL(info, 1);
  586. VERIFY_IS_EQUAL(lm.nfev, 7);
  587. VERIFY_IS_EQUAL(lm.njev, 6);
  588. // check norm^2
  589. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.1304802941E+02);
  590. // check x
  591. VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
  592. VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
  593. VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
  594. }
  595. struct misra1a_functor : Functor<double>
  596. {
  597. misra1a_functor(void) : Functor<double>(2,14) {}
  598. static const double m_x[14];
  599. static const double m_y[14];
  600. int operator()(const VectorXd &b, VectorXd &fvec)
  601. {
  602. assert(b.size()==2);
  603. assert(fvec.size()==14);
  604. for(int i=0; i<14; i++) {
  605. fvec[i] = b[0]*(1.-exp(-b[1]*m_x[i])) - m_y[i] ;
  606. }
  607. return 0;
  608. }
  609. int df(const VectorXd &b, MatrixXd &fjac)
  610. {
  611. assert(b.size()==2);
  612. assert(fjac.rows()==14);
  613. assert(fjac.cols()==2);
  614. for(int i=0; i<14; i++) {
  615. fjac(i,0) = (1.-exp(-b[1]*m_x[i]));
  616. fjac(i,1) = (b[0]*m_x[i]*exp(-b[1]*m_x[i]));
  617. }
  618. return 0;
  619. }
  620. };
  621. 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};
  622. 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};
  623. // http://www.itl.nist.gov/div898/strd/nls/data/misra1a.shtml
  624. void testNistMisra1a(void)
  625. {
  626. const int n=2;
  627. int info;
  628. VectorXd x(n);
  629. /*
  630. * First try
  631. */
  632. x<< 500., 0.0001;
  633. // do the computation
  634. misra1a_functor functor;
  635. LevenbergMarquardt<misra1a_functor> lm(functor);
  636. info = lm.minimize(x);
  637. // check return value
  638. VERIFY_IS_EQUAL(info, 1);
  639. VERIFY_IS_EQUAL(lm.nfev, 19);
  640. VERIFY_IS_EQUAL(lm.njev, 15);
  641. // check norm^2
  642. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.2455138894E-01);
  643. // check x
  644. VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
  645. VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
  646. /*
  647. * Second try
  648. */
  649. x<< 250., 0.0005;
  650. // do the computation
  651. info = lm.minimize(x);
  652. // check return value
  653. VERIFY_IS_EQUAL(info, 1);
  654. VERIFY_IS_EQUAL(lm.nfev, 5);
  655. VERIFY_IS_EQUAL(lm.njev, 4);
  656. // check norm^2
  657. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.2455138894E-01);
  658. // check x
  659. VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
  660. VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
  661. }
  662. struct hahn1_functor : Functor<double>
  663. {
  664. hahn1_functor(void) : Functor<double>(7,236) {}
  665. static const double m_x[236];
  666. int operator()(const VectorXd &b, VectorXd &fvec)
  667. {
  668. 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 ,
  669. 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 ,
  670. 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 };
  671. // int called=0; printf("call hahn1_functor with iflag=%d, called=%d\n", iflag, called); if (iflag==1) called++;
  672. assert(b.size()==7);
  673. assert(fvec.size()==236);
  674. for(int i=0; i<236; i++) {
  675. double x=m_x[i], xx=x*x, xxx=xx*x;
  676. 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];
  677. }
  678. return 0;
  679. }
  680. int df(const VectorXd &b, MatrixXd &fjac)
  681. {
  682. assert(b.size()==7);
  683. assert(fjac.rows()==236);
  684. assert(fjac.cols()==7);
  685. for(int i=0; i<236; i++) {
  686. double x=m_x[i], xx=x*x, xxx=xx*x;
  687. double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
  688. fjac(i,0) = 1.*fact;
  689. fjac(i,1) = x*fact;
  690. fjac(i,2) = xx*fact;
  691. fjac(i,3) = xxx*fact;
  692. fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
  693. fjac(i,4) = x*fact;
  694. fjac(i,5) = xx*fact;
  695. fjac(i,6) = xxx*fact;
  696. }
  697. return 0;
  698. }
  699. };
  700. 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 ,
  701. 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 ,
  702. 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};
  703. // http://www.itl.nist.gov/div898/strd/nls/data/hahn1.shtml
  704. void testNistHahn1(void)
  705. {
  706. const int n=7;
  707. int info;
  708. VectorXd x(n);
  709. /*
  710. * First try
  711. */
  712. x<< 10., -1., .05, -.00001, -.05, .001, -.000001;
  713. // do the computation
  714. hahn1_functor functor;
  715. LevenbergMarquardt<hahn1_functor> lm(functor);
  716. info = lm.minimize(x);
  717. // check return value
  718. VERIFY_IS_EQUAL(info, 1);
  719. VERIFY_IS_EQUAL(lm.nfev, 11);
  720. VERIFY_IS_EQUAL(lm.njev, 10);
  721. // check norm^2
  722. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.5324382854E+00);
  723. // check x
  724. VERIFY_IS_APPROX(x[0], 1.0776351733E+00);
  725. VERIFY_IS_APPROX(x[1],-1.2269296921E-01);
  726. VERIFY_IS_APPROX(x[2], 4.0863750610E-03);
  727. VERIFY_IS_APPROX(x[3],-1.426264e-06); // shoulde be : -1.4262662514E-06
  728. VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
  729. VERIFY_IS_APPROX(x[5], 2.4053735503E-04);
  730. VERIFY_IS_APPROX(x[6],-1.2314450199E-07);
  731. /*
  732. * Second try
  733. */
  734. x<< .1, -.1, .005, -.000001, -.005, .0001, -.0000001;
  735. // do the computation
  736. info = lm.minimize(x);
  737. // check return value
  738. VERIFY_IS_EQUAL(info, 1);
  739. VERIFY_IS_EQUAL(lm.nfev, 11);
  740. VERIFY_IS_EQUAL(lm.njev, 10);
  741. // check norm^2
  742. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.5324382854E+00);
  743. // check x
  744. VERIFY_IS_APPROX(x[0], 1.077640); // should be : 1.0776351733E+00
  745. VERIFY_IS_APPROX(x[1], -0.1226933); // should be : -1.2269296921E-01
  746. VERIFY_IS_APPROX(x[2], 0.004086383); // should be : 4.0863750610E-03
  747. VERIFY_IS_APPROX(x[3], -1.426277e-06); // shoulde be : -1.4262662514E-06
  748. VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
  749. VERIFY_IS_APPROX(x[5], 0.00024053772); // should be : 2.4053735503E-04
  750. VERIFY_IS_APPROX(x[6], -1.231450e-07); // should be : -1.2314450199E-07
  751. }
  752. struct misra1d_functor : Functor<double>
  753. {
  754. misra1d_functor(void) : Functor<double>(2,14) {}
  755. static const double x[14];
  756. static const double y[14];
  757. int operator()(const VectorXd &b, VectorXd &fvec)
  758. {
  759. assert(b.size()==2);
  760. assert(fvec.size()==14);
  761. for(int i=0; i<14; i++) {
  762. fvec[i] = b[0]*b[1]*x[i]/(1.+b[1]*x[i]) - y[i];
  763. }
  764. return 0;
  765. }
  766. int df(const VectorXd &b, MatrixXd &fjac)
  767. {
  768. assert(b.size()==2);
  769. assert(fjac.rows()==14);
  770. assert(fjac.cols()==2);
  771. for(int i=0; i<14; i++) {
  772. double den = 1.+b[1]*x[i];
  773. fjac(i,0) = b[1]*x[i] / den;
  774. fjac(i,1) = b[0]*x[i]*(den-b[1]*x[i])/den/den;
  775. }
  776. return 0;
  777. }
  778. };
  779. 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};
  780. 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};
  781. // http://www.itl.nist.gov/div898/strd/nls/data/misra1d.shtml
  782. void testNistMisra1d(void)
  783. {
  784. const int n=2;
  785. int info;
  786. VectorXd x(n);
  787. /*
  788. * First try
  789. */
  790. x<< 500., 0.0001;
  791. // do the computation
  792. misra1d_functor functor;
  793. LevenbergMarquardt<misra1d_functor> lm(functor);
  794. info = lm.minimize(x);
  795. // check return value
  796. VERIFY_IS_EQUAL(info, 3);
  797. VERIFY_IS_EQUAL(lm.nfev, 9);
  798. VERIFY_IS_EQUAL(lm.njev, 7);
  799. // check norm^2
  800. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6419295283E-02);
  801. // check x
  802. VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
  803. VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
  804. /*
  805. * Second try
  806. */
  807. x<< 450., 0.0003;
  808. // do the computation
  809. info = lm.minimize(x);
  810. // check return value
  811. VERIFY_IS_EQUAL(info, 1);
  812. VERIFY_IS_EQUAL(lm.nfev, 4);
  813. VERIFY_IS_EQUAL(lm.njev, 3);
  814. // check norm^2
  815. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6419295283E-02);
  816. // check x
  817. VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
  818. VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
  819. }
  820. struct lanczos1_functor : Functor<double>
  821. {
  822. lanczos1_functor(void) : Functor<double>(6,24) {}
  823. static const double x[24];
  824. static const double y[24];
  825. int operator()(const VectorXd &b, VectorXd &fvec)
  826. {
  827. assert(b.size()==6);
  828. assert(fvec.size()==24);
  829. for(int i=0; i<24; i++)
  830. 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];
  831. return 0;
  832. }
  833. int df(const VectorXd &b, MatrixXd &fjac)
  834. {
  835. assert(b.size()==6);
  836. assert(fjac.rows()==24);
  837. assert(fjac.cols()==6);
  838. for(int i=0; i<24; i++) {
  839. fjac(i,0) = exp(-b[1]*x[i]);
  840. fjac(i,1) = -b[0]*x[i]*exp(-b[1]*x[i]);
  841. fjac(i,2) = exp(-b[3]*x[i]);
  842. fjac(i,3) = -b[2]*x[i]*exp(-b[3]*x[i]);
  843. fjac(i,4) = exp(-b[5]*x[i]);
  844. fjac(i,5) = -b[4]*x[i]*exp(-b[5]*x[i]);
  845. }
  846. return 0;
  847. }
  848. };
  849. 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 };
  850. 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 };
  851. // http://www.itl.nist.gov/div898/strd/nls/data/lanczos1.shtml
  852. void testNistLanczos1(void)
  853. {
  854. const int n=6;
  855. int info;
  856. VectorXd x(n);
  857. /*
  858. * First try
  859. */
  860. x<< 1.2, 0.3, 5.6, 5.5, 6.5, 7.6;
  861. // do the computation
  862. lanczos1_functor functor;
  863. LevenbergMarquardt<lanczos1_functor> lm(functor);
  864. info = lm.minimize(x);
  865. // check return value
  866. VERIFY_IS_EQUAL(info, 2);
  867. VERIFY_IS_EQUAL(lm.nfev, 79);
  868. VERIFY_IS_EQUAL(lm.njev, 72);
  869. // check norm^2
  870. std::cout.precision(30);
  871. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.4290986055242372e-25); // should be 1.4307867721E-25, but nist results are on 128-bit floats
  872. // check x
  873. VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
  874. VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
  875. VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
  876. VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
  877. VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
  878. VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
  879. /*
  880. * Second try
  881. */
  882. x<< 0.5, 0.7, 3.6, 4.2, 4., 6.3;
  883. // do the computation
  884. info = lm.minimize(x);
  885. // check return value
  886. VERIFY_IS_EQUAL(info, 2);
  887. VERIFY_IS_EQUAL(lm.nfev, 9);
  888. VERIFY_IS_EQUAL(lm.njev, 8);
  889. // check norm^2
  890. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.430571737783119393e-25); // should be 1.4307867721E-25, but nist results are on 128-bit floats
  891. // check x
  892. VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
  893. VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
  894. VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
  895. VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
  896. VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
  897. VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
  898. }
  899. struct rat42_functor : Functor<double>
  900. {
  901. rat42_functor(void) : Functor<double>(3,9) {}
  902. static const double x[9];
  903. static const double y[9];
  904. int operator()(const VectorXd &b, VectorXd &fvec)
  905. {
  906. assert(b.size()==3);
  907. assert(fvec.size()==9);
  908. for(int i=0; i<9; i++) {
  909. fvec[i] = b[0] / (1.+exp(b[1]-b[2]*x[i])) - y[i];
  910. }
  911. return 0;
  912. }
  913. int df(const VectorXd &b, MatrixXd &fjac)
  914. {
  915. assert(b.size()==3);
  916. assert(fjac.rows()==9);
  917. assert(fjac.cols()==3);
  918. for(int i=0; i<9; i++) {
  919. double e = exp(b[1]-b[2]*x[i]);
  920. fjac(i,0) = 1./(1.+e);
  921. fjac(i,1) = -b[0]*e/(1.+e)/(1.+e);
  922. fjac(i,2) = +b[0]*e*x[i]/(1.+e)/(1.+e);
  923. }
  924. return 0;
  925. }
  926. };
  927. const double rat42_functor::x[9] = { 9.000E0, 14.000E0, 21.000E0, 28.000E0, 42.000E0, 57.000E0, 63.000E0, 70.000E0, 79.000E0 };
  928. const double rat42_functor::y[9] = { 8.930E0 ,10.800E0 ,18.590E0 ,22.330E0 ,39.350E0 ,56.110E0 ,61.730E0 ,64.620E0 ,67.080E0 };
  929. // http://www.itl.nist.gov/div898/strd/nls/data/ratkowsky2.shtml
  930. void testNistRat42(void)
  931. {
  932. const int n=3;
  933. int info;
  934. VectorXd x(n);
  935. /*
  936. * First try
  937. */
  938. x<< 100., 1., 0.1;
  939. // do the computation
  940. rat42_functor functor;
  941. LevenbergMarquardt<rat42_functor> lm(functor);
  942. info = lm.minimize(x);
  943. // check return value
  944. VERIFY_IS_EQUAL(info, 1);
  945. VERIFY_IS_EQUAL(lm.nfev, 10);
  946. VERIFY_IS_EQUAL(lm.njev, 8);
  947. // check norm^2
  948. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.0565229338E+00);
  949. // check x
  950. VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
  951. VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
  952. VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
  953. /*
  954. * Second try
  955. */
  956. x<< 75., 2.5, 0.07;
  957. // do the computation
  958. info = lm.minimize(x);
  959. // check return value
  960. VERIFY_IS_EQUAL(info, 1);
  961. VERIFY_IS_EQUAL(lm.nfev, 6);
  962. VERIFY_IS_EQUAL(lm.njev, 5);
  963. // check norm^2
  964. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.0565229338E+00);
  965. // check x
  966. VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
  967. VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
  968. VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
  969. }
  970. struct MGH10_functor : Functor<double>
  971. {
  972. MGH10_functor(void) : Functor<double>(3,16) {}
  973. static const double x[16];
  974. static const double y[16];
  975. int operator()(const VectorXd &b, VectorXd &fvec)
  976. {
  977. assert(b.size()==3);
  978. assert(fvec.size()==16);
  979. for(int i=0; i<16; i++)
  980. fvec[i] = b[0] * exp(b[1]/(x[i]+b[2])) - y[i];
  981. return 0;
  982. }
  983. int df(const VectorXd &b, MatrixXd &fjac)
  984. {
  985. assert(b.size()==3);
  986. assert(fjac.rows()==16);
  987. assert(fjac.cols()==3);
  988. for(int i=0; i<16; i++) {
  989. double factor = 1./(x[i]+b[2]);
  990. double e = exp(b[1]*factor);
  991. fjac(i,0) = e;
  992. fjac(i,1) = b[0]*factor*e;
  993. fjac(i,2) = -b[1]*b[0]*factor*factor*e;
  994. }
  995. return 0;
  996. }
  997. };
  998. 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 };
  999. 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 };
  1000. // http://www.itl.nist.gov/div898/strd/nls/data/mgh10.shtml
  1001. void testNistMGH10(void)
  1002. {
  1003. const int n=3;
  1004. int info;
  1005. VectorXd x(n);
  1006. /*
  1007. * First try
  1008. */
  1009. x<< 2., 400000., 25000.;
  1010. // do the computation
  1011. MGH10_functor functor;
  1012. LevenbergMarquardt<MGH10_functor> lm(functor);
  1013. info = lm.minimize(x);
  1014. // check return value
  1015. VERIFY_IS_EQUAL(info, 2);
  1016. VERIFY_IS_EQUAL(lm.nfev, 284 );
  1017. VERIFY_IS_EQUAL(lm.njev, 249 );
  1018. // check norm^2
  1019. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7945855171E+01);
  1020. // check x
  1021. VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
  1022. VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
  1023. VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
  1024. /*
  1025. * Second try
  1026. */
  1027. x<< 0.02, 4000., 250.;
  1028. // do the computation
  1029. info = lm.minimize(x);
  1030. // check return value
  1031. VERIFY_IS_EQUAL(info, 3);
  1032. VERIFY_IS_EQUAL(lm.nfev, 126);
  1033. VERIFY_IS_EQUAL(lm.njev, 116);
  1034. // check norm^2
  1035. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7945855171E+01);
  1036. // check x
  1037. VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
  1038. VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
  1039. VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
  1040. }
  1041. struct BoxBOD_functor : Functor<double>
  1042. {
  1043. BoxBOD_functor(void) : Functor<double>(2,6) {}
  1044. static const double x[6];
  1045. int operator()(const VectorXd &b, VectorXd &fvec)
  1046. {
  1047. static const double y[6] = { 109., 149., 149., 191., 213., 224. };
  1048. assert(b.size()==2);
  1049. assert(fvec.size()==6);
  1050. for(int i=0; i<6; i++)
  1051. fvec[i] = b[0]*(1.-exp(-b[1]*x[i])) - y[i];
  1052. return 0;
  1053. }
  1054. int df(const VectorXd &b, MatrixXd &fjac)
  1055. {
  1056. assert(b.size()==2);
  1057. assert(fjac.rows()==6);
  1058. assert(fjac.cols()==2);
  1059. for(int i=0; i<6; i++) {
  1060. double e = exp(-b[1]*x[i]);
  1061. fjac(i,0) = 1.-e;
  1062. fjac(i,1) = b[0]*x[i]*e;
  1063. }
  1064. return 0;
  1065. }
  1066. };
  1067. const double BoxBOD_functor::x[6] = { 1., 2., 3., 5., 7., 10. };
  1068. // http://www.itl.nist.gov/div898/strd/nls/data/boxbod.shtml
  1069. void testNistBoxBOD(void)
  1070. {
  1071. const int n=2;
  1072. int info;
  1073. VectorXd x(n);
  1074. /*
  1075. * First try
  1076. */
  1077. x<< 1., 1.;
  1078. // do the computation
  1079. BoxBOD_functor functor;
  1080. LevenbergMarquardt<BoxBOD_functor> lm(functor);
  1081. lm.parameters.ftol = 1.E6*NumTraits<double>::epsilon();
  1082. lm.parameters.xtol = 1.E6*NumTraits<double>::epsilon();
  1083. lm.parameters.factor = 10.;
  1084. info = lm.minimize(x);
  1085. // check return value
  1086. VERIFY_IS_EQUAL(info, 1);
  1087. VERIFY(lm.nfev < 31); // 31
  1088. VERIFY(lm.njev < 25); // 25
  1089. // check norm^2
  1090. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.1680088766E+03);
  1091. // check x
  1092. VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
  1093. VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
  1094. /*
  1095. * Second try
  1096. */
  1097. x<< 100., 0.75;
  1098. // do the computation
  1099. lm.resetParameters();
  1100. lm.parameters.ftol = NumTraits<double>::epsilon();
  1101. lm.parameters.xtol = NumTraits<double>::epsilon();
  1102. info = lm.minimize(x);
  1103. // check return value
  1104. VERIFY_IS_EQUAL(info, 1);
  1105. VERIFY_IS_EQUAL(lm.nfev, 15 );
  1106. VERIFY_IS_EQUAL(lm.njev, 14 );
  1107. // check norm^2
  1108. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.1680088766E+03);
  1109. // check x
  1110. VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
  1111. VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
  1112. }
  1113. struct MGH17_functor : Functor<double>
  1114. {
  1115. MGH17_functor(void) : Functor<double>(5,33) {}
  1116. static const double x[33];
  1117. static const double y[33];
  1118. int operator()(const VectorXd &b, VectorXd &fvec)
  1119. {
  1120. assert(b.size()==5);
  1121. assert(fvec.size()==33);
  1122. for(int i=0; i<33; i++)
  1123. fvec[i] = b[0] + b[1]*exp(-b[3]*x[i]) + b[2]*exp(-b[4]*x[i]) - y[i];
  1124. return 0;
  1125. }
  1126. int df(const VectorXd &b, MatrixXd &fjac)
  1127. {
  1128. assert(b.size()==5);
  1129. assert(fjac.rows()==33);
  1130. assert(fjac.cols()==5);
  1131. for(int i=0; i<33; i++) {
  1132. fjac(i,0) = 1.;
  1133. fjac(i,1) = exp(-b[3]*x[i]);
  1134. fjac(i,2) = exp(-b[4]*x[i]);
  1135. fjac(i,3) = -x[i]*b[1]*exp(-b[3]*x[i]);
  1136. fjac(i,4) = -x[i]*b[2]*exp(-b[4]*x[i]);
  1137. }
  1138. return 0;
  1139. }
  1140. };
  1141. 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 };
  1142. 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 };
  1143. // http://www.itl.nist.gov/div898/strd/nls/data/mgh17.shtml
  1144. void testNistMGH17(void)
  1145. {
  1146. const int n=5;
  1147. int info;
  1148. VectorXd x(n);
  1149. /*
  1150. * First try
  1151. */
  1152. x<< 50., 150., -100., 1., 2.;
  1153. // do the computation
  1154. MGH17_functor functor;
  1155. LevenbergMarquardt<MGH17_functor> lm(functor);
  1156. lm.parameters.ftol = NumTraits<double>::epsilon();
  1157. lm.parameters.xtol = NumTraits<double>::epsilon();
  1158. lm.parameters.maxfev = 1000;
  1159. info = lm.minimize(x);
  1160. // check norm^2
  1161. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.4648946975E-05);
  1162. // check x
  1163. VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
  1164. VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
  1165. VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
  1166. VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
  1167. VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
  1168. // check return value
  1169. VERIFY_IS_EQUAL(info, 2);
  1170. VERIFY(lm.nfev < 650); // 602
  1171. VERIFY(lm.njev < 600); // 545
  1172. /*
  1173. * Second try
  1174. */
  1175. x<< 0.5 ,1.5 ,-1 ,0.01 ,0.02;
  1176. // do the computation
  1177. lm.resetParameters();
  1178. info = lm.minimize(x);
  1179. // check return value
  1180. VERIFY_IS_EQUAL(info, 1);
  1181. VERIFY_IS_EQUAL(lm.nfev, 18);
  1182. VERIFY_IS_EQUAL(lm.njev, 15);
  1183. // check norm^2
  1184. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.4648946975E-05);
  1185. // check x
  1186. VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
  1187. VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
  1188. VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
  1189. VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
  1190. VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
  1191. }
  1192. struct MGH09_functor : Functor<double>
  1193. {
  1194. MGH09_functor(void) : Functor<double>(4,11) {}
  1195. static const double _x[11];
  1196. static const double y[11];
  1197. int operator()(const VectorXd &b, VectorXd &fvec)
  1198. {
  1199. assert(b.size()==4);
  1200. assert(fvec.size()==11);
  1201. for(int i=0; i<11; i++) {
  1202. double x = _x[i], xx=x*x;
  1203. fvec[i] = b[0]*(xx+x*b[1])/(xx+x*b[2]+b[3]) - y[i];
  1204. }
  1205. return 0;
  1206. }
  1207. int df(const VectorXd &b, MatrixXd &fjac)
  1208. {
  1209. assert(b.size()==4);
  1210. assert(fjac.rows()==11);
  1211. assert(fjac.cols()==4);
  1212. for(int i=0; i<11; i++) {
  1213. double x = _x[i], xx=x*x;
  1214. double factor = 1./(xx+x*b[2]+b[3]);
  1215. fjac(i,0) = (xx+x*b[1]) * factor;
  1216. fjac(i,1) = b[0]*x* factor;
  1217. fjac(i,2) = - b[0]*(xx+x*b[1]) * x * factor * factor;
  1218. fjac(i,3) = - b[0]*(xx+x*b[1]) * factor * factor;
  1219. }
  1220. return 0;
  1221. }
  1222. };
  1223. 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 };
  1224. 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 };
  1225. // http://www.itl.nist.gov/div898/strd/nls/data/mgh09.shtml
  1226. void testNistMGH09(void)
  1227. {
  1228. const int n=4;
  1229. int info;
  1230. VectorXd x(n);
  1231. /*
  1232. * First try
  1233. */
  1234. x<< 25., 39, 41.5, 39.;
  1235. // do the computation
  1236. MGH09_functor functor;
  1237. LevenbergMarquardt<MGH09_functor> lm(functor);
  1238. lm.parameters.maxfev = 1000;
  1239. info = lm.minimize(x);
  1240. // check return value
  1241. VERIFY_IS_EQUAL(info, 1);
  1242. VERIFY_IS_EQUAL(lm.nfev, 490 );
  1243. VERIFY_IS_EQUAL(lm.njev, 376 );
  1244. // check norm^2
  1245. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 3.0750560385E-04);
  1246. // check x
  1247. VERIFY_IS_APPROX(x[0], 0.1928077089); // should be 1.9280693458E-01
  1248. VERIFY_IS_APPROX(x[1], 0.19126423573); // should be 1.9128232873E-01
  1249. VERIFY_IS_APPROX(x[2], 0.12305309914); // should be 1.2305650693E-01
  1250. VERIFY_IS_APPROX(x[3], 0.13605395375); // should be 1.3606233068E-01
  1251. /*
  1252. * Second try
  1253. */
  1254. x<< 0.25, 0.39, 0.415, 0.39;
  1255. // do the computation
  1256. lm.resetParameters();
  1257. info = lm.minimize(x);
  1258. // check return value
  1259. VERIFY_IS_EQUAL(info, 1);
  1260. VERIFY_IS_EQUAL(lm.nfev, 18);
  1261. VERIFY_IS_EQUAL(lm.njev, 16);
  1262. // check norm^2
  1263. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 3.0750560385E-04);
  1264. // check x
  1265. VERIFY_IS_APPROX(x[0], 0.19280781); // should be 1.9280693458E-01
  1266. VERIFY_IS_APPROX(x[1], 0.19126265); // should be 1.9128232873E-01
  1267. VERIFY_IS_APPROX(x[2], 0.12305280); // should be 1.2305650693E-01
  1268. VERIFY_IS_APPROX(x[3], 0.13605322); // should be 1.3606233068E-01
  1269. }
  1270. struct Bennett5_functor : Functor<double>
  1271. {
  1272. Bennett5_functor(void) : Functor<double>(3,154) {}
  1273. static const double x[154];
  1274. static const double y[154];
  1275. int operator()(const VectorXd &b, VectorXd &fvec)
  1276. {
  1277. assert(b.size()==3);
  1278. assert(fvec.size()==154);
  1279. for(int i=0; i<154; i++)
  1280. fvec[i] = b[0]* pow(b[1]+x[i],-1./b[2]) - y[i];
  1281. return 0;
  1282. }
  1283. int df(const VectorXd &b, MatrixXd &fjac)
  1284. {
  1285. assert(b.size()==3);
  1286. assert(fjac.rows()==154);
  1287. assert(fjac.cols()==3);
  1288. for(int i=0; i<154; i++) {
  1289. double e = pow(b[1]+x[i],-1./b[2]);
  1290. fjac(i,0) = e;
  1291. fjac(i,1) = - b[0]*e/b[2]/(b[1]+x[i]);
  1292. fjac(i,2) = b[0]*e*log(b[1]+x[i])/b[2]/b[2];
  1293. }
  1294. return 0;
  1295. }
  1296. };
  1297. 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,
  1298. 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 };
  1299. 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
  1300. ,-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 ,
  1301. -31.827299E0 ,-31.821600E0 ,-31.818701E0 ,-31.812901E0 ,-31.809999E0 ,-31.807100E0 ,-31.801300E0 ,-31.798401E0 ,-31.795500E0 ,-31.789700E0 ,-31.786800E0 };
  1302. // http://www.itl.nist.gov/div898/strd/nls/data/bennett5.shtml
  1303. void testNistBennett5(void)
  1304. {
  1305. const int n=3;
  1306. int info;
  1307. VectorXd x(n);
  1308. /*
  1309. * First try
  1310. */
  1311. x<< -2000., 50., 0.8;
  1312. // do the computation
  1313. Bennett5_functor functor;
  1314. LevenbergMarquardt<Bennett5_functor> lm(functor);
  1315. lm.parameters.maxfev = 1000;
  1316. info = lm.minimize(x);
  1317. // check return value
  1318. VERIFY_IS_EQUAL(info, 1);
  1319. VERIFY_IS_EQUAL(lm.nfev, 758);
  1320. VERIFY_IS_EQUAL(lm.njev, 744);
  1321. // check norm^2
  1322. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.2404744073E-04);
  1323. // check x
  1324. VERIFY_IS_APPROX(x[0], -2.5235058043E+03);
  1325. VERIFY_IS_APPROX(x[1], 4.6736564644E+01);
  1326. VERIFY_IS_APPROX(x[2], 9.3218483193E-01);
  1327. /*
  1328. * Second try
  1329. */
  1330. x<< -1500., 45., 0.85;
  1331. // do the computation
  1332. lm.resetParameters();
  1333. info = lm.minimize(x);
  1334. // check return value
  1335. VERIFY_IS_EQUAL(info, 1);
  1336. VERIFY_IS_EQUAL(lm.nfev, 203);
  1337. VERIFY_IS_EQUAL(lm.njev, 192);
  1338. // check norm^2
  1339. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.2404744073E-04);
  1340. // check x
  1341. VERIFY_IS_APPROX(x[0], -2523.3007865); // should be -2.5235058043E+03
  1342. VERIFY_IS_APPROX(x[1], 46.735705771); // should be 4.6736564644E+01);
  1343. VERIFY_IS_APPROX(x[2], 0.93219881891); // should be 9.3218483193E-01);
  1344. }
  1345. struct thurber_functor : Functor<double>
  1346. {
  1347. thurber_functor(void) : Functor<double>(7,37) {}
  1348. static const double _x[37];
  1349. static const double _y[37];
  1350. int operator()(const VectorXd &b, VectorXd &fvec)
  1351. {
  1352. // int called=0; printf("call hahn1_functor with iflag=%d, called=%d\n", iflag, called); if (iflag==1) called++;
  1353. assert(b.size()==7);
  1354. assert(fvec.size()==37);
  1355. for(int i=0; i<37; i++) {
  1356. double x=_x[i], xx=x*x, xxx=xx*x;
  1357. 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];
  1358. }
  1359. return 0;
  1360. }
  1361. int df(const VectorXd &b, MatrixXd &fjac)
  1362. {
  1363. assert(b.size()==7);
  1364. assert(fjac.rows()==37);
  1365. assert(fjac.cols()==7);
  1366. for(int i=0; i<37; i++) {
  1367. double x=_x[i], xx=x*x, xxx=xx*x;
  1368. double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
  1369. fjac(i,0) = 1.*fact;
  1370. fjac(i,1) = x*fact;
  1371. fjac(i,2) = xx*fact;
  1372. fjac(i,3) = xxx*fact;
  1373. fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
  1374. fjac(i,4) = x*fact;
  1375. fjac(i,5) = xx*fact;
  1376. fjac(i,6) = xxx*fact;
  1377. }
  1378. return 0;
  1379. }
  1380. };
  1381. 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 };
  1382. 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};
  1383. // http://www.itl.nist.gov/div898/strd/nls/data/thurber.shtml
  1384. void testNistThurber(void)
  1385. {
  1386. const int n=7;
  1387. int info;
  1388. VectorXd x(n);
  1389. /*
  1390. * First try
  1391. */
  1392. x<< 1000 ,1000 ,400 ,40 ,0.7,0.3,0.0 ;
  1393. // do the computation
  1394. thurber_functor functor;
  1395. LevenbergMarquardt<thurber_functor> lm(functor);
  1396. lm.parameters.ftol = 1.E4*NumTraits<double>::epsilon();
  1397. lm.parameters.xtol = 1.E4*NumTraits<double>::epsilon();
  1398. info = lm.minimize(x);
  1399. // check return value
  1400. VERIFY_IS_EQUAL(info, 1);
  1401. VERIFY_IS_EQUAL(lm.nfev, 39);
  1402. VERIFY_IS_EQUAL(lm.njev, 36);
  1403. // check norm^2
  1404. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6427082397E+03);
  1405. // check x
  1406. VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
  1407. VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
  1408. VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
  1409. VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
  1410. VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
  1411. VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
  1412. VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
  1413. /*
  1414. * Second try
  1415. */
  1416. x<< 1300 ,1500 ,500 ,75 ,1 ,0.4 ,0.05 ;
  1417. // do the computation
  1418. lm.resetParameters();
  1419. lm.parameters.ftol = 1.E4*NumTraits<double>::epsilon();
  1420. lm.parameters.xtol = 1.E4*NumTraits<double>::epsilon();
  1421. info = lm.minimize(x);
  1422. // check return value
  1423. VERIFY_IS_EQUAL(info, 1);
  1424. VERIFY_IS_EQUAL(lm.nfev, 29);
  1425. VERIFY_IS_EQUAL(lm.njev, 28);
  1426. // check norm^2
  1427. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 5.6427082397E+03);
  1428. // check x
  1429. VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
  1430. VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
  1431. VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
  1432. VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
  1433. VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
  1434. VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
  1435. VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
  1436. }
  1437. struct rat43_functor : Functor<double>
  1438. {
  1439. rat43_functor(void) : Functor<double>(4,15) {}
  1440. static const double x[15];
  1441. static const double y[15];
  1442. int operator()(const VectorXd &b, VectorXd &fvec)
  1443. {
  1444. assert(b.size()==4);
  1445. assert(fvec.size()==15);
  1446. for(int i=0; i<15; i++)
  1447. fvec[i] = b[0] * pow(1.+exp(b[1]-b[2]*x[i]),-1./b[3]) - y[i];
  1448. return 0;
  1449. }
  1450. int df(const VectorXd &b, MatrixXd &fjac)
  1451. {
  1452. assert(b.size()==4);
  1453. assert(fjac.rows()==15);
  1454. assert(fjac.cols()==4);
  1455. for(int i=0; i<15; i++) {
  1456. double e = exp(b[1]-b[2]*x[i]);
  1457. double power = -1./b[3];
  1458. fjac(i,0) = pow(1.+e, power);
  1459. fjac(i,1) = power*b[0]*e*pow(1.+e, power-1.);
  1460. fjac(i,2) = -power*b[0]*e*x[i]*pow(1.+e, power-1.);
  1461. fjac(i,3) = b[0]*power*power*log(1.+e)*pow(1.+e, power);
  1462. }
  1463. return 0;
  1464. }
  1465. };
  1466. const double rat43_functor::x[15] = { 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14., 15. };
  1467. 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 };
  1468. // http://www.itl.nist.gov/div898/strd/nls/data/ratkowsky3.shtml
  1469. void testNistRat43(void)
  1470. {
  1471. const int n=4;
  1472. int info;
  1473. VectorXd x(n);
  1474. /*
  1475. * First try
  1476. */
  1477. x<< 100., 10., 1., 1.;
  1478. // do the computation
  1479. rat43_functor functor;
  1480. LevenbergMarquardt<rat43_functor> lm(functor);
  1481. lm.parameters.ftol = 1.E6*NumTraits<double>::epsilon();
  1482. lm.parameters.xtol = 1.E6*NumTraits<double>::epsilon();
  1483. info = lm.minimize(x);
  1484. // check return value
  1485. VERIFY_IS_EQUAL(info, 1);
  1486. VERIFY_IS_EQUAL(lm.nfev, 27);
  1487. VERIFY_IS_EQUAL(lm.njev, 20);
  1488. // check norm^2
  1489. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7864049080E+03);
  1490. // check x
  1491. VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
  1492. VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
  1493. VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
  1494. VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
  1495. /*
  1496. * Second try
  1497. */
  1498. x<< 700., 5., 0.75, 1.3;
  1499. // do the computation
  1500. lm.resetParameters();
  1501. lm.parameters.ftol = 1.E5*NumTraits<double>::epsilon();
  1502. lm.parameters.xtol = 1.E5*NumTraits<double>::epsilon();
  1503. info = lm.minimize(x);
  1504. // check return value
  1505. VERIFY_IS_EQUAL(info, 1);
  1506. VERIFY_IS_EQUAL(lm.nfev, 9);
  1507. VERIFY_IS_EQUAL(lm.njev, 8);
  1508. // check norm^2
  1509. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 8.7864049080E+03);
  1510. // check x
  1511. VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
  1512. VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
  1513. VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
  1514. VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
  1515. }
  1516. struct eckerle4_functor : Functor<double>
  1517. {
  1518. eckerle4_functor(void) : Functor<double>(3,35) {}
  1519. static const double x[35];
  1520. static const double y[35];
  1521. int operator()(const VectorXd &b, VectorXd &fvec)
  1522. {
  1523. assert(b.size()==3);
  1524. assert(fvec.size()==35);
  1525. for(int i=0; i<35; i++)
  1526. fvec[i] = b[0]/b[1] * exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/(b[1]*b[1])) - y[i];
  1527. return 0;
  1528. }
  1529. int df(const VectorXd &b, MatrixXd &fjac)
  1530. {
  1531. assert(b.size()==3);
  1532. assert(fjac.rows()==35);
  1533. assert(fjac.cols()==3);
  1534. for(int i=0; i<35; i++) {
  1535. double b12 = b[1]*b[1];
  1536. double e = exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/b12);
  1537. fjac(i,0) = e / b[1];
  1538. fjac(i,1) = ((x[i]-b[2])*(x[i]-b[2])/b12-1.) * b[0]*e/b12;
  1539. fjac(i,2) = (x[i]-b[2])*e*b[0]/b[1]/b12;
  1540. }
  1541. return 0;
  1542. }
  1543. };
  1544. 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};
  1545. 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 };
  1546. // http://www.itl.nist.gov/div898/strd/nls/data/eckerle4.shtml
  1547. void testNistEckerle4(void)
  1548. {
  1549. const int n=3;
  1550. int info;
  1551. VectorXd x(n);
  1552. /*
  1553. * First try
  1554. */
  1555. x<< 1., 10., 500.;
  1556. // do the computation
  1557. eckerle4_functor functor;
  1558. LevenbergMarquardt<eckerle4_functor> lm(functor);
  1559. info = lm.minimize(x);
  1560. // check return value
  1561. VERIFY_IS_EQUAL(info, 1);
  1562. VERIFY_IS_EQUAL(lm.nfev, 18);
  1563. VERIFY_IS_EQUAL(lm.njev, 15);
  1564. // check norm^2
  1565. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.4635887487E-03);
  1566. // check x
  1567. VERIFY_IS_APPROX(x[0], 1.5543827178);
  1568. VERIFY_IS_APPROX(x[1], 4.0888321754);
  1569. VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
  1570. /*
  1571. * Second try
  1572. */
  1573. x<< 1.5, 5., 450.;
  1574. // do the computation
  1575. info = lm.minimize(x);
  1576. // check return value
  1577. VERIFY_IS_EQUAL(info, 1);
  1578. VERIFY_IS_EQUAL(lm.nfev, 7);
  1579. VERIFY_IS_EQUAL(lm.njev, 6);
  1580. // check norm^2
  1581. VERIFY_IS_APPROX(lm.fvec.squaredNorm(), 1.4635887487E-03);
  1582. // check x
  1583. VERIFY_IS_APPROX(x[0], 1.5543827178);
  1584. VERIFY_IS_APPROX(x[1], 4.0888321754);
  1585. VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
  1586. }
  1587. void test_NonLinearOptimization()
  1588. {
  1589. // Tests using the examples provided by (c)minpack
  1590. CALL_SUBTEST/*_1*/(testChkder());
  1591. CALL_SUBTEST/*_1*/(testLmder1());
  1592. CALL_SUBTEST/*_1*/(testLmder());
  1593. CALL_SUBTEST/*_2*/(testHybrj1());
  1594. CALL_SUBTEST/*_2*/(testHybrj());
  1595. CALL_SUBTEST/*_2*/(testHybrd1());
  1596. CALL_SUBTEST/*_2*/(testHybrd());
  1597. CALL_SUBTEST/*_3*/(testLmstr1());
  1598. CALL_SUBTEST/*_3*/(testLmstr());
  1599. CALL_SUBTEST/*_3*/(testLmdif1());
  1600. CALL_SUBTEST/*_3*/(testLmdif());
  1601. // NIST tests, level of difficulty = "Lower"
  1602. CALL_SUBTEST/*_4*/(testNistMisra1a());
  1603. CALL_SUBTEST/*_4*/(testNistChwirut2());
  1604. // NIST tests, level of difficulty = "Average"
  1605. CALL_SUBTEST/*_5*/(testNistHahn1());
  1606. CALL_SUBTEST/*_6*/(testNistMisra1d());
  1607. CALL_SUBTEST/*_7*/(testNistMGH17());
  1608. CALL_SUBTEST/*_8*/(testNistLanczos1());
  1609. // // NIST tests, level of difficulty = "Higher"
  1610. CALL_SUBTEST/*_9*/(testNistRat42());
  1611. // CALL_SUBTEST/*_10*/(testNistMGH10());
  1612. CALL_SUBTEST/*_11*/(testNistBoxBOD());
  1613. // CALL_SUBTEST/*_12*/(testNistMGH09());
  1614. CALL_SUBTEST/*_13*/(testNistBennett5());
  1615. CALL_SUBTEST/*_14*/(testNistThurber());
  1616. CALL_SUBTEST/*_15*/(testNistRat43());
  1617. CALL_SUBTEST/*_16*/(testNistEckerle4());
  1618. }
  1619. /*
  1620. * Can be useful for debugging...
  1621. printf("info, nfev : %d, %d\n", info, lm.nfev);
  1622. printf("info, nfev, njev : %d, %d, %d\n", info, solver.nfev, solver.njev);
  1623. printf("info, nfev : %d, %d\n", info, solver.nfev);
  1624. printf("x[0] : %.32g\n", x[0]);
  1625. printf("x[1] : %.32g\n", x[1]);
  1626. printf("x[2] : %.32g\n", x[2]);
  1627. printf("x[3] : %.32g\n", x[3]);
  1628. printf("fvec.blueNorm() : %.32g\n", solver.fvec.blueNorm());
  1629. printf("fvec.blueNorm() : %.32g\n", lm.fvec.blueNorm());
  1630. printf("info, nfev, njev : %d, %d, %d\n", info, lm.nfev, lm.njev);
  1631. printf("fvec.squaredNorm() : %.13g\n", lm.fvec.squaredNorm());
  1632. std::cout << x << std::endl;
  1633. std::cout.precision(9);
  1634. std::cout << x[0] << std::endl;
  1635. std::cout << x[1] << std::endl;
  1636. std::cout << x[2] << std::endl;
  1637. std::cout << x[3] << std::endl;
  1638. */