Logo ROOT   6.07/09
Reference Guide
rf307_fullpereventerrors.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #307
5 ///
6 /// Complete example with use of full p.d.f. with per-event errors
7 ///
8 ///
9 ///
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooGaussModel.h"
18 #include "RooConstVar.h"
19 #include "RooDecay.h"
20 #include "RooLandau.h"
21 #include "RooProdPdf.h"
22 #include "RooHistPdf.h"
23 #include "RooPlot.h"
24 #include "TCanvas.h"
25 #include "TAxis.h"
26 #include "TH1.h"
27 using namespace RooFit ;
28 
29 
30 void rf307_fullpereventerrors()
31 {
32  // 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
33  // ----------------------------------------------------------------------------------------------
34 
35  // Observables
36  RooRealVar dt("dt","dt",-10,10) ;
37  RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
38 
39  // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
40  RooRealVar bias("bias","bias",0,-10,10) ;
41  RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
42  RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
43 
44  // Construct decay(dt) (x) gauss1(dt|dterr)
45  RooRealVar tau("tau","tau",1.548) ;
46  RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
47 
48 
49 
50  // 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
51  // -----------------------------------------------------------------
52 
53  // Use landau p.d.f to get empirical distribution with long tail
54  RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
55  RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
56 
57  // Construct a histogram pdf to describe the shape of the dtErr distribution
58  RooDataHist* expHistDterr = expDataDterr->binnedClone() ;
59  RooHistPdf pdfErr("pdfErr","pdfErr",dterr,*expHistDterr) ;
60 
61 
62  // 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 )
63  // ----------------------------------------------------------------------------------------------------------------------
64 
65  // Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
66  RooProdPdf model("model","model",pdfErr,Conditional(decay_gm,dt)) ;
67 
68  // (Alternatively you could also use the landau shape pdfDtErr)
69  //RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;
70 
71 
72 
73  // 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
74  // ------------------------------------------------------------------
75 
76  // Specify external dataset with dterr values to use model_dm as conditional p.d.f.
77  RooDataSet* data = model.generate(RooArgSet(dt,dterr),10000) ;
78 
79 
80 
81  // 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 )
82  // ---------------------------------------------------------------------
83 
84  // Specify dterr as conditional observable
85  model.fitTo(*data) ;
86 
87 
88 
89  // 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 )
90  // ---------------------------------------------------------------------
91 
92 
93  // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
94  TH1* hh_model = model.createHistogram("hh_model",dt,Binning(50),YVar(dterr,Binning(50))) ;
95  hh_model->SetLineColor(kBlue) ;
96 
97 
98  // Make projection of data an dt
99  RooPlot* frame = dt.frame(Title("Projection of model(dt|dterr) on dt")) ;
100  data->plotOn(frame) ;
101  model.plotOn(frame) ;
102 
103 
104  // Draw all frames on canvas
105  TCanvas* c = new TCanvas("rf307_fullpereventerrors","rf307_fullperventerrors",800, 400);
106  c->Divide(2) ;
107  c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_model->GetZaxis()->SetTitleOffset(2.5) ; hh_model->Draw("surf") ;
108  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
109 
110 
111 
112 }
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:262
return c
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
RooCmdArg Conditional(const RooArgSet &pdfSet, const RooArgSet &depSet, Bool_t depsAreCond=kFALSE)
Landau Distribution p.d.f.
Definition: RooLandau.h:24
RooCmdArg Title(const char *name)
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
const Double_t sigma
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooHistPdf implements a probablity density function sampled from a multidimensional histogram...
Definition: RooHistPdf.h:28
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2853
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
The Canvas class.
Definition: TCanvas.h:41
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
The TH1 histogram class.
Definition: TH1.h:80
RooDataHist * binnedClone(const char *newName=0, const char *newTitle=0) const
Return binned clone of this dataset.
Definition: RooDataSet.cxx:981
TAxis * GetZaxis()
Definition: TH1.h:326
RooConstVar & RooConst(Double_t val)
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1089
#define gPad
Definition: TVirtualPad.h:289
Definition: Rtypes.h:61
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26
RooCmdArg Binning(const RooAbsBinning &binning)