combinedFit.C

Go to the documentation of this file.
00001 //+ Combined (simultaneous) fit of two histogram with separate functions 
00002 //  and some common parameters
00003 //
00004 // See http://root.cern.ch/phpBB3//viewtopic.php?f=3&t=11740#p50908
00005 // for a modified version working with Fumili or GSLMultiFit 
00006 //
00007 // N.B. this macro must be compiled with ACliC 
00008 //
00009 //Author: L. Moneta - Dec 2010
00010 
00011 #include "Fit/Fitter.h"
00012 #include "Fit/BinData.h"
00013 #include "Fit/Chi2FCN.h"
00014 #include "TH1.h"
00015 #include "TList.h"
00016 #include "Math/WrappedMultiTF1.h"
00017 #include "HFitInterface.h"
00018 #include "TCanvas.h"
00019 #include "TStyle.h"
00020 
00021 
00022 // definition of shared parameter
00023 // background function 
00024 int iparB[2] = { 0,      // exp amplitude in B histo
00025                  2    // exp common parameter 
00026 };
00027 
00028 // signal + background function 
00029 int iparSB[5] = { 1, // exp amplitude in S+B histo
00030                   2, // exp common parameter
00031                   3, // gaussian amplitude
00032                   4, // gaussian mean
00033                   5  // gaussian sigma
00034 };
00035 
00036 struct GlobalChi2 { 
00037    GlobalChi2(  ROOT::Math::IMultiGenFunction & f1,  
00038                 ROOT::Math::IMultiGenFunction & f2) : 
00039       fChi2_1(&f1), fChi2_2(&f2) {}
00040 
00041    // parameter vector is first background (in common 1 and 2) 
00042    // and then is signal (only in 2)
00043    double operator() (const double *par) const {
00044       double p1[2];
00045       for (int i = 0; i < 2; ++i) p1[i] = par[iparB[i] ];
00046 
00047       double p2[5]; 
00048       for (int i = 0; i < 5; ++i) p2[i] = par[iparSB[i] ];
00049 
00050       return (*fChi2_1)(p1) + (*fChi2_2)(p2);
00051    } 
00052 
00053    const  ROOT::Math::IMultiGenFunction * fChi2_1;
00054    const  ROOT::Math::IMultiGenFunction * fChi2_2;
00055 };
00056 
00057 void combinedFit() { 
00058 
00059   TH1D * hB = new TH1D("hB","histo B",100,0,100);
00060   TH1D * hSB = new TH1D("hSB","histo S+B",100, 0,100);
00061 
00062   TF1 * fB = new TF1("fB","expo",0,100);
00063   fB->SetParameters(1,-0.05);
00064   hB->FillRandom("fB");
00065 
00066   TF1 * fS = new TF1("fS","gaus",0,100);
00067   fS->SetParameters(1,30,5);
00068 
00069   hSB->FillRandom("fB",2000);
00070   hSB->FillRandom("fS",1000);
00071 
00072   // perform now global fit 
00073 
00074   TF1 * fSB = new TF1("fSB","expo + gaus(2)",0,100);
00075 
00076   ROOT::Math::WrappedMultiTF1 wfB(*fB,1);
00077   ROOT::Math::WrappedMultiTF1 wfSB(*fSB,1);
00078 
00079   ROOT::Fit::DataOptions opt; 
00080   ROOT::Fit::DataRange rangeB; 
00081   // set the data range
00082   rangeB.SetRange(10,90);
00083   ROOT::Fit::BinData dataB(opt,rangeB); 
00084   ROOT::Fit::FillData(dataB, hB);
00085 
00086   ROOT::Fit::DataRange rangeSB; 
00087   rangeSB.SetRange(10,50);
00088   ROOT::Fit::BinData dataSB(opt,rangeSB); 
00089   ROOT::Fit::FillData(dataSB, hSB);
00090 
00091   ROOT::Fit::Chi2Function chi2_B(dataB, wfB);
00092   ROOT::Fit::Chi2Function chi2_SB(dataSB, wfSB);
00093 
00094   GlobalChi2 globalChi2(chi2_B, chi2_SB);
00095 
00096   ROOT::Fit::Fitter fitter;
00097 
00098   const int Npar = 6; 
00099   double par0[Npar] = { 5,5,-0.1,100, 30,10};
00100 
00101   // create before the parameter settings in order to fix or set range on them
00102   fitter.Config().SetParamsSettings(6,par0);
00103   // fix 5-th parameter  
00104   fitter.Config().ParSettings(4).Fix();
00105   // set limits on the third and 4-th parameter
00106   fitter.Config().ParSettings(2).SetLimits(-10,-1.E-4);
00107   fitter.Config().ParSettings(3).SetLimits(0,10000);
00108   fitter.Config().ParSettings(3).SetStepSize(5);
00109 
00110   fitter.Config().MinimizerOptions().SetPrintLevel(0);
00111   fitter.Config().SetMinimizer("Minuit2","Migrad"); 
00112 
00113   // fit FCN function directly 
00114   // (specify optionally data size and flag to indicate that is a chi2 fit)
00115   fitter.FitFCN(6,globalChi2,0,dataB.Size()+dataSB.Size(),true);
00116   ROOT::Fit::FitResult result = fitter.Result();
00117   result.Print(std::cout);
00118 
00119   TCanvas * c1 = new TCanvas("Simfit","Simultaneous fit of two histograms",
00120                              10,10,700,700);
00121   c1->Divide(1,2);
00122   c1->cd(1);
00123   gStyle->SetOptFit(1111);
00124 
00125   fB->SetFitResult( result, iparB);
00126   fB->SetRange(rangeB().first, rangeB().second);   
00127   fB->SetLineColor(kBlue);
00128   hB->GetListOfFunctions()->Add(fB);
00129   hB->Draw(); 
00130 
00131   c1->cd(2);
00132   fSB->SetFitResult( result, iparSB);
00133   fSB->SetRange(rangeSB().first, rangeSB().second);   
00134   fSB->SetLineColor(kRed);
00135   hSB->GetListOfFunctions()->Add(fSB);
00136   hSB->Draw(); 
00137     
00138 
00139 }

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