ROOT
v6-32
Reference Guide
Loading...
Searching...
No Matches
rf602_chi2fit.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
##
5
## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #602
6
##
7
## Setting up a chi^2 fit to a binned dataset
8
##
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \authors Clemens Lange, Wouter Verkerke (C version)
14
15
from
__future__
import
print_function
16
import
ROOT
17
18
19
# Set up model
20
# ---------------------
21
22
# Declare observable x
23
x =
ROOT.RooRealVar
(
"x"
,
"x"
, 0, 10)
24
25
# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
26
# their parameters
27
mean =
ROOT.RooRealVar
(
"mean"
,
"mean of gaussians"
, 5)
28
sigma1 =
ROOT.RooRealVar
(
"sigma1"
,
"width of gaussians"
, 0.5)
29
sigma2 =
ROOT.RooRealVar
(
"sigma2"
,
"width of gaussians"
, 1)
30
31
sig1 =
ROOT.RooGaussian
(
"sig1"
,
"Signal component 1"
, x, mean, sigma1)
32
sig2 =
ROOT.RooGaussian
(
"sig2"
,
"Signal component 2"
, x, mean, sigma2)
33
34
# Build Chebychev polynomial p.d.f.
35
a0 =
ROOT.RooRealVar
(
"a0"
,
"a0"
, 0.5, 0.0, 1.0)
36
a1 =
ROOT.RooRealVar
(
"a1"
,
"a1"
, 0.2, 0.0, 1.0)
37
bkg =
ROOT.RooChebychev
(
"bkg"
,
"Background"
, x, [a0, a1])
38
39
# Sum the signal components into a composite signal p.d.f.
40
sig1frac =
ROOT.RooRealVar
(
"sig1frac"
,
"fraction of component 1 in signal"
, 0.8, 0.0, 1.0)
41
sig =
ROOT.RooAddPdf
(
"sig"
,
"Signal"
, [sig1, sig2], [sig1frac])
42
43
# Sum the composite signal and background
44
bkgfrac =
ROOT.RooRealVar
(
"bkgfrac"
,
"fraction of background"
, 0.5, 0.0, 1.0)
45
model =
ROOT.RooAddPdf
(
"model"
,
"g1+g2+a"
, [bkg, sig], [bkgfrac])
46
47
# Create biuned dataset
48
# -----------------------------------------
49
50
d =
model.generate
({x}, 10000)
51
dh =
d.binnedClone
()
52
53
# Construct a chi^2 of the data and the model.
54
# When a p.d.f. is used in a chi^2 fit, probability density scaled
55
# by the number of events in the dataset to obtain the fit function
56
# If model is an extended p.d.f, expected number events is used
57
# instead of the observed number of events.
58
ll =
ROOT.RooLinkedList
()
59
model.chi2FitTo
(dh, ll)
60
61
# NB: It is also possible to fit a ROOT.RooAbsReal function to a ROOT.RooDataHist
62
# using chi2FitTo().
63
64
# Note that entries with zero bins are _not_ allowed
65
# for a proper chi^2 calculation and will give error
66
# messages
67
dsmall =
d.reduce
(
ROOT.RooFit.EventRange
(1, 100))
68
dhsmall =
dsmall.binnedClone
()
69
chi2_lowstat =
model.createChi2
(dhsmall)
70
print(
chi2_lowstat.getVal
())
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
rf602_chi2fit.py
ROOT v6-32 - Reference Guide Generated on Thu Mar 13 2025 14:14:36 (GVA Time) using Doxygen 1.10.0