This is a standard demo that can be used with any ROOT file prepared in the standard way. You specify:
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.
[#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)
[#1] INFO:Fitting -- RooAbsPdf::fitTo(simPdf) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using generic CPU library compiled with no vectorizations
[#1] INFO:Fitting -- Creation of NLL object took 9.40621 ms
Metropolis-Hastings progress: ....................................................................................................
[#1] INFO:Eval -- Proposal acceptance rate: 27.202%
[#1] INFO:Eval -- Number of steps in chain: 27202
>>>> RESULT : 95% interval on SigXsecOverSM is : [0, 2.16421]
struct BayesianMCMCOptions {
double confLevel = 0.95;
int intervalType = 2;
double maxPOI = -999;
double minPOI = -999;
int numIters = 100000;
int numBurnInSteps = 100;
};
BayesianMCMCOptions optMCMC;
void StandardBayesianMCMCDemo(const char *infile = "", const char *workspaceName = "combined",
const char *modelConfigName = "ModelConfig", const char *dataName = "obsData")
{
const char *filename = "";
if (!strcmp(infile, "")) {
filename = "results/example_combined_GaussExample_model.root";
bool fileExist = !
gSystem->AccessPathName(filename);
if (!fileExist) {
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;
if (!file) {
cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
return;
}
if (!w) {
cout << "workspace not found" << endl;
return;
}
if (!data || !mc) {
w->Print();
cout << "data or ModelConfig was not found" << endl;
return;
}
mcmc.SetConfidenceLevel(optMCMC.confLevel);
mcmc.SetProposalFunction(sp);
mcmc.SetNumIters(optMCMC.numIters);
mcmc.SetNumBurnInSteps(optMCMC.numBurnInSteps);
if (optMCMC.intervalType == 0)
if (optMCMC.intervalType == 1)
mcmc.SetLeftSideTailFraction(0.5);
if (optMCMC.intervalType == 2)
mcmc.SetLeftSideTailFraction(0.);
if (optMCMC.minPOI != -999)
firstPOI->
setMin(optMCMC.minPOI);
if (optMCMC.maxPOI != -999)
firstPOI->
setMax(optMCMC.maxPOI);
plot.Draw();
if (list->getSize() > 1) {
double n = list->getSize();
}
int iPad = 1;
plot.DrawChainScatter(*firstPOI, *nuis);
}
cout <<
"\n>>>> RESULT : " << optMCMC.confLevel * 100 <<
"% interval on " << firstPOI->
GetName() <<
" is : ["
}
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
RooAbsArg * first() const
Abstract base class for binned and unbinned datasets.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Variable that can be changed from the outside.
void setMin(const char *name, double value, bool shared=true)
Set minimum of name range to given value.
void setMax(const char *name, double value, bool shared=true)
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 ¶m)
get the highest value of param that is within the confidence interval
virtual double LowerLimit(RooRealVar ¶m)
get the lowest value of param that is within the confidence interval
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.
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A file, usually with extension .root, that stores data and code in the form of serialized objects in ...
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.
const char * GetName() const override
Returns name of object.
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...
RooStats::ModelConfig ModelConfig
Namespace for the RooStats classes.
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).