Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
StandardProfileLikelihoodDemo.py File Reference

Detailed Description

View in nbviewer Open in SWAN
Standard demo of the Profile Likelihood calculator StandardProfileLikelihoodDemo

This is a standard demo that can be used with any ROOT file prepared in the standard way. You specify:

  • name for input ROOT file
  • name of workspace inside ROOT file that holds model and data
  • name of ModelConfig that specifies details for calculator tools
  • name of dataset With the values provided below the macro will attempt to run the standard hist2workspace example and read the ROOT file that it produces.

The actual heart of the demo is only about 10 lines long.

The ProfileLikelihoodCalculator is based on Wilks's theorem and the asymptotic properties of the profile likelihood ratio (eg. that it is chi-square distributed for the true value).

[#1] INFO:InputArguments -- The deprecated RooFit::CloneData(1) option passed to createNLL() is ignored.
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minimization -- Including the following constraint terms in minimization: (alpha_syst2Constraint,alpha_syst3Constraint,gamma_stat_channel1_bin_0_constraint,gamma_stat_channel1_bin_1_constraint)
[#1] INFO:Minimization -- The following global observables have been defined and their values are taken from the model: (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
[#1] INFO:Fitting -- RooAbsPdf::fitTo(simPdf) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#0] PROGRESS:Minimization -- ProfileLikelihoodCalcultor::DoGLobalFit - find MLE
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_simPdf_obsData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#0] PROGRESS:Minimization -- ProfileLikelihoodCalcultor::DoMinimizeNLL - using Minuit2 / with strategy 1
[#1] INFO:Minimization --
RooFitResult: minimized FCN value: 15.5775, estimated distance to minimum: 1.48403e-11
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
SigXsecOverSM 1.1154e+00 +/- 5.87e-01
alpha_syst2 -8.9189e-03 +/- 9.83e-01
alpha_syst3 1.7896e-02 +/- 9.48e-01
gamma_stat_channel1_bin_0 9.9955e-01 +/- 4.93e-02
gamma_stat_channel1_bin_1 1.0036e+00 +/- 8.01e-02
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[SigXsecOverSM]) Creating instance of MINUIT
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_simPdf_obsData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[SigXsecOverSM]) determining minimum likelihood for current configurations w.r.t all observable
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[SigXsecOverSM]) minimum found at (SigXsecOverSM=1.11588)
.
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[SigXsecOverSM]) Creating instance of MINUIT
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_simPdf_obsData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[SigXsecOverSM]) determining minimum likelihood for current configurations w.r.t all observable
[#0] ERROR:InputArguments -- RooArgSet::checkForDup: ERROR argument with name SigXsecOverSM is already in this set
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[SigXsecOverSM]) minimum found at (SigXsecOverSM=1.1155)
......................................................................................................
>>>> RESULT : 95.0% interval on SigXsecOverSM is : [0.0, 2.327436551649386]
making a plot of the profile likelihood function ....(if it is taking a lot of time use less points or the TF1 drawing option)
import ROOT
workspaceName = "combined"
modelConfigName = "ModelConfig"
dataName = "obsData"
confLevel = 0.95
nScanPoints = 50
plotAsTF1 = False
poiXMin = 1
poiXMax = 0
doHypoTest = False
nullParamValue = 0
filename = "results/example_combined_GaussExample_model.root"
# if file does not exists generate with histfactory
if ROOT.gSystem.AccessPathName(filename):
# Normally this would be run on the command line
print("will run standard hist2workspace example")
ROOT.gROOT.ProcessLine(".! prepareHistFactory .")
ROOT.gROOT.ProcessLine(".! hist2workspace config/example.xml")
print("\n\n---------------------")
print("Done creating example input")
print("---------------------\n\n")
file = ROOT.TFile.Open(filename)
# -------------------------------------------------------
# Tutorial starts here
# -------------------------------------------------------
# get the workspace out of the file
w = file.Get(workspaceName)
# get the modelConfig out of the file
mc = w[modelConfigName]
# get the modelConfig out of the file
data = w[dataName]
# ---------------------------------------------
# create and use the ProfileLikelihoodCalculator
# to find and plot the 95% confidence interval
# on the parameter of interest as specified
# in the model config
pl = ROOT.RooStats.ProfileLikelihoodCalculator(data, mc)
pl.SetConfidenceLevel(confLevel)
interval = pl.GetInterval()
# print out the interval on the first Parameter of Interest
firstPOI = mc.GetParametersOfInterest().first()
limit_lower, limit_upper = interval.LowerLimit(firstPOI), interval.UpperLimit(firstPOI)
print(f"\n>>>> RESULT : {confLevel * 100}% interval on {firstPOI.GetName()} is : [{limit_lower}, {limit_upper}]\n")
# make a plot
print(
"making a plot of the profile likelihood function ....(if it is taking a lot of time use less points or the "
"TF1 drawing option)\n"
)
plot = ROOT.RooStats.LikelihoodIntervalPlot(interval)
plot.SetNPoints(nScanPoints) # do not use too many points, it could become very slow for some models
if poiXMin < poiXMax:
plot.SetRange(poiXMin, poiXMax)
opt = ""
if plotAsTF1:
opt += "tf1"
plot.Draw(opt) # use option TF1 if too slow (plot.Draw("tf1")
# if requested perform also an hypothesis test for the significance
if doHypoTest:
nullparams = ROOT.RooArgSet("nullparams")
nullparams.addClone(firstPOI)
nullparams.setRealValue(firstPOI.GetName(), nullParamValue)
pl.SetNullParameters(nullparams)
print("Perform Test of Hypothesis : null Hypothesis is " + firstPOI.GetName() + str(nullParamValue))
result = pl.GetHypoTest()
print("\n>>>> Hypotheis Test Result ")
result.Print()
Authors
Akeem Hart, Kyle Cranmer (C++ Version)

Definition in file StandardProfileLikelihoodDemo.py.