ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf306_condpereventerrors.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #306
4 ///
5 /// Complete example with use of conditional p.d.f. with per-event errors
6 ///
7 ///
8 ///
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooConstVar.h"
20 #include "RooGaussModel.h"
21 #include "RooDecay.h"
22 #include "RooLandau.h"
23 #include "RooPlot.h"
24 #include "TCanvas.h"
25 #include "TAxis.h"
26 #include "TH2D.h"
27 using namespace RooFit ;
28 
29 
30 
31 void rf306_condpereventerrors()
32 {
33  // 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
34  // ----------------------------------------------------------------------------------------------
35 
36  // Observables
37  RooRealVar dt("dt","dt",-10,10) ;
38  RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
39 
40  // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
41  RooRealVar bias("bias","bias",0,-10,10) ;
42  RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
43  RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
44 
45  // Construct decay(dt) (x) gauss1(dt|dterr)
46  RooRealVar tau("tau","tau",1.548) ;
47  RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
48 
49 
50 
51  // 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
52  // ------------------------------------------------------------------------------------------------------
53 
54  // Use landau p.d.f to get somewhat realistic distribution with long tail
55  RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
56  RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
57 
58 
59 
60  // 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 )
61  // ---------------------------------------------------------------------------------------------
62 
63  // Specify external dataset with dterr values to use decay_dm as conditional p.d.f.
64  RooDataSet* data = decay_gm.generate(dt,ProtoData(*expDataDterr)) ;
65 
66 
67 
68  // 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 )
69  // ---------------------------------------------------------------------
70 
71  // Specify dterr as conditional observable
72  decay_gm.fitTo(*data,ConditionalObservables(dterr)) ;
73 
74 
75 
76  // 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 )
77  // ---------------------------------------------------------------------
78 
79 
80  // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
81  TH1* hh_decay = decay_gm.createHistogram("hh_decay",dt,Binning(50),YVar(dterr,Binning(50))) ;
82  hh_decay->SetLineColor(kBlue) ;
83 
84 
85  // Plot decay_gm(dt|dterr) at various values of dterr
86  RooPlot* frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr")) ;
87  for (Int_t ibin=0 ; ibin<100 ; ibin+=20) {
88  dterr.setBin(ibin) ;
89  decay_gm.plotOn(frame,Normalization(5.)) ;
90  }
91 
92 
93  // Make projection of data an dt
94  RooPlot* frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt")) ;
95  data->plotOn(frame2) ;
96 
97  // Make projection of decay(dt|dterr) on dt.
98  //
99  // Instead of integrating out dterr, make a weighted average of curves
100  // at values dterr_i as given in the external dataset.
101  // (The kTRUE argument bins the data before projection to speed up the process)
102  decay_gm.plotOn(frame2,ProjWData(*expDataDterr,kTRUE)) ;
103 
104 
105 
106  // Draw all frames on canvas
107  TCanvas* c = new TCanvas("rf306_condpereventerrors","rf306_condperventerrors",1200, 400);
108  c->Divide(3) ;
109  c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_decay->GetZaxis()->SetTitleOffset(2.5) ; hh_decay->Draw("surf") ;
110  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
111  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
112 
113 
114 }
115 
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:245
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
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)
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:626
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
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
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:321
RooConstVar & RooConst(Double_t val)
#define gPad
Definition: TVirtualPad.h:288
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)