Logo ROOT   6.10/09
Reference Guide
rf306_condpereventerrors.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
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 // 07/2008 - Wouter Verkerke
10 //
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 "RooGaussModel.h"
20 #include "RooDecay.h"
21 #include "RooLandau.h"
22 #include "RooPlot.h"
23 #include "TCanvas.h"
24 #include "TH2D.h"
25 using namespace RooFit ;
26 
27 
28 
29 class TestBasic306 : public RooFitTestUnit
30 {
31 public:
32  TestBasic306(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Conditional use of per-event error p.d.f. F(t|dt)",refFile,writeRef,verbose) {} ;
33  Bool_t testCode() {
34 
35  // 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
36  // ----------------------------------------------------------------------------------------------
37 
38  // Observables
39  RooRealVar dt("dt","dt",-10,10) ;
40  RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
41 
42  // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
43  RooRealVar bias("bias","bias",0,-10,10) ;
44  RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
45  RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
46 
47  // Construct decay(dt) (x) gauss1(dt|dterr)
48  RooRealVar tau("tau","tau",1.548) ;
49  RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
50 
51 
52 
53  // 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
54  // ------------------------------------------------------------------------------------------------------
55 
56  // Use landau p.d.f to get somewhat realistic distribution with long tail
57  RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
58  RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
59 
60 
61 
62  // 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 )
63  // ---------------------------------------------------------------------------------------------
64 
65  // Specify external dataset with dterr values to use decay_dm as conditional p.d.f.
66  RooDataSet* data = decay_gm.generate(dt,ProtoData(*expDataDterr)) ;
67 
68 
69 
70  // 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 )
71  // ---------------------------------------------------------------------
72 
73  // Specify dterr as conditional observable
74  decay_gm.fitTo(*data,ConditionalObservables(dterr)) ;
75 
76 
77 
78  // 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 )
79  // ---------------------------------------------------------------------
80 
81 
82  // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
83  TH1* hh_decay = decay_gm.createHistogram("hh_decay",dt,Binning(50),YVar(dterr,Binning(50))) ;
84  hh_decay->SetLineColor(kBlue) ;
85 
86 
87  // Plot decay_gm(dt|dterr) at various values of dterr
88  RooPlot* frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr")) ;
89  for (Int_t ibin=0 ; ibin<100 ; ibin+=20) {
90  dterr.setBin(ibin) ;
91  decay_gm.plotOn(frame,Normalization(5.),Name(Form("curve_slice_%d",ibin))) ;
92  }
93 
94 
95  // Make projection of data an dt
96  RooPlot* frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt")) ;
97  data->plotOn(frame2) ;
98 
99  // Make projection of decay(dt|dterr) on dt.
100  //
101  // Instead of integrating out dterr, make a weighted average of curves
102  // at values dterr_i as given in the external dataset.
103  // (The kTRUE argument bins the data before projection to speed up the process)
104  decay_gm.plotOn(frame2,ProjWData(*expDataDterr,kTRUE)) ;
105 
106 
107  regTH(hh_decay,"rf306_model2d") ;
108  regPlot(frame,"rf306_plot1") ;
109  regPlot(frame2,"rf306_plot2") ;
110 
111  delete expDataDterr ;
112  delete data ;
113 
114  return kTRUE ;
115  }
116 } ;
117 
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
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
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:36
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
RooCmdArg Name(const char *name)
bool verbose
char * Form(const char *fmt,...)
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:56
RooConstVar & RooConst(Double_t val)
Definition: Rtypes.h:56
RooCmdArg ConditionalObservables(const RooArgSet &set)
const Bool_t kTRUE
Definition: RtypesCore.h:91
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26
RooCmdArg Binning(const RooAbsBinning &binning)