Logo ROOT  
Reference Guide
df104_HiggsToTwoPhotons.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_dataframe
3## \notebook -draw
4## This tutorial is the Higgs to two photons analysis from the ATLAS Open Data release in 2020
5## (http://opendata.atlas.cern/release/2020/documentation/). The data was taken with the ATLAS detector
6## during 2016 at a center-of-mass energy of 13 TeV. Although the Higgs to two photons decay is very rare,
7## the contribution of the Higgs can be seen as a narrow peak around 125 GeV because of the excellent
8## reconstruction and identification efficiency of photons at the ATLAS experiment.
9##
10## The analysis is translated to a RDataFrame workflow processing 1.7 GB of simulated events and data.
11##
12## \macro_image
13## \macro_code
14## \macro_output
15##
16## \date February 2020
17## \author Stefan Wunsch (KIT, CERN)
18
19import ROOT
20import os
21
22# Enable multi-threading
23ROOT.ROOT.EnableImplicitMT()
24
25# Create a ROOT dataframe for each dataset
26path = "root://eospublic.cern.ch//eos/opendata/atlas/OutreachDatasets/2020-01-22"
27df = {}
28df["data"] = ROOT.RDataFrame("mini", (os.path.join(path, "GamGam/Data/data_{}.GamGam.root".format(x)) for x in ("A", "B", "C", "D")))
29df["ggH"] = ROOT.RDataFrame("mini", os.path.join(path, "GamGam/MC/mc_343981.ggH125_gamgam.GamGam.root"))
30df["VBF"] = ROOT.RDataFrame("mini", os.path.join(path, "GamGam/MC/mc_345041.VBFH125_gamgam.GamGam.root"))
31processes = list(df.keys())
32
33# Apply scale factors and MC weight for simulated events and a weight of 1 for the data
34for p in ["ggH", "VBF"]:
35 df[p] = df[p].Define("weight",
36 "scaleFactor_PHOTON * scaleFactor_PhotonTRIGGER * scaleFactor_PILEUP * mcWeight");
37df["data"] = df["data"].Define("weight", "1.0")
38
39# Select the events for the analysis
40for p in processes:
41 # Apply preselection cut on photon trigger
42 df[p] = df[p].Filter("trigP")
43
44 # Find two good muons with tight ID, pt > 25 GeV and not in the transition region between barrel and encap
45 df[p] = df[p].Define("goodphotons", "photon_isTightID && (photon_pt > 25000) && (abs(photon_eta) < 2.37) && ((abs(photon_eta) < 1.37) || (abs(photon_eta) > 1.52))")\
46 .Filter("Sum(goodphotons) == 2")
47
48 # Take only isolated photons
49 df[p] = df[p].Filter("Sum(photon_ptcone30[goodphotons] / photon_pt[goodphotons] < 0.065) == 2")\
50 .Filter("Sum(photon_etcone20[goodphotons] / photon_pt[goodphotons] < 0.065) == 2")
51
52# Compile a function to compute the invariant mass of the diphoton system
53ROOT.gInterpreter.Declare(
54"""
55using Vec_t = const ROOT::VecOps::RVec<float>;
56float ComputeInvariantMass(Vec_t& pt, Vec_t& eta, Vec_t& phi, Vec_t& e) {
57 ROOT::Math::PtEtaPhiEVector p1(pt[0], eta[0], phi[0], e[0]);
58 ROOT::Math::PtEtaPhiEVector p2(pt[1], eta[1], phi[1], e[1]);
59 return (p1 + p2).mass() / 1000.0;
60}
61""")
62
63# Define a new column with the invariant mass and perform final event selection
64hists = {}
65for p in processes:
66 # Make four vectors and compute invariant mass
67 df[p] = df[p].Define("m_yy", "ComputeInvariantMass(photon_pt[goodphotons], photon_eta[goodphotons], photon_phi[goodphotons], photon_E[goodphotons])")
68
69 # Make additional kinematic cuts and select mass window
70 df[p] = df[p].Filter("photon_pt[goodphotons][0] / 1000.0 / m_yy > 0.35")\
71 .Filter("photon_pt[goodphotons][1] / 1000.0 / m_yy > 0.25")\
72 .Filter("m_yy > 105 && m_yy < 160")
73
74 # Book histogram of the invariant mass with this selection
75 hists[p] = df[p].Histo1D(
76 ROOT.RDF.TH1DModel(p, "Diphoton invariant mass; m_{#gamma#gamma} [GeV];Events", 30, 105, 160),
77 "m_yy", "weight")
78
79# Run the event loop
80ggh = hists["ggH"].GetValue()
81vbf = hists["VBF"].GetValue()
82data = hists["data"].GetValue()
83
84# Create the plot
85
86# Set styles
87ROOT.gROOT.SetStyle("ATLAS")
88
89# Create canvas with pads for main plot and data/MC ratio
90c = ROOT.TCanvas("c", "", 700, 750)
91
92upper_pad = ROOT.TPad("upper_pad", "", 0, 0.35, 1, 1)
93lower_pad = ROOT.TPad("lower_pad", "", 0, 0, 1, 0.35)
94for p in [upper_pad, lower_pad]:
95 p.SetLeftMargin(0.14)
96 p.SetRightMargin(0.05)
97 p.SetTickx(False)
98 p.SetTicky(False)
99upper_pad.SetBottomMargin(0)
100lower_pad.SetTopMargin(0)
101lower_pad.SetBottomMargin(0.3)
102
103upper_pad.Draw()
104lower_pad.Draw()
105
106# Fit signal + background model to data
107upper_pad.cd()
108fit = ROOT.TF1("fit", "([0]+[1]*x+[2]*x^2+[3]*x^3)+[4]*exp(-0.5*((x-[5])/[6])^2)", 105, 160)
109fit.FixParameter(5, 125.0)
110fit.FixParameter(4, 119.1)
111fit.FixParameter(6, 2.39)
112fit.SetLineColor(2)
113fit.SetLineStyle(1)
114fit.SetLineWidth(2)
115data.Fit("fit", "", "E SAME", 105, 160)
116fit.Draw("SAME")
117
118# Draw background
119bkg = ROOT.TF1("bkg", "([0]+[1]*x+[2]*x^2+[3]*x^3)", 105, 160)
120for i in range(4):
121 bkg.SetParameter(i, fit.GetParameter(i))
122bkg.SetLineColor(4)
123bkg.SetLineStyle(2)
124bkg.SetLineWidth(2)
125bkg.Draw("SAME")
126
127# Draw data
128data.SetMarkerStyle(20)
129data.SetMarkerSize(1.2)
130data.SetLineWidth(2)
131data.SetLineColor(ROOT.kBlack)
132data.Draw("E SAME")
133data.SetMinimum(1e-3)
134data.SetMaximum(8e3)
135data.GetYaxis().SetLabelSize(0.045)
136data.GetYaxis().SetTitleSize(0.05)
137data.SetStats(0)
138data.SetTitle("")
139
140# Scale simulated events with luminosity * cross-section / sum of weights
141# and merge to single Higgs signal
142lumi = 10064.0
143ggh.Scale(lumi * 0.102 / 55922617.6297)
144vbf.Scale(lumi * 0.008518764 / 3441426.13711)
145higgs = ggh.Clone()
146higgs.Add(vbf)
147higgs.Draw("HIST SAME")
148
149# Draw ratio
150lower_pad.cd()
151
152ratiobkg = ROOT.TF1("zero", "0", 105, 160)
153ratiobkg.SetLineColor(4)
154ratiobkg.SetLineStyle(2)
155ratiobkg.SetLineWidth(2)
156ratiobkg.SetMinimum(-125)
157ratiobkg.SetMaximum(250)
158ratiobkg.GetXaxis().SetLabelSize(0.08)
159ratiobkg.GetXaxis().SetTitleSize(0.12)
160ratiobkg.GetXaxis().SetTitleOffset(1.0)
161ratiobkg.GetYaxis().SetLabelSize(0.08)
162ratiobkg.GetYaxis().SetTitleSize(0.09)
163ratiobkg.GetYaxis().SetTitle("Data - Bkg.")
164ratiobkg.GetYaxis().CenterTitle()
165ratiobkg.GetYaxis().SetTitleOffset(0.7)
166ratiobkg.GetYaxis().SetNdivisions(503, False)
167ratiobkg.GetYaxis().ChangeLabel(-1, -1, 0)
168ratiobkg.GetXaxis().SetTitle("m_{#gamma#gamma} [GeV]")
169ratiobkg.Draw()
170
171ratiosig = ROOT.TH1F("ratiosig", "ratiosig", 5500, 105, 160)
172ratiosig.Eval(fit)
173ratiosig.SetLineColor(2)
174ratiosig.SetLineStyle(1)
175ratiosig.SetLineWidth(2)
176ratiosig.Add(bkg, -1)
177ratiosig.Draw("SAME")
178
179ratiodata = data.Clone()
180ratiodata.Add(bkg, -1)
181ratiodata.Draw("E SAME")
182for i in range(1, data.GetNbinsX()):
183 ratiodata.SetBinError(i, data.GetBinError(i))
184
185# Add legend
186upper_pad.cd()
187legend = ROOT.TLegend(0.55, 0.55, 0.89, 0.85)
188legend.SetTextFont(42)
189legend.SetFillStyle(0)
190legend.SetBorderSize(0)
191legend.SetTextSize(0.05)
192legend.SetTextAlign(32)
193legend.AddEntry(data, "Data" ,"lep")
194legend.AddEntry(bkg, "Background", "l")
195legend.AddEntry(fit, "Signal + Bkg.", "l")
196legend.AddEntry(higgs, "Signal", "l")
197legend.Draw("SAME")
198
199# Add ATLAS label
200text = ROOT.TLatex()
201text.SetNDC()
202text.SetTextFont(72)
203text.SetTextSize(0.05)
204text.DrawLatex(0.18, 0.84, "ATLAS")
205text.SetTextFont(42)
206text.DrawLatex(0.18 + 0.13, 0.84, "Open Data")
207text.SetTextSize(0.04)
208text.DrawLatex(0.18, 0.78, "#sqrt{s} = 13 TeV, 10 fb^{-1}");
209
210# Save the plot
211c.SaveAs("HiggsToTwoPhotons.pdf");
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:42
RVec< T > Filter(const RVec< T > &v, F &&f)
Create a new collection with the elements passing the filter expressed by the predicate.
Definition: RVec.hxx:939
A struct which stores the parameters of a TH1D.
Definition: HistoModels.hxx:27