Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
df002_dataModel.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_dataframe
3 ## \notebook -draw
4 ## This tutorial shows the possibility to use data models which are more
5 ## complex than flat ntuples with RDataFrame
6 ##
7 ## \macro_code
8 ## \macro_image
9 ##
10 ## \date May 2017
11 ## \author Danilo Piparo
12 
13 import ROOT
14 
15 # A simple helper function to fill a test tree: this makes the example stand-alone.
16 fill_tree_code = '''
17 using FourVector = ROOT::Math::XYZTVector;
18 using FourVectorVec = std::vector<FourVector>;
19 using CylFourVector = ROOT::Math::RhoEtaPhiVector;
20 
21 // A simple helper function to fill a test tree: this makes the example
22 // stand-alone.
23 void fill_tree(const char *filename, const char *treeName)
24 {
25  const double M = 0.13957; // set pi+ mass
26  TRandom3 R(1);
27 
28  auto genTracks = [&](){
29  FourVectorVec tracks;
30  const auto nPart = R.Poisson(15);
31  tracks.reserve(nPart);
32  for (int j = 0; j < nPart; ++j) {
33  const auto px = R.Gaus(0, 10);
34  const auto py = R.Gaus(0, 10);
35  const auto pt = sqrt(px * px + py * py);
36  const auto eta = R.Uniform(-3, 3);
37  const auto phi = R.Uniform(0.0, 2 * TMath::Pi());
38  CylFourVector vcyl(pt, eta, phi);
39  // set energy
40  auto E = sqrt(vcyl.R() * vcyl.R() + M * M);
41  // fill track vector
42  tracks.emplace_back(vcyl.X(), vcyl.Y(), vcyl.Z(), E);
43  }
44  return tracks;
45  };
46 
47  ROOT::RDataFrame d(64);
48  d.Define("tracks", genTracks).Snapshot<FourVectorVec>(treeName, filename, {"tracks"});
49 }
50 '''
51 
52 # We prepare an input tree to run on
53 fileName = "df002_dataModel_py.root"
54 treeName = "myTree"
55 ROOT.gInterpreter.Declare(fill_tree_code)
56 ROOT.fill_tree(fileName, treeName)
57 
58 # We read the tree from the file and create a RDataFrame, a class that
59 # allows us to interact with the data contained in the tree.
60 RDF = ROOT.ROOT.RDataFrame
61 d = RDF(treeName, fileName)
62 
63 # Operating on branches which are collection of objects
64 # Here we deal with the simplest of the cuts: we decide to accept the event
65 # only if the number of tracks is greater than 5.
66 n_cut = 'tracks.size() > 8'
67 nentries = d.Filter(n_cut).Count();
68 
69 print("%s passed all filters" %nentries.GetValue())
70 
71 # Another possibility consists in creating a new column containing the
72 # quantity we are interested in.
73 # In this example, we will cut on the number of tracks and plot their
74 # transverse momentum.
75 
76 getPt_code ='''
77 using namespace ROOT::VecOps;
78 RVec<double> getPt(const RVec<FourVector> &tracks)
79 {
80  auto pt = [](const FourVector &v) { return v.pt(); };
81  return Map(tracks, pt);
82 }
83 '''
84 ROOT.gInterpreter.Declare(getPt_code)
85 
86 getPtWeights_code ='''
87 using namespace ROOT::VecOps;
88 RVec<double> getPtWeights(const RVec<FourVector> &tracks)
89 {
90  auto ptWeight = [](const FourVector &v) { return 1. / v.Pt(); };
91  return Map(tracks, ptWeight);
92 };
93 '''
94 ROOT.gInterpreter.Declare(getPtWeights_code)
95 
96 augmented_d = d.Define('tracks_n', '(int)tracks.size()') \
97  .Filter('tracks_n > 2') \
98  .Define('tracks_pts', 'getPt( tracks )') \
99  .Define("tracks_pts_weights", 'getPtWeights( tracks )' )
100 
101 # The histogram is initialised with a tuple containing the parameters of the
102 # histogram
103 trN = augmented_d.Histo1D(("", "", 40, -.5, 39.5), "tracks_n")
104 trPts = augmented_d.Histo1D("tracks_pts")
105 trWPts = augmented_d.Histo1D("tracks_pts", "tracks_pts_weights")
106 
107 c1 = ROOT.TCanvas()
108 trN.Draw()
109 
110 c2 = ROOT.TCanvas()
111 trPts.Draw()
112 
113 c3 = ROOT.TCanvas()
114 trWPts.Draw()