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##
9## \macro_code
10## \macro_image
11##
12## \date March 2018
13## \authors Danilo Piparo (CERN), Andre Vieira Silva
14
15import ROOT
16
17filename = ROOT.gROOT.GetTutorialDir().Data() + "/dataframe/df017_vecOpsHEP.root"
18treename = "myDataset"
19
20def WithPyROOT(filename):
21 from math import sqrt
22 f = ROOT.TFile(filename)
23 h = ROOT.TH1F("pt", "With PyROOT", 16, 0, 4)
24 for event in f.myDataset:
25 for E, px, py in zip(event.E, event.px, event.py):
26 if (E > 100):
27 h.Fill(sqrt(px*px + py*py))
28 h.DrawCopy()
29
30def WithRDataFrameVecOpsJit(treename, filename):
31 f = ROOT.RDataFrame(treename, filename)
32 h = f.Define("good_pt", "sqrt(px*px + py*py)[E>100]")\
33 .Histo1D(("pt", "With RDataFrame and RVec", 16, 0, 4), "good_pt")
34 h.DrawCopy()
35
36## We plot twice the same quantity, the key is to look into the implementation
37## of the functions above
38c = ROOT.TCanvas()
39c.Divide(2,1)
40c.cd(1)
41WithPyROOT(filename)
42c.cd(2)
43WithRDataFrameVecOpsJit(treename, filename)
44c.SaveAs("df017_vecOpsHEP.png")
45
46print("Saved figure to df017_vecOpsHEP.png")
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTree,...