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