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