Logo ROOT  
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
22import ROOT
23import 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.
27fileNameUrl = "http://root.cern.ch/files/tutorials/df014_CsvDataSource_MuRun2010B.csv"
28fileName = "df014_CsvDataSource_MuRun2010B_py.csv"
29if not os.path.isfile(fileName):
30 ROOT.TFile.Cp(fileNameUrl, fileName)
31
32df = 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.
38filteredEvents = 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.
42invMass = filteredEvents.Histo1D(("invMass", "CMS Opendata: #mu#mu mass;#mu#mu mass [GeV];Events", 512, 2, 110), "m")
43
44c = ROOT.TCanvas()
45c.SetLogx()
46c.SetLogy()
47invMass.Draw()
48c.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.
54fullSpectrum = 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.
58jpsiLow = 2.95
59jpsiHigh = 3.25
60jpsiCut = 'm < %s && m > %s' % (jpsiHigh, jpsiLow)
61jpsi = 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.
65dualCanvas = ROOT.TCanvas("DualCanvas", "DualCanvas", 800, 512)
66dualCanvas.Divide(2, 1)
67leftPad = dualCanvas.cd(1)
68leftPad.SetLogx()
69leftPad.SetLogy()
70fullSpectrum.Draw("Hist")
71dualCanvas.cd(2)
72jpsi.Draw("HistP")
73dualCanvas.SaveAs("df014_jpsi.png")
74
75print("Saved figures to df014_*.png")
RDataFrame FromCSV(std::string_view fileName, bool readHeaders=true, char delimiter=',', Long64_t linesChunkSize=-1LL, std::unordered_map< std::string, char > &&colTypes={})
Factory method to create a CSV RDataFrame.
Definition RCsvDS.cxx:558