Logo ROOT  
Reference Guide
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## \author 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., 1.)
29model = ROOT.RooAddPdf(
30 "model", "model", ROOT.RooArgList(gx, px), ROOT.RooArgList(f))
31
32# Generated 10000 events in (x,y) from p.d.f. model
33modelData = model.generate(ROOT.RooArgSet(x), 10000)
34
35# Fit full range
36# ---------------------------
37
38# Fit p.d.f to all data
39r_full = model.fitTo(modelData, ROOT.RooFit.Save(ROOT.kTRUE))
40
41# Fit partial range
42# ----------------------------------
43
44# Define "signal" range in x as [-3,3]
45x.setRange("signal", -3, 3)
46
47# Fit p.d.f only to data in "signal" range
48r_sig = model.fitTo(modelData, ROOT.RooFit.Save(
49 ROOT.kTRUE), ROOT.RooFit.Range("signal"))
50
51# Plot/print results
52# ---------------------------------------
53
54# Make plot frame in x and add data and fitted model
55frame = x.frame(ROOT.RooFit.Title("Fitting a sub range"))
56modelData.plotOn(frame)
57model.plotOn(
58 frame, ROOT.RooFit.Range("Full"), ROOT.RooFit.LineStyle(
59 ROOT.kDashed), ROOT.RooFit.LineColor(
60 ROOT.kRed)) # Add shape in full ranged dashed
61model.plotOn(frame) # By default only fitted range is shown
62
63# Print fit results
64print("result of fit on all data ")
65r_full.Print()
66print("result of fit in in signal region (note increased error on signal fraction)")
67r_sig.Print()
68
69# Draw frame on canvas
70c = ROOT.TCanvas("rf203_ranges", "rf203_ranges", 600, 600)
71ROOT.gPad.SetLeftMargin(0.15)
72frame.GetYaxis().SetTitleOffset(1.4)
73frame.Draw()
74
75c.SaveAs("rf203_ranges.png")