This tutorial is the Higgs to four lepton analysis from the ATLAS Open Data release in 2020 (http://opendata.atlas.cern/release/2020/documentation/). The data was taken with the ATLAS detector during 2016 at a center-of-mass energy of 13 TeV. The decay of the Standard Model Higgs boson to two Z bosons and subsequently to four leptons is called the "golden channel". The selection leads to a narrow invariant mass peak on top a relatively smooth and small background, revealing the Higgs at 125 GeV.
The analysis is translated to a RDataFrame workflow processing about 300 MB of simulated events and data.
import ROOT
import os
my_spec = os.path.join(ROOT.gROOT.GetTutorialsDir(), "dataframe", "df106_HiggsToFourLeptons_spec.json")
df = df.DefinePerSample("xsecs", 'rdfsampleinfo_.GetD("xsecs")')
df = df.DefinePerSample("lumi", 'rdfsampleinfo_.GetD("lumi")')
df = df.DefinePerSample("sumws", 'rdfsampleinfo_.GetD("sumws")')
df = df.DefinePerSample("sample_category", 'rdfsampleinfo_.GetS("sample_category")')
ROOT.gInterpreter.Declare("""
using cRVecF = const ROOT::RVecF &;
bool GoodElectronsAndMuons(const ROOT::RVecI & type, cRVecF pt, cRVecF eta, cRVecF phi, cRVecF e, cRVecF trackd0pv, cRVecF tracksigd0pv, cRVecF z0)
{
for (size_t i = 0; i < type.size(); i++) {
ROOT::Math::PtEtaPhiEVector p(pt[i] / 1000.0, eta[i], phi[i], e[i] / 1000.0);
if (type[i] == 11) {
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;
} else {
if (abs(trackd0pv[i] / tracksigd0pv[i]) > 5 || abs(z0[i] * sin(p.Theta())) > 0.5) return false;
}
}
return true;
}
""")
df = df.Filter("trigE || trigM")
df = df.Define("good_lep", "abs(lep_eta) < 2.5 && lep_pt > 5000 && lep_ptcone30 / lep_pt < 0.3 && lep_etcone20 / lep_pt < 0.3")\
.Filter("Sum(good_lep) == 4")\
.Filter("Sum(lep_charge[good_lep]) == 0")\
.Define("goodlep_sumtypes", "Sum(lep_type[good_lep])")\
.Filter("goodlep_sumtypes == 44 || goodlep_sumtypes == 52 || goodlep_sumtypes == 48")
df = df.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])")
df = df.Define("goodlep_pt", "lep_pt[good_lep]")\
.Define("goodlep_eta", "lep_eta[good_lep]")\
.Define("goodlep_phi", "lep_phi[good_lep]")\
.Define("goodlep_E", "lep_E[good_lep]")
df = df.Filter("goodlep_pt[0] > 25000 && goodlep_pt[1] > 15000 && goodlep_pt[2] > 10000")
ROOT.gInterpreter.Declare("""
double weights(float scaleFactor_1, float scaleFactor_2, float scaleFactor_3, float scaleFactor_4, float mcWeight_, double xsecs_, double sumws_, double lumi_)
{
double weight_ = scaleFactor_1 * scaleFactor_2 * scaleFactor_3 * scaleFactor_4 * mcWeight_ * xsecs_/sumws_ * lumi_;
return (weight_);
}
""")
df = df.DefinePerSample("reweighting", 'rdfsampleinfo_.Contains("mc")')
df = df.Define("weight", 'double x; if (reweighting) {return weights(scaleFactor_ELE, scaleFactor_MUON, scaleFactor_LepTRIGGER, scaleFactor_PILEUP, mcWeight, xsecs, sumws, lumi);} else {return 1.;}')
ROOT.gInterpreter.Declare("""
float ComputeInvariantMass(cRVecF pt, cRVecF eta, cRVecF phi, cRVecF e)
{
ROOT::Math::PtEtaPhiEVector p1(pt[0], eta[0], phi[0], e[0]);
ROOT::Math::PtEtaPhiEVector p2(pt[1], eta[1], phi[1], e[1]);
ROOT::Math::PtEtaPhiEVector p3(pt[2], eta[2], phi[2], e[2]);
ROOT::Math::PtEtaPhiEVector p4(pt[3], eta[3], phi[3], e[3]);
return (p1 + p2 + p3 + p4).M() / 1000;
}
""")
df = df.Define("m4l", "ComputeInvariantMass(goodlep_pt, goodlep_eta, goodlep_phi, goodlep_E)")
histos = []
for sample_category in ["data", "higgs", "zz", "other"]:
histos.append(df.Filter(f
"sample_category == \"{sample_category}\"").Histo1D(
ROOT.RDF.TH1DModel(f
"{sample_category}",
"m4l", 24, 80, 170),
"m4l",
"weight"))
h_data = histos[0].GetValue()
h_higgs = histos[1].GetValue()
h_zz = histos[2].GetValue()
h_other = histos[3].GetValue()
h_zz.Scale(1.3)
ROOT.gROOT.SetStyle("ATLAS")
c = ROOT.TCanvas("c", "", 600, 600)
pad = ROOT.TPad("upper_pad", "", 0, 0, 1, 1)
pad.SetTickx(False)
pad.SetTicky(False)
pad.Draw()
pad.cd()
stack = ROOT.THStack()
for h, color in zip([h_other, h_zz, h_higgs], [(155, 152, 204), (100, 192, 232), (191, 34, 41)]):
h.SetLineWidth(1)
h.SetLineColor(1)
h.SetFillColor(ROOT.TColor.GetColor(*color))
stack.Add(h)
stack.Draw("HIST")
stack.GetXaxis().SetLabelSize(0.04)
stack.GetXaxis().SetTitleSize(0.045)
stack.GetXaxis().SetTitleOffset(1.3)
stack.GetXaxis().SetTitle("m_{4l}^{H#rightarrow ZZ} [GeV]")
stack.GetYaxis().SetTitle("Events")
stack.GetYaxis().SetLabelSize(0.04)
stack.GetYaxis().SetTitleSize(0.045)
stack.SetMaximum(33)
stack.GetYaxis().ChangeLabel(1, -1, 0)
h_data.SetMarkerStyle(20)
h_data.SetMarkerSize(1.2)
h_data.SetLineWidth(2)
h_data.SetLineColor(ROOT.kBlack)
h_data.Draw("E SAME")
legend = ROOT.TLegend(0.60, 0.65, 0.92, 0.92)
legend.SetTextFont(42)
legend.SetFillStyle(0)
legend.SetBorderSize(0)
legend.SetTextSize(0.04)
legend.SetTextAlign(32)
legend.AddEntry(h_data, "Data" ,"lep")
legend.AddEntry(h_higgs, "Higgs", "f")
legend.AddEntry(h_zz, "ZZ", "f")
legend.AddEntry(h_other, "Other", "f")
legend.Draw("SAME")
text = ROOT.TLatex()
text.SetNDC()
text.SetTextFont(72)
text.SetTextSize(0.045)
text.DrawLatex(0.21, 0.86, "ATLAS")
text.SetTextFont(42)
text.DrawLatex(0.21 + 0.16, 0.86, "Open Data")
text.SetTextSize(0.04)
text.DrawLatex(0.21, 0.80, "#sqrt{s} = 13 TeV, 10 fb^{-1}")
c.SaveAs("df106_HiggsToFourLeptons.png")
print("Saved figure to df106_HiggsToFourLeptons.png")
ROOT::RDataFrame FromSpec(const std::string &jsonFile)
Factory method to create an RDataFrame from a JSON specification file.
void AddProgressbar(ROOT::RDF::RNode df)
A struct which stores the parameters of a TH1D.
|=> | [Elapsed time: 0:15m processing file: 5 / 9 processed evts: 1000 / 165548 6.56e+01 evt/s 41:46m remaining time (per file)]
|=====================================================================================================================================================================================================================> | [Elapsed time: 0:16m processing file: 5 / 9 processed evts: 77000 / 165548 3.65e+04 evt/s 0:02m remaining time (per file)]
|============================================================================================================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:17m processing file: 5 / 9 processed evts: 154000 / 165548 5.00e+04 evt/s 0:00m remaining time (per file)]
|======================================================================================================================================================================================================================================================> | [Elapsed time: 0:18m processing file: 6 / 9 processed evts: 191000 / 356674 4.67e+04 evt/s 0:03m remaining time (per file)]
|========================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:19m processing file: 6 / 9 processed evts: 267000 / 356674 5.24e+04 evt/s 0:01m remaining time (per file)]
|==================================================================================================================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:20m processing file: 6 / 9 processed evts: 336000 / 356674 5.51e+04 evt/s 0:00m remaining time (per file)]
|===================================================================================================================================================================================================> | [Elapsed time: 0:21m processing file: 7 / 9 processed evts: 387000 / 910953 5.44e+04 evt/s 0:09m remaining time (per file)]
|=========================================================================================================================================================================================================================================> | [Elapsed time: 0:22m processing file: 7 / 9 processed evts: 463000 / 910953 5.71e+04 evt/s 0:07m remaining time (per file)]
|============================================================================================================================================================================================================================================================================> | [Elapsed time: 0:23m processing file: 7 / 9 processed evts: 532000 / 910953 5.84e+04 evt/s 0:06m remaining time (per file)]
|===========================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:24m processing file: 7 / 9 processed evts: 593000 / 910953 5.86e+04 evt/s 0:05m remaining time (per file)]
|==================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:25m processing file: 7 / 9 processed evts: 669000 / 910953 6.01e+04 evt/s 0:04m remaining time (per file)]
|=====================================================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:26m processing file: 7 / 9 processed evts: 739000 / 910953 6.09e+04 evt/s 0:02m remaining time (per file)]
|======================================================================================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:27m processing file: 7 / 9 processed evts: 803000 / 910953 6.11e+04 evt/s 0:01m remaining time (per file)]
|============================================================================================================================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:28m processing file: 7 / 9 processed evts: 879000 / 910953 6.21e+04 evt/s 0:00m remaining time (per file)]
|============================================================================================================================================================================================================================================================================================================================================================================================================================================================================> | [Elapsed time: 0:28m processing file: 9 / 9 processed evts: 912000 / 912535 6.29e+04 evt/s 0:00m remaining time (per file)]
Saved figure to df106_HiggsToFourLeptons.png