Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df102_NanoAODDimuonAnalysis.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_dataframe
3## \notebook -js
4## Show how NanoAOD files can be processed with RDataFrame.
5##
6## This tutorial illustrates how NanoAOD files can be processed with ROOT
7## dataframes. The NanoAOD-like input files are filled with 66 mio. events
8## from CMS OpenData containing muon candidates part of 2012 dataset
9## ([DOI: 10.7483/OPENDATA.CMS.YLIC.86ZZ](http://opendata.cern.ch/record/6004)
10## and [DOI: 10.7483/OPENDATA.CMS.M5AD.Y3V3](http://opendata.cern.ch/record/6030)).
11## The macro matches muon pairs and produces an histogram of the dimuon mass
12## spectrum showing resonances up to the Z mass.
13## Note that the bump at 30 GeV is not a resonance but a trigger effect.
14##
15## More details about the dataset can be found on [the CERN Open Data portal](http://opendata.web.cern.ch/record/12341).
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
30df = ROOT.RDataFrame("Events", "root://eospublic.cern.ch//eos/opendata/cms/derived-data/AOD2NanoAODOutreachTool/Run2012BC_DoubleMuParked_Muons.root")
31
32# For simplicity, select only events with exactly two muons and require opposite charge
33df_2mu = df.Filter("nMuon == 2", "Events with exactly two muons")
34df_os = df_2mu.Filter("Muon_charge[0] != Muon_charge[1]", "Muons with opposite charge")
35
36# Compute invariant mass of the dimuon system
37df_mass = df_os.Define("Dimuon_mass", "InvariantMass(Muon_pt, Muon_eta, Muon_phi, Muon_mass)")
38
39# Make histogram of dimuon mass spectrum. Note how we can set titles and axis labels in one go.
40h = df_mass.Histo1D(("Dimuon_mass", "Dimuon mass;m_{#mu#mu} (GeV);N_{Events}", 30000, 0.25, 300), "Dimuon_mass")
41
42# Request cut-flow report
43report = df_mass.Report()
44
45# Produce plot
46ROOT.gStyle.SetOptStat(0); ROOT.gStyle.SetTextFont(42)
47c = ROOT.TCanvas("c", "", 800, 700)
48c.SetLogx(); c.SetLogy()
49
50h.SetTitle("")
51h.GetXaxis().SetTitleSize(0.04)
52h.GetYaxis().SetTitleSize(0.04)
53h.Draw()
54
55label = ROOT.TLatex(); label.SetNDC(True)
56label.DrawLatex(0.175, 0.740, "#eta")
57label.DrawLatex(0.205, 0.775, "#rho,#omega")
58label.DrawLatex(0.270, 0.740, "#phi")
59label.DrawLatex(0.400, 0.800, "J/#psi")
60label.DrawLatex(0.415, 0.670, "#psi'")
61label.DrawLatex(0.485, 0.700, "Y(1,2,3S)")
62label.DrawLatex(0.755, 0.680, "Z")
63label.SetTextSize(0.040); label.DrawLatex(0.100, 0.920, "#bf{CMS Open Data}")
64label.SetTextSize(0.030); label.DrawLatex(0.630, 0.920, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}")
65
66c.SaveAs("dimuon_spectrum.pdf")
67
68# Print cut-flow report
69report.Print()
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...