ROOT logo
//////////////////////////////////////////////////////////////////////////
//
// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #602
// 
// Setting up a chi^2 fit to a binned dataset
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooChebychev.h"
#include "RooAddPdf.h"
#include "RooChi2Var.h"
#include "RooMinuit.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf602_chi2fit()
{

  // S e t u p   m o d e l
  // ---------------------

  // Declare observable x
  RooRealVar x("x","x",0,10) ;

  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their paramaters
  RooRealVar mean("mean","mean of gaussians",5) ;
  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
  RooRealVar sigma2("sigma2","width of gaussians",1) ;

  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;  
  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;  
  
  // Build Chebychev polynomial p.d.f.  
  RooRealVar a0("a0","a0",0.5,0.,1.) ;
  RooRealVar a1("a1","a1",-0.2,0.,1.) ;
  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;

  // Sum the signal components into a composite signal p.d.f.
  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;

  // Sum the composite signal and background 
  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
  RooAddPdf  model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;


  // C r e a t e   b i n n e d   d a t a s e t
  // -----------------------------------------

  RooDataSet* d = model.generate(x,10000) ;
  RooDataHist* dh = d->binnedClone() ;

  // Construct a chi^2 of the data and the model.
  // When a p.d.f. is used in a chi^2 fit, the probability density scaled
  // by the number of events in the dataset to obtain the fit function
  // If model is an extended p.d.f, the expected number events is used
  // instead of the observed number of events.
  model.chi2FitTo(*dh) ;

  // NB: It is also possible to fit a RooAbsReal function to a RooDataHist
  // using chi2FitTo(). 

  // Note that entries with zero bins are _not_ allowed 
  // for a proper chi^2 calculation and will give error
  // messages  
  RooDataSet* dsmall = (RooDataSet*) d->reduce(EventRange(1,100)) ;
  RooDataHist* dhsmall = dsmall->binnedClone() ;
  RooChi2Var chi2_lowstat("chi2_lowstat","chi2",model,*dhsmall) ;
  cout << chi2_lowstat.getVal() << endl ;


}
 rf602_chi2fit.C:1
 rf602_chi2fit.C:2
 rf602_chi2fit.C:3
 rf602_chi2fit.C:4
 rf602_chi2fit.C:5
 rf602_chi2fit.C:6
 rf602_chi2fit.C:7
 rf602_chi2fit.C:8
 rf602_chi2fit.C:9
 rf602_chi2fit.C:10
 rf602_chi2fit.C:11
 rf602_chi2fit.C:12
 rf602_chi2fit.C:13
 rf602_chi2fit.C:14
 rf602_chi2fit.C:15
 rf602_chi2fit.C:16
 rf602_chi2fit.C:17
 rf602_chi2fit.C:18
 rf602_chi2fit.C:19
 rf602_chi2fit.C:20
 rf602_chi2fit.C:21
 rf602_chi2fit.C:22
 rf602_chi2fit.C:23
 rf602_chi2fit.C:24
 rf602_chi2fit.C:25
 rf602_chi2fit.C:26
 rf602_chi2fit.C:27
 rf602_chi2fit.C:28
 rf602_chi2fit.C:29
 rf602_chi2fit.C:30
 rf602_chi2fit.C:31
 rf602_chi2fit.C:32
 rf602_chi2fit.C:33
 rf602_chi2fit.C:34
 rf602_chi2fit.C:35
 rf602_chi2fit.C:36
 rf602_chi2fit.C:37
 rf602_chi2fit.C:38
 rf602_chi2fit.C:39
 rf602_chi2fit.C:40
 rf602_chi2fit.C:41
 rf602_chi2fit.C:42
 rf602_chi2fit.C:43
 rf602_chi2fit.C:44
 rf602_chi2fit.C:45
 rf602_chi2fit.C:46
 rf602_chi2fit.C:47
 rf602_chi2fit.C:48
 rf602_chi2fit.C:49
 rf602_chi2fit.C:50
 rf602_chi2fit.C:51
 rf602_chi2fit.C:52
 rf602_chi2fit.C:53
 rf602_chi2fit.C:54
 rf602_chi2fit.C:55
 rf602_chi2fit.C:56
 rf602_chi2fit.C:57
 rf602_chi2fit.C:58
 rf602_chi2fit.C:59
 rf602_chi2fit.C:60
 rf602_chi2fit.C:61
 rf602_chi2fit.C:62
 rf602_chi2fit.C:63
 rf602_chi2fit.C:64
 rf602_chi2fit.C:65
 rf602_chi2fit.C:66
 rf602_chi2fit.C:67
 rf602_chi2fit.C:68
 rf602_chi2fit.C:69
 rf602_chi2fit.C:70
 rf602_chi2fit.C:71
 rf602_chi2fit.C:72
 rf602_chi2fit.C:73
 rf602_chi2fit.C:74
 rf602_chi2fit.C:75
 rf602_chi2fit.C:76
 rf602_chi2fit.C:77
 rf602_chi2fit.C:78
 rf602_chi2fit.C:79
 rf602_chi2fit.C:80
 rf602_chi2fit.C:81
 rf602_chi2fit.C:82
 rf602_chi2fit.C:83
 rf602_chi2fit.C:84
 rf602_chi2fit.C:85
 rf602_chi2fit.C:86
 rf602_chi2fit.C:87