ROOT
Version v6.34
master
v6.32
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
rf314_paramfitrange.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## Multidimensional models: working with parameterized ranges in a fit.
5
## This an example of a fit with an acceptance that changes per-event
6
##
7
## `pdf = exp(-t/tau)` with `t[tmin,5]`
8
##
9
## where `t` and `tmin` are both observables in the dataset
10
##
11
## \macro_image
12
## \macro_code
13
## \macro_output
14
##
15
## \date February 2018
16
## \authors Clemens Lange, Wouter Verkerke (C++ version)
17
18
import
ROOT
19
20
21
# Define observables and decay pdf
22
# ---------------------------------------------------------------
23
24
# Declare observables
25
t =
ROOT.RooRealVar
(
"t"
,
"t"
, 0, 5)
26
tmin =
ROOT.RooRealVar
(
"tmin"
,
"tmin"
, 0, 0, 5)
27
28
# Make parameterized range in t : [tmin,5]
29
t.setRange
(tmin,
ROOT.RooFit.RooConst
(
t.getMax
()))
30
31
# Make pdf
32
tau =
ROOT.RooRealVar
(
"tau"
,
"tau"
, -1.54, -10, -0.1)
33
model =
ROOT.RooExponential
(
"model"
,
"model"
, t, tau)
34
35
# Create input data
36
# ------------------------------------
37
38
# Generate complete dataset without acceptance cuts (for reference)
39
dall =
model.generate
({t}, 10000)
40
41
# Generate a (fake) prototype dataset for acceptance limit values
42
tmp =
ROOT.RooGaussian
(
"gmin"
,
"gmin"
, tmin, 0.0, 0.5).generate({tmin}, 5000)
43
44
# Generate dataset with t values that observe (t>tmin)
45
dacc =
model.generate
({t}, ProtoData=tmp)
46
47
# Fit pdf to data in acceptance region
48
# -----------------------------------------------------------------------
49
50
r =
model.fitTo
(dacc, Save=
True
, PrintLevel=-1)
51
52
# Plot fitted pdf on full and accepted data
53
# ---------------------------------------------------------------------------------
54
55
# Make plot frame, datasets and overlay model
56
frame =
t.frame
(Title=
"Fit to data with per-event acceptance"
)
57
dall.plotOn
(frame, MarkerColor=
"r"
, LineColor=
"r"
)
58
model.plotOn
(frame)
59
dacc.plotOn
(frame)
60
61
# Print fit results to demonstrate absence of bias
62
r.Print
(
"v"
)
63
64
c =
ROOT.TCanvas
(
"rf314_paramranges"
,
"rf314_paramranges"
, 600, 600)
65
ROOT.gPad.SetLeftMargin
(0.15)
66
frame.GetYaxis
().SetTitleOffset(1.6)
67
frame.Draw
()
68
69
c.SaveAs
(
"rf314_paramranges.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
rf314_paramfitrange.py
ROOT tags/6-34-04 - Reference Guide Generated on Fri Mar 21 2025 04:40:18 (GVA Time) using Doxygen 1.10.0