ROOT
v6-34
Reference Guide
Loading...
Searching...
No Matches
df102_NanoAODDimuonAnalysis.C
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 August 2018
22
/// \author Stefan Wunsch (KIT, CERN)
23
24
#include "
ROOT/RDataFrame.hxx
"
25
#include "
ROOT/RDFHelpers.hxx
"
26
#include "
ROOT/RVec.hxx
"
27
#include "
TCanvas.h
"
28
#include "
TH1D.h
"
29
#include "
TLatex.h
"
30
#include "
TStyle.h
"
31
32
using namespace
ROOT::VecOps
;
33
34
void
df102_NanoAODDimuonAnalysis
()
35
{
36
// Enable multi-threading
37
ROOT::EnableImplicitMT
();
38
39
// Create dataframe from NanoAOD files
40
ROOT::RDataFrame
df(
"Events"
,
"root://eospublic.cern.ch//eos/opendata/cms/derived-data/AOD2NanoAODOutreachTool/"
41
"Run2012BC_DoubleMuParked_Muons.root"
);
42
43
// Add ProgressBar
44
ROOT::RDF::Experimental::AddProgressBar
(df);
45
46
// For simplicity, select only events with exactly two muons and require opposite charge
47
auto
df_2mu
= df.Filter(
"nMuon == 2"
,
"Events with exactly two muons"
);
48
auto
df_os
=
df_2mu
.Filter(
"Muon_charge[0] != Muon_charge[1]"
,
"Muons with opposite charge"
);
49
50
// Compute invariant mass of the dimuon system
51
auto
df_mass
=
df_os
.Define(
"Dimuon_mass"
,
InvariantMass<float>
, {
"Muon_pt"
,
"Muon_eta"
,
"Muon_phi"
,
"Muon_mass"
});
52
53
// Make histogram of dimuon mass spectrum. Note how we can set title and axis labels in one go
54
auto
h
=
df_mass
.Histo1D({
"Dimuon_mass"
,
"Dimuon mass;m_{#mu#mu} (GeV);N_{Events}"
, 30000, 0.25, 300},
"Dimuon_mass"
);
55
56
// Request cut-flow report
57
auto
report
= df.Report();
58
59
// Produce plot
60
gStyle
->
SetOptStat
(0);
gStyle
->
SetTextFont
(42);
61
auto
c
=
new
TCanvas
(
"c"
,
""
, 800, 700);
62
c
->SetLogx();
c
->SetLogy();
63
64
h
->GetXaxis()->SetTitleSize(0.04);
65
h
->GetYaxis()->SetTitleSize(0.04);
66
h
->DrawClone();
67
68
TLatex
label
;
label
.SetNDC(
true
);
69
label
.DrawLatex(0.175, 0.740,
"#eta"
);
70
label
.DrawLatex(0.205, 0.775,
"#rho,#omega"
);
71
label
.DrawLatex(0.270, 0.740,
"#phi"
);
72
label
.DrawLatex(0.400, 0.800,
"J/#psi"
);
73
label
.DrawLatex(0.415, 0.670,
"#psi'"
);
74
label
.DrawLatex(0.485, 0.700,
"Y(1,2,3S)"
);
75
label
.DrawLatex(0.755, 0.680,
"Z"
);
76
label
.SetTextSize(0.040);
label
.DrawLatex(0.100, 0.920,
"#bf{CMS Open Data}"
);
77
label
.SetTextSize(0.030);
label
.DrawLatex(0.630, 0.920,
"#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}"
);
78
79
c
->SaveAs(
"dimuon_spectrum.pdf"
);
80
81
// Print cut-flow report
82
report
->Print();
83
}
84
85
int
main
()
86
{
87
df102_NanoAODDimuonAnalysis
();
88
}
main
int main()
Definition
Prototype.cxx:12
RDFHelpers.hxx
RDataFrame.hxx
c
#define c(i)
Definition
RSha256.hxx:101
h
#define h(i)
Definition
RSha256.hxx:106
RVec.hxx
TCanvas.h
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
TH1D.h
TLatex.h
TStyle.h
gStyle
R__EXTERN TStyle * gStyle
Definition
TStyle.h:436
ROOT::Detail::TRangeCast
Definition
TCollection.h:311
ROOT::RDataFrame
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
Definition
RDataFrame.hxx:41
TAttText::SetTextFont
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition
TAttText.h:46
TCanvas
The Canvas class.
Definition
TCanvas.h:23
TLatex
To draw Mathematical Formula.
Definition
TLatex.h:18
TStyle::SetOptStat
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition
TStyle.cxx:1640
ROOT::RDF::Experimental::AddProgressBar
void AddProgressBar(ROOT::RDF::RNode df)
Add ProgressBar to a ROOT::RDF::RNode.
Definition
RDFHelpers.cxx:373
ROOT::VecOps
Definition
TCollectionProxyInfo.h:42
ROOT::EnableImplicitMT
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Definition
TROOT.cxx:539
df102_NanoAODDimuonAnalysis
Definition
df102_NanoAODDimuonAnalysis.py:1
tutorials
dataframe
df102_NanoAODDimuonAnalysis.C
ROOT v6-34 - Reference Guide Generated on Fri Jan 24 2025 14:44:17 (GVA Time) using Doxygen 1.10.0