gsl_multiroots__broyden.c 9.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457
  1. /* multiroots/broyden.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_multiroots.h"
  28. #include "gsl_linalg.h"
  29. #include "gsl_multiroots__enorm.c"
  30. /* Broyden's method. It is not an efficient or modern algorithm but
  31. gives an example of a rank-1 update.
  32. C.G. Broyden, "A Class of Methods for Solving Nonlinear
  33. Simultaneous Equations", Mathematics of Computation, vol 19 (1965),
  34. p 577-593
  35. */
  36. typedef struct
  37. {
  38. gsl_matrix *H;
  39. gsl_matrix *lu;
  40. gsl_permutation *permutation;
  41. gsl_vector *v;
  42. gsl_vector *w;
  43. gsl_vector *y;
  44. gsl_vector *p;
  45. gsl_vector *fnew;
  46. gsl_vector *x_trial;
  47. double phi;
  48. }
  49. broyden_state_t;
  50. static int broyden_alloc (void *vstate, size_t n);
  51. static int broyden_set (void *vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  52. static int broyden_iterate (void *vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  53. static void broyden_free (void *vstate);
  54. static int
  55. broyden_alloc (void *vstate, size_t n)
  56. {
  57. broyden_state_t *state = (broyden_state_t *) vstate;
  58. gsl_vector *v, *w, *y, *fnew, *x_trial, *p;
  59. gsl_permutation *perm;
  60. gsl_matrix *m, *H;
  61. m = gsl_matrix_calloc (n, n);
  62. if (m == 0)
  63. {
  64. GSL_ERROR ("failed to allocate space for lu", GSL_ENOMEM);
  65. }
  66. state->lu = m;
  67. perm = gsl_permutation_calloc (n);
  68. if (perm == 0)
  69. {
  70. gsl_matrix_free (m);
  71. GSL_ERROR ("failed to allocate space for permutation", GSL_ENOMEM);
  72. }
  73. state->permutation = perm;
  74. H = gsl_matrix_calloc (n, n);
  75. if (H == 0)
  76. {
  77. gsl_permutation_free (perm);
  78. gsl_matrix_free (m);
  79. GSL_ERROR ("failed to allocate space for d", GSL_ENOMEM);
  80. }
  81. state->H = H;
  82. v = gsl_vector_calloc (n);
  83. if (v == 0)
  84. {
  85. gsl_matrix_free (H);
  86. gsl_permutation_free (perm);
  87. gsl_matrix_free (m);
  88. GSL_ERROR ("failed to allocate space for v", GSL_ENOMEM);
  89. }
  90. state->v = v;
  91. w = gsl_vector_calloc (n);
  92. if (w == 0)
  93. {
  94. gsl_vector_free (v);
  95. gsl_matrix_free (H);
  96. gsl_permutation_free (perm);
  97. gsl_matrix_free (m);
  98. GSL_ERROR ("failed to allocate space for w", GSL_ENOMEM);
  99. }
  100. state->w = w;
  101. y = gsl_vector_calloc (n);
  102. if (y == 0)
  103. {
  104. gsl_vector_free (w);
  105. gsl_vector_free (v);
  106. gsl_matrix_free (H);
  107. gsl_permutation_free (perm);
  108. gsl_matrix_free (m);
  109. GSL_ERROR ("failed to allocate space for y", GSL_ENOMEM);
  110. }
  111. state->y = y;
  112. fnew = gsl_vector_calloc (n);
  113. if (fnew == 0)
  114. {
  115. gsl_vector_free (y);
  116. gsl_vector_free (w);
  117. gsl_vector_free (v);
  118. gsl_matrix_free (H);
  119. gsl_permutation_free (perm);
  120. gsl_matrix_free (m);
  121. GSL_ERROR ("failed to allocate space for fnew", GSL_ENOMEM);
  122. }
  123. state->fnew = fnew;
  124. x_trial = gsl_vector_calloc (n);
  125. if (x_trial == 0)
  126. {
  127. gsl_vector_free (fnew);
  128. gsl_vector_free (y);
  129. gsl_vector_free (w);
  130. gsl_vector_free (v);
  131. gsl_matrix_free (H);
  132. gsl_permutation_free (perm);
  133. gsl_matrix_free (m);
  134. GSL_ERROR ("failed to allocate space for x_trial", GSL_ENOMEM);
  135. }
  136. state->x_trial = x_trial;
  137. p = gsl_vector_calloc (n);
  138. if (p == 0)
  139. {
  140. gsl_vector_free (x_trial);
  141. gsl_vector_free (fnew);
  142. gsl_vector_free (y);
  143. gsl_vector_free (w);
  144. gsl_vector_free (v);
  145. gsl_matrix_free (H);
  146. gsl_permutation_free (perm);
  147. gsl_matrix_free (m);
  148. GSL_ERROR ("failed to allocate space for p", GSL_ENOMEM);
  149. }
  150. state->p = p;
  151. return GSL_SUCCESS;
  152. }
  153. static int
  154. broyden_set (void *vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx)
  155. {
  156. broyden_state_t *state = (broyden_state_t *) vstate;
  157. size_t i, j, n = function->n;
  158. int signum = 0;
  159. GSL_MULTIROOT_FN_EVAL (function, x, f);
  160. gsl_multiroot_fdjacobian (function, x, f, GSL_SQRT_DBL_EPSILON, state->lu);
  161. gsl_linalg_LU_decomp (state->lu, state->permutation, &signum);
  162. gsl_linalg_LU_invert (state->lu, state->permutation, state->H);
  163. for (i = 0; i < n; i++)
  164. for (j = 0; j < n; j++)
  165. gsl_matrix_set(state->H,i,j,-gsl_matrix_get(state->H,i,j));
  166. for (i = 0; i < n; i++)
  167. {
  168. gsl_vector_set (dx, i, 0.0);
  169. }
  170. state->phi = enorm (f);
  171. return GSL_SUCCESS;
  172. }
  173. static int
  174. broyden_iterate (void *vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx)
  175. {
  176. broyden_state_t *state = (broyden_state_t *) vstate;
  177. double phi0, phi1, t, lambda;
  178. gsl_matrix *H = state->H;
  179. gsl_vector *p = state->p;
  180. gsl_vector *y = state->y;
  181. gsl_vector *v = state->v;
  182. gsl_vector *w = state->w;
  183. gsl_vector *fnew = state->fnew;
  184. gsl_vector *x_trial = state->x_trial;
  185. gsl_matrix *lu = state->lu;
  186. gsl_permutation *perm = state->permutation;
  187. size_t i, j, iter;
  188. size_t n = function->n;
  189. /* p = H f */
  190. for (i = 0; i < n; i++)
  191. {
  192. double sum = 0;
  193. for (j = 0; j < n; j++)
  194. {
  195. sum += gsl_matrix_get (H, i, j) * gsl_vector_get (f, j);
  196. }
  197. gsl_vector_set (p, i, sum);
  198. }
  199. t = 1;
  200. iter = 0;
  201. phi0 = state->phi;
  202. new_step:
  203. for (i = 0; i < n; i++)
  204. {
  205. double pi = gsl_vector_get (p, i);
  206. double xi = gsl_vector_get (x, i);
  207. gsl_vector_set (x_trial, i, xi + t * pi);
  208. }
  209. {
  210. int status = GSL_MULTIROOT_FN_EVAL (function, x_trial, fnew);
  211. if (status != GSL_SUCCESS)
  212. {
  213. return GSL_EBADFUNC;
  214. }
  215. }
  216. phi1 = enorm (fnew);
  217. iter++ ;
  218. if (phi1 > phi0 && iter < 10 && t > 0.1)
  219. {
  220. /* full step goes uphill, take a reduced step instead */
  221. double theta = phi1 / phi0;
  222. t *= (sqrt (1.0 + 6.0 * theta) - 1.0) / (3.0 * theta);
  223. goto new_step;
  224. }
  225. if (phi1 > phi0)
  226. {
  227. /* need to recompute Jacobian */
  228. int signum = 0;
  229. gsl_multiroot_fdjacobian (function, x, f, GSL_SQRT_DBL_EPSILON, lu);
  230. for (i = 0; i < n; i++)
  231. for (j = 0; j < n; j++)
  232. gsl_matrix_set(lu,i,j,-gsl_matrix_get(lu,i,j));
  233. gsl_linalg_LU_decomp (lu, perm, &signum);
  234. gsl_linalg_LU_invert (lu, perm, H);
  235. gsl_linalg_LU_solve (lu, perm, f, p);
  236. t = 1;
  237. for (i = 0; i < n; i++)
  238. {
  239. double pi = gsl_vector_get (p, i);
  240. double xi = gsl_vector_get (x, i);
  241. gsl_vector_set (x_trial, i, xi + t * pi);
  242. }
  243. {
  244. int status = GSL_MULTIROOT_FN_EVAL (function, x_trial, fnew);
  245. if (status != GSL_SUCCESS)
  246. {
  247. return GSL_EBADFUNC;
  248. }
  249. }
  250. phi1 = enorm (fnew);
  251. }
  252. /* y = f' - f */
  253. for (i = 0; i < n; i++)
  254. {
  255. double yi = gsl_vector_get (fnew, i) - gsl_vector_get (f, i);
  256. gsl_vector_set (y, i, yi);
  257. }
  258. /* v = H y */
  259. for (i = 0; i < n; i++)
  260. {
  261. double sum = 0;
  262. for (j = 0; j < n; j++)
  263. {
  264. sum += gsl_matrix_get (H, i, j) * gsl_vector_get (y, j);
  265. }
  266. gsl_vector_set (v, i, sum);
  267. }
  268. /* lambda = p . v */
  269. lambda = 0;
  270. for (i = 0; i < n; i++)
  271. {
  272. lambda += gsl_vector_get (p, i) * gsl_vector_get (v, i);
  273. }
  274. if (lambda == 0)
  275. {
  276. GSL_ERROR ("approximation to Jacobian has collapsed", GSL_EZERODIV) ;
  277. }
  278. /* v' = v + t * p */
  279. for (i = 0; i < n; i++)
  280. {
  281. double vi = gsl_vector_get (v, i) + t * gsl_vector_get (p, i);
  282. gsl_vector_set (v, i, vi);
  283. }
  284. /* w^T = p^T H */
  285. for (i = 0; i < n; i++)
  286. {
  287. double sum = 0;
  288. for (j = 0; j < n; j++)
  289. {
  290. sum += gsl_matrix_get (H, j, i) * gsl_vector_get (p, j);
  291. }
  292. gsl_vector_set (w, i, sum);
  293. }
  294. /* Hij -> Hij - (vi wj / lambda) */
  295. for (i = 0; i < n; i++)
  296. {
  297. double vi = gsl_vector_get (v, i);
  298. for (j = 0; j < n; j++)
  299. {
  300. double wj = gsl_vector_get (w, j);
  301. double Hij = gsl_matrix_get (H, i, j) - vi * wj / lambda;
  302. gsl_matrix_set (H, i, j, Hij);
  303. }
  304. }
  305. /* copy fnew into f */
  306. gsl_vector_memcpy (f, fnew);
  307. /* copy x_trial into x */
  308. gsl_vector_memcpy (x, x_trial);
  309. for (i = 0; i < n; i++)
  310. {
  311. double pi = gsl_vector_get (p, i);
  312. gsl_vector_set (dx, i, t * pi);
  313. }
  314. state->phi = phi1;
  315. return GSL_SUCCESS;
  316. }
  317. static void
  318. broyden_free (void *vstate)
  319. {
  320. broyden_state_t *state = (broyden_state_t *) vstate;
  321. gsl_matrix_free (state->H);
  322. gsl_matrix_free (state->lu);
  323. gsl_permutation_free (state->permutation);
  324. gsl_vector_free (state->v);
  325. gsl_vector_free (state->w);
  326. gsl_vector_free (state->y);
  327. gsl_vector_free (state->p);
  328. gsl_vector_free (state->fnew);
  329. gsl_vector_free (state->x_trial);
  330. }
  331. static const gsl_multiroot_fsolver_type broyden_type =
  332. {"broyden", /* name */
  333. sizeof (broyden_state_t),
  334. &broyden_alloc,
  335. &broyden_set,
  336. &broyden_iterate,
  337. &broyden_free};
  338. const gsl_multiroot_fsolver_type *gsl_multiroot_fsolver_broyden = &broyden_type;