Blame siman/test.c

Packit 67cb25
/* siman/test.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000 Mark Galassi
Packit 67cb25
 * 
Packit 67cb25
 * This program is free software; you can redistribute it and/or modify
Packit 67cb25
 * it under the terms of the GNU General Public License as published by
Packit 67cb25
 * the Free Software Foundation; either version 3 of the License, or (at
Packit 67cb25
 * your option) any later version.
Packit 67cb25
 * 
Packit 67cb25
 * This program is distributed in the hope that it will be useful, but
Packit 67cb25
 * WITHOUT ANY WARRANTY; without even the implied warranty of
Packit 67cb25
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Packit 67cb25
 * General Public License for more details.
Packit 67cb25
 * 
Packit 67cb25
 * You should have received a copy of the GNU General Public License
Packit 67cb25
 * along with this program; if not, write to the Free Software
Packit 67cb25
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <stdlib.h>
Packit 67cb25
#include <string.h>
Packit 67cb25
#include <math.h>
Packit 67cb25
#include <gsl/gsl_test.h>
Packit 67cb25
#include <gsl/gsl_rng.h>
Packit 67cb25
#include <gsl/gsl_siman.h>
Packit 67cb25
#include <gsl/gsl_ieee_utils.h>
Packit 67cb25
#include <stdio.h>
Packit 67cb25
Packit 67cb25
/* set up parameters for this simulated annealing run */
Packit 67cb25
#define N_TRIES 200             /* how many points do we try before stepping */
Packit 67cb25
#define ITERS_FIXED_T 1000      /* how many iterations for each T? */
Packit 67cb25
#define STEP_SIZE 1.0           /* max step size in random walk */
Packit 67cb25
#define K 1.0                   /* Boltzmann constant */
Packit 67cb25
#define T_INITIAL 0.008         /* initial temperature */
Packit 67cb25
#define MU_T 1.003              /* damping factor for temperature */
Packit 67cb25
#define T_MIN 2.0e-6
Packit 67cb25
Packit 67cb25
gsl_siman_params_t params = {N_TRIES, ITERS_FIXED_T, STEP_SIZE,
Packit 67cb25
                             K, T_INITIAL, MU_T, T_MIN};
