Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
20 int df014_CSVDataSource()
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 fileNameUrl = "http://root.cern.ch/files/tutorials/df014_CsvDataSource_MuRun2010B.csv";
25  auto fileName = "df014_CsvDataSource_MuRun2010B_cpp.csv";
26  if(gSystem->AccessPathName(fileName))
27  TFile::Cp(fileNameUrl, fileName);
28  auto tdf = ROOT::RDF::MakeCsvDataFrame(fileName);
29 
30  // Now we will apply a first filter based on two columns of the CSV,
31  // and we will define a new column that will contain the invariant mass.
32  // Note how the new invariant mass column is defined from several other
33  // columns that already existed in the CSV file.
34  auto filteredEvents =
35  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.
39  auto invMass =
40  filteredEvents.Histo1D({"invMass", "CMS Opendata: #mu#mu mass;#mu#mu mass [GeV];Events", 512, 2, 110}, "m");
41 
42  auto c = new TCanvas();
43  c->SetLogx();
44  c->SetLogy();
45  invMass->DrawClone();
46 
47  // We will now produce a plot also for the J/Psi particle. We will plot
48  // on the same canvas the full spectrum and the zoom in the J/psi particle.
49  // First we will create the full spectrum histogram from the invariant mass
50  // column, using a different histogram model than before.
51  auto fullSpectrum =
52  filteredEvents.Histo1D({"Spectrum", "Subset of CMS Run 2010B;#mu#mu mass [GeV];Events", 1024, 2, 110}, "m");
53 
54  // Next we will create the histogram for the J/psi particle, applying first
55  // the corresponding cut.
56  double jpsiLow = 2.95;
57  double jpsiHigh = 3.25;
58  auto jpsiCut = [jpsiLow, jpsiHigh](double m) { return m < jpsiHigh && m > jpsiLow; };
59  auto jpsi =
60  filteredEvents.Filter(jpsiCut, {"m"})
61  .Histo1D({"jpsi", "Subset of CMS Run 2010B: J/#psi window;#mu#mu mass [GeV];Events", 128, jpsiLow, jpsiHigh},
62  "m");
63 
64  // Finally we draw the two histograms side by side.
65  auto dualCanvas = new TCanvas("DualCanvas", "DualCanvas", 800, 512);
66  dualCanvas->Divide(2, 1);
67  auto leftPad = dualCanvas->cd(1);
68  leftPad->SetLogx();
69  leftPad->SetLogy();
70  fullSpectrum->DrawClone("Hist");
71  dualCanvas->cd(2);
72  jpsi->DrawClone("HistP");
73 
74  return 0;
75 }