8.0 How to Fit a Distribution

Introduction

In this section we will do simple fitting to the Z mass peak stored in the root file called myZPeakCRAB.root. This root file is available at /afs/cern.ch/cms/Tutorials/myZPeakCRAB.root. We will fit a Gaussian, Breit-Wigner and Convolution of relativistic Breit-Wigner plus interference term with a Gaussian. We would also use CMS.RooFit tool to fit the distribution.

Fitting the Z mass peak

The main intention of fitting the Z mass peak is to show how to fit a distribution. To do this, will need the myZPeakCRAB.root which is available as above.The different distribution that we would fit to the Z mass peak are:

  • Gaussian

   \begin{displaymath} 	G(x; \mu, \sigma) = \frac{1}{\sqrt{2\pi}\sigma}\exp\left[-\frac{(x-\mu)^2}{2\sigma^2}\right]   \end{displaymath}

  • Relativistic Breit-Wigner

   \begin{displaymath} 	B(m; M, \Gamma) = N \cdot \frac{2}{\pi} \cdot \frac{\Gamma^{2}M^{2}}{(m^{2}-M^{2})^{2}+m^{4}(\Gamma^{2}/M^{2})}   \end{displaymath}

  • Convolution of relativistic Breit-Wigner plus interference term with a Gaussian

   \begin{displaymath} 	P(m) = \int B(m'; M, \Gamma)\cdot G(m-m'; \mu, \sigma) dm'     \end{displaymath}


Some general remarks about fitting a Z peak:

To fit a generator-level Z peak a Breit-Wigner fit makes sense. However, reconstructed-level Z peaks have many detector resolutions that smear the Z mass peak. If the detector resolution is relatively poor, then it is usually good enough to fit a gaussian (since the gaussian detector resolution will overwhelm the inherent Briet-Wigner shape of the peak). If the detector resolution is fairly good, then another option is to fit a Breit-Wigner (for the inherent shape) convoluted with a gaussian (to describe the detector effects).This is in the "no-background" case. If you have backgrounds in your sample (Drell-Yan, cosmics, etc...), and you want to do the fit over a large mass range, then another function needs to be included to take care of this - an exponential is commonly used.

Fitting a Gaussian

Before going further, create a rootlogon.C in your home area on lxplus or cmslpc, for example /afs/cern.ch/user/x/xyz or /uscms/home/malik/. Paste the lines HERE in the rootlogon.C file. If you already have a rootlogon.C file, then make sure that these lines are included in it. You can obtain the input root file from /afs/cern.ch/cms/Tutorials/myZPeakCRAB.root.

There are several options to fit a Gaussian

Using the inbuilt Gaussian in ROOT

Login to ROOT as follows
root -l
and execute the following commands
TFile f("myZPeakCRAB.root");
gStyle->SetOptFit(111111);   
Zmass->Fit("gaus"); 

This will pop up the following histogram. Save this histogram as pdf or postscript or eps file using the menu (File->SaveAs) of the histogram window. myZmass_Gauss_default_nofitparam.png
Quit ROOT as follows

.q

The line gStyle->SetOptFit(111111); enables all the histogram statistics to be displayed. For more options and other information please refer to ROOT documentation.

Using a macro of your own in ROOT

The macro to run is HERE and is called FitZPeak.C. This macro calls another macro HERE named BW.C. So please copy, paste and save them with the corresponding names in YOURWORKINGAREA/CMSSW_10_2_18/src. You can also do this using the command:

wget -O FitZPeak.C https://twiki.cern.ch/twiki/pub/CMSPublic/WorkBookHowToFit/FitZPeak.C.txt
wget -O BW.C https://twiki.cern.ch/twiki/pub/CMSPublic/WorkBookHowToFit/BW.C.txt

Note that now the myZPeakCRAB.root file is opened by executing the macro itself, in addition to fitting the Z mass peak.

To run this macro execute the following command from the area YOURWORKINGAREA/CMSSW_10_2_18/src.

root -l FitZPeak.C

