stressFit.cxx

Go to the documentation of this file.
00001 // @(#)root/test:$name:  $:$id: stressFit.cxx,v 1.15 2002/10/25 10:47:51 rdm exp $
00002 // Authors: Rene Brun, Eddy Offermann  April 2006
00003    
00004 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*//
00005 //                                                                           //
00006 // Function Minimization Examples, Fred James                                //
00007 //                                                                           //
00008 // from the                                                                  //
00009 //   Proceedings of the 1972 CERN Computing and Data Processing School       // 
00010 //   Pertisau, Austria, 10-24 September, 1972 (CERN 72-21)                   //
00011 //                                                                           //
00012 // Here a collection of test problems is assembled which were found to be    //
00013 // useful in verifying and comparing minimization routines. Many of these    //
00014 // are standard functions upon which it has become conventional to try all   //
00015 // new methods, quoting the performance in the publication of the algorithm  //
00016 //                                                                           //
00017 // Each test will produce one line (Test OK or Test FAILED) . At the end of  //
00018 // the test a table is printed showing the global results Real Time and      //
00019 // Cpu Time. One single number (ROOTMARKS) is also calculated showing the    //
00020 // relative performance of your machine compared to a reference machine      //
00021 // a Pentium IV 2.4 Ghz) with 512 MBytes of memory and 120 GBytes IDE disk.  //
00022 //                                                                           //
00023 // In the main routine the fitter can be chosen through TVirtualFitter :     //
00024 //   - Minuit                                                                //
00025 //   - Minuit2                                                               //
00026 //   - Fumili                                                                //
00027 //
00028 //  To run the test, do, eg
00029 // root -b -q stressFit.cxx
00030 // root -b -q "stressFit.cxx(\"Minuit2\")"
00031 // root -b -q "stressFit.cxx+(\"Minuit2\")"
00032 //                                                                           //
00033 // The verbosity can be set through the global parameter gVerbose :          //
00034 //   -1: off  1: on                                                          //
00035 // The tolerance on the parameter deviation from the minimum can be set      //
00036 // through gAbsTolerance .                                                   //
00037 //                                                                           //
00038 // An example of output when all the tests run OK is shown below:            //
00039 // *******************************************************************       //
00040 // *  Minimization - S T R E S S suite                               *       //
00041 // *******************************************************************       //
00042 // *******************************************************************       //
00043 // *  Starting  S T R E S S                                          *       //
00044 // *******************************************************************       //
00045 // Test  1 : Wood.................................................. OK       //
00046 // Test  2 : RosenBrock............................................ OK       //
00047 // Test  3 : Powell................................................ OK       //
00048 // Test  4 : Fletcher.............................................. OK       //
00049 // Test  5 : GoldStein1............................................ OK       //
00050 // Test  6 : GoldStein2............................................ OK       //
00051 // Test  7 : TrigoFletcher......................................... OK       //
00052 // *******************************************************************       //
00053 //                                                                           //
00054 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*_*//
00055 
00056 #include <stdlib.h>
00057 #include "TVirtualFitter.h"
00058 #include "TSystem.h"
00059 #include "TROOT.h"
00060 #include "TBenchmark.h"
00061 #include "TMath.h"
00062 #include "TStopwatch.h"
00063 #include "Riostream.h"
00064 #include "TVectorD.h"
00065 #include "TMatrixD.h"
00066 
00067 Int_t stressFit(const char *theFitter="Minuit", Int_t N=2000);
00068 Int_t    gVerbose      = -1;
00069 Double_t gAbsTolerance = 0.005;
00070 
00071 //------------------------------------------------------------------------
00072 void StatusPrint(Int_t id,const TString &title,Bool_t status)
00073 {
00074   // Print test program number and its title
00075   const Int_t kMAX = 65;
00076   Char_t number[4];
00077   sprintf(number,"%2d",id);
00078   TString header = TString("Test ")+number+" : "+title;
00079   const Int_t nch = header.Length();
00080   for (Int_t i = nch; i < kMAX; i++) header += '.';
00081   cout << header << (status ? "OK" : "FAILED") << endl;
00082 }
00083 
00084 //______________________________________________________________________________
00085 void RosenBrock(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00086 {
00087   const Double_t x = par[0];
00088   const Double_t y = par[1];
00089   const Double_t tmp1 = y-x*x;
00090   const Double_t tmp2 = 1-x;
00091   f = 100*tmp1*tmp1+tmp2*tmp2;
00092 }
00093 
00094 //______________________________________________________________________________
00095 Bool_t RunRosenBrock()
00096 {
00097 //
00098 // F(x,y) = 100 (y-x^2)^2 + (1-x)^2
00099 //
00100 //   start point: F(-1.2,1.0) = 24.20
00101 //   minimum    : F(1.0,1.0)  = 0.
00102 //
00103 // This narrow, parabolic valley is probably the best known of all test cases. The floor
00104 // of the valley follows approximately the parabola y = x^2+1/200 .
00105 // There is a region where the covariance matrix is not positive-definite and even a path
00106 // where it is singular . Stepping methods tend to perform at least as well as gradient
00107 //  method for this function .
00108 // [Reference: Comput. J. 3,175 (1960).]
00109 
00110   Bool_t ok = kTRUE;
00111   TVirtualFitter *min = TVirtualFitter::Fitter(0,2);
00112   min->SetFCN(RosenBrock);
00113 
00114   Double_t arglist[100];
00115   arglist[0] = gVerbose;
00116   min->ExecuteCommand("SET PRINT",arglist,1);
00117 
00118   min->SetParameter(0,"x",-1.2,0.01,0,0);
00119   min->SetParameter(1,"y", 1.0,0.01,0,0);
00120 
00121   arglist[0] = 0;
00122   min->ExecuteCommand("SET NOW",arglist,0);
00123   arglist[0] = 1000; // number of function calls 
00124   arglist[1] = 0.001; // tolerance 
00125   min->ExecuteCommand("MIGRAD",arglist,0);
00126   min->ExecuteCommand("MIGRAD",arglist,2);
00127   min->ExecuteCommand("MINOS",arglist,0);
00128 
00129   Double_t parx,pary;
00130   Double_t we,al,bl;
00131   Char_t parName[32];
00132   min->GetParameter(0,parName,parx,we,al,bl);
00133   min->GetParameter(1,parName,pary,we,al,bl);
00134   
00135   ok = ( TMath::Abs(parx-1.) < gAbsTolerance &&
00136          TMath::Abs(pary-1.) < gAbsTolerance );
00137 
00138   delete min;
00139 
00140   return ok;
00141 }
00142 
00143 //______________________________________________________________________________
00144 void Wood4(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00145 {
00146   const Double_t w = par[0]; 
00147   const Double_t x = par[1]; 
00148   const Double_t y = par[2]; 
00149   const Double_t z = par[3]; 
00150 
00151   const Double_t w1 = w-1;
00152   const Double_t x1 = x-1;
00153   const Double_t y1 = y-1;
00154   const Double_t z1 = z-1;
00155   const Double_t tmp1 = x-w*w;
00156   const Double_t tmp2 = z-y*y;
00157 
00158   f = 100*tmp1*tmp1+w1*w1+90*tmp2*tmp2+y1*y1+10.1*(x1*x1+z1*z1)+19.8*x1*z1;
00159 }
00160 
00161 //______________________________________________________________________________
00162 Bool_t RunWood4()
00163 {
00164 //
00165 // F(w,x,y,z) = 100 (y-w^2)^2 + (w-1)^2 + 90 (z-y^2)^2
00166 //              + (1-y)^2 + 10.1 [(x-1)^2 + (z-1)^2]
00167 //              + 19.8 (x-1)(z-1)
00168 //
00169 //   start point: F(-3,-1,-3,-1) = 19192
00170 //   minimum    : F(1,1,1,1)  =   0.
00171 //
00172 // This is a fourth-degree polynomial which is reasonably well-behaved near the minimum,
00173 // but in order to get there one must cross a rather flat, four-dimensional "plateau"
00174 // which often causes minimization algorithm to get "stuck" far from the minimum. As
00175 // such it is a particularly good test of convergence criteria and simulates quite well a
00176 // feature of many physical problems in many variables where no good starting
00177 // approximation is known .
00178 // [Reference: Unpublished. See IBM Technical Report No. 320-2949.]
00179 
00180   Bool_t ok = kTRUE;
00181   TVirtualFitter *min = TVirtualFitter::Fitter(0,4);
00182   min->SetFCN(Wood4);
00183 
00184   Double_t arglist[100];
00185   arglist[0] = gVerbose;
00186   min->ExecuteCommand("SET PRINT",arglist,1);
00187 
00188   min->SetParameter(0,"w",-3.0,0.01,0,0);
00189   min->SetParameter(1,"x",-1.0,0.01,0,0);
00190   min->SetParameter(2,"y",-3.0,0.01,0,0);
00191   min->SetParameter(3,"z",-1.0,0.01,0,0);
00192 
00193   arglist[0] = 0;
00194   min->ExecuteCommand("SET NOW",arglist,0);
00195   arglist[0] = 1000; // number of function calls 
00196   arglist[1] = 0.001; // tolerance 
00197   min->ExecuteCommand("MIGRAD",arglist,0);
00198   min->ExecuteCommand("MIGRAD",arglist,2);
00199   min->ExecuteCommand("MINOS",arglist,0);
00200 
00201   Double_t parw,parx,pary,parz;
00202   Double_t we,al,bl;
00203   Char_t parName[32];
00204   min->GetParameter(1,parName,parw,we,al,bl);
00205   min->GetParameter(0,parName,parx,we,al,bl);
00206   min->GetParameter(1,parName,pary,we,al,bl);
00207   min->GetParameter(1,parName,parz,we,al,bl);
00208   
00209   ok = ( TMath::Abs(parw-1.) < gAbsTolerance &&
00210          TMath::Abs(parx-1.) < gAbsTolerance &&
00211          TMath::Abs(pary-1.) < gAbsTolerance &&
00212          TMath::Abs(parz-1.) < gAbsTolerance );
00213 
00214   delete min;
00215 
00216   return ok;
00217 }
00218 
00219 //______________________________________________________________________________
00220 void Powell(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00221 {
00222   const Double_t w = par[0]; 
00223   const Double_t x = par[1]; 
00224   const Double_t y = par[2]; 
00225   const Double_t z = par[3]; 
00226 
00227   const Double_t tmp1 = w+10*x;
00228   const Double_t tmp2 = y-z;
00229   const Double_t tmp3 = x-2*y;
00230   const Double_t tmp4 = w-z;
00231 
00232   f = tmp1*tmp1+5*tmp2*tmp2+tmp3*tmp3*tmp3*tmp3+10*tmp4*tmp4*tmp4*tmp4;
00233 }
00234 
00235 //______________________________________________________________________________
00236 Bool_t RunPowell()
00237 {
00238 //
00239 // F(w,x,y,z) = (w+10x)^2 + 5(y-z)^2 + (x-2y)^4+ 10 (w-z)^4
00240 //
00241 //   start point: F(-3,-1,0,1) = 215
00242 //   minimum    : F(0,0,0,0)  =   0.
00243 //
00244 // This function is difficult because its matrix of second derivatives becomes singular
00245 //  at the minimum. Near the minimum the function is given by (w + 10x)^2 + 5 (y-5)^2
00246 // which does not determine the minimum uniquely.
00247 // [Reference: Comput. J. 5, 147 (1962).]
00248 
00249   Bool_t ok = kTRUE;
00250   TVirtualFitter *min = TVirtualFitter::Fitter(0,4);
00251   min->SetFCN(Powell);
00252   
00253   Double_t arglist[100];
00254   arglist[0] = gVerbose;
00255   min->ExecuteCommand("SET PRINT",arglist,1);
00256   
00257   min->SetParameter(0,"w",+3.0,0.01,0,0);
00258   min->SetParameter(1,"x",-1.0,0.01,0,0);
00259   min->SetParameter(2,"y", 0.0,0.01,0,0);
00260   min->SetParameter(3,"z",+1.0,0.01,0,0);
00261     
00262   arglist[0] = 0;
00263   min->ExecuteCommand("SET NOW",arglist,0);
00264   arglist[0] = 1000; // number of function calls 
00265   arglist[1] = 0.001; // tolerance 
00266   min->ExecuteCommand("MIGRAD",arglist,0);
00267   min->ExecuteCommand("MIGRAD",arglist,2);
00268   min->ExecuteCommand("MINOS",arglist,0);
00269     
00270   Double_t parw,parx,pary,parz;
00271   Double_t we,al,bl;
00272   Char_t parName[32];
00273   min->GetParameter(1,parName,parw,we,al,bl);
00274   min->GetParameter(0,parName,parx,we,al,bl);
00275   min->GetParameter(1,parName,pary,we,al,bl);
00276   min->GetParameter(1,parName,parz,we,al,bl);
00277   
00278   ok = ( TMath::Abs(parw-0.) < gAbsTolerance &&
00279          TMath::Abs(parx-0.) < 10.*gAbsTolerance &&
00280          TMath::Abs(pary-0.) < gAbsTolerance &&
00281          TMath::Abs(parz-0.) < gAbsTolerance );
00282 
00283   delete min;
00284 
00285   return ok;
00286 }
00287 
00288 //______________________________________________________________________________
00289 void Fletcher(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00290 {
00291   const Double_t x = par[0];
00292   const Double_t y = par[1];
00293   const Double_t z = par[2];
00294 
00295   Double_t psi;
00296   if (x > 0)
00297     psi = TMath::ATan(y/x)/2/TMath::Pi();
00298   else if (x < 0)
00299     psi = 0.5+TMath::ATan(y/x)/2/TMath::Pi();
00300   else
00301     psi = 0.0;
00302 
00303   const Double_t tmp1 = z-10*psi;
00304   const Double_t tmp2 = TMath::Sqrt(x*x+y*y)-1;
00305 
00306   f = 100*(tmp1*tmp1+tmp2*tmp2)+z*z;
00307 }
00308 
00309 //______________________________________________________________________________
00310 Bool_t RunFletcher()
00311 {
00312 //
00313 // F(x,y,z) = 100 {[z - 10 G(x,y)]^2 + ( (x^2+y^2)^1/2 - 1 )^2} + z^2
00314 //
00315 //                     | arctan(y/x)        for x > 0
00316 // where 2 pi G(x,y) = |
00317 //                     | pi + arctan(y/x)   for x < 0
00318 //
00319 //   start point: F(-1,0,0) = 2500
00320 //   minimum    : F(1,0,0)  =   0.
00321 //
00322 // F is defined only for -0.25 < G(x,y) < 0.75
00323 //
00324 // This is a curved valley problem, similar to Rosenbrock's, but in three dimensions .
00325 // [Reference: Comput. J. 6, 163 (1963).]
00326 
00327   Bool_t ok = kTRUE;
00328   TVirtualFitter *min = TVirtualFitter::Fitter(0,3);
00329   min->SetFCN(Fletcher);
00330 
00331   Double_t arglist[100];
00332   arglist[0] = gVerbose;
00333   min->ExecuteCommand("SET PRINT",arglist,1);
00334 
00335   min->SetParameter(0,"x",-1.0,0.01,0,0);
00336   min->SetParameter(1,"y", 0.0,0.01,0,0);
00337   min->SetParameter(2,"z", 0.0,0.01,0,0);
00338 
00339   arglist[0] = 0;
00340   min->ExecuteCommand("SET NOW",arglist,0);
00341   arglist[0] = 1000; // number of function calls 
00342   arglist[1] = 0.001; // tolerance 
00343   min->ExecuteCommand("MIGRAD",arglist,0);
00344   min->ExecuteCommand("MIGRAD",arglist,2);
00345   min->ExecuteCommand("MINOS",arglist,0);
00346 
00347   Double_t parx,pary,parz;
00348   Double_t we,al,bl;
00349   Char_t parName[32];
00350   min->GetParameter(0,parName,parx,we,al,bl);
00351   min->GetParameter(1,parName,pary,we,al,bl);
00352   min->GetParameter(1,parName,parz,we,al,bl);
00353   
00354   ok = ( TMath::Abs(parx-1.) < gAbsTolerance &&
00355          TMath::Abs(pary-0.) < gAbsTolerance &&
00356          TMath::Abs(parz-0.) < gAbsTolerance );
00357 
00358   delete min;
00359 
00360   return ok;
00361 }
00362 
00363 //______________________________________________________________________________
00364 void GoldStein1(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00365 {
00366   const Double_t x = par[0];
00367   const Double_t y = par[1];
00368 
00369   const Double_t tmp1 = x+y+1;
00370   const Double_t tmp2 = 19-14*x+3*x*x-14*y+6*x*y+3*y*y;
00371   const Double_t tmp3 = 2*x-3*y;
00372   const Double_t tmp4 = 18-32*x+12*x*x+48*y-36*x*y+27*y*y;
00373 
00374   f = (1+tmp1*tmp1*tmp2)*(30+tmp3*tmp3*tmp4);
00375 }
00376 
00377 //______________________________________________________________________________
00378 Bool_t RunGoldStein1()
00379 {
00380 //
00381 // F(x,y) = (1 + (x+y+1)^2 * (19-14x+3x^2-14y+6xy+3y^2))
00382 //           * (30 + (2x-3y)^2 * (18-32x+12x^2+48y-36xy+27y^2))
00383 // 
00384 //   start point     : F(-0.4,-0,6) = 35
00385 //   local  minima   : F(1.2,0.8)   = 840
00386 //                     F(1.8,0.2)   = 84
00387 //                     F(-0.6,-0.4) = 30
00388 //   global minimum  : F(0.0,-1.0)  = 3
00389 //   
00390 // This is an eighth-order polynomial in two variables which is well behaved near each
00391 // minimum, but has four local minima and is of course non-positive-definite in many
00392 // regions. The saddle point between the two lowest minima occurs at F(-0.4,-0.6)=35
00393 // making this an interesting start point .
00394 // [Reference: Math. Comp. 25, 571 (1971).]
00395 
00396   Bool_t ok = kTRUE;
00397   TVirtualFitter *min = TVirtualFitter::Fitter(0,2);
00398   min->SetFCN(GoldStein1);
00399 
00400   Double_t arglist[100];
00401   arglist[0] = gVerbose;
00402   min->ExecuteCommand("SET PRINT",arglist,1);
00403   
00404   min->SetParameter(0,"x",-0.3999,0.01,-2.0,+2.0);
00405   min->SetParameter(1,"y",-0.6,0.01,-2.0,+2.0);
00406   
00407   arglist[0] = 0;
00408   min->ExecuteCommand("SET NOW",arglist,0);
00409   arglist[0] = 1000; // number of function calls 
00410   arglist[1] = 0.001; // tolerance 
00411   min->ExecuteCommand("MIGRAD",arglist,0);
00412   min->ExecuteCommand("MIGRAD",arglist,2);
00413   min->ExecuteCommand("MIGRAD",arglist,2);
00414   min->ExecuteCommand("MINOS",arglist,0);
00415 
00416   Double_t parx,pary;
00417   Double_t we,al,bl;
00418   Char_t parName[32];
00419   min->GetParameter(0,parName,parx,we,al,bl);
00420   min->GetParameter(1,parName,pary,we,al,bl);
00421 
00422   ok = ( TMath::Abs(parx-0.) < gAbsTolerance &&
00423          TMath::Abs(pary+1.) < gAbsTolerance );
00424   
00425   delete min;
00426 
00427   return ok;
00428 }
00429 
00430 //______________________________________________________________________________
00431 void GoldStein2(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00432 {
00433   const Double_t x = par[0];
00434   const Double_t y = par[1];
00435 
00436   const Double_t tmp1 = x*x+y*y-25;
00437   const Double_t tmp2 = TMath::Sin(4*x-3*y);
00438   const Double_t tmp3 = 2*x+y-10;
00439 
00440   f = TMath::Exp(0.5*tmp1*tmp1)+tmp2*tmp2*tmp2*tmp2+0.5*tmp3*tmp3;
00441 }
00442 
00443 //______________________________________________________________________________
00444 Bool_t RunGoldStein2()
00445 {
00446 //
00447 // F(x,y) = (1 + (x+y+1)^2 * (19-14x+3x^2-14y+6xy+3y^2))
00448 //           * (30 + (2x-3y)^2 * (18-32x+12x^2+48y-36xy+27y^2))
00449 // 
00450 //   start point     : F(1.6,3.4) =
00451 //   global minimum  : F(3,4)     = 1
00452 //   
00453 // This function has many local minima .
00454 // [Reference: Math. Comp. 25, 571 (1971).]
00455 
00456   Bool_t ok = kTRUE;
00457   TVirtualFitter *min = TVirtualFitter::Fitter(0,2);
00458   min->SetFCN(GoldStein2);
00459 
00460   Double_t arglist[100];
00461   arglist[0] = gVerbose;
00462   min->ExecuteCommand("SET PRINT",arglist,1);
00463 
00464   min->SetParameter(0,"x",+1.0,0.01,-5.0,+5.0);
00465   min->SetParameter(1,"y",+3.2,0.01,-5.0,+5.0);
00466 
00467   arglist[0] = 0;
00468   min->ExecuteCommand("SET NOW",arglist,0);
00469   arglist[0] = 1000; // number of function calls 
00470   arglist[1] = 0.01; // tolerance 
00471   min->ExecuteCommand("MIGRAD",arglist,0);
00472   min->ExecuteCommand("MIGRAD",arglist,2);
00473   min->ExecuteCommand("MINOS",arglist,0);
00474 
00475   Double_t parx,pary;
00476   Double_t we,al,bl;
00477   Char_t parName[32];
00478   min->GetParameter(0,parName,parx,we,al,bl);
00479   min->GetParameter(1,parName,pary,we,al,bl);
00480 
00481   ok = ( TMath::Abs(parx-3.) < gAbsTolerance &&
00482          TMath::Abs(pary-4.) < gAbsTolerance );
00483 
00484   delete min;
00485 
00486   return ok;
00487 }
00488 
00489 Double_t seed = 3;
00490 Int_t  nf;
00491 TMatrixD A;
00492 TMatrixD B;
00493 TVectorD x0;
00494 TVectorD sx0;
00495 TVectorD cx0;
00496 TVectorD sx;
00497 TVectorD cx;
00498 TVectorD v0;
00499 TVectorD v;
00500 TVectorD r;
00501 
00502 //______________________________________________________________________________
00503 void TrigoFletcher(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t /*iflag*/)
00504 {
00505   Int_t i;
00506   for (i = 0; i < nf ; i++) {
00507     cx0[i] = TMath::Cos(x0[i]);
00508     sx0[i] = TMath::Sin(x0[i]);
00509     cx [i] = TMath::Cos(par[i]);
00510     sx [i] = TMath::Sin(par[i]);
00511   }
00512 
00513   v0 = A*sx0+B*cx0;
00514   v  = A*sx +B*cx;
00515   r  = v0-v;
00516  
00517   f = r * r;
00518 }
00519 
00520 //______________________________________________________________________________
00521 Bool_t RunTrigoFletcher()
00522 {
00523 //
00524 // F(\vec{x}) = \sum_{i=1}^n ( E_i - \sum_{j=1}^n (A_{ij} \sin x_j + B_{ij} \cos x_j) )^2
00525 // 
00526 //   where E_i = \sum_{j=1}^n ( A_{ij} \sin x_{0j} + B_{ij} \cos x_{0j} )
00527 //
00528 //   B_{ij} and A_{ij} are random matrices composed of integers between -100 and 100;
00529 //   for j = 1,...,n: x_{0j} are any random numbers, -\pi < x_{0j} < \pi;
00530 //
00531 //   start point : x_j = x_{0j} + 0.1 \delta_j,  -\pi < \delta_j < \pi
00532 //   minimum     : F(\vec{x} = \vec{x}_0) = 0
00533 //   
00534 // This is a set of functions of any number of variables n, where the minimum is always
00535 // known in advance, but where the problem can be changed by choosing different
00536 // (random) values of the constants A_{ij}, B_{ij}, and x_{0j} . The difficulty can be
00537 // varied by choosing larger starting deviations \delta_j . In practice, most methods
00538 // find the "right" minimum, corresponding to \vec{x} = \vec{x}_0, but there are usually
00539 // many subsidiary minima.
00540 // [Reference: Comput. J. 6 163 (1963).]
00541 
00542  
00543   const Double_t pi = TMath::Pi();
00544   Bool_t ok = kTRUE;
00545   Double_t delta = 0.1;
00546   
00547   for (nf = 5; nf<32;nf +=5) {
00548      TVirtualFitter *min = TVirtualFitter::Fitter(0,nf);
00549      min->SetFCN(TrigoFletcher);
00550      
00551      Double_t arglist[100];
00552      arglist[0] = gVerbose;
00553      min->ExecuteCommand("SET PRINT",arglist,1);
00554      A.ResizeTo(nf,nf);
00555      B.ResizeTo(nf,nf);
00556      x0.ResizeTo(nf);
00557      sx0.ResizeTo(nf);
00558      cx0.ResizeTo(nf);
00559      sx.ResizeTo(nf);
00560      cx.ResizeTo(nf);
00561      v0.ResizeTo(nf);
00562      v.ResizeTo(nf);
00563      r.ResizeTo(nf);
00564      A.Randomize(-100.,100,seed);
00565      B.Randomize(-100.,100,seed);
00566      for (Int_t i = 0; i < nf; i++) {
00567        for (Int_t j = 0; j < nf; j++) {
00568          A(i,j) = Int_t(A(i,j));
00569          B(i,j) = Int_t(B(i,j));
00570        }
00571      }
00572 
00573      x0.Randomize(-pi,pi,seed);
00574      TVectorD x1(nf); x1.Randomize(-delta*pi,delta*pi,seed);
00575      x1+= x0;
00576 
00577      for (Int_t i = 0; i < nf; i++)
00578        min->SetParameter(i,Form("x_%d",i),x1[i],0.01,-pi*(1+delta),+pi*(1+delta));
00579 
00580      arglist[0] = 0;
00581      min->ExecuteCommand("SET NOW",arglist,0);
00582      arglist[0] = 1000; // number of function calls 
00583      arglist[1] = 0.01; // tolerance 
00584      min->ExecuteCommand("MIGRAD",arglist,0);
00585      min->ExecuteCommand("MIGRAD",arglist,2);
00586      min->ExecuteCommand("MINOS",arglist,0);
00587 
00588      Double_t par,we,al,bl;
00589      Char_t parName[32];
00590      for (Int_t i = 0; i < nf; i++) {
00591        min->GetParameter(i,parName,par,we,al,bl);
00592        ok = ok && ( TMath::Abs(par) -TMath::Abs(x0[i]) < gAbsTolerance );
00593        if (!ok) printf("nf=%d, i=%d, par=%g, x0=%g\n",nf,i,par,x0[i]);
00594      }
00595      delete min;
00596   }  
00597 
00598   return ok;
00599 }
00600 
00601 //______________________________________________________________________________
00602 Int_t stressFit(const char *theFitter, Int_t N)
00603 {
00604   TVirtualFitter::SetDefaultFitter(theFitter);
00605   
00606   cout << "******************************************************************" <<endl;
00607   cout << "*  Minimization - S T R E S S suite                              *" <<endl;
00608   cout << "******************************************************************" <<endl;
00609   cout << "******************************************************************" <<endl;
00610 
00611    TStopwatch timer;
00612    timer.Start();
00613 
00614   cout << "*  Starting  S T R E S S  with fitter : "<<TVirtualFitter::GetDefaultFitter() <<endl;
00615   cout << "******************************************************************" <<endl;
00616 
00617   gBenchmark->Start("stressFit");
00618 
00619   Bool_t okRosenBrock    = kTRUE;
00620   Bool_t okWood          = kTRUE;
00621   Bool_t okPowell        = kTRUE;
00622   Bool_t okFletcher      = kTRUE;
00623   Bool_t okGoldStein1    = kTRUE;
00624   Bool_t okGoldStein2    = kTRUE;
00625   Bool_t okTrigoFletcher = kTRUE;
00626   Int_t i;
00627   for (i = 0; i < N; i++)  okWood          = RunWood4();
00628   StatusPrint(1,"Wood",okWood);
00629   for (i = 0; i < N; i++) okRosenBrock    = RunRosenBrock();
00630   StatusPrint(2,"RosenBrock",okRosenBrock);
00631   for (i = 0; i < N; i++) okPowell        = RunPowell();
00632   StatusPrint(3,"Powell",okPowell);
00633   for (i = 0; i < N; i++) okFletcher      = RunFletcher();
00634   StatusPrint(4,"Fletcher",okFletcher);
00635   for (i = 0; i < N; i++) okGoldStein1    = RunGoldStein1();
00636   StatusPrint(5,"GoldStein1",okGoldStein1);
00637   for (i = 0; i < N; i++) okGoldStein2    = RunGoldStein2();
00638   StatusPrint(6,"GoldStein2",okGoldStein2);
00639   okTrigoFletcher = RunTrigoFletcher();
00640   StatusPrint(7,"TrigoFletcher",okTrigoFletcher);
00641 
00642   gBenchmark->Stop("stressFit");
00643 
00644 
00645   //Print table with results
00646   Bool_t UNIX = strcmp(gSystem->GetName(), "Unix") == 0;
00647   printf("******************************************************************\n");
00648   if (UNIX) {
00649      TString sp = gSystem->GetFromPipe("uname -a");
00650      sp.Resize(60);
00651      printf("*  SYS: %s\n",sp.Data());
00652      if (strstr(gSystem->GetBuildNode(),"Linux")) {
00653         sp = gSystem->GetFromPipe("lsb_release -d -s");
00654         printf("*  SYS: %s\n",sp.Data());
00655      }
00656      if (strstr(gSystem->GetBuildNode(),"Darwin")) {
00657         sp  = gSystem->GetFromPipe("sw_vers -productVersion");
00658         sp += " Mac OS X ";
00659         printf("*  SYS: %s\n",sp.Data());
00660      }
00661   } else {
00662     const Char_t *os = gSystem->Getenv("OS");
00663     if (!os) printf("*  SYS: Windows 95\n");
00664     else     printf("*  SYS: %s %s \n",os,gSystem->Getenv("PROCESSOR_IDENTIFIER"));
00665   }
00666   
00667   printf("******************************************************************\n");
00668   gBenchmark->Print("stressFit");
00669 #ifdef __CINT__
00670   Double_t reftime = 86.34; //macbrun interpreted
00671 #else
00672   Double_t reftime = 12.07; //macbrun compiled
00673 #endif
00674   const Double_t rootmarks = 800.*reftime/gBenchmark->GetCpuTime("stressFit");
00675   
00676   printf("******************************************************************\n");
00677   printf("*  ROOTMARKS =%6.1f   *  Root%-8s  %d/%d\n",rootmarks,gROOT->GetVersion(),
00678          gROOT->GetVersionDate(),gROOT->GetVersionTime());
00679   printf("******************************************************************\n");
00680 
00681    return 0;
00682 }
00683 
00684 //_____________________________batch only_____________________
00685 #ifndef __CINT__
00686 
00687 int main(int argc,const char *argv[]) 
00688 {
00689   gBenchmark = new TBenchmark();
00690   const char *fitter = "Minuit";
00691   if (argc > 1)  fitter = argv[1];
00692   if (strcmp(fitter,"Minuit") && strcmp(fitter,"Minuit2") && strcmp(fitter,"Fumili")) {
00693      printf("stressFit illegal option %s, using Minuit instead\n",fitter);
00694      fitter = "Minuit";
00695   }
00696   Int_t N = 2000;
00697   if (argc > 2) N = atoi(argv[2]);
00698   stressFit(fitter,N);  //default is Minuit
00699   return 0;
00700 }
00701 
00702 #endif

Generated on Tue Jul 5 15:15:13 2011 for ROOT_528-00b_version by  doxygen 1.5.1