Logo ROOT   6.07/09
Reference Guide
rf316_llratioplot.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #316
5 ///
6 /// Using the likelihood ratio technique to construct a signal enhanced
7 /// one-dimensional projection of a multi-dimensional p.d.f.
8 ///
9 /// \macro_image
10 /// \macro_output
11 /// \macro_code
12 /// \author 07/2008 - Wouter Verkerke
13 
14 
15 #include "RooRealVar.h"
16 #include "RooDataSet.h"
17 #include "RooGaussian.h"
18 #include "RooConstVar.h"
19 #include "RooPolynomial.h"
20 #include "RooAddPdf.h"
21 #include "RooProdPdf.h"
22 #include "TCanvas.h"
23 #include "TAxis.h"
24 #include "RooPlot.h"
25 using namespace RooFit ;
26 
27 
28 void rf316_llratioplot()
29 {
30 
31  // C r e a t e 3 D p d f a n d d a t a
32  // -------------------------------------------
33 
34  // Create observables
35  RooRealVar x("x","x",-5,5) ;
36  RooRealVar y("y","y",-5,5) ;
37  RooRealVar z("z","z",-5,5) ;
38 
39  // Create signal pdf gauss(x)*gauss(y)*gauss(z)
40  RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
41  RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
42  RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
43  RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;
44 
45  // Create background pdf poly(x)*poly(y)*poly(z)
46  RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
47  RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
48  RooPolynomial pz("pz","pz",z) ;
49  RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;
50 
51  // Create composite pdf sig+bkg
52  RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
53  RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;
54 
55  RooDataSet* data = model.generate(RooArgSet(x,y,z),20000) ;
56 
57 
58 
59  // P r o j e c t p d f a n d d a t a o n x
60  // -------------------------------------------------
61 
62  // Make plain projection of data and pdf on x observable
63  RooPlot* frame = x.frame(Title("Projection of 3D data and pdf on X"),Bins(40)) ;
64  data->plotOn(frame) ;
65  model.plotOn(frame) ;
66 
67 
68 
69  // D e f i n e p r o j e c t e d s i g n a l l i k e l i h o o d r a t i o
70  // ----------------------------------------------------------------------------------
71 
72  // Calculate projection of signal and total likelihood on (y,z) observables
73  // i.e. integrate signal and composite model over x
74  RooAbsPdf* sigyz = sig.createProjection(x) ;
75  RooAbsPdf* totyz = model.createProjection(x) ;
76 
77  // Construct the log of the signal / signal+background probability
78  RooFormulaVar llratio_func("llratio","log10(@0)-log10(@1)",RooArgList(*sigyz,*totyz)) ;
79 
80 
81 
82  // P l o t d a t a w i t h a L L r a t i o c u t
83  // -------------------------------------------------------
84 
85  // Calculate the llratio value for each event in the dataset
86  data->addColumn(llratio_func) ;
87 
88  // Extract the subset of data with large signal likelihood
89  RooDataSet* dataSel = (RooDataSet*) data->reduce(Cut("llratio>0.7")) ;
90 
91  // Make plot frame
92  RooPlot* frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"),Bins(40)) ;
93 
94  // Plot select data on frame
95  dataSel->plotOn(frame2) ;
96 
97 
98 
99  // M a k e M C p r o j e c t i o n o f p d f w i t h s a m e L L r a t i o c u t
100  // ---------------------------------------------------------------------------------------------
101 
102  // Generate large number of events for MC integration of pdf projection
103  RooDataSet* mcprojData = model.generate(RooArgSet(x,y,z),10000) ;
104 
105  // Calculate LL ratio for each generated event and select MC events with llratio)0.7
106  mcprojData->addColumn(llratio_func) ;
107  RooDataSet* mcprojDataSel = (RooDataSet*) mcprojData->reduce(Cut("llratio>0.7")) ;
108 
109  // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
110  // on set of events with the same llratio cut as was applied to data
111  model.plotOn(frame2,ProjWData(*mcprojDataSel)) ;
112 
113 
114 
115  TCanvas* c = new TCanvas("rf316_llratioplot","rf316_llratioplot",800,400) ;
116  c->Divide(2) ;
117  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
118  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
119 
120 
121 
122 }
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
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooCmdArg Cut(const char *cutSpec)
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
RooAbsData * reduce(const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg())
Create a reduced copy of this dataset.
Definition: RooAbsData.cxx:343
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
virtual RooAbsArg * addColumn(RooAbsArg &var, Bool_t adjustRange=kTRUE)
Add a column with the values of the given (function) argument to this dataset.
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
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
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Definition: RooAbsPdf.cxx:2989
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
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
Double_t y[n]
Definition: legend1.C:17
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
RooConstVar & RooConst(Double_t val)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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
RooCmdArg Bins(Int_t nbin)
#define gPad
Definition: TVirtualPad.h:289
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559