ROOT
v6-32
Reference Guide
Loading...
Searching...
No Matches
df014_CSVDataSource.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_dataframe
3
## \notebook -draw
4
## Process a CSV file with RDataFrame and the CSV data source.
5
##
6
## This tutorial illustrates how use the RDataFrame in combination with a
7
## RDataSource. In this case we use a TCsvDS. This data source allows to read
8
## a CSV file from a RDataFrame.
9
## As a result of running this tutorial, we will produce plots of the dimuon
10
## spectrum starting from a subset of the CMS collision events of Run2010B.
11
## Dataset Reference:
12
## McCauley, T. (2014). Dimuon event information derived from the Run2010B
13
## public Mu dataset. CERN Open Data Portal.
14
## DOI: [10.7483/OPENDATA.CMS.CB8H.MFFA](http://opendata.cern.ch/record/700).
15
##
16
## \macro_code
17
## \macro_image
18
##
19
## \date October 2017
20
## \author Enric Tejedor (CERN)
21
22
import
ROOT
23
import
os
24
25
# Let's first create a RDF that will read from the CSV file.
26
# The types of the columns will be automatically inferred.
27
fileNameUrl =
"http://root.cern/files/tutorials/df014_CsvDataSource_MuRun2010B.csv"
28
fileName =
"df014_CsvDataSource_MuRun2010B_py.csv"
29
if
not
os.path.isfile
(fileName):
30
ROOT.TFile.Cp
(fileNameUrl, fileName)
31
32
df =
ROOT.RDF.FromCSV
(fileName)
33
34
# Now we will apply a first filter based on two columns of the CSV,
35
# and we will define a new column that will contain the invariant mass.
36
# Note how the new invariant mass column is defined from several other
37
# columns that already existed in the CSV file.
38
filteredEvents =
df.Filter
(
"Q1 * Q2 == -1"
) \
39
.Define(
"m"
,
"sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))"
)
40
41
# Next we create a histogram to hold the invariant mass values and we draw it.
42
invMass =
filteredEvents.Histo1D
((
"invMass"
,
"CMS Opendata: #mu#mu mass;#mu#mu mass [GeV];Events"
, 512, 2, 110),
"m"
)
43
44
c =
ROOT.TCanvas
()
45
c.SetLogx
()
46
c.SetLogy
()
47
invMass.Draw
()
48
c.SaveAs
(
"df014_invMass.png"
)
49
50
# We will now produce a plot also for the J/Psi particle. We will plot
51
# on the same canvas the full spectrum and the zoom in on the J/psi particle.
52
# First we will create the full spectrum histogram from the invariant mass
53
# column, using a different histogram model than before.
54
fullSpectrum =
filteredEvents.Histo1D
((
"Spectrum"
,
"Subset of CMS Run 2010B;#mu#mu mass [GeV];Events"
, 1024, 2, 110),
"m"
)
55
56
# Next we will create the histogram for the J/psi particle, applying first
57
# the corresponding cut.
58
jpsiLow = 2.95
59
jpsiHigh = 3.25
60
jpsiCut =
'm < %s && m > %s'
% (jpsiHigh, jpsiLow)
61
jpsi =
filteredEvents.Filter
(jpsiCut) \
62
.Histo1D((
"jpsi"
,
"Subset of CMS Run 2010B: J/#psi window;#mu#mu mass [GeV];Events"
, 128, jpsiLow, jpsiHigh),
"m"
)
63
64
# Finally we draw the two histograms side by side.
65
dualCanvas =
ROOT.TCanvas
(
"DualCanvas"
,
"DualCanvas"
, 800, 512)
66
dualCanvas.Divide
(2, 1)
67
leftPad =
dualCanvas.cd
(1)
68
leftPad.SetLogx
()
69
leftPad.SetLogy
()
70
fullSpectrum.Draw
(
"Hist"
)
71
dualCanvas.cd
(2)
72
jpsi.SetMarkerStyle
(20)
73
jpsi.Draw
(
"HistP"
)
74
dualCanvas.SaveAs
(
"df014_jpsi.png"
)
75
76
print(
"Saved figures to df014_*.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
dataframe
df014_CSVDataSource.py
ROOT v6-32 - Reference Guide Generated on Fri Jan 24 2025 14:10:44 (GVA Time) using Doxygen 1.10.0