Packit 67cb25
Packit 67cb25
double square (double x) ;
Packit 67cb25
double square (double x) { return x * x ; } 
Packit 67cb25
Packit 67cb25
double E1(void *xp);
Packit 67cb25
double M1(void *xp, void *yp);
Packit 67cb25
void S1(const gsl_rng * r, void *xp, double step_size);
Packit 67cb25
void P1(void *xp);
Packit 67cb25
Packit 67cb25
/* now some functions to test in one dimension */
Packit 67cb25
double E1(void *xp)
Packit 67cb25
{
Packit 67cb25
  double x = * ((double *) xp);
Packit 67cb25
Packit 67cb25
  return exp(-square(x-1))*sin(8*x) - exp(-square(x-1000))*0.89;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double M1(void *xp, void *yp)
Packit 67cb25
{
Packit 67cb25
  double x = *((double *) xp);
Packit 67cb25
  double y = *((double *) yp);
Packit 67cb25
Packit 67cb25
  return fabs(x - y);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
void S1(const gsl_rng * r, void *xp, double step_size)
Packit 67cb25
{
Packit 67cb25
  double old_x = *((double *) xp);
Packit 67cb25
  double new_x;
Packit 67cb25
Packit 67cb25
  new_x = gsl_rng_uniform(r)*2*step_size - step_size + old_x;
Packit 67cb25
Packit 67cb25
  memcpy(xp, &new_x, sizeof(new_x));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
void P1(void *xp)
Packit 67cb25
{
Packit 67cb25
  printf(" %12g ", *((double *) xp));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int main(void)
Packit 67cb25
{
Packit 67cb25
  double x_min = 1.36312999455315182 ;
Packit 67cb25
  double x ;
Packit 67cb25
Packit 67cb25
  gsl_rng * r = gsl_rng_alloc (gsl_rng_env_setup()) ;
Packit 67cb25
Packit 67cb25
  gsl_ieee_env_setup ();
Packit 67cb25
Packit 67cb25
  /* The function tested here has multiple mimima. 
Packit 67cb25
     The global minimum is at    x = 1.36312999, (f = -0.87287)
Packit 67cb25
     There is a local minimum at x = 0.60146196, (f = -0.84893) */
Packit 67cb25
Packit 67cb25
  x = -10.0 ;
Packit 67cb25
  gsl_siman_solve(r, &x, E1, S1, M1, NULL, NULL, NULL, NULL,
Packit 67cb25
                  sizeof(double), params);
Packit 67cb25
  gsl_test_rel(x, x_min, 1e-3, "f(x)= exp(-(x-1)^2) sin(8x), x0=-10") ;
Packit 67cb25
Packit 67cb25
  x = +10.0 ;
Packit 67cb25
  gsl_siman_solve(r, &x, E1, S1, M1, NULL, NULL, NULL, NULL,
Packit 67cb25
                  sizeof(double), params);
Packit 67cb25
  gsl_test_rel(x, x_min, 1e-3, "f(x)= exp(-(x-1)^2) sin(8x), x0=10") ;
Packit 67cb25
Packit 67cb25
  /* Start at the false minimum */
Packit 67cb25
Packit 67cb25
  x = +0.6 ; 
Packit 67cb25
  gsl_siman_solve(r, &x, E1, S1, M1, NULL, NULL, NULL, NULL,
Packit 67cb25
                  sizeof(double), params);
Packit 67cb25
  gsl_test_rel(x, x_min, 1e-3, "f(x)= exp(-(x-1)^2) sin(8x), x0=0.6") ;
Packit 67cb25
Packit 67cb25
  x = +0.5 ; 
Packit 67cb25
  gsl_siman_solve(r, &x, E1, S1, M1, NULL, NULL, NULL, NULL,
Packit 67cb25
                  sizeof(double), params);
Packit 67cb25
  gsl_test_rel(x, x_min, 1e-3, "f(x)= exp(-(x-1)^2) sin(8x), x0=0.5") ;
Packit 67cb25
Packit 67cb25
  x = +0.4 ; 
Packit 67cb25
  gsl_siman_solve(r, &x, E1, S1, M1, NULL, NULL, NULL, NULL,
Packit 67cb25
                  sizeof(double), params);
Packit 67cb25
  gsl_test_rel(x, x_min, 1e-3, "f(x)= exp(-(x-1)^2) sin(8x), x0=0.4") ;
Packit 67cb25
Packit 67cb25
  gsl_rng_free(r);
Packit 67cb25
  exit (gsl_test_summary ());
Packit 67cb25
Packit 67cb25
#ifdef JUNK 
Packit 67cb25
  x0.D1 = 12.0;
Packit 67cb25
  printf("#one dimensional problem, x0 = %f\n", x0.D1);
Packit 67cb25
  gsl_siman_Usolve(r, &x0, test_E_1D, test_step_1D, distance_1D,
Packit 67cb25
                   print_pos_1D, params);
Packit 67cb25
Packit 67cb25
Packit 67cb25
  x0.D2[0] = 12.0;
Packit 67cb25
  x0.D2[1] = 5.5;
Packit 67cb25
  printf("#two dimensional problem, (x0,y0) = (%f,%f)\n",
Packit 67cb25
         x0.D2[0], x0.D2[1]);
Packit 67cb25
  gsl_siman_Usolve(r, &x0, test_E_2D, test_step_2D, distance_2D,
Packit 67cb25
                   print_pos_2D, params); 
Packit 67cb25
Packit 67cb25
  x0.D3[0] = 12.2;
Packit 67cb25
  x0.D3[1] = 5.5;
Packit 67cb25
  x0.D3[2] = -15.5;
Packit 67cb25
  printf("#three dimensional problem, (x0,y0,z0) = (%f,%f,%f)\n",
Packit 67cb25
         x0.D3[0], x0.D3[1], x0.D3[2]);
Packit 67cb25
  gsl_siman_Usolve(r, &x0, test_E_3D, test_step_3D, distance_3D, 
Packit 67cb25
                   print_pos_3D, params); 
Packit 67cb25
Packit 67cb25
  x0.D2[0] = 12.2;
Packit 67cb25
  x0.D2[1] = 5.5;
Packit 67cb25
Packit 67cb25
  gsl_siman_solve(r, &x0, test_E_2D, test_step_2D, distance_2D, print_pos_2D, params);
Packit 67cb25
  
Packit 67cb25
  x0.D3[0] = 12.2;
Packit 67cb25
  x0.D3[1] = 5.5;
Packit 67cb25
  x0.D3[2] = -15.5;
Packit 67cb25
Packit 67cb25
  gsl_siman_solve(r, &x0, test_E_3D, test_step_3D, distance_3D, print_pos_3D, params);
Packit 67cb25
Packit 67cb25
  return 0;
Packit 67cb25
#endif
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25