Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df106_HiggsToFourLeptons.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_dataframe
3## \notebook -draw
4## The Higgs to four lepton analysis from the ATLAS Open Data release of 2020, with RDataFrame.
5##
6## This tutorial is the Higgs to four lepton analysis from the ATLAS Open Data release in 2020
7## (http://opendata.atlas.cern/release/2020/documentation/). The data was taken with the ATLAS detector
8## during 2016 at a center-of-mass energy of 13 TeV. The decay of the Standard Model Higgs boson
9## to two Z bosons and subsequently to four leptons is called the "golden channel". The selection leads
10## to a narrow invariant mass peak on top a relatively smooth and small background, revealing
11## the Higgs at 125 GeV.
12##
13## The analysis is translated to a RDataFrame workflow processing about 300 MB of simulated events and data.
14##
15## \macro_image
16## \macro_code
17## \macro_output
18##
19## \date March 2020
20## \author Stefan Wunsch (KIT, CERN)
21
22import ROOT
23import json
24import os
25
26# Create a ROOT dataframe for each dataset
27# Note that we load the filenames from the external json file placed in the same folder than this script.
28path = "root://eospublic.cern.ch//eos/opendata/atlas/OutreachDatasets/2020-01-22"
29files = json.load(open(os.path.join(os.path.dirname(os.path.abspath(__file__)), "df106_HiggsToFourLeptons.json")))
30processes = files.keys()
31df = {}
32xsecs = {}
33sumws = {}
34samples = []
35for p in processes:
36 for d in files[p]:
37 # Construct the dataframes
38 folder = d[0] # Folder name
39 sample = d[1] # Sample name
40 xsecs[sample] = d[2] # Cross-section
41 sumws[sample] = d[3] # Sum of weights
42 samples.append(sample)
43 df[sample] = ROOT.RDataFrame("mini", "{}/4lep/{}/{}.4lep.root".format(path, folder, sample))
44
45# Select events for the analysis
46ROOT.gInterpreter.Declare("""
47using cRVecF = const ROOT::RVecF &;
48bool GoodElectronsAndMuons(const ROOT::RVecI & type, cRVecF pt, cRVecF eta, cRVecF phi, cRVecF e, cRVecF trackd0pv, cRVecF tracksigd0pv, cRVecF z0)
49{
50 for (size_t i = 0; i < type.size(); i++) {
51 ROOT::Math::PtEtaPhiEVector p(pt[i] / 1000.0, eta[i], phi[i], e[i] / 1000.0);
52 if (type[i] == 11) {
53 if (pt[i] < 7000 || abs(eta[i]) > 2.47 || abs(trackd0pv[i] / tracksigd0pv[i]) > 5 || abs(z0[i] * sin(p.Theta())) > 0.5) return false;
54 } else {
55 if (abs(trackd0pv[i] / tracksigd0pv[i]) > 5 || abs(z0[i] * sin(p.Theta())) > 0.5) return false;
56 }
57 }
58 return true;
59}
60""")
61
62for s in samples:
63 # Select electron or muon trigger
64 df[s] = df[s].Filter("trigE || trigM")
65
66 # Select events with exactly four good leptons conserving charge and lepton numbers
67 # Note that all collections are RVecs and good_lep is the mask for the good leptons.
68 # The lepton types are PDG numbers and set to 11 or 13 for an electron or muon
69 # irrespective of the charge.
70 df[s] = df[s].Define("good_lep", "abs(lep_eta) < 2.5 && lep_pt > 5000 && lep_ptcone30 / lep_pt < 0.3 && lep_etcone20 / lep_pt < 0.3")\
71 .Filter("Sum(good_lep) == 4")\
72 .Filter("Sum(lep_charge[good_lep]) == 0")\
73 .Define("goodlep_sumtypes", "Sum(lep_type[good_lep])")\
74 .Filter("goodlep_sumtypes == 44 || goodlep_sumtypes == 52 || goodlep_sumtypes == 48")
75
76 # Apply additional cuts depending on lepton flavour
77 df[s] = df[s].Filter("GoodElectronsAndMuons(lep_type[good_lep], lep_pt[good_lep], lep_eta[good_lep], lep_phi[good_lep], lep_E[good_lep], lep_trackd0pvunbiased[good_lep], lep_tracksigd0pvunbiased[good_lep], lep_z0[good_lep])")
78
79 # Create new columns with the kinematics of good leptons
80 df[s] = df[s].Define("goodlep_pt", "lep_pt[good_lep]")\
81 .Define("goodlep_eta", "lep_eta[good_lep]")\
82 .Define("goodlep_phi", "lep_phi[good_lep]")\
83 .Define("goodlep_E", "lep_E[good_lep]")
84
85 # Select leptons with high transverse momentum
86 df[s] = df[s].Filter("goodlep_pt[0] > 25000 && goodlep_pt[1] > 15000 && goodlep_pt[2] > 10000")
87
88# Apply luminosity, scale factors and MC weights for simulated events
89lumi = 10064.0
90for s in samples:
91 if "data" in s:
92 df[s] = df[s].Define("weight", "1.0")
93 else:
94 df[s] = df[s].Define("weight", "scaleFactor_ELE * scaleFactor_MUON * scaleFactor_LepTRIGGER * scaleFactor_PILEUP * mcWeight * {} / {} * {}".format(xsecs[s], sumws[s], lumi))
95
96# Compute invariant mass of the four lepton system and make a histogram
97ROOT.gInterpreter.Declare("""
98float ComputeInvariantMass(cRVecF pt, cRVecF eta, cRVecF phi, cRVecF e)
99{
100 ROOT::Math::PtEtaPhiEVector p1(pt[0], eta[0], phi[0], e[0]);
101 ROOT::Math::PtEtaPhiEVector p2(pt[1], eta[1], phi[1], e[1]);
102 ROOT::Math::PtEtaPhiEVector p3(pt[2], eta[2], phi[2], e[2]);
103 ROOT::Math::PtEtaPhiEVector p4(pt[3], eta[3], phi[3], e[3]);
104 return (p1 + p2 + p3 + p4).M() / 1000;
105}
106""")
107
108histos = {}
109for s in samples:
110 df[s] = df[s].Define("m4l", "ComputeInvariantMass(goodlep_pt, goodlep_eta, goodlep_phi, goodlep_E)")
111 histos[s] = df[s].Histo1D(ROOT.RDF.TH1DModel(s, "m4l", 24, 80, 170), "m4l", "weight")
112
113# Run the event loop and merge histograms of the respective processes
114
115# RunGraphs allows to run the event loops of the separate RDataFrame graphs
116# concurrently. This results in an improved usage of the available resources
117# if each separate RDataFrame can not utilize all available resources, e.g.,
118# because not enough data is available.
119ROOT.RDF.RunGraphs([histos[s] for s in samples])
120
121def merge_histos(label):
122 h = None
123 for i, d in enumerate(files[label]):
124 t = histos[d[1]].GetValue()
125 if i == 0: h = t.Clone()
126 else: h.Add(t)
127 h.SetNameTitle(label, label)
128 return h
129
130data = merge_histos("data")
131higgs = merge_histos("higgs")
132zz = merge_histos("zz")
133other = merge_histos("other")
134
135# Apply MC correction for ZZ due to missing gg->ZZ process
136zz.Scale(1.3)
137
138# Create the plot
139
140# Set styles
141ROOT.gROOT.SetStyle("ATLAS")
142
143# Create canvas with pad
144c = ROOT.TCanvas("c", "", 600, 600)
145pad = ROOT.TPad("upper_pad", "", 0, 0, 1, 1)
146pad.SetTickx(False)
147pad.SetTicky(False)
148pad.Draw()
149pad.cd()
150
151# Draw stack with MC contributions
152stack = ROOT.THStack()
153for h, color in zip([other, zz, higgs], [(155, 152, 204), (100, 192, 232), (191, 34, 41)]):
154 h.SetLineWidth(1)
155 h.SetLineColor(1)
156 h.SetFillColor(ROOT.TColor.GetColor(*color))
157 stack.Add(h)
158stack.Draw("HIST")
159stack.GetXaxis().SetLabelSize(0.04)
160stack.GetXaxis().SetTitleSize(0.045)
161stack.GetXaxis().SetTitleOffset(1.3)
162stack.GetXaxis().SetTitle("m_{4l}^{H#rightarrow ZZ} [GeV]")
163stack.GetYaxis().SetTitle("Events")
164stack.GetYaxis().SetLabelSize(0.04)
165stack.GetYaxis().SetTitleSize(0.045)
166stack.SetMaximum(33)
167stack.GetYaxis().ChangeLabel(1, -1, 0)
168
169# Draw data
170data.SetMarkerStyle(20)
171data.SetMarkerSize(1.2)
172data.SetLineWidth(2)
173data.SetLineColor(ROOT.kBlack)
174data.Draw("E SAME")
175
176# Add legend
177legend = ROOT.TLegend(0.60, 0.65, 0.92, 0.92)
178legend.SetTextFont(42)
179legend.SetFillStyle(0)
180legend.SetBorderSize(0)
181legend.SetTextSize(0.04)
182legend.SetTextAlign(32)
183legend.AddEntry(data, "Data" ,"lep")
184legend.AddEntry(higgs, "Higgs", "f")
185legend.AddEntry(zz, "ZZ", "f")
186legend.AddEntry(other, "Other", "f")
187legend.Draw("SAME")
188
189# Add ATLAS label
190text = ROOT.TLatex()
191text.SetNDC()
192text.SetTextFont(72)
193text.SetTextSize(0.045)
194text.DrawLatex(0.21, 0.86, "ATLAS")
195text.SetTextFont(42)
196text.DrawLatex(0.21 + 0.16, 0.86, "Open Data")
197text.SetTextSize(0.04)
198text.DrawLatex(0.21, 0.80, "#sqrt{s} = 13 TeV, 10 fb^{-1}")
199
200# Save the plot
201c.SaveAs("df106_HiggsToFourLeptons.png")
202print("Saved figure to df106_HiggsToFourLeptons.png")
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTree,...
void RunGraphs(std::vector< RResultHandle > handles)
Trigger the event loop of multiple RDataFrames concurrently.
A struct which stores the parameters of a TH1D.