gsl_multiroots__hybrid.c 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662
  1. /* multiroots/hybrid.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__dogleg.c"
  30. typedef struct
  31. {
  32. size_t iter;
  33. size_t ncfail;
  34. size_t ncsuc;
  35. size_t nslow1;
  36. size_t nslow2;
  37. double fnorm;
  38. double delta;
  39. gsl_matrix *J;
  40. gsl_matrix *q;
  41. gsl_matrix *r;
  42. gsl_vector *tau;
  43. gsl_vector *diag;
  44. gsl_vector *qtf;
  45. gsl_vector *newton;
  46. gsl_vector *gradient;
  47. gsl_vector *x_trial;
  48. gsl_vector *f_trial;
  49. gsl_vector *df;
  50. gsl_vector *qtdf;
  51. gsl_vector *rdx;
  52. gsl_vector *w;
  53. gsl_vector *v;
  54. }
  55. hybrid_state_t;
  56. static int hybrid_alloc (void *vstate, size_t n);
  57. static int hybrid_set (void *vstate, gsl_multiroot_function * func,
  58. gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  59. static int hybrids_set (void *vstate, gsl_multiroot_function * func,
  60. gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  61. static int set (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  62. gsl_vector * f, gsl_vector * dx, int scale);
  63. static int hybrid_iterate (void *vstate, gsl_multiroot_function * func,
  64. gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  65. static void hybrid_free (void *vstate);
  66. static int iterate (void *vstate, gsl_multiroot_function * func,
  67. gsl_vector * x, gsl_vector * f, gsl_vector * dx,
  68. int scale);
  69. static int
  70. hybrid_alloc (void *vstate, size_t n)
  71. {
  72. hybrid_state_t *state = (hybrid_state_t *) vstate;
  73. gsl_matrix *J, *q, *r;
  74. gsl_vector *tau, *diag, *qtf, *newton, *gradient, *x_trial, *f_trial,
  75. *df, *qtdf, *rdx, *w, *v;
  76. J = gsl_matrix_calloc (n, n);
  77. if (J == 0)
  78. {
  79. GSL_ERROR ("failed to allocate space for J", GSL_ENOMEM);
  80. }
  81. state->J = J;
  82. q = gsl_matrix_calloc (n, n);
  83. if (q == 0)
  84. {
  85. gsl_matrix_free (J);
  86. GSL_ERROR ("failed to allocate space for q", GSL_ENOMEM);
  87. }
  88. state->q = q;
  89. r = gsl_matrix_calloc (n, n);
  90. if (r == 0)
  91. {
  92. gsl_matrix_free (J);
  93. gsl_matrix_free (q);
  94. GSL_ERROR ("failed to allocate space for r", GSL_ENOMEM);
  95. }
  96. state->r = r;
  97. tau = gsl_vector_calloc (n);
  98. if (tau == 0)
  99. {
  100. gsl_matrix_free (J);
  101. gsl_matrix_free (q);
  102. gsl_matrix_free (r);
  103. GSL_ERROR ("failed to allocate space for tau", GSL_ENOMEM);
  104. }
  105. state->tau = tau;
  106. diag = gsl_vector_calloc (n);
  107. if (diag == 0)
  108. {
  109. gsl_matrix_free (J);
  110. gsl_matrix_free (q);
  111. gsl_matrix_free (r);
  112. gsl_vector_free (tau);
  113. GSL_ERROR ("failed to allocate space for diag", GSL_ENOMEM);
  114. }
  115. state->diag = diag;
  116. qtf = gsl_vector_calloc (n);
  117. if (qtf == 0)
  118. {
  119. gsl_matrix_free (J);
  120. gsl_matrix_free (q);
  121. gsl_matrix_free (r);
  122. gsl_vector_free (tau);
  123. gsl_vector_free (diag);
  124. GSL_ERROR ("failed to allocate space for qtf", GSL_ENOMEM);
  125. }
  126. state->qtf = qtf;
  127. newton = gsl_vector_calloc (n);
  128. if (newton == 0)
  129. {
  130. gsl_matrix_free (J);
  131. gsl_matrix_free (q);
  132. gsl_matrix_free (r);
  133. gsl_vector_free (tau);
  134. gsl_vector_free (diag);
  135. gsl_vector_free (qtf);
  136. GSL_ERROR ("failed to allocate space for newton", GSL_ENOMEM);
  137. }
  138. state->newton = newton;
  139. gradient = gsl_vector_calloc (n);
  140. if (gradient == 0)
  141. {
  142. gsl_matrix_free (J);
  143. gsl_matrix_free (q);
  144. gsl_matrix_free (r);
  145. gsl_vector_free (tau);
  146. gsl_vector_free (diag);
  147. gsl_vector_free (qtf);
  148. gsl_vector_free (newton);
  149. GSL_ERROR ("failed to allocate space for gradient", GSL_ENOMEM);
  150. }
  151. state->gradient = gradient;
  152. x_trial = gsl_vector_calloc (n);
  153. if (x_trial == 0)
  154. {
  155. gsl_matrix_free (J);
  156. gsl_matrix_free (q);
  157. gsl_matrix_free (r);
  158. gsl_vector_free (tau);
  159. gsl_vector_free (diag);
  160. gsl_vector_free (qtf);
  161. gsl_vector_free (newton);
  162. gsl_vector_free (gradient);
  163. GSL_ERROR ("failed to allocate space for x_trial", GSL_ENOMEM);
  164. }
  165. state->x_trial = x_trial;
  166. f_trial = gsl_vector_calloc (n);
  167. if (f_trial == 0)
  168. {
  169. gsl_matrix_free (J);
  170. gsl_matrix_free (q);
  171. gsl_matrix_free (r);
  172. gsl_vector_free (tau);
  173. gsl_vector_free (diag);
  174. gsl_vector_free (qtf);
  175. gsl_vector_free (newton);
  176. gsl_vector_free (gradient);
  177. gsl_vector_free (x_trial);
  178. GSL_ERROR ("failed to allocate space for f_trial", GSL_ENOMEM);
  179. }
  180. state->f_trial = f_trial;
  181. df = gsl_vector_calloc (n);
  182. if (df == 0)
  183. {
  184. gsl_matrix_free (J);
  185. gsl_matrix_free (q);
  186. gsl_matrix_free (r);
  187. gsl_vector_free (tau);
  188. gsl_vector_free (diag);
  189. gsl_vector_free (qtf);
  190. gsl_vector_free (newton);
  191. gsl_vector_free (gradient);
  192. gsl_vector_free (x_trial);
  193. gsl_vector_free (f_trial);
  194. GSL_ERROR ("failed to allocate space for df", GSL_ENOMEM);
  195. }
  196. state->df = df;
  197. qtdf = gsl_vector_calloc (n);
  198. if (qtdf == 0)
  199. {
  200. gsl_matrix_free (J);
  201. gsl_matrix_free (q);
  202. gsl_matrix_free (r);
  203. gsl_vector_free (tau);
  204. gsl_vector_free (diag);
  205. gsl_vector_free (qtf);
  206. gsl_vector_free (newton);
  207. gsl_vector_free (gradient);
  208. gsl_vector_free (x_trial);
  209. gsl_vector_free (f_trial);
  210. gsl_vector_free (df);
  211. GSL_ERROR ("failed to allocate space for qtdf", GSL_ENOMEM);
  212. }
  213. state->qtdf = qtdf;
  214. rdx = gsl_vector_calloc (n);
  215. if (rdx == 0)
  216. {
  217. gsl_matrix_free (J);
  218. gsl_matrix_free (q);
  219. gsl_matrix_free (r);
  220. gsl_vector_free (tau);
  221. gsl_vector_free (diag);
  222. gsl_vector_free (qtf);
  223. gsl_vector_free (newton);
  224. gsl_vector_free (gradient);
  225. gsl_vector_free (x_trial);
  226. gsl_vector_free (f_trial);
  227. gsl_vector_free (df);
  228. gsl_vector_free (qtdf);
  229. GSL_ERROR ("failed to allocate space for rdx", GSL_ENOMEM);
  230. }
  231. state->rdx = rdx;
  232. w = gsl_vector_calloc (n);
  233. if (w == 0)
  234. {
  235. gsl_matrix_free (J);
  236. gsl_matrix_free (q);
  237. gsl_matrix_free (r);
  238. gsl_vector_free (tau);
  239. gsl_vector_free (diag);
  240. gsl_vector_free (qtf);
  241. gsl_vector_free (newton);
  242. gsl_vector_free (gradient);
  243. gsl_vector_free (x_trial);
  244. gsl_vector_free (f_trial);
  245. gsl_vector_free (df);
  246. gsl_vector_free (qtdf);
  247. gsl_vector_free (rdx);
  248. GSL_ERROR ("failed to allocate space for w", GSL_ENOMEM);
  249. }
  250. state->w = w;
  251. v = gsl_vector_calloc (n);
  252. if (v == 0)
  253. {
  254. gsl_matrix_free (J);
  255. gsl_matrix_free (q);
  256. gsl_matrix_free (r);
  257. gsl_vector_free (tau);
  258. gsl_vector_free (diag);
  259. gsl_vector_free (qtf);
  260. gsl_vector_free (newton);
  261. gsl_vector_free (gradient);
  262. gsl_vector_free (x_trial);
  263. gsl_vector_free (f_trial);
  264. gsl_vector_free (df);
  265. gsl_vector_free (qtdf);
  266. gsl_vector_free (rdx);
  267. gsl_vector_free (w);
  268. GSL_ERROR ("failed to allocate space for v", GSL_ENOMEM);
  269. }
  270. state->v = v;
  271. return GSL_SUCCESS;
  272. }
  273. static int
  274. hybrid_set (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  275. gsl_vector * f, gsl_vector * dx)
  276. {
  277. int status = set (vstate, func, x, f, dx, 0);
  278. return status;
  279. }
  280. static int
  281. hybrids_set (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  282. gsl_vector * f, gsl_vector * dx)
  283. {
  284. int status = set (vstate, func, x, f, dx, 1);
  285. return status;
  286. }
  287. static int
  288. set (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  289. gsl_vector * f, gsl_vector * dx, int scale)
  290. {
  291. hybrid_state_t *state = (hybrid_state_t *) vstate;
  292. gsl_matrix *J = state->J;
  293. gsl_matrix *q = state->q;
  294. gsl_matrix *r = state->r;
  295. gsl_vector *tau = state->tau;
  296. gsl_vector *diag = state->diag;
  297. int status;
  298. status = GSL_MULTIROOT_FN_EVAL (func, x, f);
  299. if (status)
  300. {
  301. return status;
  302. }
  303. status = gsl_multiroot_fdjacobian (func, x, f, GSL_SQRT_DBL_EPSILON, J);
  304. if (status)
  305. {
  306. return status;
  307. }
  308. state->iter = 1;
  309. state->fnorm = enorm (f);
  310. state->ncfail = 0;
  311. state->ncsuc = 0;
  312. state->nslow1 = 0;
  313. state->nslow2 = 0;
  314. gsl_vector_set_all (dx, 0.0);
  315. /* Store column norms in diag */
  316. if (scale)
  317. compute_diag (J, diag);
  318. else
  319. gsl_vector_set_all (diag, 1.0);
  320. /* Set delta to factor |D x| or to factor if |D x| is zero */
  321. state->delta = compute_delta (diag, x);
  322. /* Factorize J into QR decomposition */
  323. status = gsl_linalg_QR_decomp (J, tau);
  324. if (status)
  325. {
  326. return status;
  327. }
  328. status = gsl_linalg_QR_unpack (J, tau, q, r);
  329. return status;
  330. }
  331. static int
  332. hybrid_iterate (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  333. gsl_vector * f, gsl_vector * dx)
  334. {
  335. int status = iterate (vstate, func, x, f, dx, 0);
  336. return status;
  337. }
  338. static int
  339. hybrids_iterate (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  340. gsl_vector * f, gsl_vector * dx)
  341. {
  342. int status = iterate (vstate, func, x, f, dx, 1);
  343. return status;
  344. }
  345. static int
  346. iterate (void *vstate, gsl_multiroot_function * func, gsl_vector * x,
  347. gsl_vector * f, gsl_vector * dx, int scale)
  348. {
  349. hybrid_state_t *state = (hybrid_state_t *) vstate;
  350. const double fnorm = state->fnorm;
  351. gsl_matrix *J = state->J;
  352. gsl_matrix *q = state->q;
  353. gsl_matrix *r = state->r;
  354. gsl_vector *tau = state->tau;
  355. gsl_vector *diag = state->diag;
  356. gsl_vector *qtf = state->qtf;
  357. gsl_vector *x_trial = state->x_trial;
  358. gsl_vector *f_trial = state->f_trial;
  359. gsl_vector *df = state->df;
  360. gsl_vector *qtdf = state->qtdf;
  361. gsl_vector *rdx = state->rdx;
  362. gsl_vector *w = state->w;
  363. gsl_vector *v = state->v;
  364. double prered, actred;
  365. double pnorm, fnorm1, fnorm1p;
  366. double ratio;
  367. double p1 = 0.1, p5 = 0.5, p001 = 0.001, p0001 = 0.0001;
  368. /* Compute qtf = Q^T f */
  369. compute_qtf (q, f, qtf);
  370. /* Compute dogleg step */
  371. dogleg (r, qtf, diag, state->delta, state->newton, state->gradient, dx);
  372. /* Take a trial step */
  373. compute_trial_step (x, dx, state->x_trial);
  374. pnorm = scaled_enorm (diag, dx);
  375. if (state->iter == 1)
  376. {
  377. if (pnorm < state->delta)
  378. {
  379. state->delta = pnorm;
  380. }
  381. }
  382. /* Evaluate function at x + p */
  383. {
  384. int status = GSL_MULTIROOT_FN_EVAL (func, x_trial, f_trial);
  385. if (status != GSL_SUCCESS)
  386. {
  387. return GSL_EBADFUNC;
  388. }
  389. }
  390. /* Set df = f_trial - f */
  391. compute_df (f_trial, f, df);
  392. /* Compute the scaled actual reduction */
  393. fnorm1 = enorm (f_trial);
  394. actred = compute_actual_reduction (fnorm, fnorm1);
  395. /* Compute rdx = R dx */
  396. compute_rdx (r, dx, rdx);
  397. /* Compute the scaled predicted reduction phi1p = |Q^T f + R dx| */
  398. fnorm1p = enorm_sum (qtf, rdx);
  399. prered = compute_predicted_reduction (fnorm, fnorm1p);
  400. /* Compute the ratio of the actual to predicted reduction */
  401. if (prered > 0)
  402. {
  403. ratio = actred / prered;
  404. }
  405. else
  406. {
  407. ratio = 0;
  408. }
  409. /* Update the step bound */
  410. if (ratio < p1)
  411. {
  412. state->ncsuc = 0;
  413. state->ncfail++;
  414. state->delta *= p5;
  415. }
  416. else
  417. {
  418. state->ncfail = 0;
  419. state->ncsuc++;
  420. if (ratio >= p5 || state->ncsuc > 1)
  421. state->delta = GSL_MAX (state->delta, pnorm / p5);
  422. if (fabs (ratio - 1) <= p1)
  423. state->delta = pnorm / p5;
  424. }
  425. /* Test for successful iteration */
  426. if (ratio >= p0001)
  427. {
  428. gsl_vector_memcpy (x, x_trial);
  429. gsl_vector_memcpy (f, f_trial);
  430. state->fnorm = fnorm1;
  431. state->iter++;
  432. }
  433. /* Determine the progress of the iteration */
  434. state->nslow1++;
  435. if (actred >= p001)
  436. state->nslow1 = 0;
  437. if (actred >= p1)
  438. state->nslow2 = 0;
  439. if (state->ncfail == 2)
  440. {
  441. gsl_multiroot_fdjacobian (func, x, f, GSL_SQRT_DBL_EPSILON, J);
  442. state->nslow2++;
  443. if (state->iter == 1)
  444. {
  445. if (scale)
  446. compute_diag (J, diag);
  447. state->delta = compute_delta (diag, x);
  448. }
  449. else
  450. {
  451. if (scale)
  452. update_diag (J, diag);
  453. }
  454. /* Factorize J into QR decomposition */
  455. gsl_linalg_QR_decomp (J, tau);
  456. gsl_linalg_QR_unpack (J, tau, q, r);
  457. return GSL_SUCCESS;
  458. }
  459. /* Compute qtdf = Q^T df, w = (Q^T df - R dx)/|dx|, v = D^2 dx/|dx| */
  460. compute_qtf (q, df, qtdf);
  461. compute_wv (qtdf, rdx, dx, diag, pnorm, w, v);
  462. /* Rank-1 update of the jacobian Q'R' = Q(R + w v^T) */
  463. gsl_linalg_QR_update (q, r, w, v);
  464. /* No progress as measured by jacobian evaluations */
  465. if (state->nslow2 == 5)
  466. {
  467. return GSL_ENOPROGJ;
  468. }
  469. /* No progress as measured by function evaluations */
  470. if (state->nslow1 == 10)
  471. {
  472. return GSL_ENOPROG;
  473. }
  474. return GSL_SUCCESS;
  475. }
  476. static void
  477. hybrid_free (void *vstate)
  478. {
  479. hybrid_state_t *state = (hybrid_state_t *) vstate;
  480. gsl_vector_free (state->v);
  481. gsl_vector_free (state->w);
  482. gsl_vector_free (state->rdx);
  483. gsl_vector_free (state->qtdf);
  484. gsl_vector_free (state->df);
  485. gsl_vector_free (state->f_trial);
  486. gsl_vector_free (state->x_trial);
  487. gsl_vector_free (state->gradient);
  488. gsl_vector_free (state->newton);
  489. gsl_vector_free (state->qtf);
  490. gsl_vector_free (state->diag);
  491. gsl_vector_free (state->tau);
  492. gsl_matrix_free (state->r);
  493. gsl_matrix_free (state->q);
  494. gsl_matrix_free (state->J);
  495. }
  496. static const gsl_multiroot_fsolver_type hybrid_type = {
  497. "hybrid", /* name */
  498. sizeof (hybrid_state_t),
  499. &hybrid_alloc,
  500. &hybrid_set,
  501. &hybrid_iterate,
  502. &hybrid_free
  503. };
  504. static const gsl_multiroot_fsolver_type hybrids_type = {
  505. "hybrids", /* name */
  506. sizeof (hybrid_state_t),
  507. &hybrid_alloc,
  508. &hybrids_set,
  509. &hybrids_iterate,
  510. &hybrid_free
  511. };
  512. const gsl_multiroot_fsolver_type *gsl_multiroot_fsolver_hybrid = &hybrid_type;
  513. const gsl_multiroot_fsolver_type *gsl_multiroot_fsolver_hybrids =
  514. &hybrids_type;