Logo ROOT  
Reference Guide
df102_NanoAODDimuonAnalysis.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_dataframe
3## \notebook -js
4## This tutorial illustrates how NanoAOD files can be processed with ROOT
5## dataframes. The NanoAOD-like input files are filled with 66 mio. events
6## from CMS OpenData containing muon candidates part of 2012 dataset
7## ([DOI: 10.7483/OPENDATA.CMS.YLIC.86ZZ](http://opendata.cern.ch/record/6004)
8## and [DOI: 10.7483/OPENDATA.CMS.M5AD.Y3V3](http://opendata.cern.ch/record/6030)).
9## The macro matches muon pairs and produces an histogram of the dimuon mass
10## spectrum showing resonances up to the Z mass.
11## Note that the bump at 30 GeV is not a resonance but a trigger effect.
12##
13## Some more details about the dataset:
14## - It contains about 66 millions events (muon and electron collections, plus some other information, e.g. about primary vertices)
15## - It spans two compressed ROOT files located on EOS for about a total size of 7.5 GB.
16##
17## \macro_image
18## \macro_code
19## \macro_output
20##
21## \date April 2019
22## \author Stefan Wunsch (KIT, CERN)
23
24import ROOT
25
26# Enable multi-threading
27ROOT.ROOT.EnableImplicitMT()
28
29# Create dataframe from NanoAOD files
30files = ROOT.std.vector("string")(2)
31files[0] = "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012B_DoubleMuParked.root"
32files[1] = "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012C_DoubleMuParked.root"
33df = ROOT.RDataFrame("Events", files)
34
35# For simplicity, select only events with exactly two muons and require opposite charge
36df_2mu = df.Filter("nMuon == 2", "Events with exactly two muons")
37df_os = df_2mu.Filter("Muon_charge[0] != Muon_charge[1]", "Muons with opposite charge")
38
39# Compute invariant mass of the dimuon system
40df_mass = df_os.Define("Dimuon_mass", "InvariantMass(Muon_pt, Muon_eta, Muon_phi, Muon_mass)")
41
42# Make histogram of dimuon mass spectrum
43h = df_mass.Histo1D(("Dimuon_mass", "Dimuon_mass", 30000, 0.25, 300), "Dimuon_mass")
44
45# Request cut-flow report
46report = df_mass.Report()
47
48# Produce plot
49ROOT.gStyle.SetOptStat(0); ROOT.gStyle.SetTextFont(42)
50c = ROOT.TCanvas("c", "", 800, 700)
51c.SetLogx(); c.SetLogy()
52
53h.SetTitle("")
54h.GetXaxis().SetTitle("m_{#mu#mu} (GeV)"); h.GetXaxis().SetTitleSize(0.04)
55h.GetYaxis().SetTitle("N_{Events}"); h.GetYaxis().SetTitleSize(0.04)
56h.Draw()
57
58label = ROOT.TLatex(); label.SetNDC(True)
59label.DrawLatex(0.175, 0.740, "#eta")
60label.DrawLatex(0.205, 0.775, "#rho,#omega")
61label.DrawLatex(0.270, 0.740, "#phi")
62label.DrawLatex(0.400, 0.800, "J/#psi")
63label.DrawLatex(0.415, 0.670, "#psi'")
64label.DrawLatex(0.485, 0.700, "Y(1,2,3S)")
65label.DrawLatex(0.755, 0.680, "Z")
66label.SetTextSize(0.040); label.DrawLatex(0.100, 0.920, "#bf{CMS Open Data}")
67label.SetTextSize(0.030); label.DrawLatex(0.630, 0.920, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}")
68
69c.SaveAs("dimuon_spectrum.pdf")
70
71# Print cut-flow report
72report.Print()
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:42