Logo ROOT   6.16/01
Reference Guide
rf203_ranges.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #203
5## Fitting and plotting in sub ranges
6##
7## \macro_code
8##
9## \date February 2018
10## \author Clemens Lange
11## \author Wouter Verkerke (C version)
12
13from __future__ import print_function
14import ROOT
15
16# Set up model
17# ---------------------
18
19# Construct observables x
20x = ROOT.RooRealVar("x", "x", -10, 10)
21
22# Construct gaussx(x,mx,1)
23mx = ROOT.RooRealVar("mx", "mx", 0, -10, 10)
24gx = ROOT.RooGaussian("gx", "gx", x, mx, ROOT.RooFit.RooConst(1))
25
26# px = 1 (flat in x)
27px = ROOT.RooPolynomial("px", "px", x)
28
29# model = f*gx + (1-f)px
30f = ROOT.RooRealVar("f", "f", 0., 1.)
31model = ROOT.RooAddPdf(
32 "model", "model", ROOT.RooArgList(gx, px), ROOT.RooArgList(f))
33
34# Generated 10000 events in (x,y) from p.d.f. model
35modelData = model.generate(ROOT.RooArgSet(x), 10000)
36
37# Fit full range
38# ---------------------------
39
40# Fit p.d.f to all data
41r_full = model.fitTo(modelData, ROOT.RooFit.Save(ROOT.kTRUE))
42
43# Fit partial range
44# ----------------------------------
45
46# Define "signal" range in x as [-3,3]
47x.setRange("signal", -3, 3)
48
49# Fit p.d.f only to data in "signal" range
50r_sig = model.fitTo(modelData, ROOT.RooFit.Save(
51 ROOT.kTRUE), ROOT.RooFit.Range("signal"))
52
53# Plot/print results
54# ---------------------------------------
55
56# Make plot frame in x and add data and fitted model
57frame = x.frame(ROOT.RooFit.Title("Fitting a sub range"))
58modelData.plotOn(frame)
59model.plotOn(frame, ROOT.RooFit.Range("Full"), ROOT.RooFit.LineStyle(
60 ROOT.kDashed), ROOT.RooFit.LineColor(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")