Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf311_rangeplot.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ##
5 ## Multidimensional models: projecting p.d.f and data ranges in continuous observables
6 ##
7 ## \macro_code
8 ##
9 ## \date February 2018
10 ## \author Clemens Lange, Wouter Verkerke (C++ version)
11 
12 import ROOT
13 
14 # Create 3D pdf and data
15 # -------------------------------------------
16 
17 # Create observables
18 x = ROOT.RooRealVar("x", "x", -5, 5)
19 y = ROOT.RooRealVar("y", "y", -5, 5)
20 z = ROOT.RooRealVar("z", "z", -5, 5)
21 
22 # Create signal pdf gauss(x)*gauss(y)*gauss(z)
23 gx = ROOT.RooGaussian(
24  "gx", "gx", x, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
25 gy = ROOT.RooGaussian(
26  "gy", "gy", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
27 gz = ROOT.RooGaussian(
28  "gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
29 sig = ROOT.RooProdPdf("sig", "sig", ROOT.RooArgList(gx, gy, gz))
30 
31 # Create background pdf poly(x)*poly(y)*poly(z)
32 px = ROOT.RooPolynomial("px", "px", x, ROOT.RooArgList(
33  ROOT.RooFit.RooConst(-0.1), ROOT.RooFit.RooConst(0.004)))
34 py = ROOT.RooPolynomial("py", "py", y, ROOT.RooArgList(
35  ROOT.RooFit.RooConst(0.1), ROOT.RooFit.RooConst(-0.004)))
36 pz = ROOT.RooPolynomial("pz", "pz", z)
37 bkg = ROOT.RooProdPdf("bkg", "bkg", ROOT.RooArgList(px, py, pz))
38 
39 # Create composite pdf sig+bkg
40 fsig = ROOT.RooRealVar("fsig", "signal fraction", 0.1, 0., 1.)
41 model = ROOT.RooAddPdf(
42  "model", "model", ROOT.RooArgList(
43  sig, bkg), ROOT.RooArgList(fsig))
44 
45 data = model.generate(ROOT.RooArgSet(x, y, z), 20000)
46 
47 # Project pdf and data on x
48 # -------------------------------------------------
49 
50 # Make plain projection of data and pdf on x observable
51 frame = x.frame(ROOT.RooFit.Title(
52  "Projection of 3D data and pdf on X"), ROOT.RooFit.Bins(40))
53 data.plotOn(frame)
54 model.plotOn(frame)
55 
56 # Project pdf and data on x in signal range
57 # ----------------------------------------------------------------------------------
58 
59 # Define signal region in y and z observables
60 y.setRange("sigRegion", -1, 1)
61 z.setRange("sigRegion", -1, 1)
62 
63 # Make plot frame
64 frame2 = x.frame(ROOT.RooFit.Title(
65  "Same projection on X in signal range of (Y,Z)"), ROOT.RooFit.Bins(40))
66 
67 # Plot subset of data in which all observables are inside "sigRegion"
68 # For observables that do not have an explicit "sigRegion" range defined (e.g. observable)
69 # an implicit definition is used that is identical to the full range (i.e.
70 # [-5,5] for x)
71 data.plotOn(frame2, ROOT.RooFit.CutRange("sigRegion"))
72 
73 # Project model on x, projected observables (y,z) only in "sigRegion"
74 model.plotOn(frame2, ROOT.RooFit.ProjectionRange("sigRegion"))
75 
76 c = ROOT.TCanvas("rf311_rangeplot", "rf310_rangeplot", 800, 400)
77 c.Divide(2)
78 c.cd(1)
79 ROOT.gPad.SetLeftMargin(0.15)
80 frame.GetYaxis().SetTitleOffset(1.4)
81 frame.Draw()
82 c.cd(2)
83 ROOT.gPad.SetLeftMargin(0.15)
84 frame2.GetYaxis().SetTitleOffset(1.4)
85 frame2.Draw()
86 
87 c.SaveAs("rf311_rangeplot.png")