Logo ROOT   6.16/01
Reference Guide
df014_CSVDataSource.C
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
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.
24 auto fileName = "df014_CsvDataSource_MuRun2010B.csv";
25 auto tdf = ROOT::RDF::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 auto filteredEvents =
32 tdf.Filter("Q1 * Q2 == -1")
33 .Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))");
34
35 // Next we create a histogram to hold the invariant mass values and we draw it.
36 auto invMass =
37 filteredEvents.Histo1D({"invMass", "CMS Opendata: #mu#mu mass;#mu#mu mass [GeV];Events", 512, 2, 110}, "m");
38
39 auto c = new TCanvas();
40 c->SetLogx();
41 c->SetLogy();
42 invMass->DrawClone();
43
44 // We will now produce a plot also for the J/Psi particle. We will plot
45 // on the same canvas the full spectrum and the zoom in the J/psi particle.
46 // First we will create the full spectrum histogram from the invariant mass
47 // column, using a different histogram model than before.
48 auto fullSpectrum =
49 filteredEvents.Histo1D({"Spectrum", "Subset of CMS Run 2010B;#mu#mu mass [GeV];Events", 1024, 2, 110}, "m");
50
51 // Next we will create the histogram for the J/psi particle, applying first
52 // the corresponding cut.
53 double jpsiLow = 2.95;
54 double jpsiHigh = 3.25;
55 auto jpsiCut = [jpsiLow, jpsiHigh](double m) { return m < jpsiHigh && m > jpsiLow; };
56 auto jpsi =
57 filteredEvents.Filter(jpsiCut, {"m"})
58 .Histo1D({"jpsi", "Subset of CMS Run 2010B: J/#psi window;#mu#mu mass [GeV];Events", 128, jpsiLow, jpsiHigh},
59 "m");
60
61 // Finally we draw the two histograms side by side.
62 auto dualCanvas = new TCanvas("DualCanvas", "DualCanvas", 800, 512);
63 dualCanvas->Divide(2, 1);
64 auto leftPad = dualCanvas->cd(1);
65 leftPad->SetLogx();
66 leftPad->SetLogy();
67 fullSpectrum->DrawClone("Hist");
68 dualCanvas->cd(2);
69 jpsi->DrawClone("HistP");
70
71 return 0;
72}
#define c(i)
Definition: RSha256.hxx:101
The Canvas class.
Definition: TCanvas.h:31
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:454
auto * m
Definition: textangle.C:8