Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df017_vecOpsHEP.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_dataframe
3## \notebook -draw
4## Use RVecs to plot the transverse momentum of selected particles.
5##
6## This tutorial shows how VecOps can be used to slim down the programming
7## model typically adopted in HEP for analysis.
8## In this case we have a dataset containing the kinematic properties of
9## particles stored in individual arrays.
10## We want to plot the transverse momentum of these particles if the energy is
11## greater than 100 MeV.
12
13## \macro_code
14## \macro_image
15##
16## \date March 2018
17## \authors Danilo Piparo (CERN), Andre Vieira Silva
18
19import ROOT
20
21filename = ROOT.gROOT.GetTutorialDir().Data() + "/dataframe/df017_vecOpsHEP.root"
22treename = "myDataset"
23
24def WithPyROOT(filename):
25 from math import sqrt
26 f = ROOT.TFile(filename)
27 h = ROOT.TH1F("pt", "With PyROOT", 16, 0, 4)
28 for event in f.myDataset:
29 for E, px, py in zip(event.E, event.px, event.py):
30 if (E > 100):
31 h.Fill(sqrt(px*px + py*py))
32 h.DrawCopy()
33
34def WithRDataFrameVecOpsJit(treename, filename):
35 f = ROOT.RDataFrame(treename, filename)
36 h = f.Define("good_pt", "sqrt(px*px + py*py)[E>100]")\
37 .Histo1D(("pt", "With RDataFrame and RVec", 16, 0, 4), "good_pt")
38 h.DrawCopy()
39
40## We plot twice the same quantity, the key is to look into the implementation
41## of the functions above.
42c = ROOT.TCanvas()
43c.Divide(2,1)
44c.cd(1)
45WithPyROOT(filename)
46c.cd(2)
47WithRDataFrameVecOpsJit(treename, filename)
48c.SaveAs("df017_vecOpsHEP.png")
49
50print("Saved figure to df017_vecOpsHEP.png")
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...