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