Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf103_interprfuncs.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ## Basic functionality: interpreted functions and p.d.f.s
5 ##
6 ## \macro_code
7 ##
8 ## \date February 2018
9 ## \author Clemens Lange, Wouter Verkerke (C++ version)
10 
11 import ROOT
12 
13 # Generic interpreted p.d.f.
14 # ------------------------------
15 
16 # Declare observable x
17 x = ROOT.RooRealVar("x", "x", -20, 20)
18 
19 # Construct generic pdf from interpreted expression
20 # ------------------------------------------------------
21 
22 # ROOT.To construct a proper p.d.f, the formula expression is explicitly normalized internally by dividing
23 # it by a numeric integral of the expresssion over x in the range [-20,20]
24 #
25 alpha = ROOT.RooRealVar("alpha", "alpha", 5, 0.1, 10)
26 genpdf = ROOT.RooGenericPdf(
27  "genpdf",
28  "genpdf",
29  "(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))",
30  ROOT.RooArgList(
31  x,
32  alpha))
33 
34 # Sample, fit and plot generic pdf
35 # ---------------------------------------------------------------
36 
37 # Generate a toy dataset from the interpreted p.d.f
38 data = genpdf.generate(ROOT.RooArgSet(x), 10000)
39 
40 # Fit the interpreted p.d.f to the generated data
41 genpdf.fitTo(data)
42 
43 # Make a plot of the data and the p.d.f overlaid
44 xframe = x.frame(ROOT.RooFit.Title("Interpreted expression pdf"))
45 data.plotOn(xframe)
46 genpdf.plotOn(xframe)
47 
48 # Standard p.d.f. adjust with interpreted helper function
49 # ------------------------------------------------------------------------------------------------------------
50 # Make a gauss(x,sqrt(mean2),sigma) from a standard ROOT.RooGaussian #
51 #
52 # Construct standard pdf with formula replacing parameter
53 # ------------------------------------------------------------------------------------------------------------
54 
55 # Construct parameter mean2 and sigma
56 mean2 = ROOT.RooRealVar("mean2", "mean^2", 10, 0, 200)
57 sigma = ROOT.RooRealVar("sigma", "sigma", 3, 0.1, 10)
58 
59 # Construct interpreted function mean = sqrt(mean^2)
60 mean = ROOT.RooFormulaVar(
61  "mean", "mean", "sqrt(mean2)", ROOT.RooArgList(mean2))
62 
63 # Construct a gaussian g2(x,sqrt(mean2),sigma)
64 g2 = ROOT.RooGaussian("g2", "h2", x, mean, sigma)
65 
66 # Generate toy data
67 # ---------------------------------
68 
69 # Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian
70 # dataset with mean 10 and width 3
71 g1 = ROOT.RooGaussian("g1", "g1", x, ROOT.RooFit.RooConst(
72  10), ROOT.RooFit.RooConst(3))
73 data2 = g1.generate(ROOT.RooArgSet(x), 1000)
74 
75 # Fit and plot tailored standard pdf
76 # -------------------------------------------------------------------
77 
78 # Fit g2 to data from g1
79 r = g2.fitTo(data2, ROOT.RooFit.Save()) # ROOT.RooFitResult
80 r.Print()
81 
82 # Plot data on frame and overlay projection of g2
83 xframe2 = x.frame(ROOT.RooFit.Title("Tailored Gaussian pdf"))
84 data2.plotOn(xframe2)
85 g2.plotOn(xframe2)
86 
87 # Draw all frames on a canvas
88 c = ROOT.TCanvas("rf103_interprfuncs", "rf103_interprfuncs", 800, 400)
89 c.Divide(2)
90 c.cd(1)
91 ROOT.gPad.SetLeftMargin(0.15)
92 xframe.GetYaxis().SetTitleOffset(1.4)
93 xframe.Draw()
94 c.cd(2)
95 ROOT.gPad.SetLeftMargin(0.15)
96 xframe2.GetYaxis().SetTitleOffset(1.4)
97 xframe2.Draw()
98 
99 c.SaveAs("rf103_interprfuncs.png")