This should pop up a histogram (shown below) and you will find yourself in a ROOT session. You can save this plot from the menu on top of the histogram (File->SaveAs) and quit ROOT session by doing the following:

.q

myZmass_Gauss_macro_nofitparam.png

Here is some explanation of the macro. We have defined the Gaussian (also defined Breit-Wigner, fitted in the next exercise) distribution that we want to fit in the file BW.C as:


Double_t mygauss(Double_t * x, Double_t * par)
{
  Double_t arg = 0;
  if (par[2]<0) par[2]=-par[2];  // par[2]: sigma
  if (par[2] != 0) arg = (x[0] - par[1])/par[2];  // par[1]: mean
 
 //return par[0]*BIN_SIZE*TMath::Exp(-0.5*arg*arg)/
  //   (TMath::Sqrt(2*TMath::Pi())*par[2]); 
   return par[0]*TMath::Exp(-0.5*arg*arg)/
     (TMath::Sqrt(2*TMath::Pi())*par[2]); // par[0] is constant
 
}
 

par[0], par[1] and par[2], respectively, are the constant, mean and the sigma parameters. Also x[0] mean the x-axis variable. BW.C is called by FitZPeak.C in the code line gROOT->LoadMacro("BW.C"); . The initial values of the three fitted are defined as follows in the FitZPeak.C

func->SetParameter(0,1.0);   func->SetParName(0,"const");  
func->SetParameter(2,5.0);   func->SetParName(2,"sigma");  
func->SetParameter(1,95.0);  func->SetParName(1,"mean");

Using a macro in RooFit

Before we start, have a look at the twiki CMS.RooFit to get a feeling for it. The macro to fit Z mass peak in RooFit is HERE and is called RooFitMacro.C. So please copy and paste the lines from HERE and save them with the name RooFitMacro.C in YOURWORKINGAREA/CMSSW_10_2_18/src. You can also do this using the command:

wget https://twiki.cern.ch/twiki/pub/CMSPublic/WorkBookHowToFit/RooFitMacro.C

Then execute the following:

root -l RooFitMacro.C

This should pop a histogram (shown below) and you will find yourself in a ROOT session. You can save this plot from the menu on top of the histogram (File->SaveAs) and quit ROOT session by doing the following:

.q

Note that if you look at the macro HERE, it has the following code lines:

 RooGaussian gauss("gauss","gauss",x,mean,sigma);
 // RooBreitWigner gauss("gauss","gauss",x,mean,sigma);
 // RooVoigtian gauss("gauss","gauss",x,mean,width,sigma);
To do a Gaussian, we simply commented out the Breit-Wigner and Voigtian (convolution of Breit-Wigner and Gaussian).

myZmass_Gauss_RooFit_nofitparam.png

Fitting a Breit-Wigner

Using a macro of your own ROOT

To make fit Breit-Wigner, we first uncomment out the Breit-Wigner and comment out the Gaussian part in FitZPeak.C as follows (using /* and */) :

////////////////
//For Gaussian//
///////////////
/*
TF1 *func = new TF1("mygauss",mygauss,massMIN, massMAX,3); 
func->SetParameter(0,1.0);   func->SetParName(0,"const");  
func->SetParameter(2,5.0);   func->SetParName(2,"sigma");  
func->SetParameter(1,95.0);     func->SetParName(1,"mean");

Z_mass->Fit("mygauss","QR");
TF1 *fit = Z_mass->GetFunction("mygauss");
*/
/////////////////////
// For Breit-Wigner//
////////////////////
TF1 *func = new TF1("mybw",mybw,massMIN, massMAX,3);
func->SetParameter(0,1.0);   func->SetParName(0,"const");
func->SetParameter(2,5.0);     func->SetParName(1,"sigma");
func->SetParameter(1,95.0);    func->SetParName(2,"mean");

Z_mass->Fit("mybw","QR");
TF1 *fit = Z_mass->GetFunction("mybw");

Then execute the following:

root -l FitZPeak.C

