27wspace = ROOT.RooWorkspace()
29 "Poisson::countingModel(obs[150,0,300], " "sum(s[50,0,120]*ratioSigEff[1.,0,3.],b[100]*ratioBkgEff[1.,0.,3.]))"
31wspace.factory(
"Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)")
32wspace.factory(
"Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)")
33wspace.factory(
"PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)")
36modelWithConstraints = wspace[
"modelWithConstraints"]
42ratioSigEff = wspace[
"ratioSigEff"]
43ratioBkgEff = wspace[
"ratioBkgEff"]
44constrainedParams = {ratioSigEff, ratioBkgEff}
46gSigEff = wspace[
"gSigEff"]
47gSigBkg = wspace[
"gSigBkg"]
53data = ROOT.RooDataSet(
"exampleData",
"exampleData", {obs})
57modelWithConstraints.fitTo(data, Constrain=constrainedParams, PrintLevel=-1)
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)
70wspace.writeToFile(
"rs101_ws.root")
73plc = ROOT.RooStats.ProfileLikelihoodCalculator(data, modelConfig)
75lrinterval = plc.GetInterval()
78dataCanvas = ROOT.TCanvas(
"dataCanvas")
79dataCanvas.Divide(2, 1)
82plotInt = ROOT.RooStats.LikelihoodIntervalPlot(lrinterval)
83plotInt.SetTitle(
"Profile Likelihood Ratio and Posterior for S")
87fc = ROOT.RooStats.FeldmanCousins(data, modelConfig)
88fc.UseAdaptiveSampling(
True)
89fc.FluctuateNumDataEntries(
False)
93fcint = fc.GetInterval()
95fit = modelWithConstraints.fitTo(data, Save=
True, PrintLevel=-1)
100ph = ROOT.RooStats.ProposalHelper()
101ph.SetVariables(fit.floatParsFinal())
102ph.SetCovMatrix(fit.covarianceMatrix())
103ph.SetUpdateProposalParameters(
True)
105pdfProp = ph.GetProposalFunction()
107mc = ROOT.RooStats.MCMCCalculator(data, modelConfig)
110mc.SetNumBurnInSteps(40)
111mc.SetProposalFunction(pdfProp)
112mc.SetLeftSideTailFraction(0.5)
113mcInt = mc.GetInterval()
116print(
"Profile lower limit on s = ", lrinterval.LowerLimit(s))
117print(
"Profile upper limit on s = ", lrinterval.UpperLimit(s))
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")
134mcPlot = ROOT.RooStats.MCMCIntervalPlot(mcInt)
135mcPlot.SetLineColor(ROOT.kMagenta)
136mcPlot.SetLineWidth(2)
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())
149chainData = mcInt.GetChainAsDataSet()
151print(
"plotting the chain data - nentries = ", chainData.numEntries())
152chain = ROOT.RooStats.GetAsTTree(
"chainTreeData",
"chainTreeData", chainData)
153chain.SetMarkerStyle(6)
154chain.SetMarkerColor(ROOT.kRed)
156chain.Draw(
"s:ratioSigEff:ratioBkgEff",
"nll_MarkovChain_local_",
"box")
159parScanData = fc.GetPointsToScan()
160print(
"plotting the scanned points used in the frequentist construction - npoints = ", parScanData.numEntries())
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)
177dataCanvas.SaveAs(
"rs101_limitexample.png")