Logo ROOT  
Reference Guide
rf316_llratioplot.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4///
5/// Multidimensional models: using the likelihood ratio technique to construct a signal
6/// enhanced one-dimensional projection of a multi-dimensional p.d.f.
7///
8/// \macro_image
9/// \macro_output
10/// \macro_code
11///
12/// \date 07/2008
13/// \author Wouter Verkerke
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"
25using namespace RooFit;
26
28{
29
30 // C r e a t e 3 D p d f a n d d a t a
31 // -------------------------------------------
32
33 // Create observables
34 RooRealVar x("x", "x", -5, 5);
35 RooRealVar y("y", "y", -5, 5);
36 RooRealVar z("z", "z", -5, 5);
37
38 // Create signal pdf gauss(x)*gauss(y)*gauss(z)
39 RooGaussian gx("gx", "gx", x, RooConst(0), RooConst(1));
40 RooGaussian gy("gy", "gy", y, RooConst(0), RooConst(1));
41 RooGaussian gz("gz", "gz", z, RooConst(0), RooConst(1));
42 RooProdPdf sig("sig", "sig", RooArgSet(gx, gy, gz));
43
44 // Create background pdf poly(x)*poly(y)*poly(z)
45 RooPolynomial px("px", "px", x, RooArgSet(RooConst(-0.1), RooConst(0.004)));
46 RooPolynomial py("py", "py", y, RooArgSet(RooConst(0.1), RooConst(-0.004)));
47 RooPolynomial pz("pz", "pz", z);
48 RooProdPdf bkg("bkg", "bkg", RooArgSet(px, py, pz));
49
50 // Create composite pdf sig+bkg
51 RooRealVar fsig("fsig", "signal fraction", 0.1, 0., 1.);
52 RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
53
54 RooDataSet *data = model.generate(RooArgSet(x, y, z), 20000);
55
56 // P r o j e c t p d f a n d d a t a o n x
57 // -------------------------------------------------
58
59 // Make plain projection of data and pdf on x observable
60 RooPlot *frame = x.frame(Title("Projection of 3D data and pdf on X"), Bins(40));
61 data->plotOn(frame);
62 model.plotOn(frame);
63
64 // 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
65 // ----------------------------------------------------------------------------------
66
67 // Calculate projection of signal and total likelihood on (y,z) observables
68 // i.e. integrate signal and composite model over x
69 RooAbsPdf *sigyz = sig.createProjection(x);
70 RooAbsPdf *totyz = model.createProjection(x);
71
72 // Construct the log of the signal / signal+background probability
73 RooFormulaVar llratio_func("llratio", "log10(@0)-log10(@1)", RooArgList(*sigyz, *totyz));
74
75 // P l o t d a t a w i t h a L L r a t i o c u t
76 // -------------------------------------------------------
77
78 // Calculate the llratio value for each event in the dataset
79 data->addColumn(llratio_func);
80
81 // Extract the subset of data with large signal likelihood
82 RooDataSet *dataSel = (RooDataSet *)data->reduce(Cut("llratio>0.7"));
83
84 // Make plot frame
85 RooPlot *frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"), Bins(40));
86
87 // Plot select data on frame
88 dataSel->plotOn(frame2);
89
90 // 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
91 // ---------------------------------------------------------------------------------------------
92
93 // Generate large number of events for MC integration of pdf projection
94 RooDataSet *mcprojData = model.generate(RooArgSet(x, y, z), 10000);
95
96 // Calculate LL ratio for each generated event and select MC events with llratio)0.7
97 mcprojData->addColumn(llratio_func);
98 RooDataSet *mcprojDataSel = (RooDataSet *)mcprojData->reduce(Cut("llratio>0.7"));
99
100 // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
101 // on set of events with the same llratio cut as was applied to data
102 model.plotOn(frame2, ProjWData(*mcprojDataSel));
103
104 TCanvas *c = new TCanvas("rf316_llratioplot", "rf316_llratioplot", 800, 400);
105 c->Divide(2);
106 c->cd(1);
107 gPad->SetLeftMargin(0.15);
108 frame->GetYaxis()->SetTitleOffset(1.4);
109 frame->Draw();
110 c->cd(2);
111 gPad->SetLeftMargin(0.15);
112 frame2->GetYaxis()->SetTitleOffset(1.4);
113 frame2->Draw();
114}
#define c(i)
Definition: RSha256.hxx:101
#define gPad
Definition: TVirtualPad.h:287
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:381
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
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:546
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:3342
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
virtual RooAbsArg * addColumn(RooAbsArg &var, Bool_t adjustRange=kTRUE)
Add a column with the values of the given (function) argument to this dataset.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Definition: RooFormulaVar.h:29
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooPolynomial implements a polynomial p.d.f of the form.
Definition: RooPolynomial.h:28
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:27
RooConstVar & RooConst(Double_t val)
RooCmdArg Bins(Int_t nbin)
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg Cut(const char *cutSpec)
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
const char * Title
Definition: TXMLSetup.cxx:67