Logo ROOT   6.10/00
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Namespaces
tdf002_dataModel.py File Reference

Namespaces

 tdf002_dataModel
 

Detailed Description

This tutorial shows the possibility to use data models which are more complex than flat ntuples with TDataFrame.

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

Definition in file tdf002_dataModel.py.