36ROOT.ROOT.EnableImplicitMT()
39higgs_header_path = os.path.join(os.sep, str(ROOT.gROOT.GetTutorialDir()) + os.sep,
"dataframe" + os.sep,
40 "df103_NanoAODHiggsAnalysis_python.h")
42ROOT.gInterpreter.Declare(
'#include "{}"'.format(higgs_header_path))
46def reco_higgs_to_2el2mu(df):
47 """Reconstruct Higgs from two electrons and two muons"""
49 df_base = selection_2el2mu(df)
51 df_z_mass = df_base.Define(
"Z_mass",
"compute_z_masses_2el2mu(Electron_pt, Electron_eta, Electron_phi,"
52 " Electron_mass, Muon_pt, Muon_eta, Muon_phi, Muon_mass)")
54 df_z_cut = filter_z_candidates(df_z_mass)
56 df_h_mass = df_z_cut.Define(
"H_mass",
"compute_higgs_mass_2el2mu(Electron_pt, Electron_eta, Electron_phi,"
57 " Electron_mass, Muon_pt, Muon_eta, Muon_phi, Muon_mass)")
62def selection_2el2mu(df):
63 """Select interesting events with two electrons and two muons"""
64 df_ge2el2mu = df.Filter(
"nElectron>=2 && nMuon>=2",
"At least two electrons and two muons")
65 df_eta = df_ge2el2mu.Filter(
"All(abs(Electron_eta)<2.5) && All(abs(Muon_eta)<2.4)",
"Eta cuts")
67 df_pt = df_eta.Filter(
"pt_cuts(Muon_pt, Electron_pt)",
"Pt cuts")
69 df_dr = df_pt.Filter(
"dr_cuts(Muon_eta, Muon_phi, Electron_eta, Electron_phi)",
"Dr cuts")
71 df_iso = df_dr.Filter(
"All(abs(Electron_pfRelIso03_all)<0.40) && All(abs(Muon_pfRelIso04_all)<0.40)",
72 "Require good isolation")
73 df_el_ip3d = df_iso.Define(
"Electron_ip3d_el",
"sqrt(Electron_dxy*Electron_dxy + Electron_dz*Electron_dz)")
74 df_el_sip3d = df_el_ip3d.Define(
"Electron_sip3d_el",
75 "Electron_ip3d_el/sqrt(Electron_dxyErr*Electron_dxyErr + "
76 "Electron_dzErr*Electron_dzErr)")
77 df_el_track = df_el_sip3d.Filter(
"All(Electron_sip3d_el<4) && All(abs(Electron_dxy)<0.5) &&"
78 " All(abs(Electron_dz)<1.0)",
79 "Electron track close to primary vertex with small uncertainty")
80 df_mu_ip3d = df_el_track.Define(
"Muon_ip3d_mu",
"sqrt(Muon_dxy*Muon_dxy + Muon_dz*Muon_dz)")
82 df_mu_sip3d = df_mu_ip3d.Define(
"Muon_sip3d_mu",
83 "Muon_ip3d_mu/sqrt(Muon_dxyErr*Muon_dxyErr + Muon_dzErr*Muon_dzErr)")
84 df_mu_track = df_mu_sip3d.Filter(
"All(Muon_sip3d_mu<4) && All(abs(Muon_dxy)<0.5) && All(abs(Muon_dz)<1.0)",
85 "Muon track close to primary vertex with small uncertainty")
86 df_2p2n = df_mu_track.Filter(
"Sum(Electron_charge)==0 && Sum(Muon_charge)==0",
87 "Two opposite charged electron and muon pairs")
92def reco_higgs_to_4mu(df):
93 """Reconstruct Higgs from four muons"""
95 df_base = selection_4mu(df)
98 df_z_idx = df_base.Define(
"Z_idx",
"reco_zz_to_4l(Muon_pt, Muon_eta, Muon_phi, Muon_mass, Muon_charge)")
101 df_z_dr = df_z_idx.Filter(
"filter_z_dr(Z_idx, Muon_eta, Muon_phi)",
"Delta R separation of muons building Z system")
104 df_z_mass = df_z_dr.Define(
"Z_mass",
"compute_z_masses_4l(Z_idx, Muon_pt, Muon_eta, Muon_phi, Muon_mass)")
107 df_z_cut = filter_z_candidates(df_z_mass)
110 df_h_mass = df_z_cut.Define(
"H_mass",
"compute_higgs_mass_4l(Z_idx, Muon_pt, Muon_eta, Muon_phi, Muon_mass)")
115def selection_4mu(df):
116 """Select interesting events with four muons"""
117 df_ge4m = df.Filter(
"nMuon>=4",
"At least four muons")
119 df_iso = df_ge4m.Filter(
"All(abs(Muon_pfRelIso04_all)<0.40)",
"Require good isolation")
120 df_kin = df_iso.Filter(
"All(Muon_pt>5) && All(abs(Muon_eta)<2.4)",
"Good muon kinematics")
121 df_ip3d = df_kin.Define(
"Muon_ip3d",
"sqrt(Muon_dxy*Muon_dxy + Muon_dz*Muon_dz)")
122 df_sip3d = df_ip3d.Define(
"Muon_sip3d",
"Muon_ip3d/sqrt(Muon_dxyErr*Muon_dxyErr + Muon_dzErr*Muon_dzErr)")
123 df_pv = df_sip3d.Filter(
"All(Muon_sip3d<4) && All(abs(Muon_dxy)<0.5) && All(abs(Muon_dz)<1.0)",
124 "Track close to primary vertex with small uncertainty")
125 df_2p2n = df_pv.Filter(
"nMuon==4 && Sum(Muon_charge==1)==2 && Sum(Muon_charge==-1)==2",
126 "Two positive and two negative muons")
131def filter_z_candidates(df):
132 """Apply selection on reconstructed Z candidates"""
133 df_z1_cut = df.Filter(
"Z_mass[0] > 40 && Z_mass[0] < 120",
"Mass of first Z candidate in [40, 120]")
134 df_z2_cut = df_z1_cut.Filter(
"Z_mass[1] > 12 && Z_mass[1] < 120",
"Mass of second Z candidate in [12, 120]")
139def reco_higgs_to_4el(df):
140 """Reconstruct Higgs from four electrons"""
142 df_base = selection_4el(df)
145 df_z_idx = df_base.Define(
"Z_idx",
146 "reco_zz_to_4l(Electron_pt, Electron_eta, Electron_phi, Electron_mass, Electron_charge)")
149 df_z_dr = df_z_idx.Filter(
"filter_z_dr(Z_idx, Electron_eta, Electron_phi)",
150 "Delta R separation of Electrons building Z system")
153 df_z_mass = df_z_dr.Define(
"Z_mass",
154 "compute_z_masses_4l(Z_idx, Electron_pt, Electron_eta, Electron_phi, Electron_mass)")
157 df_z_cut = filter_z_candidates(df_z_mass)
160 df_h_mass = df_z_cut.Define(
"H_mass",
161 "compute_higgs_mass_4l(Z_idx, Electron_pt, Electron_eta, Electron_phi, Electron_mass)")
166def selection_4el(df):
167 """Select interesting events with four electrons"""
168 df_ge4el = df.Filter(
"nElectron>=4",
"At least our electrons")
169 df_iso = df_ge4el.Filter(
"All(abs(Electron_pfRelIso03_all)<0.40)",
"Require good isolation")
170 df_kin = df_iso.Filter(
"All(Electron_pt>7) && All(abs(Electron_eta)<2.5)",
"Good Electron kinematics")
171 df_ip3d = df_kin.Define(
"Electron_ip3d",
"sqrt(Electron_dxy*Electron_dxy + Electron_dz*Electron_dz)")
172 df_sip3d = df_ip3d.Define(
"Electron_sip3d",
173 "Electron_ip3d/sqrt(Electron_dxyErr*Electron_dxyErr + Electron_dzErr*Electron_dzErr)")
174 df_pv = df_sip3d.Filter(
"All(Electron_sip3d<4) && All(abs(Electron_dxy)<0.5) && All(abs(Electron_dz)<1.0)",
175 "Track close to primary vertex with small uncertainty")
176 df_2p2n = df_pv.Filter(
"nElectron==4 && Sum(Electron_charge==1)==2 && Sum(Electron_charge==-1)==2",
177 "Two positive and two negative electrons")
182def plot(sig, bkg, data, x_label, filename):
184 Plot invariant mass for signal
and background processes
from simulated
185 events overlay the measured data.
188 ROOT.gStyle.SetOptStat(0)
189 ROOT.gStyle.SetTextFont(42)
190 d = ROOT.TCanvas(
"d",
"", 800, 700)
192 ROOT.SetOwnership(d,
False)
193 d.SetLeftMargin(0.15)
202 h_cmb.GetXaxis().SetTitle(x_label)
203 h_cmb.GetXaxis().SetTitleSize(0.04)
204 h_cmb.GetYaxis().SetTitle(
"N_{Events}")
205 h_cmb.GetYaxis().SetTitleSize(0.04)
206 h_cmb.SetLineColor(ROOT.kRed)
207 h_cmb.SetLineWidth(2)
209 h_bkg.SetLineWidth(2)
210 h_bkg.SetFillStyle(1001)
211 h_bkg.SetLineColor(ROOT.kBlack)
212 h_bkg.SetFillColor(ROOT.kAzure - 9)
216 h_data.SetLineWidth(1)
217 h_data.SetMarkerStyle(20)
218 h_data.SetMarkerSize(1.0)
219 h_data.SetMarkerColor(ROOT.kBlack)
220 h_data.SetLineColor(ROOT.kBlack)
223 h_cmb.DrawCopy(
"HIST")
224 h_bkg.DrawCopy(
"HIST SAME")
225 h_data.DrawCopy(
"PE1 SAME")
228 legend = ROOT.TLegend(0.62, 0.70, 0.82, 0.88)
229 legend.SetFillColor(0)
230 legend.SetBorderSize(0)
231 legend.SetTextSize(0.03)
232 legend.AddEntry(h_data,
"Data",
"PE1")
233 legend.AddEntry(h_bkg,
"ZZ",
"f")
234 legend.AddEntry(h_cmb,
"m_{H} = 125 GeV",
"f")
238 cms_label = ROOT.TLatex()
239 cms_label.SetTextSize(0.04)
240 cms_label.DrawLatexNDC(0.16, 0.92,
"#bf{CMS Open Data}")
241 header = ROOT.TLatex()
242 header.SetTextSize(0.03)
243 header.DrawLatexNDC(0.63, 0.92,
"#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}")
254 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/SMHiggsToZZTo4L.root")
261 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/ZZTo4mu.root")
264 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/ZZTo4e.root")
267 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/ZZTo2e2mu.root")
270 doublemu_files = ROOT.std.vector(
"string")(2)
271 doublemu_files[0] =
"root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012B_DoubleMuParked.root"
272 doublemu_files[1] =
"root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012C_DoubleMuParked.root"
276 doubleel_files = ROOT.std.vector(
"string")(2)
277 doubleel_files[0] =
"root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012B_DoubleElectron.root"
278 doubleel_files[1] =
"root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012C_DoubleElectron.root"
289 nevt_ZZTo4mu = 1499064.0
292 nevt_ZZTo4el = 1499093.0
294 xsec_ZZTo2el2mu = 0.18
295 nevt_ZZTo2el2mu = 1497445.0
297 xsec_SMHiggsToZZTo4L = 0.0065
298 nevt_SMHiggsToZZTo4L = 299973.0
301 weight_sig_4mu = luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L
302 weight_bkg_4mu = luminosity * xsec_ZZTo4mu * scale_ZZTo4l / nevt_ZZTo4mu
304 weight_sig_4el = luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L
305 weight_bkg_4el = luminosity * xsec_ZZTo4el * scale_ZZTo4l / nevt_ZZTo4el
307 weight_sig_2el2mu = luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L
308 weight_bkg_2el2mu = luminosity * xsec_ZZTo2el2mu * scale_ZZTo4l / nevt_ZZTo2el2mu
311 df_sig_4mu_reco = reco_higgs_to_4mu(df_sig_4l)
313 df_h_sig_4mu = df_sig_4mu_reco.Define(
"weight",
"{}".format(weight_sig_4mu))\
314 .Histo1D((
"h_sig_4mu",
"", nbins, 70, 180),
"H_mass",
"weight")
316 df_bkg_4mu_reco = reco_higgs_to_4mu(df_bkg_4mu)
318 df_h_bkg_4mu = df_bkg_4mu_reco.Define(
"weight",
"{}".format(weight_bkg_4mu))\
319 .Histo1D((
"h_bkg_4mu",
"", nbins, 70, 180),
"H_mass",
"weight")
321 df_data_4mu_reco = reco_higgs_to_4mu(df_data_doublemu)
323 df_h_data_4mu = df_data_4mu_reco.Define(
"weight",
"1.0")\
324 .Histo1D((
"h_data_4mu",
"", nbins, 70, 180),
"H_mass",
"weight")
327 df_sig_4el_reco = reco_higgs_to_4el(df_sig_4l)
329 df_h_sig_4el = df_sig_4el_reco.Define(
"weight",
"{}".format(weight_sig_4el))\
330 .Histo1D((
"h_sig_4el",
"", nbins, 70, 180),
"H_mass",
"weight")
332 df_bkg_4el_reco = reco_higgs_to_4el(df_bkg_4el)
334 df_h_bkg_4el = df_bkg_4el_reco.Define(
"weight",
"{}".format(weight_bkg_4el))\
335 .Histo1D((
"h_bkg_4el",
"", nbins, 70, 180),
"H_mass",
"weight")
337 df_data_4el_reco = reco_higgs_to_4el(df_data_doubleel)
339 df_h_data_4el = df_data_4el_reco.Define(
"weight",
"1.0")\
340 .Histo1D((
"h_data_4el",
"", nbins, 70, 180),
"H_mass",
"weight")
343 df_sig_2el2mu_reco = reco_higgs_to_2el2mu(df_sig_4l)
345 df_h_sig_2el2mu = df_sig_2el2mu_reco.Define(
"weight",
"{}".format(weight_sig_2el2mu))\
346 .Histo1D((
"h_sig_2el2mu",
"", nbins, 70, 180),
"H_mass",
"weight")
348 df_bkg_2el2mu_reco = reco_higgs_to_2el2mu(df_bkg_2el2mu)
350 df_h_bkg_2el2mu = df_bkg_2el2mu_reco.Define(
"weight",
"{}".format(weight_bkg_2el2mu))\
351 .Histo1D((
"h_bkg_2el2mu",
"", nbins, 70, 180),
"H_mass",
"weight")
353 df_data_2el2mu_reco = reco_higgs_to_2el2mu(df_data_doublemu)
355 df_h_data_2el2mu = df_data_2el2mu_reco.Define(
"weight",
"1.0")\
356 .Histo1D((
"h_data_2el2mu_doublemu",
"", nbins, 70, 180),
"H_mass",
"weight")
359 signal_4mu = df_h_sig_4mu.GetValue()
360 background_4mu = df_h_bkg_4mu.GetValue()
361 data_4mu = df_h_data_4mu.GetValue()
363 signal_4el = df_h_sig_4el.GetValue()
364 background_4el = df_h_bkg_4el.GetValue()
365 data_4el = df_h_data_4el.GetValue()
367 signal_2el2mu = df_h_sig_2el2mu.GetValue()
368 background_2el2mu = df_h_bkg_2el2mu.GetValue()
369 data_2el2mu = df_h_data_2el2mu.GetValue()
372 plot(signal_4mu, background_4mu, data_4mu,
"m_{4#mu} (GeV)",
"higgs_4mu.pdf")
373 plot(signal_4el, background_4el, data_4el,
"m_{4e} (GeV)",
"higgs_4el.pdf")
374 plot(signal_2el2mu, background_2el2mu, data_2el2mu,
"m_{2e2#mu} (GeV)",
"higgs_2el2mu.pdf")
381 h_sig_4l = signal_4mu
382 h_sig_4l.Add(signal_4el)
383 h_sig_4l.Add(signal_2el2mu)
385 h_bkg_4l = background_4mu
386 h_bkg_4l.Add(background_4el)
387 h_bkg_4l.Add(background_2el2mu)
390 h_data_4l.Add(data_4el)
391 h_data_4l.Add(data_2el2mu)
394 plot(h_sig_4l, h_bkg_4l, h_data_4l,
"m_{4l} (GeV)",
"higgs_4l.pdf")
397if __name__ ==
"__main__":
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...