Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf401_importttreethx.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'DATA AND CATEGORIES' RooFit tutorial macro #401
5##
6## Overview of advanced option for importing data from ROOT ROOT.TTree and ROOT.THx histograms
7## Basic import options are demonstrated in rf102_dataimport.py
8##
9## \macro_code
10## \macro_output
11##
12## \date February 2018
13## \authors Clemens Lange, Wouter Verkerke (C version)
14
15import ROOT
16from array import array
17
18
19def makeTH1(name, mean, sigma):
20 """Create ROOT TH1 filled with a Gaussian distribution."""
21
22 hh = ROOT.TH1D(name, name, 100, -10, 10)
23 for i in range(1000):
24 hh.Fill(ROOT.gRandom.Gaus(mean, sigma))
25
26 return hh
27
28
29def makeTTree():
30 """Create ROOT ROOT.TTree filled with a Gaussian distribution in x and a uniform distribution in y."""
31
32 tree = ROOT.TTree("tree", "tree")
33 px = array("d", [0])
34 py = array("d", [0])
35 pz = array("d", [0])
36 pi = array("i", [0])
37 tree.Branch("x", px, "x/D")
38 tree.Branch("y", py, "y/D")
39 tree.Branch("z", pz, "z/D")
40 tree.Branch("i", pi, "i/I")
41 for i in range(100):
42 px[0] = ROOT.gRandom.Gaus(0, 3)
43 py[0] = ROOT.gRandom.Uniform() * 30 - 15
44 pz[0] = ROOT.gRandom.Gaus(0, 5)
45 pi[0] = i % 3
46 tree.Fill()
47
48 return tree
49
50
51# Import multiple TH1 into a RooDataHist
52# ----------------------------------------------------------
53
54# Create thee ROOT ROOT.TH1 histograms
55hh_1 = makeTH1("hh1", 0, 3)
56hh_2 = makeTH1("hh2", -3, 1)
57hh_3 = makeTH1("hh3", +3, 4)
58
59# Declare observable x
60x = ROOT.RooRealVar("x", "x", -10, 10)
61
62# Create category observable c that serves as index for the ROOT histograms
63c = ROOT.RooCategory("c", "c")
64c.defineType("SampleA")
65c.defineType("SampleB")
66c.defineType("SampleC")
67
68# Create a binned dataset that imports contents of all ROOT.TH1 mapped by
69# index category c
70dh = ROOT.RooDataHist("dh", "dh", [x], Index=c, Import={"SampleA": hh_1, "SampleB": hh_2, "SampleC": hh_3})
71dh.Print()
72
73dh2 = ROOT.RooDataHist("dh", "dh", [x], Index=c, Import={"SampleA": hh_1, "SampleB": hh_2, "SampleC": hh_3})
74dh2.Print()
75
76# Importing a ROOT TTree into a RooDataSet with cuts
77# --------------------------------------------------------------------
78
79tree = makeTTree()
80
81# Define observables y,z
82y = ROOT.RooRealVar("y", "y", -10, 10)
83z = ROOT.RooRealVar("z", "z", -10, 10)
84
85# Import only observables (y,z)
86ds = ROOT.RooDataSet("ds", "ds", {x, y}, Import=tree)
87ds.Print()
88
89# Import observables (x,y,z) but only event for which (y+z<0) is ROOT.True
90# Import observables (x,y,z) but only event for which (y+z<0) is ROOT.True
91ds2 = ROOT.RooDataSet("ds2", "ds2", {x, y, z}, Import=tree, Cut="y+z<0")
92ds2.Print()
93
94# Importing integer ROOT TTree branches
95# ---------------------------------------------------------------
96
97# Import integer tree branch as ROOT.RooRealVar
98i = ROOT.RooRealVar("i", "i", 0, 5)
99ds3 = ROOT.RooDataSet("ds3", "ds3", {i, x}, Import=tree)
100ds3.Print()
101
102# Define category i
103icat = ROOT.RooCategory("i", "i", {"State0": 0, "State1": 1})
104
105# Import integer tree branch as ROOT.RooCategory (only events with i==0 and i==1
106# will be imported as those are the only defined states)
107ds4 = ROOT.RooDataSet("ds4", "ds4", {icat, x}, Import=tree)
108ds4.Print()
109
110# Import multiple RooDataSets into a RooDataSet
111# ----------------------------------------------------------------------------------------
112
113# Create three ROOT.RooDataSets in (y,z)
114dsA = ds2.reduce({x, y}, "z<-5")
115dsB = ds2.reduce({x, y}, "abs(z)<5")
116dsC = ds2.reduce({x, y}, "z>5")
117
118# Create a dataset that imports contents of all the above datasets mapped
119# by index category c
120dsABC = ROOT.RooDataSet("dsABC", "dsABC", {x, y}, Index=c, Import={"SampleA": dsA, "SampleB": dsB, "SampleC": dsC})
121
122dsABC.Print()