Logo ROOT   6.08/07
Reference Guide
rf306_condpereventerrors.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #306
5 ///
6 /// Complete example with use of conditional p.d.f. with per-event errors
7 ///
8 /// \macro_image
9 /// \macro_output
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 "RooConstVar.h"
18 #include "RooGaussModel.h"
19 #include "RooDecay.h"
20 #include "RooLandau.h"
21 #include "RooPlot.h"
22 #include "TCanvas.h"
23 #include "TAxis.h"
24 #include "TH2D.h"
25 using namespace RooFit ;
26 
27 
28 
29 void rf306_condpereventerrors()
30 {
31  // 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
32  // ----------------------------------------------------------------------------------------------
33 
34  // Observables
35  RooRealVar dt("dt","dt",-10,10) ;
36  RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
37 
38  // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
39  RooRealVar bias("bias","bias",0,-10,10) ;
40  RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
41  RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
42 
43  // Construct decay(dt) (x) gauss1(dt|dterr)
44  RooRealVar tau("tau","tau",1.548) ;
45  RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
46 
47 
48 
49  // C o n s t r u c t f a k e ' e x t e r n a l ' d a t a w i t h p e r - e v e n t e r r o r
50  // ------------------------------------------------------------------------------------------------------
51 
52  // Use landau p.d.f to get somewhat realistic distribution with long tail
53  RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
54  RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
55 
56 
57 
58  // S a m p l e d a t a f r o m c o n d i t i o n a l d e c a y _ g m ( d t | d t e r r )
59  // ---------------------------------------------------------------------------------------------
60 
61  // Specify external dataset with dterr values to use decay_dm as conditional p.d.f.
62  RooDataSet* data = decay_gm.generate(dt,ProtoData(*expDataDterr)) ;
63 
64 
65 
66  // 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 )
67  // ---------------------------------------------------------------------
68 
69  // Specify dterr as conditional observable
70  decay_gm.fitTo(*data,ConditionalObservables(dterr)) ;
71 
72 
73 
74  // 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 )
75  // ---------------------------------------------------------------------
76 
77 
78  // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
79  TH1* hh_decay = decay_gm.createHistogram("hh_decay",dt,Binning(50),YVar(dterr,Binning(50))) ;
80  hh_decay->SetLineColor(kBlue) ;
81 
82 
83  // Plot decay_gm(dt|dterr) at various values of dterr
84  RooPlot* frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr")) ;
85  for (Int_t ibin=0 ; ibin<100 ; ibin+=20) {
86  dterr.setBin(ibin) ;
87  decay_gm.plotOn(frame,Normalization(5.)) ;
88  }
89 
90 
91  // Make projection of data an dt
92  RooPlot* frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt")) ;
93  data->plotOn(frame2) ;
94 
95  // Make projection of decay(dt|dterr) on dt.
96  //
97  // Instead of integrating out dterr, make a weighted average of curves
98  // at values dterr_i as given in the external dataset.
99  // (The kTRUE argument bins the data before projection to speed up the process)
100  decay_gm.plotOn(frame2,ProjWData(*expDataDterr,kTRUE)) ;
101 
102 
103 
104  // Draw all frames on canvas
105  TCanvas* c = new TCanvas("rf306_condpereventerrors","rf306_condperventerrors",1200, 400);
106  c->Divide(3) ;
107  c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_decay->GetZaxis()->SetTitleOffset(2.5) ; hh_decay->Draw("surf") ;
108  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
109  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
110 
111 
112 }
113 
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
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
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
return c
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
int Int_t
Definition: RtypesCore.h:41
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
Landau Distribution p.d.f.
Definition: RooLandau.h:24
RooCmdArg Title(const char *name)
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
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
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:2851
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 ProtoData(const RooDataSet &protoData, Bool_t randomizeOrder=kFALSE, Bool_t resample=kFALSE)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg Normalization(Double_t scaleFactor)
The TH1 histogram class.
Definition: TH1.h:80
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
RooCmdArg ConditionalObservables(const RooArgSet &set)
const Bool_t kTRUE
Definition: Rtypes.h:91
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)