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