Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf408_RDataFrameToRooFit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Fill RooDataSet/RooDataHist in RDataFrame.
5##
6## This tutorial shows how to fill RooFit data classes directly from RDataFrame.
7## Using two small helpers, we tell RDataFrame where the data has to go.
8##
9## \macro_code
10## \macro_output
11##
12## \date July 2021
13## \author Harshal Shende, Stephan Hageboeck (C++ version)
14
15import ROOT
16import math
17
18
19# Set up
20# ------------------------
21
22# We create an RDataFrame with two columns filled with 2 million random numbers.
23d = ROOT.RDataFrame(2000000)
24dd = d.Define("x", "gRandom->Uniform(-5., 5.)").Define("y", "gRandom->Gaus(1., 3.)")
25
26
27# We create RooFit variables that will represent the dataset.
28x = ROOT.RooRealVar("x", "x", -5.0, 5.0)
29y = ROOT.RooRealVar("y", "y", -50.0, 50.0)
30x.setBins(10)
31y.setBins(20)
32
33
34# Booking the creation of RooDataSet / RooDataHist in RDataFrame
35# ----------------------------------------------------------------
36
37# Method 1:
38# ---------
39# We directly book the RooDataSetHelper action.
40# We need to pass
41# - the RDataFrame column types as template parameters
42# - the constructor arguments for RooDataSet (they follow the same syntax as the usual RooDataSet constructors)
43# - the column names that RDataFrame should fill into the dataset
44#
45# NOTE: RDataFrame columns are matched to RooFit variables by position, *not by name*!
46rooDataSet = dd.Book(
47 ROOT.std.move(ROOT.RooDataSetHelper("dataset", "Title of dataset", ROOT.RooArgSet(x, y))), ("x", "y")
48)
49
50
51# Method 2:
52# ---------
53# We first declare the RooDataHistHelper
54rdhMaker = ROOT.RooDataHistHelper("dataset", "Title of dataset", ROOT.RooArgSet(x, y))
55
56# Then, we move it into an RDataFrame action:
57rooDataHist = dd.Book(ROOT.std.move(rdhMaker), ("x", "y"))
58
59
60# Run it and inspect the results
61# -------------------------------
62
63# Let's inspect the dataset / datahist.
64# Note that the first time we touch one of those objects, the RDataFrame event loop will run.
65
66
67def printData(data):
68 print("")
69 data.Print()
70 for i in range(min(data.numEntries(), 20)):
71 print(
72 "("
73 + ", ".join(["{0:8.3f}".format(var.getVal()) for var in data.get(i)])
74 + ", ) weight={0:10.3f}".format(data.weight())
75 )
76
77 print("mean(x) = {0:.3f}".format(data.mean(x)) + "\tsigma(x) = {0:.3f}".format(math.sqrt(data.moment(x, 2.0))))
78 print("mean(y) = {0:.3f}".format(data.mean(y)) + "\tsigma(y) = {0:.3f}\n".format(math.sqrt(data.moment(y, 2.0))))
79
80
81printData(rooDataSet)
82printData(rooDataHist)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...