ROOT logo

From $ROOTSYS/tutorials/roofit/rf307_fullpereventerrors.C

//////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #307
// 
// Complete example with use of full p.d.f. with per-event errors
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooGaussModel.h"
#include "RooConstVar.h"
#include "RooDecay.h"
#include "RooLandau.h"
#include "RooProdPdf.h"
#include "RooHistPdf.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
using namespace RooFit ;


void rf307_fullpereventerrors()
{
  // B - p h y s i c s   p d f   w i t h   p e r - e v e n t  G a u s s i a n   r e s o l u t i o n
  // ----------------------------------------------------------------------------------------------

  // Observables
  RooRealVar dt("dt","dt",-10,10) ;
  RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;

  // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
  RooRealVar bias("bias","bias",0,-10,10) ;
  RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
  RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;

  // Construct decay(dt) (x) gauss1(dt|dterr)
  RooRealVar tau("tau","tau",1.548) ;
  RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;



  // C o n s t r u c t   e m p i r i c a l   p d f   f o r   p e r - e v e n t   e r r o r
  // -----------------------------------------------------------------

  // Use landau p.d.f to get empirical distribution with long tail
  RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
  RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;

  // Construct a histogram pdf to describe the shape of the dtErr distribution
  RooDataHist* expHistDterr = expDataDterr->binnedClone() ;
  RooHistPdf pdfErr("pdfErr","pdfErr",dterr,*expHistDterr) ;


  // C o n s t r u c t   c o n d i t i o n a l   p r o d u c t   d e c a y _ d m ( d t | d t e r r ) * p d f ( d t e r r )
  // ----------------------------------------------------------------------------------------------------------------------

  // Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
  RooProdPdf model("model","model",pdfErr,Conditional(decay_gm,dt)) ;

  // (Alternatively you could also use the landau shape pdfDtErr)
  //RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;

  

  // S a m p l e,   f i t   a n d   p l o t   p r o d u c t   m o d e l 
  // ------------------------------------------------------------------

  // Specify external dataset with dterr values to use model_dm as conditional p.d.f.
  RooDataSet* data = model.generate(RooArgSet(dt,dterr),10000) ;

  

  // F i t   c o n d i t i o n a l   d e c a y _ d m ( d t | d t e r r )
  // ---------------------------------------------------------------------

  // Specify dterr as conditional observable
  model.fitTo(*data) ;


  
  // P l o t   c o n d i t i o n a l   d e c a y _ d m ( d t | d t e r r )
  // ---------------------------------------------------------------------


  // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
  TH1* hh_model = model.createHistogram("hh_model",dt,Binning(50),YVar(dterr,Binning(50))) ;
  hh_model->SetLineColor(kBlue) ;


  // Make projection of data an dt
  RooPlot* frame = dt.frame(Title("Projection of model(dt|dterr) on dt")) ;
  data->plotOn(frame) ;
  model.plotOn(frame) ;


  // Draw all frames on canvas
  TCanvas* c = new TCanvas("rf307_fullpereventerrors","rf307_fullperventerrors",800, 400);
  c->Divide(2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_model->GetZaxis()->SetTitleOffset(2.5) ; hh_model->Draw("surf") ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;



}
 rf307_fullpereventerrors.C:1
 rf307_fullpereventerrors.C:2
 rf307_fullpereventerrors.C:3
 rf307_fullpereventerrors.C:4
 rf307_fullpereventerrors.C:5
 rf307_fullpereventerrors.C:6
 rf307_fullpereventerrors.C:7
 rf307_fullpereventerrors.C:8
 rf307_fullpereventerrors.C:9
 rf307_fullpereventerrors.C:10
 rf307_fullpereventerrors.C:11
 rf307_fullpereventerrors.C:12
 rf307_fullpereventerrors.C:13
 rf307_fullpereventerrors.C:14
 rf307_fullpereventerrors.C:15
 rf307_fullpereventerrors.C:16
 rf307_fullpereventerrors.C:17
 rf307_fullpereventerrors.C:18
 rf307_fullpereventerrors.C:19
 rf307_fullpereventerrors.C:20
 rf307_fullpereventerrors.C:21
 rf307_fullpereventerrors.C:22
 rf307_fullpereventerrors.C:23
 rf307_fullpereventerrors.C:24
 rf307_fullpereventerrors.C:25
 rf307_fullpereventerrors.C:26
 rf307_fullpereventerrors.C:27
 rf307_fullpereventerrors.C:28
 rf307_fullpereventerrors.C:29
 rf307_fullpereventerrors.C:30
 rf307_fullpereventerrors.C:31
 rf307_fullpereventerrors.C:32
 rf307_fullpereventerrors.C:33
 rf307_fullpereventerrors.C:34
 rf307_fullpereventerrors.C:35
 rf307_fullpereventerrors.C:36
 rf307_fullpereventerrors.C:37
 rf307_fullpereventerrors.C:38
 rf307_fullpereventerrors.C:39
 rf307_fullpereventerrors.C:40
 rf307_fullpereventerrors.C:41
 rf307_fullpereventerrors.C:42
 rf307_fullpereventerrors.C:43
 rf307_fullpereventerrors.C:44
 rf307_fullpereventerrors.C:45
 rf307_fullpereventerrors.C:46
 rf307_fullpereventerrors.C:47
 rf307_fullpereventerrors.C:48
 rf307_fullpereventerrors.C:49
 rf307_fullpereventerrors.C:50
 rf307_fullpereventerrors.C:51
 rf307_fullpereventerrors.C:52
 rf307_fullpereventerrors.C:53
 rf307_fullpereventerrors.C:54
 rf307_fullpereventerrors.C:55
 rf307_fullpereventerrors.C:56
 rf307_fullpereventerrors.C:57
 rf307_fullpereventerrors.C:58
 rf307_fullpereventerrors.C:59
 rf307_fullpereventerrors.C:60
 rf307_fullpereventerrors.C:61
 rf307_fullpereventerrors.C:62
 rf307_fullpereventerrors.C:63
 rf307_fullpereventerrors.C:64
 rf307_fullpereventerrors.C:65
 rf307_fullpereventerrors.C:66
 rf307_fullpereventerrors.C:67
 rf307_fullpereventerrors.C:68
 rf307_fullpereventerrors.C:69
 rf307_fullpereventerrors.C:70
 rf307_fullpereventerrors.C:71
 rf307_fullpereventerrors.C:72
 rf307_fullpereventerrors.C:73
 rf307_fullpereventerrors.C:74
 rf307_fullpereventerrors.C:75
 rf307_fullpereventerrors.C:76
 rf307_fullpereventerrors.C:77
 rf307_fullpereventerrors.C:78
 rf307_fullpereventerrors.C:79
 rf307_fullpereventerrors.C:80
 rf307_fullpereventerrors.C:81
 rf307_fullpereventerrors.C:82
 rf307_fullpereventerrors.C:83
 rf307_fullpereventerrors.C:84
 rf307_fullpereventerrors.C:85
 rf307_fullpereventerrors.C:86
 rf307_fullpereventerrors.C:87
 rf307_fullpereventerrors.C:88
 rf307_fullpereventerrors.C:89
 rf307_fullpereventerrors.C:90
 rf307_fullpereventerrors.C:91
 rf307_fullpereventerrors.C:92
 rf307_fullpereventerrors.C:93
 rf307_fullpereventerrors.C:94
 rf307_fullpereventerrors.C:95
 rf307_fullpereventerrors.C:96
 rf307_fullpereventerrors.C:97
 rf307_fullpereventerrors.C:98
 rf307_fullpereventerrors.C:99
 rf307_fullpereventerrors.C:100
 rf307_fullpereventerrors.C:101
 rf307_fullpereventerrors.C:102
 rf307_fullpereventerrors.C:103
 rf307_fullpereventerrors.C:104
 rf307_fullpereventerrors.C:105
 rf307_fullpereventerrors.C:106
 rf307_fullpereventerrors.C:107
 rf307_fullpereventerrors.C:108
 rf307_fullpereventerrors.C:109
 rf307_fullpereventerrors.C:110
 rf307_fullpereventerrors.C:111
 rf307_fullpereventerrors.C:112
 rf307_fullpereventerrors.C:113
 rf307_fullpereventerrors.C:114
 rf307_fullpereventerrors.C:115