Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf203_ranges.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Addition and convolution: fitting and plotting in sub ranges
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11from __future__ import print_function
12import ROOT
13
14# Set up model
15# ---------------------
16
17# Construct observables x
18x = ROOT.RooRealVar("x", "x", -10, 10)
19
20# Construct gaussx(x,mx,1)
21mx = ROOT.RooRealVar("mx", "mx", 0, -10, 10)
22gx = ROOT.RooGaussian("gx", "gx", x, mx, ROOT.RooFit.RooConst(1))
23
24# px = 1 (flat in x)
25px = ROOT.RooPolynomial("px", "px", x)
26
27# model = f*gx + (1-f)px
28f = ROOT.RooRealVar("f", "f", 0.0, 1.0)
29model = ROOT.RooAddPdf("model", "model", [gx, px], [f])
30
31# Generated 10000 events in (x,y) from pdf model
32modelData = model.generate({x}, 10000)
33
34# Fit full range
35# ---------------------------
36
37# Fit pdf to all data
38r_full = model.fitTo(modelData, Save=True)
39
40# Fit partial range
41# ----------------------------------
42
43# Define "signal" range in x as [-3,3]
44x.setRange("signal", -3, 3)
45
46# Fit pdf only to data in "signal" range
47r_sig = model.fitTo(modelData, Save=True, Range="signal")
48
49# Plot/print results
50# ---------------------------------------
51
52# Make plot frame in x and add data and fitted model
53frame = x.frame(Title="Fitting a sub range")
54modelData.plotOn(frame)
55model.plotOn(frame, Range="Full", LineColor="r", LineStyle="--") # Add shape in full ranged dashed
56model.plotOn(frame) # By default only fitted range is shown
57
58# Print fit results
59print("result of fit on all data ")
60r_full.Print()
61print("result of fit in in signal region (note increased error on signal fraction)")
62r_sig.Print()
63
64# Draw frame on canvas
65c = ROOT.TCanvas("rf203_ranges", "rf203_ranges", 600, 600)
66ROOT.gPad.SetLeftMargin(0.15)
67frame.GetYaxis().SetTitleOffset(1.4)
68frame.Draw()
69
70c.SaveAs("rf203_ranges.png")