This should pop a histogram (shown below) and you will find yourself in ROOT session. You can save this plot from the menu on top of the histogram (File->SaveAs) and quit ROOT session by doing the following:

.q

myZmass_BW_macro_nofitparam.png

Using a macro in RooFit

Before we proceed we need to uncomment and comment out few lines in RooFitMacro.C to have them look as follows:

 //RooGaussian gauss("gauss","gauss",x,mean,sigma);
 RooBreitWigner gauss("gauss","gauss",x,mean,sigma);
 // RooVoigtian gauss("gauss","gauss",x,mean,width,sigma);

Then execute

root -l RooFitMacro.C

This should pop a histogram and you will find yourself in ROOT session. You can save this plot from the menu on top of the histogram (File->SaveAs) and quit ROOT session by doing the following:

.q

myZmass_BW_RooFit_nofitparam.png

Fitting a Convolution of Gaussian and Breit-Wigner

Using a macro in RooFit

Before we proceed we need to uncomment and comment out few lines in RooFitMacro.C to have them look as follows:

 //RooGaussian gauss("gauss","gauss",x,mean,sigma);
 // RooBreitWigner gauss("gauss","gauss",x,mean,sigma);
  RooVoigtian gauss("gauss","gauss",x,mean,width,sigma);

Then execute

root -l RooFitMacro.C

This should pop a histogram and you will find yourself in ROOT session. You can save this plot from the menu on top of the histogram (File->SaveAs) and quit ROOT session by doing the following:

.q

myZmass_ConvBWGauss_RooFit_nofitparam.png

Review status

Reviewer/Editor and Date (copy from screen) Comments
Samir Guragain - 27 Apr 2012 Changed CMSSW_3_3_0 to 424 and correctly added/modified RooFitMacro

Responsible: SudhirMalik
Last reviewed by: SamirGuragain - 27-Apr-2012

Topic attachments
I Attachment History Action Size Date Who Comment
Texttxt BW.C.txt r1 manage 0.8 K 2010-03-29 - 05:23 SudhirMalik  
Texttxt FitZPeak.C.txt r1 manage 1.4 K 2010-03-29 - 05:22 SudhirMalik  
C source code filec RooFitMacro.C r1 manage 1.6 K 2012-04-27 - 21:42 UnknownUser  
PNGpng myZmass_BW_RooFit_nofitparam.png r1 manage 28.2 K 2010-03-29 - 05:26 SudhirMalik  
PNGpng myZmass_BW_macro_nofitparam.png r1 manage 21.7 K 2010-03-29 - 05:25 SudhirMalik  
PNGpng myZmass_ConvBWGauss_RooFit_nofitparam.png r1 manage 21.9 K 2010-03-29 - 05:27 SudhirMalik  
PNGpng myZmass_Gauss_RooFit_nofitparam.png r1 manage 21.0 K 2010-03-29 - 05:41 SudhirMalik  
PNGpng myZmass_Gauss_default_nofitparam.png r1 manage 18.8 K 2010-03-29 - 05:31 SudhirMalik  
PNGpng myZmass_Gauss_macro_nofitparam.png r1 manage 21.3 K 2010-03-29 - 05:25 SudhirMalik  
Texttxt rootlogon.C.txt r2 r1 manage 0.2 K 2010-03-29 - 19:44 SudhirMalik  
Edit | Attach | Watch | Print version | History: r10 < r9 < r8 < r7 < r6 | Backlinks | Raw View | WYSIWYG | More topic actions
Topic revision: r10 - 2020-10-13 - NitishDhingra


ESSENTIALS

ADVANCED TOPICS


 
    • Cern Search Icon Cern Search
    • TWiki Search Icon TWiki Search
    • Google Search Icon Google Search

    CMSPublic All webs login

This site is powered by the TWiki collaboration platform Powered by PerlCopyright &© 2008-2023 by the contributing authors. All material on this collaboration platform is the property of the contributing authors.
or Ideas, requests, problems regarding TWiki? use Discourse or Send feedback