rgsl/types/siman.rs
1//
2// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
3//
4
5/*!
6# 25 Simulated Annealing
7
8Stochastic search techniques are used when the structure of a space is not well understood or
9is not smooth, so that techniques like Newton’s method (which requires calculating Jacobian
10derivative matrices) cannot be used. In particular, these techniques are frequently used to
11solve combinatorial optimization problems, such as the traveling salesman problem.
12
13The goal is to find a point in the space at which a real valued energy function (or cost
14function) is minimized. Simulated annealing is a minimization technique which has given
15good results in avoiding local minima; it is based on the idea of taking a random walk
16through the space at successively lower temperatures, where the probability of taking a
17step is given by a Boltzmann distribution.
18
19The functions described in this chapter are declared in the header file gsl_siman.h.
20
21## Simulated Annealing algorithm
22
23The simulated annealing algorithm takes random walks through the problem space, looking
24for points with low energies; in these random walks, the probability of taking a step is
25determined by the Boltzmann distribution,
26
27> p = e −(E i+1 −E i )/(kT )
28> if E i+1 > E i , and p = 1 when E i+1 ≤ E i
29
30In other words, a step will occur if the new energy is lower. If the new energy is higher,
31the transition can still occur, and its likelihood is proportional to the temperature T and
32inversely proportional to the energy difference E i+1 − E i .
33
34The temperature T is initially set to a high value, and a random walk is carried out
35at that temperature. Then the temperature is lowered very slightly according to a cooling
36schedule, for example: T → T /μ T where μ T is slightly greater than 1.
37
38The slight probability of taking a step that gives higher energy is what allows simulated
39annealing to frequently get out of local minima.
40!*/
41
42const GSL_LOG_DBL_MIN: f64 = -7.0839641853226408e+02;
43
44pub struct SimAnnealing<T: Clone> {
45 x0_p: T,
46 params: SimAnnealingParams,
47 Efunc_t: gsl_siman_Efunc_t<T>,
48 step_t: gsl_siman_step_t<T>,
49 metric_t: gsl_siman_metric_t<T>,
50 print_t: Option<gsl_siman_print_t<T>>,
51}
52
53type gsl_siman_Efunc_t<T> = fn(&T) -> f64;
54type gsl_siman_step_t<T> = fn(&mut crate::Rng, &mut T, f64);
55type gsl_siman_metric_t<T> = fn(&T, &T) -> f64;
56type gsl_siman_print_t<T> = fn(&T);
57
58fn boltzmann(E: f64, new_E: f64, T: f64, params: &SimAnnealingParams) -> f64 {
59 let x = -(new_E - E) / (params.k * T);
60 // avoid underflow errors for large uphill steps
61 if x < GSL_LOG_DBL_MIN {
62 0.0
63 } else {
64 x.exp()
65 }
66}
67
68impl<T> SimAnnealing<T>
69where
70 T: Clone,
71{
72 pub fn new(
73 x0_p: T,
74 ef: gsl_siman_Efunc_t<T>,
75 take_step: gsl_siman_step_t<T>,
76 distance: gsl_siman_metric_t<T>,
77 print_pos: Option<gsl_siman_print_t<T>>,
78 params: SimAnnealingParams,
79 ) -> SimAnnealing<T> {
80 SimAnnealing {
81 x0_p,
82 params,
83 Efunc_t: ef,
84 step_t: take_step,
85 metric_t: distance,
86 print_t: print_pos,
87 }
88 }
89
90 /*
91 /* implementation of a basic simulated annealing algorithm */
92
93 void
94 gsl_siman_solve (const gsl_rng * r, void *x0_p, gsl_siman_Efunc_t Ef,
95 gsl_siman_step_t take_step,
96 gsl_siman_metric_t distance,
97 gsl_siman_print_t print_position,
98 gsl_siman_copy_t copyfunc,
99 gsl_siman_copy_construct_t copy_constructor,
100 gsl_siman_destroy_t destructor,
101 size_t element_size,
102 gsl_siman_params_t params)
103 {
104 void *x, *new_x, *best_x;
105 double E, new_E, best_E;
106 int i;
107 double T, T_factor;
108 int n_evals = 1, n_iter = 0, n_accepts, n_rejects, n_eless;
109
110 /* this function requires that either the dynamic functions (copy,
111 copy_constructor and destrcutor) are passed, or that an element
112 size is given */
113 assert((copyfunc != NULL && copy_constructor != NULL && destructor != NULL)
114 || (element_size != 0));
115
116 distance = 0 ; /* This parameter is not currently used */
117 E = Ef(x0_p);
118
119 if (copyfunc) {
120 x = copy_constructor(x0_p);
121 new_x = copy_constructor(x0_p);
122 best_x = copy_constructor(x0_p);
123 } else {
124 x = (void *) malloc (element_size);
125 memcpy (x, x0_p, element_size);
126 new_x = (void *) malloc (element_size);
127 best_x = (void *) malloc (element_size);
128 memcpy (best_x, x0_p, element_size);
129 }
130
131 best_E = E;
132
133 T = params.t_initial;
134 T_factor = 1.0 / params.mu_t;
135
136 if (print_position) {
137 printf ("#-iter #-evals temperature position energy\n");
138 }
139
140 while (1) {
141
142 n_accepts = 0;
143 n_rejects = 0;
144 n_eless = 0;
145
146 for (i = 0; i < params.iters_fixed_T; ++i) {
147
148 copy_state(x, new_x, element_size, copyfunc);
149
150 take_step (r, new_x, params.step_size);
151 new_E = Ef (new_x);
152
153 if(new_E <= best_E){
154 if (copyfunc) {
155 copyfunc(new_x,best_x);
156 } else {
157 memcpy (best_x, new_x, element_size);
158 }
159 best_E=new_E;
160 }
161
162 ++n_evals; /* keep track of Ef() evaluations */
163 /* now take the crucial step: see if the new point is accepted
164 or not, as determined by the boltzmann probability */
165 if (new_E < E) {
166 if (new_E < best_E) {
167 copy_state(new_x, best_x, element_size, copyfunc);
168 best_E = new_E;
169 }
170
171 /* yay! take a step */
172 copy_state(new_x, x, element_size, copyfunc);
173 E = new_E;
174 ++n_eless;
175
176 } else if (gsl_rng_uniform(r) < boltzmann(E, new_E, T, ¶ms)) {
177 /* yay! take a step */
178 copy_state(new_x, x, element_size, copyfunc);
179 E = new_E;
180 ++n_accepts;
181
182 } else {
183 ++n_rejects;
184 }
185 }
186
187 if (print_position) {
188 /* see if we need to print stuff as we go */
189 /* printf("%5d %12g %5d %3d %3d %3d", n_iter, T, n_evals, */
190 /* 100*n_eless/n_steps, 100*n_accepts/n_steps, */
191 /* 100*n_rejects/n_steps); */
192 printf ("%5d %7d %12g", n_iter, n_evals, T);
193 print_position (x);
194 printf (" %12g %12g\n", E, best_E);
195 }
196
197 /* apply the cooling schedule to the temperature */
198 /* FIXME: I should also introduce a cooling schedule for the iters */
199 T *= T_factor;
200 ++n_iter;
201 if (T < params.t_min) {
202 break;
203 }
204 }
205
206 /* at the end, copy the result onto the initial point, so we pass it
207 back to the caller */
208 copy_state(best_x, x0_p, element_size, copyfunc);
209
210 if (copyfunc) {
211 destructor(x);
212 destructor(new_x);
213 destructor(best_x);
214 } else {
215 free (x);
216 free (new_x);
217 free (best_x);
218 }
219 }
220 */
221 /// This function performs a simulated annealing search through a given space. The space
222 /// is specified by providing the functions Ef and distance. The simulated annealing steps
223 /// are generated using the random number generator `rng` and the function `take_step`.
224 ///
225 /// The starting configuration of the system should be given by x0\_p.
226 ///
227 /// The params structure (described below) controls the run by providing the temperature
228 /// schedule and other tunable parameters to the algorithm.
229 ///
230 /// On exit the best result achieved during the search is returned. If the annealing
231 /// process has been successful this should be a good approximation to the optimal point
232 /// in the space.
233 ///
234 /// If the argument `print_pos` is not None, a debugging log will be printed to
235 /// stdout with the following columns: ```#-iter #-evals temperature position energy best_energy```
236 /// and the output of the function print position itself.
237 pub fn solve(&self, rng: &mut crate::Rng) -> T {
238 let mut x = self.x0_p.clone();
239 let mut new_x = self.x0_p.clone();
240 let mut best_x = self.x0_p.clone();
241
242 let mut n_evals = 0_usize;
243
244 let mut E = (self.Efunc_t)(&self.x0_p);
245 let mut best_E = E;
246
247 let mut Temp = self.params.t_initial;
248 let Temp_factor = 1.0 / self.params.mu_t;
249
250 if self.print_t.is_some() {
251 println!(
252 "{i:^6} | {e:^7} | {t:^12} | {p:^15} | {E:^13}",
253 i = "#-iter",
254 e = "#-evals",
255 t = "temperature",
256 p = "position",
257 E = "energy"
258 );
259 }
260
261 let mut iter = 0;
262 loop {
263 //let mut n_accepts = 0;
264 //let mut n_rejects = 0;
265 //let mut n_eless = 0;
266
267 for _ in 0..self.params.iters_fixed_T {
268 x = new_x.clone();
269
270 (self.step_t)(rng, &mut new_x, self.params.step_size);
271 let new_E = (self.Efunc_t)(&new_x);
272 n_evals += 1; // keep track of Ef() evaluations
273
274 if new_E <= best_E {
275 best_x = new_x.clone();
276 best_E = new_E;
277 }
278
279 // now take the crucial step: see if the new point is accepted
280 // or not, as determined by the boltzmann probability
281 if new_E < E {
282 if new_E < best_E {
283 best_x = new_x.clone();
284 best_E = new_E;
285 }
286 // yay! take a step
287 x = new_x.clone();
288 E = new_E;
289 //n_eless += 1;
290 } else if rng.uniform() < boltzmann(E, new_E, Temp, &self.params) {
291 // yay! take a step
292 x = new_x.clone();
293 E = new_E;
294 //n_accepts += 1;
295 }
296 }
297
298 if let Some(ref printf) = self.print_t {
299 print!("{:>06} | {:>07} | {:>12.10} | ", iter, n_evals, Temp);
300 printf(&x);
301 println!(" | {:+>13.12}", E);
302 }
303
304 Temp *= Temp_factor;
305 iter += 1;
306 if Temp < self.params.t_min {
307 break;
308 }
309 }
310
311 best_x
312 }
313
314 /*
315 /* implementation of a simulated annealing algorithm with many tries */
316
317 void
318 gsl_siman_solve_many (const gsl_rng * r, void *x0_p, gsl_siman_Efunc_t Ef,
319 gsl_siman_step_t take_step,
320 gsl_siman_metric_t distance,
321 gsl_siman_print_t print_position,
322 size_t element_size,
323 gsl_siman_params_t params)
324 {
325 /* the new set of trial points, and their energies and probabilities */
326 void *x, *new_x;
327 double *energies, *probs, *sum_probs;
328 double Ex; /* energy of the chosen point */
329 double T, T_factor; /* the temperature and a step multiplier */
330 int i;
331 double u; /* throw the die to choose a new "x" */
332 int n_iter;
333
334 if (print_position) {
335 printf ("#-iter temperature position");
336 printf (" delta_pos energy\n");
337 }
338
339 x = (void *) malloc (params.n_tries * element_size);
340 new_x = (void *) malloc (params.n_tries * element_size);
341 energies = (double *) malloc (params.n_tries * sizeof (double));
342 probs = (double *) malloc (params.n_tries * sizeof (double));
343 sum_probs = (double *) malloc (params.n_tries * sizeof (double));
344
345 T = params.t_initial;
346 T_factor = 1.0 / params.mu_t;
347
348 memcpy (x, x0_p, element_size);
349
350 n_iter = 0;
351 while (1)
352 {
353 Ex = Ef (x);
354 for (i = 0; i < params.n_tries - 1; ++i)
355 { /* only go to N_TRIES-2 */
356 /* center the new_x[] around x, then pass it to take_step() */
357 sum_probs[i] = 0;
358 memcpy ((char *)new_x + i * element_size, x, element_size);
359 take_step (r, (char *)new_x + i * element_size, params.step_size);
360 energies[i] = Ef ((char *)new_x + i * element_size);
361 probs[i] = boltzmann(Ex, energies[i], T, ¶ms);
362 }
363 /* now add in the old value of "x", so it is a contendor */
364 memcpy ((char *)new_x + (params.n_tries - 1) * element_size, x, element_size);
365 energies[params.n_tries - 1] = Ex;
366 probs[params.n_tries - 1] = boltzmann(Ex, energies[i], T, ¶ms);
367
368 /* now throw biased die to see which new_x[i] we choose */
369 sum_probs[0] = probs[0];
370 for (i = 1; i < params.n_tries; ++i)
371 {
372 sum_probs[i] = sum_probs[i - 1] + probs[i];
373 }
374 u = gsl_rng_uniform (r) * sum_probs[params.n_tries - 1];
375 for (i = 0; i < params.n_tries; ++i)
376 {
377 if (u < sum_probs[i])
378 {
379 memcpy (x, (char *) new_x + i * element_size, element_size);
380 break;
381 }
382 }
383 if (print_position)
384 {
385 printf ("%5d\t%12g\t", n_iter, T);
386 print_position (x);
387 printf ("\t%12g\t%12g\n", distance (x, x0_p), Ex);
388 }
389 T *= T_factor;
390 ++n_iter;
391 if (T < params.t_min)
392 {
393 break;
394 }
395 }
396
397 /* now return the value via x0_p */
398 memcpy (x0_p, x, element_size);
399
400 /* printf("the result is: %g (E=%g)\n", x, Ex); */
401
402 free (x);
403 free (new_x);
404 free (energies);
405 free (probs);
406 free (sum_probs);
407 }
408 */
409 /// Like the function solve, but performs multiple runs and returns the best result.
410 pub fn solve_many(&self, rng: &mut crate::Rng) -> T {
411 let mut x = self.x0_p.clone();
412 let mut new_x = Vec::with_capacity(self.params.n_tries);
413
414 let mut energies = Vec::with_capacity(self.params.n_tries);
415 let mut probs = Vec::with_capacity(self.params.n_tries);
416 let mut sum_probs = Vec::with_capacity(self.params.n_tries);
417
418 let mut Temp = self.params.t_initial;
419 let Temp_factor = 1.0 / self.params.mu_t;
420
421 if self.print_t.is_some() {
422 println!(
423 "{i:^6} | {t:^12} | {p:^15} | {d:^15} | {E:^13}",
424 i = "#-iter",
425 t = "temperature",
426 p = "position",
427 d = "delta_pos",
428 E = "energy"
429 );
430 }
431
432 let mut iter = 0;
433 loop {
434 let Ex = (self.Efunc_t)(&x);
435 for i in 0..self.params.n_tries - 1 {
436 // only go to N_TRIES-2
437 sum_probs.push(0.0);
438 new_x.push(x.clone());
439
440 (self.step_t)(rng, &mut new_x[i], self.params.step_size);
441 energies.push((self.Efunc_t)(&new_x[i]));
442 probs.push(boltzmann(Ex, energies[i], Temp, &self.params));
443 }
444
445 // now add in the old value of "x", so it is a contendor
446 new_x.push(x.clone());
447 energies.push(Ex);
448 probs.push(boltzmann(
449 Ex,
450 energies[self.params.n_tries - 2],
451 Temp,
452 &self.params,
453 ));
454 sum_probs.push(0.0);
455
456 // now throw biased die to see which new_x[i] we choose
457 sum_probs[0] = probs[0];
458 for i in 1..self.params.n_tries {
459 sum_probs[i] = sum_probs[i - 1] + probs[i];
460 }
461 let u = rng.uniform() * *sum_probs.last().unwrap();
462 for i in 0..self.params.n_tries {
463 if u < sum_probs[i] {
464 x = new_x[i].clone();
465 break;
466 }
467 }
468
469 if let Some(ref printf) = self.print_t {
470 print!("{:>06} | {:>12.10} | ", iter, Temp);
471 printf(&x);
472 println!(
473 " | {: >15.12} | {: >13.12}",
474 (self.metric_t)(&x, &self.x0_p),
475 Ex
476 );
477 }
478
479 Temp *= Temp_factor;
480 iter += 1;
481 if Temp < self.params.t_min {
482 break;
483 }
484 }
485 x
486 }
487}
488
489pub struct SimAnnealingParams {
490 n_tries: usize,
491 iters_fixed_T: usize,
492 step_size: f64,
493 k: f64,
494 t_initial: f64,
495 mu_t: f64,
496 t_min: f64,
497}
498
499impl SimAnnealingParams {
500 /// These are the parameters that control a run of the simulated annealing algorithm.
501 /// This structure contains all the information needed to control the search,
502 /// beyond the energy function, the step function and the initial guess.
503 ///
504 /// - n_tries: The number of points to try for each step.
505 /// - iters: The number of iterations at each temperature.
506 /// - step_size: The maximum step size in the random walk.
507 /// - k, t_initial, mu_t, t_min: The parameters of the Boltzmann distribution and
508 /// cooling schedule.
509 pub fn new(
510 n_tries: usize,
511 iters: usize,
512 step_size: f64,
513 k: f64,
514 t_initial: f64,
515 mu_t: f64,
516 t_min: f64,
517 ) -> SimAnnealingParams {
518 SimAnnealingParams {
519 n_tries,
520 iters_fixed_T: iters,
521 step_size,
522 k,
523 t_initial,
524 mu_t,
525 t_min,
526 }
527 }
528}