Logo ROOT   6.14/05
Reference Guide
df014_CSVDataSource.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_dataframe
3 ## \notebook
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 ##
16 ## \date October 2017
17 ## \author Enric Tejedor
18 
19 import ROOT
20 
21 # Let's first create a RDF that will read from the CSV file.
22 # The types of the columns will be automatically inferred.
23 fileName = "df014_CsvDataSource_MuRun2010B.csv"
24 MakeCsvDataFrame = ROOT.ROOT.RDF.MakeCsvDataFrame
25 tdf = MakeCsvDataFrame(fileName)
26 
27 # Now we will apply a first filter based on two columns of the CSV,
28 # and we will define a new column that will contain the invariant mass.
29 # Note how the new invariant mass column is defined from several other
30 # columns that already existed in the CSV file.
31 filteredEvents = tdf.Filter("Q1 * Q2 == -1") \
32  .Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))")
33 
34 # Next we create a histogram to hold the invariant mass values and we draw it.
35 invMass = filteredEvents.Histo1D(("invMass", "CMS Opendata: #mu#mu mass;#mu#mu mass [GeV];Events", 512, 2, 110), "m")
36 
37 c = ROOT.TCanvas()
38 c.SetLogx()
39 c.SetLogy()
40 invMass.Draw()
41 
42 # We will now produce a plot also for the J/Psi particle. We will plot
43 # on the same canvas the full spectrum and the zoom in the J/psi particle.
44 # First we will create the full spectrum histogram from the invariant mass
45 # column, using a different histogram model than before.
46 fullSpectrum = filteredEvents.Histo1D(("Spectrum", "Subset of CMS Run 2010B;#mu#mu mass [GeV];Events", 1024, 2, 110), "m")
47 
48 # Next we will create the histogram for the J/psi particle, applying first
49 # the corresponding cut.
50 jpsiLow = 2.95
51 jpsiHigh = 3.25
52 jpsiCut = 'm < %s && m > %s' % (jpsiHigh, jpsiLow)
53 jpsi = filteredEvents.Filter(jpsiCut) \
54  .Histo1D(("jpsi", "Subset of CMS Run 2010B: J/#psi window;#mu#mu mass [GeV];Events", 128, jpsiLow, jpsiHigh), "m")
55 
56 # Finally we draw the two histograms side by side.
57 dualCanvas = ROOT.TCanvas("DualCanvas", "DualCanvas", 800, 512)
58 dualCanvas.Divide(2, 1)
59 leftPad = dualCanvas.cd(1)
60 leftPad.SetLogx()
61 leftPad.SetLogy()
62 fullSpectrum.Draw("Hist")
63 dualCanvas.cd(2)
64 jpsi.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:439