gsl_siman__siman.c 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274
  1. /* siman/siman.c
  2. *
  3. * Copyright (C) 2007 Brian Gough
  4. * Copyright (C) 1996, 1997, 1998, 1999, 2000 Mark Galassi
  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. #include "gsl__config.h"
  21. #include <stdio.h>
  22. #include <math.h>
  23. #include <stdlib.h>
  24. #include <string.h>
  25. #include <assert.h>
  26. #include "gsl_machine.h"
  27. #include "gsl_rng.h"
  28. #include "gsl_siman.h"
  29. static inline double
  30. boltzmann(double E, double new_E, double T, gsl_siman_params_t *params)
  31. {
  32. double x = -(new_E - E) / (params->k * T);
  33. /* avoid underflow errors for large uphill steps */
  34. return (x < GSL_LOG_DBL_MIN) ? 0.0 : exp(x);
  35. }
  36. static inline void
  37. copy_state(void *src, void *dst, size_t size, gsl_siman_copy_t copyfunc)
  38. {
  39. if (copyfunc) {
  40. copyfunc(src, dst);
  41. } else {
  42. memcpy(dst, src, size);
  43. }
  44. }
  45. /* implementation of a basic simulated annealing algorithm */
  46. void
  47. gsl_siman_solve (const gsl_rng * r, void *x0_p, gsl_siman_Efunc_t Ef,
  48. gsl_siman_step_t take_step,
  49. gsl_siman_metric_t distance,
  50. gsl_siman_print_t print_position,
  51. gsl_siman_copy_t copyfunc,
  52. gsl_siman_copy_construct_t copy_constructor,
  53. gsl_siman_destroy_t destructor,
  54. size_t element_size,
  55. gsl_siman_params_t params)
  56. {
  57. void *x, *new_x, *best_x;
  58. double E, new_E, best_E;
  59. int i;
  60. double T, T_factor;
  61. int n_evals = 1, n_iter = 0, n_accepts, n_rejects, n_eless;
  62. /* this function requires that either the dynamic functions (copy,
  63. copy_constructor and destrcutor) are passed, or that an element
  64. size is given */
  65. assert((copyfunc != NULL && copy_constructor != NULL && destructor != NULL)
  66. || (element_size != 0));
  67. distance = 0 ; /* This parameter is not currently used */
  68. E = Ef(x0_p);
  69. if (copyfunc) {
  70. x = copy_constructor(x0_p);
  71. new_x = copy_constructor(x0_p);
  72. best_x = copy_constructor(x0_p);
  73. } else {
  74. x = (void *) malloc (element_size);
  75. memcpy (x, x0_p, element_size);
  76. new_x = (void *) malloc (element_size);
  77. best_x = (void *) malloc (element_size);
  78. memcpy (best_x, x0_p, element_size);
  79. }
  80. best_E = E;
  81. T = params.t_initial;
  82. T_factor = 1.0 / params.mu_t;
  83. if (print_position) {
  84. printf ("#-iter #-evals temperature position energy\n");
  85. }
  86. while (1) {
  87. n_accepts = 0;
  88. n_rejects = 0;
  89. n_eless = 0;
  90. for (i = 0; i < params.iters_fixed_T; ++i) {
  91. copy_state(x, new_x, element_size, copyfunc);
  92. take_step (r, new_x, params.step_size);
  93. new_E = Ef (new_x);
  94. if(new_E <= best_E){
  95. if (copyfunc) {
  96. copyfunc(new_x,best_x);
  97. } else {
  98. memcpy (best_x, new_x, element_size);
  99. }
  100. best_E=new_E;
  101. }
  102. ++n_evals; /* keep track of Ef() evaluations */
  103. /* now take the crucial step: see if the new point is accepted
  104. or not, as determined by the boltzmann probability */
  105. if (new_E < E) {
  106. if (new_E < best_E) {
  107. copy_state(new_x, best_x, element_size, copyfunc);
  108. best_E = new_E;
  109. }
  110. /* yay! take a step */
  111. copy_state(new_x, x, element_size, copyfunc);
  112. E = new_E;
  113. ++n_eless;
  114. } else if (gsl_rng_uniform(r) < boltzmann(E, new_E, T, &params)) {
  115. /* yay! take a step */
  116. copy_state(new_x, x, element_size, copyfunc);
  117. E = new_E;
  118. ++n_accepts;
  119. } else {
  120. ++n_rejects;
  121. }
  122. }
  123. if (print_position) {
  124. /* see if we need to print stuff as we go */
  125. /* printf("%5d %12g %5d %3d %3d %3d", n_iter, T, n_evals, */
  126. /* 100*n_eless/n_steps, 100*n_accepts/n_steps, */
  127. /* 100*n_rejects/n_steps); */
  128. printf ("%5d %7d %12g", n_iter, n_evals, T);
  129. print_position (x);
  130. printf (" %12g %12g\n", E, best_E);
  131. }
  132. /* apply the cooling schedule to the temperature */
  133. /* FIXME: I should also introduce a cooling schedule for the iters */
  134. T *= T_factor;
  135. ++n_iter;
  136. if (T < params.t_min) {
  137. break;
  138. }
  139. }
  140. /* at the end, copy the result onto the initial point, so we pass it
  141. back to the caller */
  142. copy_state(best_x, x0_p, element_size, copyfunc);
  143. if (copyfunc) {
  144. destructor(x);
  145. destructor(new_x);
  146. destructor(best_x);
  147. } else {
  148. free (x);
  149. free (new_x);
  150. free (best_x);
  151. }
  152. }
  153. /* implementation of a simulated annealing algorithm with many tries */
  154. void
  155. gsl_siman_solve_many (const gsl_rng * r, void *x0_p, gsl_siman_Efunc_t Ef,
  156. gsl_siman_step_t take_step,
  157. gsl_siman_metric_t distance,
  158. gsl_siman_print_t print_position,
  159. size_t element_size,
  160. gsl_siman_params_t params)
  161. {
  162. /* the new set of trial points, and their energies and probabilities */
  163. void *x, *new_x;
  164. double *energies, *probs, *sum_probs;
  165. double Ex; /* energy of the chosen point */
  166. double T, T_factor; /* the temperature and a step multiplier */
  167. int i;
  168. double u; /* throw the die to choose a new "x" */
  169. int n_iter;
  170. if (print_position) {
  171. printf ("#-iter temperature position");
  172. printf (" delta_pos energy\n");
  173. }
  174. x = (void *) malloc (params.n_tries * element_size);
  175. new_x = (void *) malloc (params.n_tries * element_size);
  176. energies = (double *) malloc (params.n_tries * sizeof (double));
  177. probs = (double *) malloc (params.n_tries * sizeof (double));
  178. sum_probs = (double *) malloc (params.n_tries * sizeof (double));
  179. T = params.t_initial;
  180. T_factor = 1.0 / params.mu_t;
  181. memcpy (x, x0_p, element_size);
  182. n_iter = 0;
  183. while (1)
  184. {
  185. Ex = Ef (x);
  186. for (i = 0; i < params.n_tries - 1; ++i)
  187. { /* only go to N_TRIES-2 */
  188. /* center the new_x[] around x, then pass it to take_step() */
  189. sum_probs[i] = 0;
  190. memcpy ((char *)new_x + i * element_size, x, element_size);
  191. take_step (r, (char *)new_x + i * element_size, params.step_size);
  192. energies[i] = Ef ((char *)new_x + i * element_size);
  193. probs[i] = boltzmann(Ex, energies[i], T, &params);
  194. }
  195. /* now add in the old value of "x", so it is a contendor */
  196. memcpy ((char *)new_x + (params.n_tries - 1) * element_size, x, element_size);
  197. energies[params.n_tries - 1] = Ex;
  198. probs[params.n_tries - 1] = boltzmann(Ex, energies[i], T, &params);
  199. /* now throw biased die to see which new_x[i] we choose */
  200. sum_probs[0] = probs[0];
  201. for (i = 1; i < params.n_tries; ++i)
  202. {
  203. sum_probs[i] = sum_probs[i - 1] + probs[i];
  204. }
  205. u = gsl_rng_uniform (r) * sum_probs[params.n_tries - 1];
  206. for (i = 0; i < params.n_tries; ++i)
  207. {
  208. if (u < sum_probs[i])
  209. {
  210. memcpy (x, (char *) new_x + i * element_size, element_size);
  211. break;
  212. }
  213. }
  214. if (print_position)
  215. {
  216. printf ("%5d\t%12g\t", n_iter, T);
  217. print_position (x);
  218. printf ("\t%12g\t%12g\n", distance (x, x0_p), Ex);
  219. }
  220. T *= T_factor;
  221. ++n_iter;
  222. if (T < params.t_min)
  223. {
  224. break;
  225. }
  226. }
  227. /* now return the value via x0_p */
  228. memcpy (x0_p, x, element_size);
  229. /* printf("the result is: %g (E=%g)\n", x, Ex); */
  230. free (x);
  231. free (new_x);
  232. free (energies);
  233. free (probs);
  234. free (sum_probs);
  235. }