00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include <cmath>
00015 #include <vector>
00016 #include <algorithm>
00017
00018 #include "Math/GSLSimAnnealing.h"
00019 #include "Math/GSLRndmEngines.h"
00020 #include "Math/SMatrix.h"
00021 #include "Math/Math.h"
00022 #include "TH1.h"
00023 #include "TGraph.h"
00024 #include "TCanvas.h"
00025 #include "TApplication.h"
00026
00027 bool showGraphics = true;
00028
00029 using namespace ROOT::Math;
00030
00031
00032
00033
00034 struct s_tsp_city {
00035 const char * name;
00036 double lat, longitude;
00037 };
00038 typedef struct s_tsp_city Stsp_city;
00039
00040 Stsp_city cities[] = {{"Santa Fe", 35.68, 105.95},
00041 {"Phoenix", 33.54, 112.07},
00042 {"Albuquerque", 35.12, 106.62},
00043 {"Clovis", 34.41, 103.20},
00044 {"Durango", 37.29, 107.87},
00045 {"Dallas", 32.79, 96.77},
00046 {"Tesuque", 35.77, 105.92},
00047 {"Grants", 35.15, 107.84},
00048 {"Los Alamos", 35.89, 106.28},
00049 {"Las Cruces", 32.34, 106.76},
00050 {"Cortez", 37.35, 108.58},
00051 {"Gallup", 35.52, 108.74}};
00052
00053 #define N_CITIES (sizeof(cities)/sizeof(Stsp_city))
00054
00055 double distance_matrix[N_CITIES][N_CITIES];
00056
00057
00058 double city_distance(Stsp_city c1, Stsp_city c2)
00059 {
00060 const double earth_radius = 6375.000;
00061
00062 double sla1 = std::sin(c1.lat*M_PI/180), cla1 = std::cos(c1.lat*M_PI/180),
00063 slo1 = std::sin(c1.longitude*M_PI/180), clo1 = std::cos(c1.longitude*M_PI/180);
00064 double sla2 = std::sin(c2.lat*M_PI/180), cla2 = std::cos(c2.lat*M_PI/180),
00065 slo2 = std::sin(c2.longitude*M_PI/180), clo2 = std::cos(c2.longitude*M_PI/180);
00066
00067 double x1 = cla1*clo1;
00068 double x2 = cla2*clo2;
00069
00070 double y1 = cla1*slo1;
00071 double y2 = cla2*slo2;
00072
00073 double z1 = sla1;
00074 double z2 = sla2;
00075
00076 double dot_product = x1*x2 + y1*y2 + z1*z2;
00077
00078 double angle = std::acos(dot_product);
00079
00080
00081 return angle*earth_radius;
00082 }
00083
00084
00085 void print_distance_matrix()
00086 {
00087 unsigned int i, j;
00088
00089 for (i = 0; i < N_CITIES; ++i) {
00090 printf("# ");
00091 for (j = 0; j < N_CITIES; ++j) {
00092 printf("%15.8f ", distance_matrix[i][j]);
00093 }
00094 printf("\n");
00095 }
00096 }
00097
00098
00099
00100
00101
00102 class MySimAnFunc : public GSLSimAnFunc {
00103
00104 public:
00105
00106 MySimAnFunc( std::vector<double> & allDist)
00107 {
00108 calculate_distance_matrix();
00109
00110 for (unsigned int i = 0; i < N_CITIES; ++i)
00111 fRoute[i] = i;
00112
00113
00114 fDist = &allDist;
00115 }
00116
00117
00118 virtual ~MySimAnFunc() {}
00119
00120 unsigned int Route(unsigned int i) const { return fRoute[i]; }
00121
00122 const unsigned int * Route() const { return fRoute; }
00123 unsigned int * Route() { return fRoute; }
00124
00125 virtual MySimAnFunc * Clone() const { return new MySimAnFunc(*this); }
00126
00127 std::vector<double> & AllDist() { return *fDist; }
00128
00129 virtual double Energy() const {
00130
00131
00132
00133 double enrg = 0;
00134 for (unsigned int i = 0; i < N_CITIES; ++i) {
00135
00136
00137 enrg += fDistanceMatrix( fRoute[i] , fRoute[ (i + 1) % N_CITIES ] );
00138 }
00139
00140
00141 return enrg;
00142 }
00143
00144 virtual double Distance(const GSLSimAnFunc & f) const {
00145 const MySimAnFunc * f2 = dynamic_cast<const MySimAnFunc *> (&f);
00146 assert (f2 != 0);
00147 double d = 0;
00148
00149 for (unsigned int i = 0; i < N_CITIES; ++i) {
00150 d += (( fRoute[i] == f2->Route(i) ) ? 0 : 1 );
00151 }
00152 return d;
00153 }
00154 virtual void Step(const GSLRandomEngine & r, double ) {
00155
00156 int x1, x2, dummy;
00157
00158
00159
00160 x1 = r.RndmInt(N_CITIES);
00161 do {
00162 x2 = r.RndmInt(N_CITIES);
00163 } while (x2 == x1);
00164
00165
00166 dummy = fRoute[x1];
00167 fRoute[x1] = fRoute[x2];
00168 fRoute[x2] = dummy;
00169
00170
00171
00172 }
00173
00174 virtual void Print() {
00175 printf(" [");
00176 for (unsigned i = 0; i < N_CITIES; ++i) {
00177 printf(" %d ", fRoute[i]);
00178 }
00179 printf("] ");
00180 fDist->push_back(Energy());
00181 }
00182
00183
00184 virtual MySimAnFunc & FastCopy(const GSLSimAnFunc & f) {
00185 const MySimAnFunc * rhs = dynamic_cast<const MySimAnFunc *>(&f);
00186 assert (rhs != 0);
00187 std::copy(rhs->fRoute, rhs->fRoute + N_CITIES, fRoute);
00188 return *this;
00189 }
00190
00191 double Distance(int i, int j) const { return fDistanceMatrix(i,j); }
00192
00193 void PrintRoute();
00194
00195 void SetRoute(unsigned int * r) { std::copy(r,r+N_CITIES,fRoute); }
00196
00197 private:
00198
00199 void calculate_distance_matrix();
00200
00201
00202 unsigned int fRoute[N_CITIES];
00203 typedef SMatrix<double,N_CITIES,N_CITIES, MatRepSym<double,N_CITIES> > Matrix;
00204 Matrix fDistanceMatrix;
00205
00206 std::vector<double> * fDist;
00207
00208 };
00209
00210
00211
00212 void MySimAnFunc::calculate_distance_matrix()
00213 {
00214 unsigned int i, j;
00215 double dist;
00216
00217 for (i = 0; i < N_CITIES; ++i) {
00218 for (j = 0; j <= i; ++j) {
00219 if (i == j) {
00220 dist = 0;
00221 } else {
00222 dist = city_distance(cities[i], cities[j]);
00223 }
00224 fDistanceMatrix(i,j) = dist;
00225 }
00226 }
00227 }
00228
00229 void MySimAnFunc::PrintRoute() {
00230
00231 double dtot = 0;
00232 for (unsigned int i = 0; i < N_CITIES; ++i) {
00233 std::cout << std::setw(20) << cities[Route(i)].name << " \t " << Route(i);
00234 int j = Route(i);
00235 int k = Route( (i+ 1) % N_CITIES );
00236 dtot += Distance(j,k);
00237 std::cout << "\tdistance [" << j << " , " << k << " ]\t= " << Distance(j,k) << "\tTotal Distance\t = " << dtot << std::endl;
00238 }
00239 std::cout << "Total Route energy is " << dtot << std::endl;
00240 }
00241
00242
00243
00244 void simanTSP(bool debug = true) {
00245
00246 std::vector<double> allDist;
00247 allDist.reserve(5000);
00248
00249
00250 MySimAnFunc f(allDist);
00251
00252 GSLSimAnnealing siman;
00253
00254 GSLSimAnParams & p = siman.Params();
00255 p.n_tries = 200;
00256 p.iters_fixed_T = 10;
00257 p.step_size = 1;
00258 p.k = 1;
00259 p.t_initial = 5000;
00260 p.mu = 1.01;
00261 p.t_min = 0.5;
00262
00263
00264
00265
00266 siman.Solve(f, debug);
00267
00268 unsigned int niter = allDist.size();
00269 std::cout << "minimum found is for distance " << f.Energy() << std::endl;
00270 if (debug) std::cout << "number of iterations is " << niter << std::endl;
00271
00272 std::cout << "Best Route is \n";
00273 f.PrintRoute();
00274
00275
00276 double x0[N_CITIES+1];
00277 double y0[N_CITIES+1];
00278 double xmin[N_CITIES+1];
00279 double ymin[N_CITIES+1];
00280
00281
00282
00283 TH1 * h1 = new TH1D("h1","Distance",niter+1,0,niter+1);
00284 for (unsigned int i = 1; i <=niter; ++i) {
00285 h1->SetBinContent(i,allDist[i]);
00286 }
00287
00288 for (unsigned int i = 0; i < N_CITIES; ++i) {
00289
00290 x0[i] = -cities[i].longitude;
00291 y0[i] = cities[i].lat;
00292 xmin[i] = -cities[f.Route(i)].longitude;
00293 ymin[i] = cities[f.Route(i)].lat;
00294 }
00295
00296 x0[N_CITIES] = x0[0];
00297 y0[N_CITIES] = y0[0];
00298 xmin[N_CITIES] = xmin[0];
00299 ymin[N_CITIES] = ymin[0];
00300
00301 if ( showGraphics )
00302 {
00303
00304 TGraph *g0 = new TGraph(N_CITIES+1,x0,y0);
00305 TGraph *gmin = new TGraph(N_CITIES+1,xmin,ymin);
00306
00307 TCanvas * c1 = new TCanvas("c","TSP",10,10,1000,800);
00308 c1->Divide(2,2);
00309 c1->cd(1);
00310 g0->Draw("alp");
00311 g0->SetMarkerStyle(20);
00312 c1->cd(2);
00313 gmin->SetMarkerStyle(20);
00314 gmin->Draw("alp");
00315 c1->cd(3);
00316 h1->Draw();
00317 }
00318
00319 }
00320
00321 unsigned int r1[N_CITIES];
00322 unsigned int r2[N_CITIES];
00323 unsigned int r3[N_CITIES];
00324 unsigned int nfiter = 0;
00325 double best_E, second_E, third_E;
00326
00327
00328
00329 void do_all_perms(MySimAnFunc & f, int offset )
00330 {
00331 unsigned int * r = f.Route();
00332
00333
00334 while (std::next_permutation(r+offset,r+N_CITIES) ) {
00335
00336 double E = f.Energy();
00337 nfiter++;
00338
00339 if (E < best_E) {
00340 third_E = second_E;
00341 std::copy(r2,r2+N_CITIES,r3);
00342 second_E = best_E;
00343 std::copy(r1,r1+N_CITIES,r2);
00344 best_E = E;
00345 std::copy(r,r+N_CITIES,r1);
00346 } else if (E < second_E) {
00347 third_E = second_E;
00348 std::copy(r2,r2+N_CITIES,r3);
00349 second_E = E;
00350 std::copy(r,r+N_CITIES,r2);
00351 } else if (E < third_E) {
00352 third_E = E;
00353 std::copy(r,r+N_CITIES,r3);
00354 }
00355 }
00356 }
00357
00358 #ifdef O
00359
00360
00361 void do_all_perms(MySimAnFunc & f, int n) {
00362
00363 if (n == (N_CITIES-1)) {
00364
00365 const unsigned int * r = f.Route();
00366
00367 if (E < best_E) {
00368 third_E = second_E;
00369 std::copy(r2,r2+N_CITIES,r3);
00370 second_E = best_E;
00371 std::copy(r1,r1+N_CITIES,r2);
00372 best_E = E;
00373 std::copy(r,r+N_CITIES,r1);
00374 } else if (E < second_E) {
00375 third_E = second_E;
00376 std::copy(r2,r2+N_CITIES,r3);
00377 second_E = E;
00378 std::copy(r,r+N_CITIES,r2);
00379 } else if (E < third_E) {
00380 third_E = E;
00381 std::copy(r,r+N_CITIES,r3);
00382 }
00383 } else {
00384 unsigned int newr[N_CITIES];
00385 unsigned int j;
00386 int swap_tmp;
00387 const unsigned int * r = f.Route();
00388 std::copy(r,r+N_CITIES,newr);
00389 for (j = n; j < N_CITIES; ++j) {
00390
00391 swap_tmp = newr[j];
00392 newr[j] = newr[n];
00393 newr[n] = swap_tmp;
00394 f.SetRoute(newr);
00395 do_all_perms(f, n+1);
00396 }
00397 }
00398 }
00399 #endif
00400
00401
00402
00403
00404
00405
00406 void FullSearch()
00407 {
00408
00409 std::vector<double> dummy;
00410
00411 MySimAnFunc f(dummy);
00412
00413
00414
00415 const unsigned int * r = f.Route();
00416 std::copy(r,r+N_CITIES,r1);
00417 std::copy(r,r+N_CITIES,r2);
00418 std::copy(r,r+N_CITIES,r3);
00419 best_E = f.Energy();
00420 second_E = 1.E100;
00421 third_E = 1.E100;
00422
00423
00424 do_all_perms(f, 1 );
00425
00426 std::cout << "number of calls " << nfiter << std::endl;
00427
00428
00429 printf("\n# exhaustive best route: ");
00430 std::cout << best_E << std::endl;
00431 f.SetRoute(r1); f.PrintRoute();
00432
00433
00434 printf("\n# exhaustive second_best route: ");
00435 std::cout << second_E << std::endl;
00436 f.SetRoute(r2); f.PrintRoute();
00437
00438 printf("\n# exhaustive third_best route: ");
00439 std::cout << third_E << std::endl;
00440 f.SetRoute(r3); f.PrintRoute();
00441
00442 }
00443
00444
00445
00446
00447 #ifndef __CINT__
00448 int main(int argc, char **argv)
00449 {
00450 using std::cout;
00451 using std::endl;
00452 using std::cerr;
00453
00454 if ( argc > 1 && argc != 2 )
00455 {
00456 cerr << "Usage: " << argv[0] << " [-ng]\n";
00457 cerr << " where:\n";
00458 cerr << " -ng : no graphics mode";
00459 cerr << endl;
00460 exit(1);
00461 }
00462
00463 if ( argc == 2 && strcmp( argv[1], "-ng") == 0 )
00464 {
00465 showGraphics = false;
00466 }
00467
00468 if ( showGraphics )
00469 {
00470 TApplication* theApp = 0;
00471 theApp = new TApplication("App",&argc,argv);
00472 simanTSP(true);
00473 theApp->Run();
00474 delete theApp;
00475 theApp = 0;
00476 }
00477 else
00478 {
00479 simanTSP(false);
00480 }
00481
00482
00483
00484 return 0;
00485 }
00486 #endif
00487