ROOT
Version master
v6.36
v6.34
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
rf803_mcstudy_addons2.py
Go to the documentation of this file.
1
## \ingroup tutorial_roofit_main
2
## \notebook
3
##
4
## 'VALIDATION AND MC STUDIES' RooFit tutorial macro #803
5
##
6
## RooMCStudy: Using the randomizer and profile likelihood add-on models
7
##
8
## \macro_image
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \author Clemens Lange
14
15
16
import
ROOT
17
18
19
# Create model
20
# -----------------------
21
22
# Simulation of signal and background of top quark decaying into
23
# 3 jets with background
24
25
# Observable
26
mjjj =
ROOT.RooRealVar
(
"mjjj"
,
"m(3jet) (GeV)"
, 100, 85.0, 350.0)
27
28
# Signal component (Gaussian)
29
mtop =
ROOT.RooRealVar
(
"mtop"
,
"m(top)"
, 162)
30
wtop =
ROOT.RooRealVar
(
"wtop"
,
"m(top) resolution"
, 15.2)
31
sig =
ROOT.RooGaussian
(
"sig"
,
"top signal"
, mjjj, mtop, wtop)
32
33
# Background component (Chebychev)
34
c0 =
ROOT.RooRealVar
(
"c0"
,
"Chebychev coefficient 0"
, -0.846, -1.0, 1.0)
35
c1 =
ROOT.RooRealVar
(
"c1"
,
"Chebychev coefficient 1"
, 0.112, -1.0, 1.0)
36
c2 =
ROOT.RooRealVar
(
"c2"
,
"Chebychev coefficient 2"
, 0.076, -1.0, 1.0)
37
bkg =
ROOT.RooChebychev
(
"bkg"
,
"combinatorial background"
, mjjj, [c0, c1, c2])
38
39
# Composite model
40
nsig =
ROOT.RooRealVar
(
"nsig"
,
"number of signal events"
, 53, 0, 1e3)
41
nbkg =
ROOT.RooRealVar
(
"nbkg"
,
"number of background events"
, 103, 0, 5e3)
42
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [sig, bkg], [nsig, nbkg])
43
44
# Create manager
45
# ---------------------------
46
47
# Configure manager to perform binned extended likelihood fits (Binned=True, Extended=True) on data generated
48
# with a Poisson fluctuation on Nobs (Extended=True)
49
mcs =
ROOT.RooMCStudy
(
50
model, {mjjj}, Binned=
True
, Silence=
True
, Extended=
True
, FitOptions={
"Extended"
:
True
,
"PrintEvalErrors"
: -1}
51
)
52
53
# Customize manager
54
# ---------------------------------
55
56
# Add module that randomizes the summed value of nsig+nbkg
57
# sampling from a uniform distribution between 0 and 1000
58
#
59
# In general one can randomize a single parameter, a
60
# sum of N parameters, either a uniform or a Gaussian
61
# distribution. Multiple randomization can be executed
62
# by a single randomizer module
63
64
randModule =
ROOT.RooRandomizeParamMCSModule
()
65
randModule.sampleSumUniform
({nsig, nbkg}, 50, 500)
66
mcs.addModule
(randModule)
67
68
# Add profile likelihood calculation of significance. Redo each
69
# fit while keeping parameter nsig fixed to zero. For each toy,
70
# the difference in -log(L) of both fits is stored, well
71
# a simple significance interpretation of the delta(-logL)
72
# Dnll = 0.5 sigma^2
73
74
sigModule =
ROOT.RooDLLSignificanceMCSModule
(nsig, 0)
75
mcs.addModule
(sigModule)
76
77
# Run manager, make plots
78
# ---------------------------------------------
79
80
# Run 1000 experiments. ROOT.This configuration will generate a fair number
81
# of (harmless) MINUIT warnings due to the instability of the Chebychev polynomial fit
82
# at low statistics.
83
mcs.generateAndFit
(500)
84
85
# Make some plots
86
binning =
ROOT.RooFit.AutoBinning
(40)
87
dll_vs_ngen =
mcs.fitParDataSet
().createHistogram(
"ngen,dll_nullhypo_nsig"
, binning, binning)
88
z_vs_ngen =
mcs.fitParDataSet
().createHistogram(
"ngen,significance_nullhypo_nsig"
, binning, binning)
89
errnsig_vs_ngen =
mcs.fitParDataSet
().createHistogram(
"ngen,nsigerr"
, binning, binning)
90
errnsig_vs_nsig =
mcs.fitParDataSet
().createHistogram(
"nsig,nsigerr"
, binning, binning)
91
92
# Draw plots on canvas
93
c =
ROOT.TCanvas
(
"rf803_mcstudy_addons2"
,
"rf802_mcstudy_addons2"
, 800, 800)
94
c.Divide
(2, 2)
95
c.cd
(1)
96
ROOT.gPad.SetLeftMargin
(0.15)
97
dll_vs_ngen.GetYaxis
().SetTitleOffset(1.6)
98
dll_vs_ngen.Draw
(
"box"
)
99
c.cd
(2)
100
ROOT.gPad.SetLeftMargin
(0.15)
101
z_vs_ngen.GetYaxis
().SetTitleOffset(1.6)
102
z_vs_ngen.Draw
(
"box"
)
103
c.cd
(3)
104
ROOT.gPad.SetLeftMargin
(0.15)
105
errnsig_vs_ngen.GetYaxis
().SetTitleOffset(1.6)
106
errnsig_vs_ngen.Draw
(
"box"
)
107
c.cd
(4)
108
ROOT.gPad.SetLeftMargin
(0.15)
109
errnsig_vs_nsig.GetYaxis
().SetTitleOffset(1.6)
110
errnsig_vs_nsig.Draw
(
"box"
)
111
112
c.SaveAs
(
"rf803_mcstudy_addons2.png"
)
113
114
# Make ROOT.RooMCStudy object available on command line after
115
# macro finishes
116
ROOT.gDirectory.Add
(mcs)
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
roofit
rf803_mcstudy_addons2.py
ROOT master - Reference Guide Generated on Sun Apr 27 2025 15:45:24 (GVA Time) using Doxygen 1.10.0