Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
StandardBayesianMCMCDemo.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Standard demo of the Bayesian MCMC calculator

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 default parameters 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 MCMCCalculator is a Bayesian tool that uses the Metropolis-Hastings algorithm to efficiently integrate in many dimensions. It is not as accurate as the BayesianCalculator for simple problems, but it scales to much more complicated cases.

[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
RooAbsTestStatistic::initSimMode: creating slave calculator #0 for state channel1 (2 dataset entries)
[#1] INFO:Fitting -- RooAbsTestStatistic::initSimMode: created 1 slave calculators.
[#1] INFO:Minimization -- Including the following constraint terms in minimization: (lumiConstraint,alpha_syst1Constraint,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)
Metropolis-Hastings progress: ....................................................................................................
[#1] INFO:Eval -- Proposal acceptance rate: 27.377%
[#1] INFO:Eval -- Number of steps in chain: 27377
>>>> RESULT : 95% interval on SigXsecOverSM is : [0, 2.17276]
#include "TFile.h"
#include "TROOT.h"
#include "TCanvas.h"
#include "TMath.h"
#include "TSystem.h"
#include "RooWorkspace.h"
#include "RooAbsData.h"
#include "RooFitResult.h"
using namespace RooFit;
using namespace RooStats;
struct BayesianMCMCOptions {
double confLevel = 0.95;
int intervalType = 2; // type of interval (0 is shortest, 1 central, 2 upper limit)
double maxPOI = -999; // force different values of POI for doing the scan (default is given value)
double minPOI = -999;
int numIters = 100000; // number of iterations
int numBurnInSteps = 100; // number of burn in steps to be ignored
};
BayesianMCMCOptions optMCMC;
void StandardBayesianMCMCDemo(const char *infile = "", const char *workspaceName = "combined",
const char *modelConfigName = "ModelConfig", const char *dataName = "obsData")
{
// -------------------------------------------------------
// First part is just to access a user-defined file
// or create the standard example file if it doesn't exist
const char *filename = "";
if (!strcmp(infile, "")) {
filename = "results/example_combined_GaussExample_model.root";
bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
// if file does not exists generate with histfactory
if (!fileExist) {
// Normally this would be run on the command line
cout << "will run standard hist2workspace example" << endl;
gROOT->ProcessLine(".! prepareHistFactory .");
gROOT->ProcessLine(".! hist2workspace config/example.xml");
cout << "\n\n---------------------" << endl;
cout << "Done creating example input" << endl;
cout << "---------------------\n\n" << endl;
}
} else
filename = infile;
// Try to open the file
// if input file was specified byt not found, quit
if (!file) {
cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
return;
}
// -------------------------------------------------------
// Tutorial starts here
// -------------------------------------------------------
// get the workspace out of the file
RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
if (!w) {
cout << "workspace not found" << endl;
return;
}
// get the modelConfig out of the file
ModelConfig *mc = (ModelConfig *)w->obj(modelConfigName);
// get the modelConfig out of the file
RooAbsData *data = w->data(dataName);
// make sure ingredients are found
if (!data || !mc) {
w->Print();
cout << "data or ModelConfig was not found" << endl;
return;
}
// Want an efficient proposal function
// default is uniform.
/*
// this one is based on the covariance matrix of fit
RooFitResult* fit = mc->GetPdf()->fitTo(*data,Save());
ProposalHelper ph;
ph.SetVariables((RooArgSet&)fit->floatParsFinal());
ph.SetCovMatrix(fit->covarianceMatrix());
ph.SetUpdateProposalParameters(kTRUE); // auto-create mean vars and add mappings
ph.SetCacheSize(100);
ProposalFunction* pf = ph.GetProposalFunction();
*/
// this proposal function seems fairly robust
// -------------------------------------------------------
// create and use the MCMCCalculator
// to find and plot the 95% credible interval
// on the parameter of interest as specified
// in the model config
MCMCCalculator mcmc(*data, *mc);
mcmc.SetConfidenceLevel(optMCMC.confLevel); // 95% interval
// mcmc.SetProposalFunction(*pf);
mcmc.SetProposalFunction(sp);
mcmc.SetNumIters(optMCMC.numIters); // Metropolis-Hastings algorithm iterations
mcmc.SetNumBurnInSteps(optMCMC.numBurnInSteps); // first N steps to be ignored as burn-in
// default is the shortest interval.
if (optMCMC.intervalType == 0)
mcmc.SetIntervalType(MCMCInterval::kShortest); // for shortest interval (not really needed)
if (optMCMC.intervalType == 1)
mcmc.SetLeftSideTailFraction(0.5); // for central interval
if (optMCMC.intervalType == 2)
mcmc.SetLeftSideTailFraction(0.); // for upper limit
if (optMCMC.minPOI != -999)
firstPOI->setMin(optMCMC.minPOI);
if (optMCMC.maxPOI != -999)
firstPOI->setMax(optMCMC.maxPOI);
MCMCInterval *interval = mcmc.GetInterval();
// make a plot
// TCanvas* c1 =
auto c1 = new TCanvas("IntervalPlot");
MCMCIntervalPlot plot(*interval);
TCanvas *c2 = new TCanvas("extraPlots");
const RooArgSet *list = mc->GetNuisanceParameters();
if (list->getSize() > 1) {
double n = list->getSize();
int ny = TMath::CeilNint(sqrt(n));
int nx = TMath::CeilNint(double(n) / ny);
c2->Divide(nx, ny);
}
// draw a scatter plot of chain results for poi vs each nuisance parameters
RooRealVar *nuis = NULL;
int iPad = 1; // iPad, that's funny
while ((nuis = (RooRealVar *)it.Next())) {
c2->cd(iPad++);
plot.DrawChainScatter(*firstPOI, *nuis);
}
// print out the interval on the first Parameter of Interest
cout << "\n>>>> RESULT : " << optMCMC.confLevel * 100 << "% interval on " << firstPOI->GetName() << " is : ["
<< interval->LowerLimit(*firstPOI) << ", " << interval->UpperLimit(*firstPOI) << "] " << endl;
gPad = c1;
}
winID h TVirtualViewer3D TVirtualGLPainter char TVirtualGLPainter plot
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
#define gROOT
Definition TROOT.h:407
R__EXTERN TSystem * gSystem
Definition TSystem.h:560
#define gPad
Int_t getSize() const
Return the number of elements in the collection.
RooAbsArg * first() const
TIterator * createIterator(bool dir=kIterForward) const
TIterator-style iteration over contained elements.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
void setMin(const char *name, double value)
Set minimum of name range to given value.
void setMax(const char *name, double value)
Set maximum of name range to given value.
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.
virtual double UpperLimit(RooRealVar &param)
get the highest value of param that is within the confidence interval
virtual double 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:35
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
Class implementing a proposal function that samples the parameter space by moving only in one coordin...
Persistable container for RooFit projects.
The Canvas class.
Definition TCanvas.h:23
A ROOT file is composed of a header, followed by consecutive data records (TKey instances) with a wel...
Definition TFile.h:53
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4075
TObject * Next()
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
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:1283
return c1
Definition legend1.C:41
const Int_t n
Definition legend1.C:16
return c2
Definition legend2.C:14
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
Namespace for the RooStats classes.
Definition Asimov.h:19
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:674
Definition file.py:1
Author
Kyle Cranmer

Definition in file StandardBayesianMCMCDemo.C.