Logo ROOT   6.18/05
Reference Guide
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 enhanced one-dimensional
5/// projection of a multi-dimensional p.d.f.
6///
7/// \macro_image
8/// \macro_output
9/// \macro_code
10/// \author 07/2008 - Wouter Verkerke
11
12#include "RooRealVar.h"
13#include "RooDataSet.h"
14#include "RooGaussian.h"
15#include "RooConstVar.h"
16#include "RooPolynomial.h"
17#include "RooAddPdf.h"
18#include "RooProdPdf.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit;
23
25{
26
27 // C r e a t e 3 D p d f a n d d a t a
28 // -------------------------------------------
29
30 // Create observables
31 RooRealVar x("x", "x", -5, 5);
32 RooRealVar y("y", "y", -5, 5);
33 RooRealVar z("z", "z", -5, 5);
34
35 // Create signal pdf gauss(x)*gauss(y)*gauss(z)
36 RooGaussian gx("gx", "gx", x, RooConst(0), RooConst(1));
37 RooGaussian gy("gy", "gy", y, RooConst(0), RooConst(1));
38 RooGaussian gz("gz", "gz", z, RooConst(0), RooConst(1));
39 RooProdPdf sig("sig", "sig", RooArgSet(gx, gy, gz));
40
41 // Create background pdf poly(x)*poly(y)*poly(z)
42 RooPolynomial px("px", "px", x, RooArgSet(RooConst(-0.1), RooConst(0.004)));
43 RooPolynomial py("py", "py", y, RooArgSet(RooConst(0.1), RooConst(-0.004)));
44 RooPolynomial pz("pz", "pz", z);
45 RooProdPdf bkg("bkg", "bkg", RooArgSet(px, py, pz));
46
47 // Create composite pdf sig+bkg
48 RooRealVar fsig("fsig", "signal fraction", 0.1, 0., 1.);
49 RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
50
51 RooDataSet *data = model.generate(RooArgSet(x, y, z), 20000);
52
53 // P r o j e c t p d f a n d d a t a o n x
54 // -------------------------------------------------
55
56 // Make plain projection of data and pdf on x observable
57 RooPlot *frame = x.frame(Title("Projection of 3D data and pdf on X"), Bins(40));
58 data->plotOn(frame);
59 model.plotOn(frame);
60
61 // 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
62 // ----------------------------------------------------------------------------------
63
64 // Calculate projection of signal and total likelihood on (y,z) observables
65 // i.e. integrate signal and composite model over x
66 RooAbsPdf *sigyz = sig.createProjection(x);
67 RooAbsPdf *totyz = model.createProjection(x);
68
69 // Construct the log of the signal / signal+background probability
70 RooFormulaVar llratio_func("llratio", "log10(@0)-log10(@1)", RooArgList(*sigyz, *totyz));
71
72 // P l o t d a t a w i t h a L L r a t i o c u t
73 // -------------------------------------------------------
74
75 // Calculate the llratio value for each event in the dataset
76 data->addColumn(llratio_func);
77
78 // Extract the subset of data with large signal likelihood
79 RooDataSet *dataSel = (RooDataSet *)data->reduce(Cut("llratio>0.7"));
80
81 // Make plot frame
82 RooPlot *frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"), Bins(40));
83
84 // Plot select data on frame
85 dataSel->plotOn(frame2);
86
87 // 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
88 // ---------------------------------------------------------------------------------------------
89
90 // Generate large number of events for MC integration of pdf projection
91 RooDataSet *mcprojData = model.generate(RooArgSet(x, y, z), 10000);
92
93 // Calculate LL ratio for each generated event and select MC events with llratio)0.7
94 mcprojData->addColumn(llratio_func);
95 RooDataSet *mcprojDataSel = (RooDataSet *)mcprojData->reduce(Cut("llratio>0.7"));
96
97 // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
98 // on set of events with the same llratio cut as was applied to data
99 model.plotOn(frame2, ProjWData(*mcprojDataSel));
100
101 TCanvas *c = new TCanvas("rf316_llratioplot", "rf316_llratioplot", 800, 400);
102 c->Divide(2);
103 c->cd(1);
104 gPad->SetLeftMargin(0.15);
105 frame->GetYaxis()->SetTitleOffset(1.4);
106 frame->Draw();
107 c->cd(2);
108 gPad->SetLeftMargin(0.15);
109 frame2->GetYaxis()->SetTitleOffset(1.4);
110 frame2->Draw();
111}
#define c(i)
Definition: RSha256.hxx:101
#define gPad
Definition: TVirtualPad.h:286
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: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:2979
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:31
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:27
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:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
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 fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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:294
The Canvas class.
Definition: TCanvas.h:31
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
Template specialisation used in RooAbsArg:
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooConstVar & RooConst(Double_t val)
RooCmdArg Cut(const char *cutSpec)
RooCmdArg Bins(Int_t nbin)
const char * Title
Definition: TXMLSetup.cxx:67