gsl_multifit__lmder.c 9.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388
  1. /* multfit/lmder.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__config.h"
  20. #include <stddef.h>
  21. #include <stdlib.h>
  22. #include <stdio.h>
  23. #include <math.h>
  24. #include <float.h>
  25. #include "gsl_math.h"
  26. #include "gsl_errno.h"
  27. #include "gsl_multifit_nlin.h"
  28. #include "gsl_blas.h"
  29. #include "gsl_linalg.h"
  30. #include "gsl_permutation.h"
  31. typedef struct
  32. {
  33. size_t iter;
  34. double xnorm;
  35. double fnorm;
  36. double delta;
  37. double par;
  38. gsl_matrix *r;
  39. gsl_vector *tau;
  40. gsl_vector *diag;
  41. gsl_vector *qtf;
  42. gsl_vector *newton;
  43. gsl_vector *gradient;
  44. gsl_vector *x_trial;
  45. gsl_vector *f_trial;
  46. gsl_vector *df;
  47. gsl_vector *sdiag;
  48. gsl_vector *rptdx;
  49. gsl_vector *w;
  50. gsl_vector *work1;
  51. gsl_permutation * perm;
  52. }
  53. lmder_state_t;
  54. static int lmder_alloc (void *vstate, size_t n, size_t p);
  55. static int lmder_set (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  56. static int lmsder_set (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  57. static int set (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx, int scale);
  58. static int lmder_iterate (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  59. static void lmder_free (void *vstate);
  60. static int iterate (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx, int scale);
  61. #include "gsl_multifit__lmutil.c"
  62. #include "gsl_multifit__lmpar.c"
  63. #include "gsl_multifit__lmset.c"
  64. #include "gsl_multifit__lmiterate.c"
  65. static int
  66. lmder_alloc (void *vstate, size_t n, size_t p)
  67. {
  68. lmder_state_t *state = (lmder_state_t *) vstate;
  69. gsl_matrix *r;
  70. gsl_vector *tau, *diag, *qtf, *newton, *gradient, *x_trial, *f_trial,
  71. *df, *sdiag, *rptdx, *w, *work1;
  72. gsl_permutation *perm;
  73. r = gsl_matrix_calloc (n, p);
  74. if (r == 0)
  75. {
  76. GSL_ERROR ("failed to allocate space for r", GSL_ENOMEM);
  77. }
  78. state->r = r;
  79. tau = gsl_vector_calloc (GSL_MIN(n, p));
  80. if (tau == 0)
  81. {
  82. gsl_matrix_free (r);
  83. GSL_ERROR ("failed to allocate space for tau", GSL_ENOMEM);
  84. }
  85. state->tau = tau;
  86. diag = gsl_vector_calloc (p);
  87. if (diag == 0)
  88. {
  89. gsl_matrix_free (r);
  90. gsl_vector_free (tau);
  91. GSL_ERROR ("failed to allocate space for diag", GSL_ENOMEM);
  92. }
  93. state->diag = diag;
  94. qtf = gsl_vector_calloc (n);
  95. if (qtf == 0)
  96. {
  97. gsl_matrix_free (r);
  98. gsl_vector_free (tau);
  99. gsl_vector_free (diag);
  100. GSL_ERROR ("failed to allocate space for qtf", GSL_ENOMEM);
  101. }
  102. state->qtf = qtf;
  103. newton = gsl_vector_calloc (p);
  104. if (newton == 0)
  105. {
  106. gsl_matrix_free (r);
  107. gsl_vector_free (tau);
  108. gsl_vector_free (diag);
  109. gsl_vector_free (qtf);
  110. GSL_ERROR ("failed to allocate space for newton", GSL_ENOMEM);
  111. }
  112. state->newton = newton;
  113. gradient = gsl_vector_calloc (p);
  114. if (gradient == 0)
  115. {
  116. gsl_matrix_free (r);
  117. gsl_vector_free (tau);
  118. gsl_vector_free (diag);
  119. gsl_vector_free (qtf);
  120. gsl_vector_free (newton);
  121. GSL_ERROR ("failed to allocate space for gradient", GSL_ENOMEM);
  122. }
  123. state->gradient = gradient;
  124. x_trial = gsl_vector_calloc (p);
  125. if (x_trial == 0)
  126. {
  127. gsl_matrix_free (r);
  128. gsl_vector_free (tau);
  129. gsl_vector_free (diag);
  130. gsl_vector_free (qtf);
  131. gsl_vector_free (newton);
  132. gsl_vector_free (gradient);
  133. GSL_ERROR ("failed to allocate space for x_trial", GSL_ENOMEM);
  134. }
  135. state->x_trial = x_trial;
  136. f_trial = gsl_vector_calloc (n);
  137. if (f_trial == 0)
  138. {
  139. gsl_matrix_free (r);
  140. gsl_vector_free (tau);
  141. gsl_vector_free (diag);
  142. gsl_vector_free (qtf);
  143. gsl_vector_free (newton);
  144. gsl_vector_free (gradient);
  145. gsl_vector_free (x_trial);
  146. GSL_ERROR ("failed to allocate space for f_trial", GSL_ENOMEM);
  147. }
  148. state->f_trial = f_trial;
  149. df = gsl_vector_calloc (n);
  150. if (df == 0)
  151. {
  152. gsl_matrix_free (r);
  153. gsl_vector_free (tau);
  154. gsl_vector_free (diag);
  155. gsl_vector_free (qtf);
  156. gsl_vector_free (newton);
  157. gsl_vector_free (gradient);
  158. gsl_vector_free (x_trial);
  159. gsl_vector_free (f_trial);
  160. GSL_ERROR ("failed to allocate space for df", GSL_ENOMEM);
  161. }
  162. state->df = df;
  163. sdiag = gsl_vector_calloc (p);
  164. if (sdiag == 0)
  165. {
  166. gsl_matrix_free (r);
  167. gsl_vector_free (tau);
  168. gsl_vector_free (diag);
  169. gsl_vector_free (qtf);
  170. gsl_vector_free (newton);
  171. gsl_vector_free (gradient);
  172. gsl_vector_free (x_trial);
  173. gsl_vector_free (f_trial);
  174. gsl_vector_free (df);
  175. GSL_ERROR ("failed to allocate space for sdiag", GSL_ENOMEM);
  176. }
  177. state->sdiag = sdiag;
  178. rptdx = gsl_vector_calloc (n);
  179. if (rptdx == 0)
  180. {
  181. gsl_matrix_free (r);
  182. gsl_vector_free (tau);
  183. gsl_vector_free (diag);
  184. gsl_vector_free (qtf);
  185. gsl_vector_free (newton);
  186. gsl_vector_free (gradient);
  187. gsl_vector_free (x_trial);
  188. gsl_vector_free (f_trial);
  189. gsl_vector_free (df);
  190. gsl_vector_free (sdiag);
  191. GSL_ERROR ("failed to allocate space for rptdx", GSL_ENOMEM);
  192. }
  193. state->rptdx = rptdx;
  194. w = gsl_vector_calloc (n);
  195. if (w == 0)
  196. {
  197. gsl_matrix_free (r);
  198. gsl_vector_free (tau);
  199. gsl_vector_free (diag);
  200. gsl_vector_free (qtf);
  201. gsl_vector_free (newton);
  202. gsl_vector_free (gradient);
  203. gsl_vector_free (x_trial);
  204. gsl_vector_free (f_trial);
  205. gsl_vector_free (df);
  206. gsl_vector_free (sdiag);
  207. gsl_vector_free (rptdx);
  208. GSL_ERROR ("failed to allocate space for w", GSL_ENOMEM);
  209. }
  210. state->w = w;
  211. work1 = gsl_vector_calloc (p);
  212. if (work1 == 0)
  213. {
  214. gsl_matrix_free (r);
  215. gsl_vector_free (tau);
  216. gsl_vector_free (diag);
  217. gsl_vector_free (qtf);
  218. gsl_vector_free (newton);
  219. gsl_vector_free (gradient);
  220. gsl_vector_free (x_trial);
  221. gsl_vector_free (f_trial);
  222. gsl_vector_free (df);
  223. gsl_vector_free (sdiag);
  224. gsl_vector_free (rptdx);
  225. gsl_vector_free (w);
  226. GSL_ERROR ("failed to allocate space for work1", GSL_ENOMEM);
  227. }
  228. state->work1 = work1;
  229. perm = gsl_permutation_calloc (p);
  230. if (perm == 0)
  231. {
  232. gsl_matrix_free (r);
  233. gsl_vector_free (tau);
  234. gsl_vector_free (diag);
  235. gsl_vector_free (qtf);
  236. gsl_vector_free (newton);
  237. gsl_vector_free (gradient);
  238. gsl_vector_free (x_trial);
  239. gsl_vector_free (f_trial);
  240. gsl_vector_free (df);
  241. gsl_vector_free (sdiag);
  242. gsl_vector_free (rptdx);
  243. gsl_vector_free (w);
  244. gsl_vector_free (work1);
  245. GSL_ERROR ("failed to allocate space for perm", GSL_ENOMEM);
  246. }
  247. state->perm = perm;
  248. return GSL_SUCCESS;
  249. }
  250. static int
  251. lmder_set (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
  252. {
  253. int status = set (vstate, fdf, x, f, J, dx, 0);
  254. return status ;
  255. }
  256. static int
  257. lmsder_set (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
  258. {
  259. int status = set (vstate, fdf, x, f, J, dx, 1);
  260. return status ;
  261. }
  262. static int
  263. lmder_iterate (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
  264. {
  265. int status = iterate (vstate, fdf, x, f, J, dx, 0);
  266. return status;
  267. }
  268. static int
  269. lmsder_iterate (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
  270. {
  271. int status = iterate (vstate, fdf, x, f, J, dx, 1);
  272. return status;
  273. }
  274. static void
  275. lmder_free (void *vstate)
  276. {
  277. lmder_state_t *state = (lmder_state_t *) vstate;
  278. gsl_permutation_free (state->perm);
  279. gsl_vector_free (state->work1);
  280. gsl_vector_free (state->w);
  281. gsl_vector_free (state->rptdx);
  282. gsl_vector_free (state->sdiag);
  283. gsl_vector_free (state->df);
  284. gsl_vector_free (state->f_trial);
  285. gsl_vector_free (state->x_trial);
  286. gsl_vector_free (state->gradient);
  287. gsl_vector_free (state->newton);
  288. gsl_vector_free (state->qtf);
  289. gsl_vector_free (state->diag);
  290. gsl_vector_free (state->tau);
  291. gsl_matrix_free (state->r);
  292. }
  293. static const gsl_multifit_fdfsolver_type lmder_type =
  294. {
  295. "lmder", /* name */
  296. sizeof (lmder_state_t),
  297. &lmder_alloc,
  298. &lmder_set,
  299. &lmder_iterate,
  300. &lmder_free
  301. };
  302. static const gsl_multifit_fdfsolver_type lmsder_type =
  303. {
  304. "lmsder", /* name */
  305. sizeof (lmder_state_t),
  306. &lmder_alloc,
  307. &lmsder_set,
  308. &lmsder_iterate,
  309. &lmder_free
  310. };
  311. const gsl_multifit_fdfsolver_type *gsl_multifit_fdfsolver_lmder = &lmder_type;
  312. const gsl_multifit_fdfsolver_type *gsl_multifit_fdfsolver_lmsder = &lmsder_type;