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)
53dataOrig = ROOT.RooDataSet("exampleData", "exampleData", {obs})
54dataOrig.add(obs)
55
56# not necessary
57modelWithConstraints.fitTo(dataOrig, 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(dataOrig)
69wspace.SetName("w")
70# wspace.writeToFile("rs101_ws.root")
71
72# Make sure we reference the data in the workspace from now on
73data = wspace[dataOrig.GetName()]
74
75# First, let's use a Calculator based on the Profile Likelihood Ratio
76plc = ROOT.RooStats.ProfileLikelihoodCalculator(data, modelConfig)
77plc.SetTestSize(0.05)
78lrinterval = plc.GetInterval()
79
80# Let's make a plot
81dataCanvas = ROOT.TCanvas("dataCanvas")
82dataCanvas.Divide(2, 1)
83dataCanvas.cd(1)
84
85plotInt = ROOT.RooStats.LikelihoodIntervalPlot(lrinterval)
86plotInt.SetTitle("Profile Likelihood Ratio and Posterior for S")
87plotInt.Draw()
88
89# Second, use a Calculator based on the Feldman Cousins technique
90fc = ROOT.RooStats.FeldmanCousins(data, modelConfig)
91fc.UseAdaptiveSampling(True)
92fc.FluctuateNumDataEntries(False) # number counting analysis: dataset always has 1 entry with N events observed
93fc.SetNBins(100) # number of points to test per parameter
94fc.SetTestSize(0.05)
95# fc.SaveBeltToFile(True) # optional
96fcint = fc.GetInterval()
97
98fit = 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
103ph = ROOT.RooStats.ProposalHelper()
104ph.SetVariables(fit.floatParsFinal())
105ph.SetCovMatrix(fit.covarianceMatrix())
106ph.SetUpdateProposalParameters(True)
107ph.SetCacheSize(100)
108pdfProp = ph.GetProposalFunction()
109
110mc = ROOT.RooStats.MCMCCalculator(data, modelConfig)
111mc.SetNumIters(20000) # steps to propose in the chain
112mc.SetTestSize(0.05) # 95% CL
113mc.SetNumBurnInSteps(40) # ignore first N steps in chain as "burn in"
114mc.SetProposalFunction(pdfProp)
115mc.SetLeftSideTailFraction(0.5) # find a "central" interval
116mcInt = mc.GetInterval()
117
118# Get Lower and Upper limits from Profile Calculator
119print("Profile lower limit on s = ", lrinterval.LowerLimit(s))
120print("Profile upper limit on s = ", lrinterval.UpperLimit(s))
121
122# Get Lower and Upper limits from FeldmanCousins with profile construction
123if 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
137mcPlot = ROOT.RooStats.MCMCIntervalPlot(mcInt)
138mcPlot.SetLineColor(ROOT.kMagenta)
139mcPlot.SetLineWidth(2)
140mcPlot.Draw("same")
141
142mcul = mcInt.UpperLimit(s)
143mcll = mcInt.LowerLimit(s)
144print("MCMC lower limit on s = ", mcll)
145print("MCMC upper limit on s = ", mcul)
146print("MCMC Actual confidence level: ", mcInt.GetActualConfidenceLevel())
147
148# 3-d plot of the parameter points
149dataCanvas.cd(2)
150
151# also plot the points in the markov chain
152chainData = mcInt.GetChainAsDataSet()
153
154print("plotting the chain data - nentries = ", chainData.numEntries())
155chain = ROOT.RooStats.GetAsTTree("chainTreeData", "chainTreeData", chainData)
156chain.SetMarkerStyle(6)
157chain.SetMarkerColor(ROOT.kRed)
158
159chain.Draw("s:ratioSigEff:ratioBkgEff", "nll_MarkovChain_local_", "box") # 3-d box proportional to posterior
160
161# the points used in the profile construction
162parScanData = fc.GetPointsToScan()
163print("plotting the scanned points used in the frequentist construction - npoints = ", parScanData.numEntries())
164
165gr = ROOT.TGraph2D(parScanData.numEntries())
166for 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
173gr.SetMarkerStyle(24)
174gr.Draw("P SAME")
175
176# print timing info
177t.Stop()
178t.Print()
179
180dataCanvas.SaveAs("rs101_limitexample.png")