gsl_multifit__lmpar.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500
  1. /* multifit/lmpar.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. #include "gsl_permute_vector_double.h"
  20. #include "gsl_multifit__qrsolv.c"
  21. static size_t
  22. count_nsing (const gsl_matrix * r)
  23. {
  24. /* Count the number of nonsingular entries. Returns the index of the
  25. first entry which is singular. */
  26. size_t n = r->size2;
  27. size_t i;
  28. for (i = 0; i < n; i++)
  29. {
  30. double rii = gsl_matrix_get (r, i, i);
  31. if (rii == 0)
  32. {
  33. break;
  34. }
  35. }
  36. return i;
  37. }
  38. static void
  39. compute_newton_direction (const gsl_matrix * r, const gsl_permutation * perm,
  40. const gsl_vector * qtf, gsl_vector * x)
  41. {
  42. /* Compute and store in x the Gauss-Newton direction. If the
  43. Jacobian is rank-deficient then obtain a least squares
  44. solution. */
  45. const size_t n = r->size2;
  46. size_t i, j, nsing;
  47. for (i = 0 ; i < n ; i++)
  48. {
  49. double qtfi = gsl_vector_get (qtf, i);
  50. gsl_vector_set (x, i, qtfi);
  51. }
  52. nsing = count_nsing (r);
  53. #ifdef DEBUG
  54. printf("nsing = %d\n", nsing);
  55. printf("r = "); gsl_matrix_fprintf(stdout, r, "%g"); printf("\n");
  56. printf("qtf = "); gsl_vector_fprintf(stdout, x, "%g"); printf("\n");
  57. #endif
  58. for (i = nsing; i < n; i++)
  59. {
  60. gsl_vector_set (x, i, 0.0);
  61. }
  62. if (nsing > 0)
  63. {
  64. for (j = nsing; j > 0 && j--;)
  65. {
  66. double rjj = gsl_matrix_get (r, j, j);
  67. double temp = gsl_vector_get (x, j) / rjj;
  68. gsl_vector_set (x, j, temp);
  69. for (i = 0; i < j; i++)
  70. {
  71. double rij = gsl_matrix_get (r, i, j);
  72. double xi = gsl_vector_get (x, i);
  73. gsl_vector_set (x, i, xi - rij * temp);
  74. }
  75. }
  76. }
  77. gsl_permute_vector_inverse (perm, x);
  78. }
  79. static void
  80. compute_newton_correction (const gsl_matrix * r, const gsl_vector * sdiag,
  81. const gsl_permutation * p, gsl_vector * x,
  82. double dxnorm,
  83. const gsl_vector * diag, gsl_vector * w)
  84. {
  85. size_t n = r->size2;
  86. size_t i, j;
  87. for (i = 0; i < n; i++)
  88. {
  89. size_t pi = gsl_permutation_get (p, i);
  90. double dpi = gsl_vector_get (diag, pi);
  91. double xpi = gsl_vector_get (x, pi);
  92. gsl_vector_set (w, i, dpi * (dpi * xpi) / dxnorm);
  93. }
  94. for (j = 0; j < n; j++)
  95. {
  96. double sj = gsl_vector_get (sdiag, j);
  97. double wj = gsl_vector_get (w, j);
  98. double tj = wj / sj;
  99. gsl_vector_set (w, j, tj);
  100. for (i = j + 1; i < n; i++)
  101. {
  102. double rij = gsl_matrix_get (r, i, j);
  103. double wi = gsl_vector_get (w, i);
  104. gsl_vector_set (w, i, wi - rij * tj);
  105. }
  106. }
  107. }
  108. static void
  109. compute_newton_bound (const gsl_matrix * r, const gsl_vector * x,
  110. double dxnorm, const gsl_permutation * perm,
  111. const gsl_vector * diag, gsl_vector * w)
  112. {
  113. /* If the jacobian is not rank-deficient then the Newton step
  114. provides a lower bound for the zero of the function. Otherwise
  115. set this bound to zero. */
  116. size_t n = r->size2;
  117. size_t i, j;
  118. size_t nsing = count_nsing (r);
  119. if (nsing < n)
  120. {
  121. gsl_vector_set_zero (w);
  122. return;
  123. }
  124. for (i = 0; i < n; i++)
  125. {
  126. size_t pi = gsl_permutation_get (perm, i);
  127. double dpi = gsl_vector_get (diag, pi);
  128. double xpi = gsl_vector_get (x, pi);
  129. gsl_vector_set (w, i, dpi * (dpi * xpi / dxnorm));
  130. }
  131. for (j = 0; j < n; j++)
  132. {
  133. double sum = 0;
  134. for (i = 0; i < j; i++)
  135. {
  136. sum += gsl_matrix_get (r, i, j) * gsl_vector_get (w, i);
  137. }
  138. {
  139. double rjj = gsl_matrix_get (r, j, j);
  140. double wj = gsl_vector_get (w, j);
  141. gsl_vector_set (w, j, (wj - sum) / rjj);
  142. }
  143. }
  144. }
  145. static void
  146. compute_gradient_direction (const gsl_matrix * r, const gsl_permutation * p,
  147. const gsl_vector * qtf, const gsl_vector * diag,
  148. gsl_vector * g)
  149. {
  150. const size_t n = r->size2;
  151. size_t i, j;
  152. for (j = 0; j < n; j++)
  153. {
  154. double sum = 0;
  155. for (i = 0; i <= j; i++)
  156. {
  157. sum += gsl_matrix_get (r, i, j) * gsl_vector_get (qtf, i);
  158. }
  159. {
  160. size_t pj = gsl_permutation_get (p, j);
  161. double dpj = gsl_vector_get (diag, pj);
  162. gsl_vector_set (g, j, sum / dpj);
  163. }
  164. }
  165. }
  166. static int
  167. lmpar (gsl_matrix * r, const gsl_permutation * perm, const gsl_vector * qtf,
  168. const gsl_vector * diag, double delta, double * par_inout,
  169. gsl_vector * newton, gsl_vector * gradient, gsl_vector * sdiag,
  170. gsl_vector * x, gsl_vector * w)
  171. {
  172. double dxnorm, gnorm, fp, fp_old, par_lower, par_upper, par_c;
  173. double par = *par_inout;
  174. size_t iter = 0;
  175. #ifdef DEBUG
  176. printf("ENTERING lmpar\n");
  177. #endif
  178. compute_newton_direction (r, perm, qtf, newton);
  179. #ifdef DEBUG
  180. printf ("newton = ");
  181. gsl_vector_fprintf (stdout, newton, "%g");
  182. printf ("\n");
  183. printf ("diag = ");
  184. gsl_vector_fprintf (stdout, diag, "%g");
  185. printf ("\n");
  186. #endif
  187. /* Evaluate the function at the origin and test for acceptance of
  188. the Gauss-Newton direction. */
  189. dxnorm = scaled_enorm (diag, newton);
  190. fp = dxnorm - delta;
  191. #ifdef DEBUG
  192. printf ("dxnorm = %g, delta = %g, fp = %g\n", dxnorm, delta, fp);
  193. #endif
  194. if (fp <= 0.1 * delta)
  195. {
  196. gsl_vector_memcpy (x, newton);
  197. #ifdef DEBUG
  198. printf ("took newton (fp = %g, delta = %g)\n", fp, delta);
  199. #endif
  200. *par_inout = 0;
  201. return GSL_SUCCESS;
  202. }
  203. #ifdef DEBUG
  204. printf ("r = ");
  205. gsl_matrix_fprintf (stdout, r, "%g");
  206. printf ("\n");
  207. printf ("newton = ");
  208. gsl_vector_fprintf (stdout, newton, "%g");
  209. printf ("\n");
  210. printf ("dxnorm = %g\n", dxnorm);
  211. #endif
  212. compute_newton_bound (r, newton, dxnorm, perm, diag, w);
  213. #ifdef DEBUG
  214. printf("perm = "); gsl_permutation_fprintf(stdout, perm, "%d");
  215. printf ("diag = ");
  216. gsl_vector_fprintf (stdout, diag, "%g");
  217. printf ("\n");
  218. printf ("w = ");
  219. gsl_vector_fprintf (stdout, w, "%g");
  220. printf ("\n");
  221. #endif
  222. {
  223. double wnorm = enorm (w);
  224. double phider = wnorm * wnorm;
  225. /* w == zero if r rank-deficient,
  226. then set lower bound to zero form MINPACK, lmder.f
  227. Hans E. Plesser 2002-02-25 (hans.plesser@itf.nlh.no) */
  228. if ( wnorm > 0 )
  229. par_lower = fp / (delta * phider);
  230. else
  231. par_lower = 0.0;
  232. }
  233. #ifdef DEBUG
  234. printf("par = %g\n", par );
  235. printf("par_lower = %g\n", par_lower);
  236. #endif
  237. compute_gradient_direction (r, perm, qtf, diag, gradient);
  238. gnorm = enorm (gradient);
  239. #ifdef DEBUG
  240. printf("gradient = "); gsl_vector_fprintf(stdout, gradient, "%g"); printf("\n");
  241. printf("gnorm = %g\n", gnorm);
  242. #endif
  243. par_upper = gnorm / delta;
  244. if (par_upper == 0)
  245. {
  246. par_upper = GSL_DBL_MIN / GSL_MIN_DBL(delta, 0.1);
  247. }
  248. #ifdef DEBUG
  249. printf("par_upper = %g\n", par_upper);
  250. #endif
  251. if (par > par_upper)
  252. {
  253. #ifdef DEBUG
  254. printf("set par to par_upper\n");
  255. #endif
  256. par = par_upper;
  257. }
  258. else if (par < par_lower)
  259. {
  260. #ifdef DEBUG
  261. printf("set par to par_lower\n");
  262. #endif
  263. par = par_lower;
  264. }
  265. if (par == 0)
  266. {
  267. par = gnorm / dxnorm;
  268. #ifdef DEBUG
  269. printf("set par to gnorm/dxnorm = %g\n", par);
  270. #endif
  271. }
  272. /* Beginning of iteration */
  273. iteration:
  274. iter++;
  275. #ifdef DEBUG
  276. printf("lmpar iteration = %d\n", iter);
  277. #endif
  278. #ifdef BRIANSFIX
  279. /* Seems like this is described in the paper but not in the MINPACK code */
  280. if (par < par_lower || par > par_upper)
  281. {
  282. par = GSL_MAX_DBL (0.001 * par_upper, sqrt(par_lower * par_upper));
  283. }
  284. #endif
  285. /* Evaluate the function at the current value of par */
  286. if (par == 0)
  287. {
  288. par = GSL_MAX_DBL (0.001 * par_upper, GSL_DBL_MIN);
  289. #ifdef DEBUG
  290. printf("par = 0, set par to = %g\n", par);
  291. #endif
  292. }
  293. /* Compute the least squares solution of [ R P x - Q^T f, sqrt(par) D x]
  294. for A = Q R P^T */
  295. #ifdef DEBUG
  296. printf ("calling qrsolv with par = %g\n", par);
  297. #endif
  298. {
  299. double sqrt_par = sqrt(par);
  300. qrsolv (r, perm, sqrt_par, diag, qtf, x, sdiag, w);
  301. }
  302. dxnorm = scaled_enorm (diag, x);
  303. fp_old = fp;
  304. fp = dxnorm - delta;
  305. #ifdef DEBUG
  306. printf ("After qrsolv dxnorm = %g, delta = %g, fp = %g\n", dxnorm, delta, fp);
  307. printf ("sdiag = ") ; gsl_vector_fprintf(stdout, sdiag, "%g"); printf("\n");
  308. printf ("x = ") ; gsl_vector_fprintf(stdout, x, "%g"); printf("\n");
  309. printf ("r = ") ; gsl_matrix_fprintf(stdout, r, "%g"); printf("\nXXX\n");
  310. #endif
  311. /* If the function is small enough, accept the current value of par */
  312. if (fabs (fp) <= 0.1 * delta)
  313. goto line220;
  314. if (par_lower == 0 && fp <= fp_old && fp_old < 0)
  315. goto line220;
  316. /* Check for maximum number of iterations */
  317. if (iter == 10)
  318. goto line220;
  319. /* Compute the Newton correction */
  320. compute_newton_correction (r, sdiag, perm, x, dxnorm, diag, w);
  321. #ifdef DEBUG
  322. printf ("newton_correction = ");
  323. gsl_vector_fprintf(stdout, w, "%g"); printf("\n");
  324. #endif
  325. {
  326. double wnorm = enorm (w);
  327. par_c = fp / (delta * wnorm * wnorm);
  328. }
  329. #ifdef DEBUG
  330. printf("fp = %g\n", fp);
  331. printf("par_lower = %g\n", par_lower);
  332. printf("par_upper = %g\n", par_upper);
  333. printf("par_c = %g\n", par_c);
  334. #endif
  335. /* Depending on the sign of the function, update par_lower or par_upper */
  336. if (fp > 0)
  337. {
  338. if (par > par_lower)
  339. {
  340. par_lower = par;
  341. #ifdef DEBUG
  342. printf("fp > 0: set par_lower = par = %g\n", par);
  343. #endif
  344. }
  345. }
  346. else if (fp < 0)
  347. {
  348. if (par < par_upper)
  349. {
  350. #ifdef DEBUG
  351. printf("fp < 0: set par_upper = par = %g\n", par);
  352. #endif
  353. par_upper = par;
  354. }
  355. }
  356. /* Compute an improved estimate for par */
  357. #ifdef DEBUG
  358. printf("improved estimate par = MAX(%g, %g) \n", par_lower, par+par_c);
  359. #endif
  360. par = GSL_MAX_DBL (par_lower, par + par_c);
  361. #ifdef DEBUG
  362. printf("improved estimate par = %g \n", par);
  363. #endif
  364. goto iteration;
  365. line220:
  366. #ifdef DEBUG
  367. printf("LEAVING lmpar, par = %g\n", par);
  368. #endif
  369. *par_inout = par;
  370. return GSL_SUCCESS;
  371. }