-
Notifications
You must be signed in to change notification settings - Fork 18
/
rf111_derivatives.py
88 lines (68 loc) · 2.73 KB
/
rf111_derivatives.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
#####################################
#
# 'BASIC FUNCTIONALITY' ROOT.RooFit tutorial macro #111
#
# Numerical 1st, and 3rd order derivatives w.r.t. observables and parameters
#
# pdf = gauss(x,m,s)
#
#
# 07/2008 - Wouter Verkerke
#
# /
import ROOT
def rf111_derivatives():
# S e t u p m o d e l
# ---------------------
# Declare variables x,mean, with associated name, title, value and allowed
# range
x = ROOT.RooRealVar("x", "x", -10, 10)
mean = ROOT.RooRealVar("mean", "mean of gaussian", 1, -10, 10)
sigma = ROOT.RooRealVar("sigma", "width of gaussian", 1, 0.1, 10)
# Build gaussian p.d.f in terms of x, and sigma
gauss = ROOT.RooGaussian("gauss", "gaussian PDF", x, mean, sigma)
# C r e a t e a n d p l o t d e r i v a t i v e s w . r . t . x
# ----------------------------------------------------------------------
# Derivative of normalized gauss(x) w.r.t. observable x
dgdx = gauss.derivative(x, 1)
# Second and third derivative of normalized gauss(x) w.r.t. observable x
d2gdx2 = gauss.derivative(x, 2)
d3gdx3 = gauss.derivative(x, 3)
# Construct plot frame in 'x'
xframe = x.frame(ROOT.RooFit.Title("d(Gauss)/dx"))
# Plot gauss in frame (i.e. in x)
gauss.plotOn(xframe)
# Plot derivatives in same frame
dgdx.plotOn(xframe, ROOT.RooFit.LineColor(ROOT.kMagenta))
d2gdx2.plotOn(xframe, ROOT.RooFit.LineColor(ROOT.kRed))
d3gdx3.plotOn(xframe, ROOT.RooFit.LineColor(ROOT.kOrange))
# C r e a t e a n d p l o t d e r i v a t i v e s w . r . t . s i g m a
# ------------------------------------------------------------------------------
# Derivative of normalized gauss(x) w.r.t. parameter sigma
dgds = gauss.derivative(sigma, 1)
# Second and third derivative of normalized gauss(x) w.r.t. parameter sigma
d2gds2 = gauss.derivative(sigma, 2)
d3gds3 = gauss.derivative(sigma, 3)
# Construct plot frame in 'sigma'
sframe = sigma.frame(ROOT.RooFit.Title(
"d(Gauss)/d(sigma)"), ROOT.RooFit.Range(0., 2.))
# Plot gauss in frame (i.e. in x)
gauss.plotOn(sframe)
# Plot derivatives in same frame
dgds.plotOn(sframe, ROOT.RooFit.LineColor(ROOT.kMagenta))
d2gds2.plotOn(sframe, ROOT.RooFit.LineColor(ROOT.kRed))
d3gds3.plotOn(sframe, ROOT.RooFit.LineColor(ROOT.kOrange))
# Draw all frames on a canvas
c = ROOT.TCanvas("rf111_derivatives", "rf111_derivatives", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.6)
xframe.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
sframe.GetYaxis().SetTitleOffset(1.6)
sframe.Draw()
c.SaveAs("rf111_derivatives.png")
if __name__ == "__main__":
rf111_derivatives()