ROOT
master
Reference Guide
Loading...
Searching...
No Matches
rf705_linearmorph.py
2
## \ingroup tutorial_roofit_main
3
## \notebook
4
## 'SPECIAL PDFS' RooFit tutorial macro #705
5
##
6
## Linear interpolation between p.d.f shapes using the 'Alex Read' algorithm
7
##
8
## \macro_image
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \authors Clemens Lange, Wouter Verkerke (C version)
14
15
16
import
ROOT
17
18
19
# Create end point pdf shapes
20
# ------------------------------------------------------
21
22
# Observable
23
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -20, 20)
24
25
# Lower end point shape: a Gaussian
26
g1mean =
ROOT.RooRealVar
(
"g1mean"
,
"g1mean"
, -10)
27
g1 =
ROOT.RooGaussian
(
"g1"
,
"g1"
, x, g1mean, 2.0)
28
29
# Upper end point shape: a Polynomial
30
g2 =
ROOT.RooPolynomial
(
"g2"
,
"g2"
, x, [-0.03, -0.001])
31
32
# Create interpolating pdf
33
# -----------------------------------------------
34
35
# Create interpolation variable
36
alpha =
ROOT.RooRealVar
(
"alpha"
,
"alpha"
, 0, 1.0)
37
38
# Specify sampling density on observable and interpolation variable
39
x.setBins
(1000,
"cache"
)
40
alpha.setBins
(50,
"cache"
)
41
42
# Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
43
# and g2(x) at a=a_max
44
lmorph =
ROOT.RooIntegralMorph
(
"lmorph"
,
"lmorph"
, g1, g2, x, alpha)
45
46
# Plot interpolating pdf aat various alphas a l p h a
47
# -----------------------------------------------------------------------------
48
49
# Show end points as blue curves
50
frame1 =
x.frame
()
51
g1.plotOn
(frame1)
52
g2.plotOn
(frame1)
53
54
# Show interpolated shapes in red
55
alpha.setVal
(0.125)
56
lmorph.plotOn
(frame1, LineColor=
"r"
)
57
alpha.setVal
(0.25)
58
lmorph.plotOn
(frame1, LineColor=
"r"
)
59
alpha.setVal
(0.375)
60
lmorph.plotOn
(frame1, LineColor=
"r"
)
61
alpha.setVal
(0.50)
62
lmorph.plotOn
(frame1, LineColor=
"r"
)
63
alpha.setVal
(0.625)
64
lmorph.plotOn
(frame1, LineColor=
"r"
)
65
alpha.setVal
(0.75)
66
lmorph.plotOn
(frame1, LineColor=
"r"
)
67
alpha.setVal
(0.875)
68
lmorph.plotOn
(frame1, LineColor=
"r"
)
69
alpha.setVal
(0.95)
70
lmorph.plotOn
(frame1, LineColor=
"r"
)
71
72
# Show 2D distribution of pdf(x,alpha)
73
# -----------------------------------------------------------------------
74
75
# Create 2D histogram
76
hh =
lmorph.createHistogram
(
"hh"
, x, Binning=40, YVar=dict(var=alpha, Binning=40))
77
hh.SetLineColor
(
"kBlue"
)
78
79
# Fit pdf to dataset with alpha=0.8
80
# -----------------------------------------------------------------
81
82
# Generate a toy dataset alpha = 0.8
83
alpha.setVal
(0.8)
84
data =
lmorph.generate
({x}, 1000)
85
86
# Fit pdf to toy data
87
lmorph.setCacheAlpha
(
True
)
88
lmorph.fitTo
(data, Verbose=
True
, PrintLevel=-1)
89
90
# Plot fitted pdf and data overlaid
91
frame2 =
x.frame
(Bins=100)
92
data.plotOn
(frame2)
93
lmorph.plotOn
(frame2)
94
95
# Scan -log(L) vs alpha
96
# -----------------------------------------
97
98
# Show scan -log(L) of dataset w.r.t alpha
99
frame3 =
alpha.frame
(Bins=100, Range=(0.1, 0.9))
100
101
# Make 2D pdf of histogram
102
nll =
lmorph.createNLL
(data)
103
nll.plotOn
(frame3, ShiftToZero=
True
)
104
105
lmorph.setCacheAlpha
(
False
)
106
107
c =
ROOT.TCanvas
(
"rf705_linearmorph"
,
"rf705_linearmorph"
, 800, 800)
108
c.Divide
(2, 2)
109
c.cd
(1)
110
ROOT.gPad.SetLeftMargin
(0.15)
111
frame1.GetYaxis
().SetTitleOffset(1.6)
112
frame1.Draw
()
113
c.cd
(2)
114
ROOT.gPad.SetLeftMargin
(0.20)
115
hh.GetZaxis
().SetTitleOffset(2.5)
116
hh.Draw
(
"surf"
)
117
c.cd
(3)
118
ROOT.gPad.SetLeftMargin
(0.15)
119
frame3.GetYaxis
().SetTitleOffset(1.4)
120
frame3.Draw
()
121
c.cd
(4)
122
ROOT.gPad.SetLeftMargin
(0.15)
123
frame2.GetYaxis
().SetTitleOffset(1.4)
124
frame2.Draw
()
125
126
c.SaveAs
(
"rf705_linearmorph.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
roofit
rf705_linearmorph.py
ROOT master - Reference Guide Generated on Mon Feb 17 2025 15:03:49 (GVA Time) using Doxygen 1.10.0