24from ROOT.Experimental import RCanvas, RText, RAttrText, RAttrFont, RAttrFill, RAttrLine, RLegend, RPadPos, RPadExtent, TObjectDrawable
27ROOT.ROOT.EnableImplicitMT()
30path =
"root://eospublic.cern.ch//eos/opendata/atlas/OutreachDatasets/2020-01-22"
32df[
"data"] =
ROOT.RDataFrame(
"mini", (os.path.join(path,
"GamGam/Data/data_{}.GamGam.root".format(x))
for x
in (
"A",
"B",
"C",
"D")))
33df[
"ggH"] =
ROOT.RDataFrame(
"mini", os.path.join(path,
"GamGam/MC/mc_343981.ggH125_gamgam.GamGam.root"))
34df[
"VBF"] =
ROOT.RDataFrame(
"mini", os.path.join(path,
"GamGam/MC/mc_345041.VBFH125_gamgam.GamGam.root"))
35processes = list(df.keys())
38for p
in [
"ggH",
"VBF"]:
39 df[p] = df[p].Define(
"weight",
40 "scaleFactor_PHOTON * scaleFactor_PhotonTRIGGER * scaleFactor_PILEUP * mcWeight");
41df[
"data"] = df[
"data"].Define(
"weight",
"1.0")
46 df[p] = df[p].Filter(
"trigP")
49 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))")\
50 .Filter(
"Sum(goodphotons) == 2")
53 df[p] = df[p].Filter(
"Sum(photon_ptcone30[goodphotons] / photon_pt[goodphotons] < 0.065) == 2")\
54 .Filter(
"Sum(photon_etcone20[goodphotons] / photon_pt[goodphotons] < 0.065) == 2")
57ROOT.gInterpreter.Declare(
59using Vec_t = const ROOT::VecOps::RVec<float>;
60float ComputeInvariantMass(Vec_t& pt, Vec_t& eta, Vec_t& phi, Vec_t& e) {
61 ROOT::Math::PtEtaPhiEVector p1(pt[0], eta[0], phi[0], e[0]);
62 ROOT::Math::PtEtaPhiEVector p2(pt[1], eta[1], phi[1], e[1]);
63 return (p1 + p2).mass() / 1000.0;
71 df[p] = df[p].Define(
"m_yy",
"ComputeInvariantMass(photon_pt[goodphotons], photon_eta[goodphotons], photon_phi[goodphotons], photon_E[goodphotons])")
74 df[p] = df[p].Filter(
"photon_pt[goodphotons][0] / 1000.0 / m_yy > 0.35")\
75 .Filter(
"photon_pt[goodphotons][1] / 1000.0 / m_yy > 0.25")\
76 .Filter(
"m_yy > 105 && m_yy < 160")
79 hists[p] = df[p].Histo1D(
80 ROOT.RDF.TH1DModel(p,
"Diphoton invariant mass; m_{#gamma#gamma} [GeV];Events", 30, 105, 160),
91ggh = hists[
"ggH"].GetValue()
92vbf = hists[
"VBF"].GetValue()
93data = hists[
"data"].GetValue()
101c = RCanvas.Create(
"df104_HiggsToTwoPhotons")
103lower_pad = c.AddPad(RPadPos(0,0.65), RPadExtent(1, 0.35))
104upper_pad = c.AddPad(RPadPos(0,0), RPadExtent(1, 0.65))
106upper_frame = upper_pad.AddFrame()
107upper_frame.margins.bottom = 0
108upper_frame.margins.left = 0.14
109upper_frame.margins.right = 0.05
110upper_frame.x.labels.hide =
True
112lower_frame = lower_pad.AddFrame()
113lower_frame.margins.top = 0
114lower_frame.margins.left = 0.14
115lower_frame.margins.right = 0.05
116lower_frame.margins.bottom = 0.3
119fit = ROOT.TF1(
"fit",
"([0]+[1]*x+[2]*x^2+[3]*x^3)+[4]*exp(-0.5*((x-[5])/[6])^2)", 105, 160)
120fit.FixParameter(5, 125.0)
121fit.FixParameter(4, 119.1)
122fit.FixParameter(6, 2.39)
127data.Fit(
"fit",
"0",
"", 105, 160)
130data.SetMarkerStyle(20)
131data.SetMarkerSize(1.2)
133data.SetLineColor(ROOT.kBlack)
136data.GetYaxis().SetLabelSize(0.045)
137data.GetYaxis().SetTitleSize(0.05)
141data_drawable = upper_pad.Add[TObjectDrawable]()
142data_drawable.Set(data,
"E")
145fit_drawable = upper_pad.Add[TObjectDrawable]()
146fit_drawable.Set(fit,
"SAME")
149bkg = ROOT.TF1(
"bkg",
"([0]+[1]*x+[2]*x^2+[3]*x^3)", 105, 160)
151 bkg.SetParameter(i, fit.GetParameter(i))
155bkg_drawable = upper_pad.Add[TObjectDrawable]()
156bkg_drawable.Set(bkg,
"SAME")
161ggh.Scale(lumi * 0.102 / 55922617.6297)
162vbf.Scale(lumi * 0.008518764 / 3441426.13711)
165higgs_drawable = upper_pad.Add[TObjectDrawable]()
166higgs_drawable.Set(higgs,
"HIST SAME")
170ratiobkg = ROOT.TH1I(
"zero",
"", 10, 105, 160)
171ratiobkg.SetLineColor(4)
172ratiobkg.SetLineStyle(2)
173ratiobkg.SetLineWidth(2)
174ratiobkg.SetMinimum(-125)
175ratiobkg.SetMaximum(250)
176ratiobkg.GetXaxis().SetLabelSize(0.08)
177ratiobkg.GetXaxis().SetTitleSize(0.12)
178ratiobkg.GetXaxis().SetTitleOffset(1.0)
179ratiobkg.GetYaxis().SetLabelSize(0.08)
180ratiobkg.GetYaxis().SetTitleSize(0.09)
181ratiobkg.GetYaxis().SetTitle(
"Data - Bkg.")
182ratiobkg.GetYaxis().CenterTitle()
183ratiobkg.GetYaxis().SetNdivisions(503,
False)
184ratiobkg.GetYaxis().ChangeLabel(-1, -1, 0)
185ratiobkg.GetXaxis().SetTitle(
"m_{#gamma#gamma} [GeV]")
186lower_pad.Add[TObjectDrawable]().Set(ratiobkg,
"AXIS")
188ratiosig = ROOT.TH1F(
"ratiosig",
"ratiosig", 5500, 105, 160)
190ratiosig.SetLineColor(2)
191ratiosig.SetLineStyle(1)
192ratiosig.SetLineWidth(2)
194lower_pad.Add[TObjectDrawable]().Set(ratiosig,
"SAME")
196ratiodata = data.Clone()
197ratiodata.Add(bkg, -1)
198for i
in range(1, data.GetNbinsX()):
199 ratiodata.SetBinError(i, data.GetBinError(i))
201lower_pad.Add[TObjectDrawable]().Set(ratiodata,
"E SAME")
204legend = upper_pad.Draw[RLegend](RPadPos(0.01, 0.01), RPadExtent(0.3, 0.4))
205legend.text.size = 0.05
206legend.text.align = RAttrText.kRightCenter
207legend.text.font = RAttrFont.kArial
208legend.border.style = RAttrLine.kNone
209legend.border.width = 0
210legend.fill.style = RAttrFill.kHollow
212legend.AddEntry(data_drawable,
"Data",
"lme")
213legend.AddEntry(bkg_drawable,
"Background",
"l")
214legend.AddEntry(fit_drawable,
"Signal + Bkg.",
"l")
215legend.AddEntry(higgs_drawable,
"Signal",
"l")
218lbl1 = upper_pad.Draw[RText](RPadPos(0.05, 0.88),
"ATLAS")
220lbl1.text.font = RAttrFont.kArialBoldOblique
222lbl1.text.align = RAttrText.kLeftBottom
223lbl2 = upper_pad.Draw[RText](RPadPos(0.05 + 0.16, 0.88),
"Open Data")
225lbl2.text.font = RAttrFont.kArial
227lbl2.text.align = RAttrText.kLeftBottom
228lbl3 = upper_pad.Draw[RText](RPadPos(0.05, 0.82),
"#sqrt{s} = 13 TeV, 10 fb^{-1}")
230lbl3.text.font = RAttrFont.kArial
232lbl3.text.align = RAttrText.kLeftBottom
240print(
"Saved figure to df104.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.