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_code
9/// \macro_output
10///
11/// \date July 2008
12/// \author Wouter Verkerke
13
14#include "RooRealVar.h"
15#include "RooDataSet.h"
16#include "RooGaussian.h"
17#include "RooPolynomial.h"
18#include "RooAddPdf.h"
19#include "RooProdPdf.h"
20#include "TCanvas.h"
21#include "TAxis.h"
22#include "RooPlot.h"
23using namespace RooFit;
24
26{
27
28 // C r e a t e 3 D p d f a n d d a t a
29 // -------------------------------------------
30
31 // Create observables
32 RooRealVar x("x", "x", -5, 5);
33 RooRealVar y("y", "y", -5, 5);
34 RooRealVar z("z", "z", -5, 5);
35
36 // Create signal pdf gauss(x)*gauss(y)*gauss(z)
37 RooGaussian gx("gx", "gx", x, 0.0, 1.0);
38 RooGaussian gy("gy", "gy", y, 0.0, 1.0);
39 RooGaussian gz("gz", "gz", z, 0.0, 1.0);
40 RooProdPdf sig("sig", "sig", RooArgSet(gx, gy, gz));
41
42 // Create background pdf poly(x)*poly(y)*poly(z)
43 RooPolynomial px("px", "px", x, RooArgSet(-0.1, 0.004));
44 RooPolynomial py("py", "py", y, RooArgSet(0.1, -0.004));
45 RooPolynomial pz("pz", "pz", z);
46 RooProdPdf bkg("bkg", "bkg", RooArgSet(px, py, pz));
47
48 // Create composite pdf sig+bkg
49 RooRealVar fsig("fsig", "signal fraction", 0.1, 0., 1.);
50 RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
51
52 std::unique_ptr<RooDataSet> data{model.generate({x, y, z}, 20000)};
53
54 // P r o j e c t p d f a n d d a t a o n x
55 // -------------------------------------------------
56
57 // Make plain projection of data and pdf on x observable
58 RooPlot *frame = x.frame(Title("Projection of 3D data and pdf on X"), Bins(40));
59 data->plotOn(frame);
60 model.plotOn(frame);
61
62 // 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
63 // ----------------------------------------------------------------------------------
64
65 // Calculate projection of signal and total likelihood on (y,z) observables
66 // i.e. integrate signal and composite model over x
67 RooAbsPdf *sigyz = sig.createProjection(x);
68 RooAbsPdf *totyz = model.createProjection(x);
69
70 // Construct the log of the signal / signal+background probability
71 RooFormulaVar llratio_func("llratio", "log10(@0)-log10(@1)", RooArgList(*sigyz, *totyz));
72
73 // P l o t d a t a w i t h a L L r a t i o c u t
74 // -------------------------------------------------------
75
76 // Calculate the llratio value for each event in the dataset
77 data->addColumn(llratio_func);
78
79 // Extract the subset of data with large signal likelihood
80 std::unique_ptr<RooAbsData> dataSel{data->reduce(Cut("llratio>0.7"))};
81
82 // Make plot frame
83 RooPlot *frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"), Bins(40));
84
85 // Plot select data on frame
86 dataSel->plotOn(frame2);
87
88 // 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
89 // ---------------------------------------------------------------------------------------------
90
91 // Generate large number of events for MC integration of pdf projection
92 std::unique_ptr<RooDataSet> mcprojData{model.generate({x, y, z}, 10000)};
93
94 // Calculate LL ratio for each generated event and select MC events with llratio)0.7
95 mcprojData->addColumn(llratio_func);
96 std::unique_ptr<RooAbsData> mcprojDataSel{mcprojData->reduce(Cut("llratio>0.7"))};
97
98 // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
99 // on set of events with the same llratio cut as was applied to data
100 model.plotOn(frame2, ProjWData(*mcprojDataSel));
101
102 TCanvas *c = new TCanvas("rf316_llratioplot", "rf316_llratioplot", 800, 400);
103 c->Divide(2);
104 c->cd(1);
105 gPad->SetLeftMargin(0.15);
106 frame->GetYaxis()->SetTitleOffset(1.4);
107 frame->Draw();
108 c->cd(2);
109 gPad->SetLeftMargin(0.15);
110 frame2->GetYaxis()->SetTitleOffset(1.4);
111 frame2->Draw();
112}
#define c(i)
Definition RSha256.hxx:101
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
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
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:185
TAxis * GetYaxis() const
Definition RooPlot.cxx:1224
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:597
RooPolynomial implements a polynomial p.d.f of the form.
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:39
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
RooCmdArg ProjWData(const RooAbsData &projData, bool binData=false)
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...
Definition CodegenImpl.h:64
const char * Title
Definition TXMLSetup.cxx:68