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
rs101_limitexample.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roostats
3
## \notebook
4
## Limits: number counting experiment with uncertainty on both the background rate and signal efficiency.
5
##
6
## The usage of a Confidence Interval Calculator to set a limit on the signal is illustrated
7
##
8
## \macro_image
9
## \macro_output
10
## \macro_code
11
##
12
## \date June 2022
13
## \authors Artem Busorgin, Kyle Cranmer (C++ version)
14
15
import
ROOT
16
17
# --------------------------------------
18
# An example of setting a limit in a number counting experiment with uncertainty on background and signal
19
20
# to time the macro
21
t =
ROOT.TStopwatch
()
22
t.Start
()
23
24
# --------------------------------------
25
# The Model building stage
26
# --------------------------------------
27
wspace =
ROOT.RooWorkspace
()
28
wspace.factory
(
29
"Poisson::countingModel(obs[150,0,300], "
"sum(s[50,0,120]*ratioSigEff[1.,0,3.],b[100]*ratioBkgEff[1.,0.,3.]))"
30
)
# counting model
31
wspace.factory
(
"Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)"
)
# 5% signal efficiency uncertainty
32
wspace.factory
(
"Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)"
)
# 10% background efficiency uncertainty
33
wspace.factory
(
"PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)"
)
# product of terms
34
wspace.Print
()
35
36
modelWithConstraints = wspace[
"modelWithConstraints"
]
# get the model
37
obs = wspace[
"obs"
]
# get the observable
38
s = wspace[
"s"
]
# get the signal we care about
39
b = wspace[
"b"
]
# get the background and set it to a constant. Uncertainty included in ratioBkgEff
40
b.setConstant
()
41
42
ratioSigEff = wspace[
"ratioSigEff"
]
# get uncertain parameter to constrain
43
ratioBkgEff = wspace[
"ratioBkgEff"
]
# get uncertain parameter to constrain
44
constrainedParams = {ratioSigEff, ratioBkgEff}
# need to constrain these in the fit (should change default behavior)
45
46
gSigEff = wspace[
"gSigEff"
]
# global observables for signal efficiency
47
gSigBkg = wspace[
"gSigBkg"
]
# global obs for background efficiency
48
gSigEff.setConstant
()
49
gSigBkg.setConstant
()
50
51
# Create an example dataset with 160 observed events
52
obs.setVal
(160.0)
53
dataOrig =
ROOT.RooDataSet
(
"exampleData"
,
"exampleData"
, {obs})
54
dataOrig.add
(obs)
55
56
# not necessary
57
modelWithConstraints.fitTo
(dataOrig, Constrain=constrainedParams, PrintLevel=-1)
58
59
# Now let's make some confidence intervals for s, our parameter of interest
60
modelConfig =
ROOT.RooStats.ModelConfig
(wspace)
61
modelConfig.SetPdf
(modelWithConstraints)
62
modelConfig.SetParametersOfInterest
({s})
63
modelConfig.SetNuisanceParameters
(constrainedParams)
64
modelConfig.SetObservables
(obs)
65
modelConfig.SetGlobalObservables
({gSigEff, gSigBkg})
66
modelConfig.SetName
(
"ModelConfig"
)
67
wspace.Import
(modelConfig)
68
wspace.Import
(dataOrig)
69
wspace.SetName
(
"w"
)
70
# wspace.writeToFile("rs101_ws.root")
71
72
# Make sure we reference the data in the workspace from now on
73
data = wspace[
dataOrig.GetName
()]
74
75
# First, let's use a Calculator based on the Profile Likelihood Ratio
76
plc =
ROOT.RooStats.ProfileLikelihoodCalculator
(data, modelConfig)
77
plc.SetTestSize
(0.05)
78
lrinterval =
plc.GetInterval
()
79
80
# Let's make a plot
81
dataCanvas =
ROOT.TCanvas
(
"dataCanvas"
)
82
dataCanvas.Divide
(2, 1)
83
dataCanvas.cd
(1)
84
85
plotInt =
ROOT.RooStats.LikelihoodIntervalPlot
(lrinterval)
86
plotInt.SetTitle
(
"Profile Likelihood Ratio and Posterior for S"
)
87
plotInt.Draw
()
88
89
# Second, use a Calculator based on the Feldman Cousins technique
90
fc =
ROOT.RooStats.FeldmanCousins
(data, modelConfig)
91
fc.UseAdaptiveSampling
(
True
)
92
fc.FluctuateNumDataEntries
(
False
)
# number counting analysis: dataset always has 1 entry with N events observed
93
fc.SetNBins
(100)
# number of points to test per parameter
94
fc.SetTestSize
(0.05)
95
# fc.SaveBeltToFile(True) # optional
96
fcint =
fc.GetInterval
()
97
98
fit =
modelWithConstraints.fitTo
(data, Save=
True
, PrintLevel=-1)
99
100
# Third, use a Calculator based on Markov Chain monte carlo
101
# Before configuring the calculator, let's make a ProposalFunction
102
# that will achieve a high acceptance rate
103
ph =
ROOT.RooStats.ProposalHelper
()
104
ph.SetVariables
(
fit.floatParsFinal
())
105
ph.SetCovMatrix
(
fit.covarianceMatrix
())
106
ph.SetUpdateProposalParameters
(
True
)
107
ph.SetCacheSize
(100)
108
pdfProp =
ph.GetProposalFunction
()
109
110
mc =
ROOT.RooStats.MCMCCalculator
(data, modelConfig)
111
mc.SetNumIters
(20000)
# steps to propose in the chain
112
mc.SetTestSize
(0.05)
# 95% CL
113
mc.SetNumBurnInSteps
(40)
# ignore first N steps in chain as "burn in"
114
mc.SetProposalFunction
(pdfProp)
115
mc.SetLeftSideTailFraction
(0.5)
# find a "central" interval
116
mcInt =
mc.GetInterval
()
117
118
# Get Lower and Upper limits from Profile Calculator
119
print(
"Profile lower limit on s = "
,
lrinterval.LowerLimit
(s))
120
print(
"Profile upper limit on s = "
,
lrinterval.UpperLimit
(s))
121
122
# Get Lower and Upper limits from FeldmanCousins with profile construction
123
if
fcint:
124
fcul =
fcint.UpperLimit
(s)
125
fcll =
fcint.LowerLimit
(s)
126
print(
"FC lower limit on s = "
, fcll)
127
print(
"FC upper limit on s = "
, fcul)
128
fcllLine =
ROOT.TLine
(fcll, 0, fcll, 1)
129
fculLine =
ROOT.TLine
(fcul, 0, fcul, 1)
130
fcllLine.SetLineColor
(
ROOT.kRed
)
131
fculLine.SetLineColor
(
ROOT.kRed
)
132
fcllLine.Draw
(
"same"
)
133
fculLine.Draw
(
"same"
)
134
dataCanvas.Update
()
135
136
# Plot MCMC interval and print some statistics
137
mcPlot =
ROOT.RooStats.MCMCIntervalPlot
(mcInt)
138
mcPlot.SetLineColor
(
ROOT.kMagenta
)
139
mcPlot.SetLineWidth
(2)
140
mcPlot.Draw
(
"same"
)
141
142
mcul =
mcInt.UpperLimit
(s)
143
mcll =
mcInt.LowerLimit
(s)
144
print(
"MCMC lower limit on s = "
, mcll)
145
print(
"MCMC upper limit on s = "
, mcul)
146
print(
"MCMC Actual confidence level: "
,
mcInt.GetActualConfidenceLevel
())
147
148
# 3-d plot of the parameter points
149
dataCanvas.cd
(2)
150
151
# also plot the points in the markov chain
152
chainData =
mcInt.GetChainAsDataSet
()
153
154
print(
"plotting the chain data - nentries = "
,
chainData.numEntries
())
155
chain =
ROOT.RooStats.GetAsTTree
(
"chainTreeData"
,
"chainTreeData"
, chainData)
156
chain.SetMarkerStyle
(6)
157
chain.SetMarkerColor
(
ROOT.kRed
)
158
159
chain.Draw
(
"s:ratioSigEff:ratioBkgEff"
,
"nll_MarkovChain_local_"
,
"box"
)
# 3-d box proportional to posterior
160
161
# the points used in the profile construction
162
parScanData =
fc.GetPointsToScan
()
163
print(
"plotting the scanned points used in the frequentist construction - npoints = "
,
parScanData.numEntries
())
164
165
gr =
ROOT.TGraph2D
(
parScanData.numEntries
())
166
for
ievt
in
range
(
parScanData.numEntries
()):
167
evt =
parScanData.get
(ievt)
168
x =
evt.getRealValue
(
"ratioBkgEff"
)
169
y =
evt.getRealValue
(
"ratioSigEff"
)
170
z =
evt.getRealValue
(
"s"
)
171
gr.SetPoint
(ievt, x, y, z)
172
173
gr.SetMarkerStyle
(24)
174
gr.Draw
(
"P SAME"
)
175
176
# print timing info
177
t.Stop
()
178
t.Print
()
179
180
dataCanvas.SaveAs
(
"rs101_limitexample.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
ROOT::Detail::TRangeCast
Definition
TCollection.h:311
tutorials
roostats
rs101_limitexample.py
ROOT v6-34 - Reference Guide Generated on Fri Mar 28 2025 14:54:41 (GVA Time) using Doxygen 1.10.0