gsl_multiroots__dogleg.c 9.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414
  1. /* multiroots/dogleg.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_multiroots__enorm.c"
  20. static void compute_diag (const gsl_matrix * J, gsl_vector * diag);
  21. static void update_diag (const gsl_matrix * J, gsl_vector * diag);
  22. static double compute_delta (gsl_vector * diag, gsl_vector * x);
  23. static void compute_df (const gsl_vector * f_trial, const gsl_vector * f, gsl_vector * df);
  24. static void compute_wv (const gsl_vector * qtdf, const gsl_vector *rdx, const gsl_vector *dx, const gsl_vector *diag, double pnorm, gsl_vector * w, gsl_vector * v);
  25. static double scaled_enorm (const gsl_vector * d, const gsl_vector * f);
  26. static double scaled_enorm (const gsl_vector * d, const gsl_vector * f) {
  27. double e2 = 0 ;
  28. size_t i, n = f->size ;
  29. for (i = 0; i < n ; i++) {
  30. double fi= gsl_vector_get(f, i);
  31. double di= gsl_vector_get(d, i);
  32. double u = di * fi;
  33. e2 += u * u ;
  34. }
  35. return sqrt(e2);
  36. }
  37. static double enorm_sum (const gsl_vector * a, const gsl_vector * b);
  38. static double enorm_sum (const gsl_vector * a, const gsl_vector * b) {
  39. double e2 = 0 ;
  40. size_t i, n = a->size ;
  41. for (i = 0; i < n ; i++) {
  42. double ai= gsl_vector_get(a, i);
  43. double bi= gsl_vector_get(b, i);
  44. double u = ai + bi;
  45. e2 += u * u ;
  46. }
  47. return sqrt(e2);
  48. }
  49. static void
  50. compute_wv (const gsl_vector * qtdf, const gsl_vector *rdx, const gsl_vector *dx, const gsl_vector *diag, double pnorm, gsl_vector * w, gsl_vector * v)
  51. {
  52. size_t i, n = qtdf->size;
  53. for (i = 0; i < n; i++)
  54. {
  55. double qtdfi = gsl_vector_get (qtdf, i);
  56. double rdxi = gsl_vector_get (rdx, i);
  57. double dxi = gsl_vector_get (dx, i);
  58. double diagi = gsl_vector_get (diag, i);
  59. gsl_vector_set (w, i, (qtdfi - rdxi) / pnorm);
  60. gsl_vector_set (v, i, diagi * diagi * dxi / pnorm);
  61. }
  62. }
  63. static void
  64. compute_df (const gsl_vector * f_trial, const gsl_vector * f, gsl_vector * df)
  65. {
  66. size_t i, n = f->size;
  67. for (i = 0; i < n; i++)
  68. {
  69. double dfi = gsl_vector_get (f_trial, i) - gsl_vector_get (f, i);
  70. gsl_vector_set (df, i, dfi);
  71. }
  72. }
  73. static void
  74. compute_diag (const gsl_matrix * J, gsl_vector * diag)
  75. {
  76. size_t i, j, n = diag->size;
  77. for (j = 0; j < n; j++)
  78. {
  79. double sum = 0;
  80. for (i = 0; i < n; i++)
  81. {
  82. double Jij = gsl_matrix_get (J, i, j);
  83. sum += Jij * Jij;
  84. }
  85. if (sum == 0)
  86. sum = 1.0;
  87. gsl_vector_set (diag, j, sqrt (sum));
  88. }
  89. }
  90. static void
  91. update_diag (const gsl_matrix * J, gsl_vector * diag)
  92. {
  93. size_t i, j, n = diag->size;
  94. for (j = 0; j < n; j++)
  95. {
  96. double cnorm, diagj, sum = 0;
  97. for (i = 0; i < n; i++)
  98. {
  99. double Jij = gsl_matrix_get (J, i, j);
  100. sum += Jij * Jij;
  101. }
  102. if (sum == 0)
  103. sum = 1.0;
  104. cnorm = sqrt (sum);
  105. diagj = gsl_vector_get (diag, j);
  106. if (cnorm > diagj)
  107. gsl_vector_set (diag, j, cnorm);
  108. }
  109. }
  110. static double
  111. compute_delta (gsl_vector * diag, gsl_vector * x)
  112. {
  113. double Dx = scaled_enorm (diag, x);
  114. double factor = 100;
  115. return (Dx > 0) ? factor * Dx : factor;
  116. }
  117. static double
  118. compute_actual_reduction (double fnorm, double fnorm1)
  119. {
  120. double actred;
  121. if (fnorm1 < fnorm)
  122. {
  123. double u = fnorm1 / fnorm;
  124. actred = 1 - u * u;
  125. }
  126. else
  127. {
  128. actred = -1;
  129. }
  130. return actred;
  131. }
  132. static double
  133. compute_predicted_reduction (double fnorm, double fnorm1)
  134. {
  135. double prered;
  136. if (fnorm1 < fnorm)
  137. {
  138. double u = fnorm1 / fnorm;
  139. prered = 1 - u * u;
  140. }
  141. else
  142. {
  143. prered = 0;
  144. }
  145. return prered;
  146. }
  147. static void
  148. compute_qtf (const gsl_matrix * q, const gsl_vector * f, gsl_vector * qtf)
  149. {
  150. size_t i, j, N = f->size ;
  151. for (j = 0; j < N; j++)
  152. {
  153. double sum = 0;
  154. for (i = 0; i < N; i++)
  155. sum += gsl_matrix_get (q, i, j) * gsl_vector_get (f, i);
  156. gsl_vector_set (qtf, j, sum);
  157. }
  158. }
  159. static void
  160. compute_rdx (const gsl_matrix * r, const gsl_vector * dx, gsl_vector * rdx)
  161. {
  162. size_t i, j, N = dx->size ;
  163. for (i = 0; i < N; i++)
  164. {
  165. double sum = 0;
  166. for (j = i; j < N; j++)
  167. {
  168. sum += gsl_matrix_get (r, i, j) * gsl_vector_get (dx, j);
  169. }
  170. gsl_vector_set (rdx, i, sum);
  171. }
  172. }
  173. static void
  174. compute_trial_step (gsl_vector *x, gsl_vector * dx, gsl_vector * x_trial)
  175. {
  176. size_t i, N = x->size;
  177. for (i = 0; i < N; i++)
  178. {
  179. double pi = gsl_vector_get (dx, i);
  180. double xi = gsl_vector_get (x, i);
  181. gsl_vector_set (x_trial, i, xi + pi);
  182. }
  183. }
  184. static int
  185. newton_direction (const gsl_matrix * r, const gsl_vector * qtf, gsl_vector * p)
  186. {
  187. const size_t N = r->size2;
  188. size_t i;
  189. int status;
  190. status = gsl_linalg_R_solve (r, qtf, p);
  191. #ifdef DEBUG
  192. printf("rsolve status = %d\n", status);
  193. #endif
  194. for (i = 0; i < N; i++)
  195. {
  196. double pi = gsl_vector_get (p, i);
  197. gsl_vector_set (p, i, -pi);
  198. }
  199. return status;
  200. }
  201. static void
  202. gradient_direction (const gsl_matrix * r, const gsl_vector * qtf,
  203. const gsl_vector * diag, gsl_vector * g)
  204. {
  205. const size_t M = r->size1;
  206. const size_t N = r->size2;
  207. size_t i, j;
  208. for (j = 0; j < M; j++)
  209. {
  210. double sum = 0;
  211. double dj;
  212. for (i = 0; i < N; i++)
  213. {
  214. sum += gsl_matrix_get (r, i, j) * gsl_vector_get (qtf, i);
  215. }
  216. dj = gsl_vector_get (diag, j);
  217. gsl_vector_set (g, j, -sum / dj);
  218. }
  219. }
  220. static void
  221. minimum_step (double gnorm, const gsl_vector * diag, gsl_vector * g)
  222. {
  223. const size_t N = g->size;
  224. size_t i;
  225. for (i = 0; i < N; i++)
  226. {
  227. double gi = gsl_vector_get (g, i);
  228. double di = gsl_vector_get (diag, i);
  229. gsl_vector_set (g, i, (gi / gnorm) / di);
  230. }
  231. }
  232. static void
  233. compute_Rg (const gsl_matrix * r, const gsl_vector * gradient, gsl_vector * Rg)
  234. {
  235. const size_t N = r->size2;
  236. size_t i, j;
  237. for (i = 0; i < N; i++)
  238. {
  239. double sum = 0;
  240. for (j = i; j < N; j++)
  241. {
  242. double gj = gsl_vector_get (gradient, j);
  243. double rij = gsl_matrix_get (r, i, j);
  244. sum += rij * gj;
  245. }
  246. gsl_vector_set (Rg, i, sum);
  247. }
  248. }
  249. static void
  250. scaled_addition (double alpha, gsl_vector * newton, double beta, gsl_vector * gradient, gsl_vector * p)
  251. {
  252. const size_t N = p->size;
  253. size_t i;
  254. for (i = 0; i < N; i++)
  255. {
  256. double ni = gsl_vector_get (newton, i);
  257. double gi = gsl_vector_get (gradient, i);
  258. gsl_vector_set (p, i, alpha * ni + beta * gi);
  259. }
  260. }
  261. static int
  262. dogleg (const gsl_matrix * r, const gsl_vector * qtf,
  263. const gsl_vector * diag, double delta,
  264. gsl_vector * newton, gsl_vector * gradient, gsl_vector * p)
  265. {
  266. double qnorm, gnorm, sgnorm, bnorm, temp;
  267. newton_direction (r, qtf, newton);
  268. #ifdef DEBUG
  269. printf("newton = "); gsl_vector_fprintf(stdout, newton, "%g"); printf("\n");
  270. #endif
  271. qnorm = scaled_enorm (diag, newton);
  272. if (qnorm <= delta)
  273. {
  274. gsl_vector_memcpy (p, newton);
  275. #ifdef DEBUG
  276. printf("took newton (qnorm = %g <= delta = %g)\n", qnorm, delta);
  277. #endif
  278. return GSL_SUCCESS;
  279. }
  280. gradient_direction (r, qtf, diag, gradient);
  281. #ifdef DEBUG
  282. printf("grad = "); gsl_vector_fprintf(stdout, gradient, "%g"); printf("\n");
  283. #endif
  284. gnorm = enorm (gradient);
  285. if (gnorm == 0)
  286. {
  287. double alpha = delta / qnorm;
  288. double beta = 0;
  289. scaled_addition (alpha, newton, beta, gradient, p);
  290. #ifdef DEBUG
  291. printf("took scaled newton because gnorm = 0\n");
  292. #endif
  293. return GSL_SUCCESS;
  294. }
  295. minimum_step (gnorm, diag, gradient);
  296. compute_Rg (r, gradient, p); /* Use p as temporary space to compute Rg */
  297. #ifdef DEBUG
  298. printf("mingrad = "); gsl_vector_fprintf(stdout, gradient, "%g"); printf("\n");
  299. printf("Rg = "); gsl_vector_fprintf(stdout, p, "%g"); printf("\n");
  300. #endif
  301. temp = enorm (p);
  302. sgnorm = (gnorm / temp) / temp;
  303. if (sgnorm > delta)
  304. {
  305. double alpha = 0;
  306. double beta = delta;
  307. scaled_addition (alpha, newton, beta, gradient, p);
  308. #ifdef DEBUG
  309. printf("took gradient\n");
  310. #endif
  311. return GSL_SUCCESS;
  312. }
  313. bnorm = enorm (qtf);
  314. {
  315. double bg = bnorm / gnorm;
  316. double bq = bnorm / qnorm;
  317. double dq = delta / qnorm;
  318. double dq2 = dq * dq;
  319. double sd = sgnorm / delta;
  320. double sd2 = sd * sd;
  321. double t1 = bg * bq * sd;
  322. double u = t1 - dq;
  323. double t2 = t1 - dq * sd2 + sqrt (u * u + (1-dq2) * (1 - sd2));
  324. double alpha = dq * (1 - sd2) / t2;
  325. double beta = (1 - alpha) * sgnorm;
  326. #ifdef DEBUG
  327. printf("bnorm = %g\n", bnorm);
  328. printf("gnorm = %g\n", gnorm);
  329. printf("qnorm = %g\n", qnorm);
  330. printf("delta = %g\n", delta);
  331. printf("alpha = %g beta = %g\n", alpha, beta);
  332. printf("took scaled combination of newton and gradient\n");
  333. #endif
  334. scaled_addition (alpha, newton, beta, gradient, p);
  335. }
  336. return GSL_SUCCESS;
  337. }