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