gsl_multimin__simplex.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464
  1. /* multimin/simplex.c
  2. *
  3. * Copyright (C) 2007 Brian Gough
  4. * Copyright (C) 2002 Tuomo Keskitalo, Ivo Alxneit
  5. *
  6. * This program is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU General Public License as published by
  8. * the Free Software Foundation; either version 3 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * This program is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  14. * General Public License for more details.
  15. *
  16. * You should have received a copy of the GNU General Public License
  17. * along with this program; if not, write to the Free Software
  18. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  19. */
  20. /*
  21. - Originally written by Tuomo Keskitalo <tuomo.keskitalo@iki.fi>
  22. - Corrections to nmsimplex_iterate and other functions
  23. by Ivo Alxneit <ivo.alxneit@psi.ch>
  24. - Additional help by Brian Gough <bjg@network-theory.co.uk>
  25. */
  26. /* The Simplex method of Nelder and Mead,
  27. also known as the polytope search alogorithm. Ref:
  28. Nelder, J.A., Mead, R., Computer Journal 7 (1965) pp. 308-313.
  29. This implementation uses n+1 corner points in the simplex.
  30. */
  31. #include "gsl__config.h"
  32. #include <stdlib.h>
  33. #include "gsl_blas.h"
  34. #include "gsl_multimin.h"
  35. typedef struct
  36. {
  37. gsl_matrix *x1; /* simplex corner points */
  38. gsl_vector *y1; /* function value at corner points */
  39. gsl_vector *ws1; /* workspace 1 for algorithm */
  40. gsl_vector *ws2; /* workspace 2 for algorithm */
  41. }
  42. nmsimplex_state_t;
  43. static double
  44. nmsimplex_move_corner (const double coeff, const nmsimplex_state_t * state,
  45. size_t corner, gsl_vector * xc,
  46. const gsl_multimin_function * f)
  47. {
  48. /* moves a simplex corner scaled by coeff (negative value represents
  49. mirroring by the middle point of the "other" corner points)
  50. and gives new corner in xc and function value at xc as a
  51. return value
  52. */
  53. gsl_matrix *x1 = state->x1;
  54. size_t i, j;
  55. double newval, mp;
  56. for (j = 0; j < x1->size2; j++)
  57. {
  58. mp = 0.0;
  59. for (i = 0; i < x1->size1; i++)
  60. {
  61. if (i != corner)
  62. {
  63. mp += (gsl_matrix_get (x1, i, j));
  64. }
  65. }
  66. mp /= (double) (x1->size1 - 1);
  67. newval = mp - coeff * (mp - gsl_matrix_get (x1, corner, j));
  68. gsl_vector_set (xc, j, newval);
  69. }
  70. newval = GSL_MULTIMIN_FN_EVAL (f, xc);
  71. return newval;
  72. }
  73. static int
  74. nmsimplex_contract_by_best (nmsimplex_state_t * state, size_t best,
  75. gsl_vector * xc, gsl_multimin_function * f)
  76. {
  77. /* Function contracts the simplex in respect to
  78. best valued corner. That is, all corners besides the
  79. best corner are moved. */
  80. /* the xc vector is simply work space here */
  81. gsl_matrix *x1 = state->x1;
  82. gsl_vector *y1 = state->y1;
  83. size_t i, j;
  84. double newval;
  85. int status = GSL_SUCCESS;
  86. for (i = 0; i < x1->size1; i++)
  87. {
  88. if (i != best)
  89. {
  90. for (j = 0; j < x1->size2; j++)
  91. {
  92. newval = 0.5 * (gsl_matrix_get (x1, i, j)
  93. + gsl_matrix_get (x1, best, j));
  94. gsl_matrix_set (x1, i, j, newval);
  95. }
  96. /* evaluate function in the new point */
  97. gsl_matrix_get_row (xc, x1, i);
  98. newval = GSL_MULTIMIN_FN_EVAL (f, xc);
  99. gsl_vector_set (y1, i, newval);
  100. /* notify caller that we found at least one bad function value.
  101. we finish the contraction (and do not abort) to allow the user
  102. to handle the situation */
  103. if(!gsl_finite(newval))
  104. {
  105. status = GSL_EBADFUNC;
  106. }
  107. }
  108. }
  109. return status;
  110. }
  111. static int
  112. nmsimplex_calc_center (const nmsimplex_state_t * state, gsl_vector * mp)
  113. {
  114. /* calculates the center of the simplex to mp */
  115. gsl_matrix *x1 = state->x1;
  116. size_t i, j;
  117. double val;
  118. for (j = 0; j < x1->size2; j++)
  119. {
  120. val = 0.0;
  121. for (i = 0; i < x1->size1; i++)
  122. {
  123. val += gsl_matrix_get (x1, i, j);
  124. }
  125. val /= x1->size1;
  126. gsl_vector_set (mp, j, val);
  127. }
  128. return GSL_SUCCESS;
  129. }
  130. static double
  131. nmsimplex_size (nmsimplex_state_t * state)
  132. {
  133. /* calculates simplex size as average sum of length of vectors
  134. from simplex center to corner points:
  135. ( sum ( || y - y_middlepoint || ) ) / n
  136. */
  137. gsl_vector *s = state->ws1;
  138. gsl_vector *mp = state->ws2;
  139. gsl_matrix *x1 = state->x1;
  140. size_t i;
  141. double ss = 0.0;
  142. /* Calculate middle point */
  143. nmsimplex_calc_center (state, mp);
  144. for (i = 0; i < x1->size1; i++)
  145. {
  146. gsl_matrix_get_row (s, x1, i);
  147. gsl_blas_daxpy (-1.0, mp, s);
  148. ss += gsl_blas_dnrm2 (s);
  149. }
  150. return ss / (double) (x1->size1);
  151. }
  152. static int
  153. nmsimplex_alloc (void *vstate, size_t n)
  154. {
  155. nmsimplex_state_t *state = (nmsimplex_state_t *) vstate;
  156. if (n == 0)
  157. {
  158. GSL_ERROR("invalid number of parameters specified", GSL_EINVAL);
  159. }
  160. state->x1 = gsl_matrix_alloc (n + 1, n);
  161. if (state->x1 == NULL)
  162. {
  163. GSL_ERROR ("failed to allocate space for x1", GSL_ENOMEM);
  164. }
  165. state->y1 = gsl_vector_alloc (n + 1);
  166. if (state->y1 == NULL)
  167. {
  168. gsl_matrix_free(state->x1);
  169. GSL_ERROR ("failed to allocate space for y", GSL_ENOMEM);
  170. }
  171. state->ws1 = gsl_vector_alloc (n);
  172. if (state->ws1 == NULL)
  173. {
  174. gsl_matrix_free(state->x1);
  175. gsl_vector_free(state->y1);
  176. GSL_ERROR ("failed to allocate space for ws1", GSL_ENOMEM);
  177. }
  178. state->ws2 = gsl_vector_alloc (n);
  179. if (state->ws2 == NULL)
  180. {
  181. gsl_matrix_free(state->x1);
  182. gsl_vector_free(state->y1);
  183. gsl_vector_free(state->ws1);
  184. GSL_ERROR ("failed to allocate space for ws2", GSL_ENOMEM);
  185. }
  186. return GSL_SUCCESS;
  187. }
  188. static int
  189. nmsimplex_set (void *vstate, gsl_multimin_function * f,
  190. const gsl_vector * x,
  191. double *size, const gsl_vector * step_size)
  192. {
  193. int status;
  194. size_t i;
  195. double val;
  196. nmsimplex_state_t *state = (nmsimplex_state_t *) vstate;
  197. gsl_vector *xtemp = state->ws1;
  198. if (xtemp->size != x->size)
  199. {
  200. GSL_ERROR("incompatible size of x", GSL_EINVAL);
  201. }
  202. if (xtemp->size != step_size->size)
  203. {
  204. GSL_ERROR("incompatible size of step_size", GSL_EINVAL);
  205. }
  206. /* first point is the original x0 */
  207. val = GSL_MULTIMIN_FN_EVAL (f, x);
  208. if (!gsl_finite(val))
  209. {
  210. GSL_ERROR("non-finite function value encountered", GSL_EBADFUNC);
  211. }
  212. gsl_matrix_set_row (state->x1, 0, x);
  213. gsl_vector_set (state->y1, 0, val);
  214. /* following points are initialized to x0 + step_size */
  215. for (i = 0; i < x->size; i++)
  216. {
  217. status = gsl_vector_memcpy (xtemp, x);
  218. if (status != 0)
  219. {
  220. GSL_ERROR ("vector memcopy failed", GSL_EFAILED);
  221. }
  222. val = gsl_vector_get (xtemp, i) + gsl_vector_get (step_size, i);
  223. gsl_vector_set (xtemp, i, val);
  224. val = GSL_MULTIMIN_FN_EVAL (f, xtemp);
  225. if (!gsl_finite(val))
  226. {
  227. GSL_ERROR("non-finite function value encountered", GSL_EBADFUNC);
  228. }
  229. gsl_matrix_set_row (state->x1, i + 1, xtemp);
  230. gsl_vector_set (state->y1, i + 1, val);
  231. }
  232. /* Initialize simplex size */
  233. *size = nmsimplex_size (state);
  234. return GSL_SUCCESS;
  235. }
  236. static void
  237. nmsimplex_free (void *vstate)
  238. {
  239. nmsimplex_state_t *state = (nmsimplex_state_t *) vstate;
  240. gsl_matrix_free (state->x1);
  241. gsl_vector_free (state->y1);
  242. gsl_vector_free (state->ws1);
  243. gsl_vector_free (state->ws2);
  244. }
  245. static int
  246. nmsimplex_iterate (void *vstate, gsl_multimin_function * f,
  247. gsl_vector * x, double *size, double *fval)
  248. {
  249. /* Simplex iteration tries to minimize function f value */
  250. /* Includes corrections from Ivo Alxneit <ivo.alxneit@psi.ch> */
  251. nmsimplex_state_t *state = (nmsimplex_state_t *) vstate;
  252. /* xc and xc2 vectors store tried corner point coordinates */
  253. gsl_vector *xc = state->ws1;
  254. gsl_vector *xc2 = state->ws2;
  255. gsl_vector *y1 = state->y1;
  256. gsl_matrix *x1 = state->x1;
  257. size_t n = y1->size;
  258. size_t i;
  259. size_t hi = 0, s_hi = 0, lo = 0;
  260. double dhi, ds_hi, dlo;
  261. int status;
  262. double val, val2;
  263. if (xc->size != x->size)
  264. {
  265. GSL_ERROR("incompatible size of x", GSL_EINVAL);
  266. }
  267. /* get index of highest, second highest and lowest point */
  268. dhi = ds_hi = dlo = gsl_vector_get (y1, 0);
  269. for (i = 1; i < n; i++)
  270. {
  271. val = (gsl_vector_get (y1, i));
  272. if (val < dlo)
  273. {
  274. dlo = val;
  275. lo = i;
  276. }
  277. else if (val > dhi)
  278. {
  279. ds_hi = dhi;
  280. s_hi = hi;
  281. dhi = val;
  282. hi = i;
  283. }
  284. else if (val > ds_hi)
  285. {
  286. ds_hi = val;
  287. s_hi = i;
  288. }
  289. }
  290. /* reflect the highest value */
  291. val = nmsimplex_move_corner (-1.0, state, hi, xc, f);
  292. if (gsl_finite(val) && val < gsl_vector_get (y1, lo))
  293. {
  294. /* reflected point becomes lowest point, try expansion */
  295. val2 = nmsimplex_move_corner (-2.0, state, hi, xc2, f);
  296. if (gsl_finite(val2) && val2 < gsl_vector_get (y1, lo))
  297. {
  298. gsl_matrix_set_row (x1, hi, xc2);
  299. gsl_vector_set (y1, hi, val2);
  300. }
  301. else
  302. {
  303. gsl_matrix_set_row (x1, hi, xc);
  304. gsl_vector_set (y1, hi, val);
  305. }
  306. }
  307. /* reflection does not improve things enough
  308. or
  309. we got a non-finite (illegal) function value */
  310. else if (!gsl_finite(val) || val > gsl_vector_get (y1, s_hi))
  311. {
  312. if (gsl_finite(val) && val <= gsl_vector_get (y1, hi))
  313. {
  314. /* if trial point is better than highest point, replace
  315. highest point */
  316. gsl_matrix_set_row (x1, hi, xc);
  317. gsl_vector_set (y1, hi, val);
  318. }
  319. /* try one dimensional contraction */
  320. val2 = nmsimplex_move_corner (0.5, state, hi, xc2, f);
  321. if (gsl_finite(val2) && val2 <= gsl_vector_get (y1, hi))
  322. {
  323. gsl_matrix_set_row (state->x1, hi, xc2);
  324. gsl_vector_set (y1, hi, val2);
  325. }
  326. else
  327. {
  328. /* contract the whole simplex in respect to the best point */
  329. status = nmsimplex_contract_by_best (state, lo, xc, f);
  330. if (status != GSL_SUCCESS)
  331. {
  332. GSL_ERROR ("nmsimplex_contract_by_best failed", GSL_EFAILED);
  333. }
  334. }
  335. }
  336. else
  337. {
  338. /* trial point is better than second highest point.
  339. Replace highest point by it */
  340. gsl_matrix_set_row (x1, hi, xc);
  341. gsl_vector_set (y1, hi, val);
  342. }
  343. /* return lowest point of simplex as x */
  344. lo = gsl_vector_min_index (y1);
  345. gsl_matrix_get_row (x, x1, lo);
  346. *fval = gsl_vector_get (y1, lo);
  347. /* Update simplex size */
  348. *size = nmsimplex_size (state);
  349. return GSL_SUCCESS;
  350. }
  351. static const gsl_multimin_fminimizer_type nmsimplex_type =
  352. { "nmsimplex", /* name */
  353. sizeof (nmsimplex_state_t),
  354. &nmsimplex_alloc,
  355. &nmsimplex_set,
  356. &nmsimplex_iterate,
  357. &nmsimplex_free
  358. };
  359. const gsl_multimin_fminimizer_type
  360. * gsl_multimin_fminimizer_nmsimplex = &nmsimplex_type;