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