Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
df015_LazyDataSource.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// \notebook -draw
4 /// This tutorial illustrates how to take advantage of a *lazy data source*
5 /// creating a data frame from columns of one or multiple parent dataframe(s),
6 /// delaying the creation of the columns to the actual usage of the daughter
7 /// data frame.
8 /// Dataset Reference:
9 /// McCauley, T. (2014). Dimuon event information derived from the Run2010B
10 /// public Mu dataset. CERN Open Data Portal.
11 /// DOI: [10.7483/OPENDATA.CMS.CB8H.MFFA](http://opendata.cern.ch/record/700).
12 /// From the ROOT website: https://root.cern.ch/files/tutorials/tdf014_CsvDataSource_MuRun2010B.csv
13 ///
14 /// \macro_code
15 /// \macro_image
16 ///
17 /// \date February 2018
18 /// \author Danilo Piparo
19 
20 int df015_LazyDataSource()
21 {
22  using namespace ROOT::RDF;
23 
24  // Let's first create a RDF that will read from the CSV file.
25  // See the tutorial relative to CSV data sources for more details!
26  auto fileNameUrl = "http://root.cern.ch/files/tutorials/df014_CsvDataSource_MuRun2010B.csv";
27  auto fileName = "df015_CsvDataSource_MuRun2010B.csv";
28  if(gSystem->AccessPathName(fileName))
29  TFile::Cp(fileNameUrl, fileName);
30 
31  auto csv_rdf = MakeCsvDataFrame(fileName);
32 
33  // Now we take out four columns: px and py of the first muon in the muon pair
34  std::string px1Name = "px1";
35  auto px1 = csv_rdf.Take<double>(px1Name);
36  std::string py1Name = "py1";
37  auto py1 = csv_rdf.Take<double>(py1Name);
38 
39  // Now we create a new dataframe built on top of the columns above. Note that up to now, no event loop
40  // has been carried out!
41  auto tdf = MakeLazyDataFrame(std::make_pair(px1Name, px1), std::make_pair(py1Name, py1));
42 
43  // We build a histogram of the transverse momentum the muons.
44  auto ptFormula = [](double px, double py) { return sqrt(px * px + py * py); };
45  auto pt_h = tdf.Define("pt", ptFormula, {"px1", "py1"})
46  .Histo1D<double>({"pt", "Muon p_{T};p_{T} [GeV/c];", 128, 0, 128}, "pt");
47 
48  auto can = new TCanvas();
49  can->SetLogy();
50  pt_h->DrawCopy();
51 
52  return 0;
53 }