Logo ROOT   6.16/01
Reference Guide
StandardBayesianMCMCDemo.C
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
29
30#include "TFile.h"
31#include "TROOT.h"
32#include "TCanvas.h"
33#include "TMath.h"
34#include "TSystem.h"
35#include "RooWorkspace.h"
36#include "RooAbsData.h"
37
45#include "RooFitResult.h"
46
47
48using namespace RooFit;
49using namespace RooStats;
50
51struct BayesianMCMCOptions {
52
53 double confLevel = 0.95;
54 int intervalType = 2; // type of interval (0 is shortest, 1 central, 2 upper limit)
55 double maxPOI = -999; // force different values of POI for doing the scan (default is given value)
56 double minPOI = -999;
57 int numIters = 100000; // number of iterations
58 int numBurnInSteps = 100; // number of burn in steps to be ignored
59};
60
61BayesianMCMCOptions optMCMC;
62
63void StandardBayesianMCMCDemo(const char* infile = "",
64 const char* workspaceName = "combined",
65 const char* modelConfigName = "ModelConfig",
66 const char* dataName = "obsData"){
67
68 // -------------------------------------------------------
69 // First part is just to access a user-defined file
70 // or create the standard example file if it doesn't exist
71
72
73
74 const char* filename = "";
75 if (!strcmp(infile,"")) {
76 filename = "results/example_combined_GaussExample_model.root";
77 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
78 // if file does not exists generate with histfactory
79 if (!fileExist) {
80#ifdef _WIN32
81 cout << "HistFactory file cannot be generated on Windows - exit" << endl;
82 return;
83#endif
84 // Normally this would be run on the command line
85 cout <<"will run standard hist2workspace example"<<endl;
86 gROOT->ProcessLine(".! prepareHistFactory .");
87 gROOT->ProcessLine(".! hist2workspace config/example.xml");
88 cout <<"\n\n---------------------"<<endl;
89 cout <<"Done creating example input"<<endl;
90 cout <<"---------------------\n\n"<<endl;
91 }
92
93 }
94 else
95 filename = infile;
96
97 // Try to open the file
98 TFile *file = TFile::Open(filename);
99
100 // if input file was specified byt not found, quit
101 if(!file ){
102 cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
103 return;
104 }
105
106
107
108 // -------------------------------------------------------
109 // Tutorial starts here
110 // -------------------------------------------------------
111
112 // get the workspace out of the file
113 RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
114 if(!w){
115 cout <<"workspace not found" << endl;
116 return;
117 }
118
119 // get the modelConfig out of the file
120 ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
121
122 // get the modelConfig out of the file
123 RooAbsData* data = w->data(dataName);
124
125 // make sure ingredients are found
126 if(!data || !mc){
127 w->Print();
128 cout << "data or ModelConfig was not found" <<endl;
129 return;
130 }
131
132 // Want an efficient proposal function
133 // default is uniform.
134
135 /*
136 // this one is based on the covariance matrix of fit
137 RooFitResult* fit = mc->GetPdf()->fitTo(*data,Save());
138 ProposalHelper ph;
139 ph.SetVariables((RooArgSet&)fit->floatParsFinal());
140 ph.SetCovMatrix(fit->covarianceMatrix());
141 ph.SetUpdateProposalParameters(kTRUE); // auto-create mean vars and add mappings
142 ph.SetCacheSize(100);
143 ProposalFunction* pf = ph.GetProposalFunction();
144 */
145
146 // this proposal function seems fairly robust
147 SequentialProposal sp(0.1);
148 // -------------------------------------------------------
149 // create and use the MCMCCalculator
150 // to find and plot the 95% credible interval
151 // on the parameter of interest as specified
152 // in the model config
153 MCMCCalculator mcmc(*data,*mc);
154 mcmc.SetConfidenceLevel(optMCMC.confLevel); // 95% interval
155 // mcmc.SetProposalFunction(*pf);
156 mcmc.SetProposalFunction(sp);
157 mcmc.SetNumIters(optMCMC.numIters); // Metropolis-Hastings algorithm iterations
158 mcmc.SetNumBurnInSteps(optMCMC.numBurnInSteps); // first N steps to be ignored as burn-in
159
160 // default is the shortest interval.
161 if (optMCMC.intervalType == 0) mcmc.SetIntervalType(MCMCInterval::kShortest); // for shortest interval (not really needed)
162 if (optMCMC.intervalType == 1) mcmc.SetLeftSideTailFraction(0.5); // for central interval
163 if (optMCMC.intervalType == 2) mcmc.SetLeftSideTailFraction(0.); // for upper limit
164
165 RooRealVar* firstPOI = (RooRealVar*) mc->GetParametersOfInterest()->first();
166 if (optMCMC.minPOI != -999)
167 firstPOI->setMin(optMCMC.minPOI);
168 if (optMCMC.maxPOI != -999)
169 firstPOI->setMax(optMCMC.maxPOI);
170
171 MCMCInterval* interval = mcmc.GetInterval();
172
173 // make a plot
174 //TCanvas* c1 =
175 auto c1 = new TCanvas("IntervalPlot");
176 MCMCIntervalPlot plot(*interval);
177 plot.Draw();
178
179 TCanvas* c2 = new TCanvas("extraPlots");
180 const RooArgSet* list = mc->GetNuisanceParameters();
181 if(list->getSize()>1){
182 double n = list->getSize();
183 int ny = TMath::CeilNint( sqrt(n) );
184 int nx = TMath::CeilNint(double(n)/ny);
185 c2->Divide( nx,ny);
186 }
187
188 // draw a scatter plot of chain results for poi vs each nuisance parameters
190 RooRealVar* nuis = NULL;
191 int iPad=1; // iPad, that's funny
192 while( (nuis = (RooRealVar*) it->Next() )){
193 c2->cd(iPad++);
194 plot.DrawChainScatter(*firstPOI,*nuis);
195 }
196
197 // print out the interval on the first Parameter of Interest
198 cout << "\n>>>> RESULT : " << optMCMC.confLevel*100 << "% interval on " <<firstPOI->GetName()<<" is : ["<<
199 interval->LowerLimit(*firstPOI) << ", "<<
200 interval->UpperLimit(*firstPOI) <<"] "<<endl;
201
202
203 gPad = c1;
204}
double sqrt(double)
#define gROOT
Definition: TROOT.h:410
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
#define gPad
Definition: TVirtualPad.h:286
Int_t getSize() const
RooAbsArg * first() const
TIterator * createIterator(Bool_t dir=kIterForward) const
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
void setMin(const char *name, Double_t value)
Set minimum of name range to given value.
Definition: RooRealVar.cxx:387
void setMax(const char *name, Double_t value)
Set maximum of name range to given value.
Definition: RooRealVar.cxx:417
Bayesian Calculator estimating an interval or a credible region using the Markov-Chain Monte Carlo me...
This class provides simple and straightforward utilities to plot a MCMCInterval object.
MCMCInterval is a concrete implementation of the RooStats::ConfInterval interface.
Definition: MCMCInterval.h:30
virtual Double_t UpperLimit(RooRealVar &param)
get the highest value of param that is within the confidence interval
virtual Double_t LowerLimit(RooRealVar &param)
get the lowest value of param that is within the confidence interval
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:231
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
Definition: ModelConfig.h:234
Class implementing a proposal function that samples the parameter space by moving only in one coordin...
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
void Print(Option_t *opts=0) const
Print contents of the workspace.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
The Canvas class.
Definition: TCanvas.h:31
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseGeneralPurpose, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3975
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1286
return c1
Definition: legend1.C:41
const Int_t n
Definition: legend1.C:16
return c2
Definition: legend2.C:14
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20
Int_t CeilNint(Double_t x)
Definition: TMath.h:687
Definition: file.py:1