Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf109_chi2residpull.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'BASIC FUNCTIONALITY' RooFit tutorial macro #109
5## Calculating chi^2 from histograms and curves in ROOT.RooPlots,
6## making histogram of residual and pull distributions
7##
8## \macro_image
9## \macro_code
10## \macro_output
11##
12## \date February 2018
13## \authors Clemens Lange, Wouter Verkerke (C version)
14
15from __future__ import print_function
16import ROOT
17
18# Set up model
19# ---------------------
20
21# Create observables
22x = ROOT.RooRealVar("x", "x", -10, 10)
23
24# Create Gaussian
25sigma = ROOT.RooRealVar("sigma", "sigma", 3, 0.1, 10)
26mean = ROOT.RooRealVar("mean", "mean", 0, -10, 10)
27gauss = ROOT.RooGaussian("gauss", "gauss", x, mean, sigma)
28
29# Generate a sample of 1000 events with sigma=3
30data = gauss.generate({x}, 10000)
31
32# Change sigma to 3.15
33sigma.setVal(3.15)
34
35# Plot data and slightly distorted model
36# ---------------------------------------------------------------------------
37
38# Overlay projection of gauss with sigma=3.15 on data with sigma=3.0
39frame1 = x.frame(Title="Data with distorted Gaussian pdf", Bins=40)
40data.plotOn(frame1, DataError="SumW2")
41gauss.plotOn(frame1)
42
43# Calculate chi^2
44# ------------------------------
45
46# Show the chi^2 of the curve w.r.t. the histogram
47# If multiple curves or datasets live in the frame you can specify
48# the name of the relevant curve and/or dataset in chiSquare()
49print("chi^2 = ", frame1.chiSquare())
50
51# Show residual and pull dists
52# -------------------------------------------------------
53
54# Construct a histogram with the residuals of the data w.r.t. the curve
55hresid = frame1.residHist()
56
57# Construct a histogram with the pulls of the data w.r.t the curve
58hpull = frame1.pullHist()
59
60# Create a frame to draw the residual distribution and add the
61# distribution to the frame
62frame2 = x.frame(Title="Residual Distribution")
63frame2.addPlotable(hresid, "P")
64
65# Create a frame to draw the pull distribution and add the distribution to
66# the frame
67frame3 = x.frame(Title="Pull Distribution")
68frame3.addPlotable(hpull, "P")
69
70c = ROOT.TCanvas("rf109_chi2residpull", "rf109_chi2residpull", 900, 300)
71c.Divide(3)
72c.cd(1)
73ROOT.gPad.SetLeftMargin(0.15)
74frame1.GetYaxis().SetTitleOffset(1.6)
75frame1.Draw()
76c.cd(2)
77ROOT.gPad.SetLeftMargin(0.15)
78frame2.GetYaxis().SetTitleOffset(1.6)
79frame2.Draw()
80c.cd(3)
81ROOT.gPad.SetLeftMargin(0.15)
82frame3.GetYaxis().SetTitleOffset(1.6)
83frame3.Draw()
84
85c.SaveAs("rf109_chi2residpull.png")