Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf101_basics.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## This tutorial illustrates the basic features of RooFit.
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Set up model
14# ---------------------
15# Declare variables x,mean,sigma with associated name, title, initial
16# value and allowed range
17x = ROOT.RooRealVar("x", "x", -10, 10)
18mean = ROOT.RooRealVar("mean", "mean of gaussian", 1, -10, 10)
19sigma = ROOT.RooRealVar("sigma", "width of gaussian", 1, 0.1, 10)
20
21# Build gaussian pdf in terms of x,mean and sigma
22gauss = ROOT.RooGaussian("gauss", "gaussian PDF", x, mean, sigma)
23
24# Construct plot frame in 'x'
25xframe = x.frame(Title="Gaussian pdf") # RooPlot
26
27# Plot model and change parameter values
28# ---------------------------------------------------------------------------
29# Plot gauss in frame (i.e. in x)
30gauss.plotOn(xframe)
31
32# Change the value of sigma to 3
33sigma.setVal(3)
34
35# Plot gauss in frame (i.e. in x) and draw frame on canvas
36gauss.plotOn(xframe, LineColor="r")
37
38# Generate events
39# -----------------------------
40# Generate a dataset of 1000 events in x from gauss
41data = gauss.generate({x}, 10000) # ROOT.RooDataSet
42
43# Make a second plot frame in x and draw both the
44# data and the pdf in the frame
45xframe2 = x.frame(Title="Gaussian pdf with data") # RooPlot
46data.plotOn(xframe2)
47gauss.plotOn(xframe2)
48
49# Fit model to data
50# -----------------------------
51# Fit pdf to data
52gauss.fitTo(data)
53
54# Print values of mean and sigma (that now reflect fitted values and
55# errors)
56mean.Print()
57sigma.Print()
58
59# Draw all frames on a canvas
60c = ROOT.TCanvas("rf101_basics", "rf101_basics", 800, 400)
61c.Divide(2)
62c.cd(1)
63ROOT.gPad.SetLeftMargin(0.15)
64xframe.GetYaxis().SetTitleOffset(1.6)
65xframe.Draw()
66c.cd(2)
67ROOT.gPad.SetLeftMargin(0.15)
68xframe2.GetYaxis().SetTitleOffset(1.6)
69xframe2.Draw()
70
71c.SaveAs("rf101_basics.png")