Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
StandardBayesianMCMCDemo.py
Go to the documentation of this file.
1# \file
2# \ingroup tutorial_roostats
3# \notebook -js
4# Standard demo of the Bayesian MCMC calculator
5#
6# This is a standard demo that can be used with any ROOT file
7# prepared in the standard way. You specify:
8# - name for input ROOT file
9# - name of workspace inside ROOT file that holds model and data
10# - name of ModelConfig that specifies details for calculator tools
11# - name of dataset
12#
13# With default parameters the macro will attempt to run the
14# standard hist2workspace example and read the ROOT file
15# that it produces.
16#
17# The actual heart of the demo is only about 10 lines long.
18#
19# The MCMCCalculator is a Bayesian tool that uses
20# the Metropolis-Hastings algorithm to efficiently integrate
21# in many dimensions. It is not as accurate as the BayesianCalculator
22# for simple problems, but it scales to much more complicated cases.
23#
24# \macro_image
25# \macro_output
26# \macro_code
27#
28# \author Kyle Cranmer (C++ version), and P. P. (Python translation)
29
30import ROOT
31
32
33# general Structure definition
34class Struct:
35 def __init__(self, **entries):
36 self.__dict__.update(entries)
37
38
39BayesianMCMCOptions = Struct(
40 confLevel=0.95, # type of interval (0 is shortest, 1 central, 2 upper limit)
41 intervalType=2, # force different values of POI for doing the scan (default is given value)
42 maxPOI=-999,
43 minPOI=-999, # number of iterations
44 numIters=100000, # number of burn in steps to be ignored
45 numBurnInSteps=100,
46)
47
48optMCMC = BayesianMCMCOptions
49
50
51def StandardBayesianMCMCDemo(infile="", workspaceName="combined", modelConfigName="ModelConfig", dataName="obsData"):
52
53 # -------------------------------------------------------
54 # First part is just to access a user-defined file
55 # or create the standard example file if it doesn't exist
56
57 filename = ""
58 if infile == "":
59 filename = "results/example_combined_GaussExample_model.root"
60 fileExist = ROOT.gSystem.AccessPathName(filename)
61 # if file does not exists generate with histfactory
62 if not fileExist:
63 # Normally this would be run on the command line
64 print(f"will run standard hist2workspace example")
65 ROOT.gROOT.ProcessLine(".! prepareHistFactory .")
66 ROOT.gROOT.ProcessLine(".! hist2workspace config/example.xml")
67 print(f"\n\n---------------------")
68 print(f"Done creating example input")
69 print(f"---------------------\n\n")
70 else:
71 filename = infile
72 # Try to open the file
73 file = ROOT.TFile.Open(filename)
74
75 # if input file was specified but not found, quit
76 if not file:
77 print(f"StandardRooStatsDemoMacro: Input file {filename} is not found")
78 return
79
80 # -------------------------------------------------------
81 # Tutorial starts here
82 # -------------------------------------------------------
83
84 # get the workspace out of the file
85 # w = (RooWorkspace )file.Get(workspaceName)
86 w = file.Get(workspaceName)
87 if not w:
88 print(f"workspace not found")
89 return
90
91 # get the modelConfig out of the file
92 mc = w.obj(modelConfigName)
93
94 # get the modelConfig out of the file
95 data = w.data(dataName)
96
97 # make sure ingredients are found
98 if not data or not mc:
99 w.Print()
100 print(f"data or ModelConfig was not found")
101 return
102
103 # Want an efficient proposal function
104 # default is uniform.
105 """
106 #
107 # this one is based on the covariance matrix of fit
108 fit = mc.GetPdf().fitTo(data,Save())
109 ProposalHelper ph
110 ph.SetVariables((RooArgSet&)fit.floatParsFinal())
111 ph.SetCovMatrix(fit.covarianceMatrix())
112 ph.SetUpdateProposalParameters(True) # auto-create mean vars and add mappings
113 ph.SetCacheSize(100)
114 pf = ph.GetProposalFunction()
115 """
116
117 # this proposal function seems fairly robust
118 sp = ROOT.RooStats.SequentialProposal(0.1)
119 # -------------------------------------------------------
120 # create and use the MCMCCalculator
121 # to find and plot the 95% credible interval
122 # on the parameter of interest as specified
123 # in the model config
124 mcmc = ROOT.RooStats.MCMCCalculator(data, mc)
125 mcmc.SetConfidenceLevel(optMCMC.confLevel) # 95% interval
126 # mcmc.SetProposalFunction(*pf);
127 mcmc.SetProposalFunction(sp)
128 mcmc.SetNumIters(optMCMC.numIters) # Metropolis-Hastings algorithm iterations
129 mcmc.SetNumBurnInSteps(optMCMC.numBurnInSteps) # first N steps to be ignored as burn-in
130
131 # default is the shortest interval.
132 (optMCMC.intervalType == 0)
133 mcmc.SetIntervalType(ROOT.RooStats.MCMCInterval.kShortest) # for shortest interval (not really needed)
134 (optMCMC.intervalType == 1)
135 mcmc.SetLeftSideTailFraction(0.5) # for central interval
136 (optMCMC.intervalType == 2)
137 mcmc.SetLeftSideTailFraction(0.0) # for upper limit
138
139 firstPOI = mc.GetParametersOfInterest().first()
140 if optMCMC.minPOI != -999:
141 firstPOI.setMin(optMCMC.minPOI)
142 if optMCMC.maxPOI != -999:
143 firstPOI.setMax(optMCMC.maxPOI)
144
145 interval = mcmc.GetInterval()
146
147 # make a plot
148 c1 = ROOT.TCanvas("IntervalPlot")
149 plot = ROOT.RooStats.MCMCIntervalPlot(interval)
150 plot.Draw()
151
152 c2 = ROOT.TCanvas("extraPlots")
153 list = mc.GetNuisanceParameters()
154 if list.getSize() > 1:
155 n = list.getSize()
156 ny = ROOT.TMath.CeilNint(ROOT.sqrt(n))
157 nx = ROOT.TMath.CeilNint(ROOT.double(n) / ny)
158 c2.Divide(nx, ny)
159
160 # draw a scatter plot of chain results for poi vs each nuisance parameters
161 nuis = ROOT.kNone
162 iPad = 1 # iPad, that's funny
163
164 for nuis in mc.GetNuisanceParameters():
165 c2.cd(iPad)
166 iPad += 1
167 plot.DrawChainScatter(firstPOI, nuis)
168
169 # print out the interval on the first Parameter of Interest
170 print("\n>>>> RESULT : ", optMCMC.confLevel * 100, "interval on ", firstPOI.GetName(), "is : [")
171 print(interval.LowerLimit(firstPOI), interval.UpperLimit(firstPOI))
172 print("] ")
173 gPad = c1
174
175 c1.SaveAs("StandardBayesianMCMCDemo.1.IntervalPlot.png")
176 c2.SaveAs("StandardBayesianMCMCDemo.2.extraPlots.png")
177
178
179StandardBayesianMCMCDemo("", "combined", "ModelConfig", "obsData")
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)