00001 #include "TH1.h"
00002 #include "TF1.h"
00003 #include "Fit/BinData.h"
00004 #include "Fit/Chi2FCN.h"
00005
00006 #include "Math/WrappedMultiTF1.h"
00007 #include "Math/Minimizer.h"
00008 #include "Math/GeneticMinimizer.h"
00009 #include "Math/Factory.h"
00010 #include "HFitInterface.h"
00011
00012 #include "TMath.h"
00013
00014 #include <TStopwatch.h>
00015
00016
00017
00018
00019 Double_t background(Double_t *x, Double_t *par) {
00020 return par[0] + par[1]*x[0];
00021 }
00022
00023
00024 Double_t gaussianPeak(Double_t *x, Double_t *par) {
00025 return par[0]*TMath::Gaus(x[0],par[1], par[2]);
00026 }
00027
00028
00029 Double_t fitFunction(Double_t *x, Double_t *par) {
00030 return background(x,par) + gaussianPeak(x,&par[2]) + gaussianPeak(x,&par[5]);
00031 }
00032
00033
00034 double par0[8] = { 1, 0.05, 10 , 2, 0.5 , 10 , 7 , 1. };
00035 const int ndata = 10000;
00036 const double gAbsTolerance = 0.1;
00037 const int gVerbose = 0;
00038
00039 using std::cout;
00040 using std::endl;
00041
00042 int GAMinimize(ROOT::Math::IMultiGenFunction& chi2Func, double& xm1, double& xm2)
00043 {
00044
00045 ROOT::Math::GeneticMinimizer* min = new ROOT::Math::GeneticMinimizer();
00046 if (min == 0) {
00047 cout << "Error creating minimizer " << endl;
00048 return -1;
00049 }
00050
00051
00052 min->SetMaxFunctionCalls(1000000);
00053 min->SetMaxIterations(100000);
00054 min->SetTolerance(gAbsTolerance);
00055 min->SetPrintLevel(gVerbose);
00056 min->SetFunction(chi2Func);
00057 min->SetParameters(100, 300);
00058
00059
00060
00061 double x0[8];
00062 std::copy(par0, par0 + 8, x0);
00063 x0[3] = xm1;
00064 x0[6] = xm2;
00065
00066 for (unsigned int i = 0; i < chi2Func.NDim(); ++i) {
00067 #ifdef DEBUG
00068 cout << "set variable " << i << " to value " << x0[i] << endl;
00069 #endif
00070 if ( i == 3 || i == 6 )
00071 min->SetLimitedVariable(i,"x" + ROOT::Math::Util::ToString(i),x0[i], 0.1,2,8);
00072 else
00073 min->SetLimitedVariable(i,"x" + ROOT::Math::Util::ToString(i),x0[i], 0.1,x0[i]-2,x0[i]+2);
00074 }
00075
00076
00077
00078 if ( !min->Minimize() ) return 1;
00079 min->MinValue();
00080
00081
00082 cout << "Min values by GeneticMinimizer: " << min->X()[3] << " " << min->X()[6] << endl;
00083 xm1 = min->X()[3];
00084 xm2 = min->X()[6];
00085
00086 return 0;
00087 }
00088
00089
00090 int GAMinTutorial()
00091 {
00092 double x1=0, x2=0;
00093
00094
00095 TF1 * fitFunc = new TF1("fitFunc",fitFunction,0,10,8);
00096 fitFunc->SetParameters(par0);
00097
00098
00099 TH1D * h1 = new TH1D("h1","h1",100, 0, 10);
00100 for (int i = 0; i < ndata; ++i) {
00101 h1->Fill(fitFunc->GetRandom() );
00102 }
00103
00104
00105 ROOT::Fit::BinData d;
00106 ROOT::Fit::FillData(d,h1);
00107 ROOT::Math::WrappedMultiTF1 f(*fitFunc);
00108
00109
00110 ROOT::Fit::Chi2FCN<ROOT::Math::IMultiGenFunction> chi2Func(d,f);
00111
00112
00113 TStopwatch t;
00114 t.Start();
00115 GAMinimize(chi2Func, x1,x2);
00116 t.Stop();
00117 cout << "Time :\t " << t.RealTime() << " " << t.CpuTime() << endl;
00118
00119 return 0;
00120 }
00121
00122 int main(int argc, char **argv)
00123 {
00124 int status = GAMinTutorial();
00125
00126 return status;
00127 }