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 BayesianCalculator is based on Bayes's theorem and performs the integration using ROOT's numeric integration utilities
␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:ObjectHandling -- RooWorkspace::import(combined) importing RooUniform::prior
[#1] INFO:Minization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minization -- createNLL: caching constraint set under name CONSTR_OF_PDF_simPdf_FOR_OBS_channelCat:obs_x_channel1 with 4 entries
[#1] INFO:Minization -- Including the following contraint terms in minimization: (alpha_syst2Constraint,alpha_syst3Constraint,gamma_stat_channel1_bin_0_constraint,gamma_stat_channel1_bin_1_constraint)
[#1] INFO:Minization -- The following global observables have been defined: (nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
RooAbsTestStatistic::initSimMode: creating slave calculator #0 for state channel1 (2 dataset entries)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:Fitting -- RooAbsTestStatistic::initSimMode: created 1 slave calculators.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:Eval -- BayesianCalculator::GetPosteriorFunction : nll value -1033.75 poi value = 1
[#1] INFO:Eval -- BayesianCalculator::GetPosteriorFunction : minimum of NLL vs POI for POI = 1.12121 min NLL = -1033.78
[#1] INFO:InputArguments -- BayesianCalculator:GetInterval Compute the interval from the posterior cdf
[#1] INFO:NumericIntegration -- PosteriorCdfFunction - integral of posterior = 0.0560849 +/- 0.000423427
[#0] WARNING:Eval -- BayesianCalculator::GetInterval : 315 errors reported in evaluating log-likelihood function
[#1] INFO:Eval -- BayesianCalculator::GetInterval - found a valid interval : [0.170773 , 2.33979 ]
>>>> RESULT : 95% interval on SigXsecOverSM is : [0.170773, 2.33979]
#include <cassert>
struct BayesianNumericalOptions {
double confLevel = 0.95;
int nToys =
10000;
bool scanPosterior =
false;
bool plotPosterior = false;
int nScanPoints = 50;
int intervalType = 1;
double maxPOI = -999;
double nSigmaNuisance = -1;
};
BayesianNumericalOptions optBayes;
void StandardBayesianNumericalDemo(const char *infile = "", const char *workspaceName = "combined",
const char *modelConfigName = "ModelConfig", const char *dataName = "obsData")
{
double confLevel = optBayes.confLevel;
TString integrationType = optBayes.integrationType;
int nToys = optBayes.nToys;
bool scanPosterior = optBayes.scanPosterior;
bool plotPosterior = optBayes.plotPosterior;
int nScanPoints = optBayes.nScanPoints;
int intervalType = optBayes.intervalType;
int maxPOI = optBayes.maxPOI;
double nSigmaNuisance = optBayes.nSigmaNuisance;
const char *filename = "";
if (!strcmp(infile, "")) {
filename = "results/example_combined_GaussExample_model.root";
if (!fileExist) {
#ifdef _WIN32
cout << "HistFactory file cannot be generated on Windows - exit" << endl;
return;
#endif
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;
cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
return;
}
if (!w) {
cout << "workspace not found" << endl;
return;
}
ModelConfig *mc = (ModelConfig *)w->
obj(modelConfigName);
if (!data || !mc) {
cout << "data or ModelConfig was not found" << endl;
return;
}
RooUniform prior(
"prior",
"", *mc->GetParametersOfInterest());
mc->SetPriorPdf(*w->
pdf(
"prior"));
if (nSigmaNuisance > 0) {
assert(pdf);
for (int i = 0; i < nuisPar.getSize(); ++i) {
v->setMin(
TMath::Max(
v->getMin(),
v->getVal() - nSigmaNuisance *
v->getError()));
v->setMax(
TMath::Min(
v->getMax(),
v->getVal() + nSigmaNuisance *
v->getError()));
std::cout <<
"setting interval for nuisance " <<
v->GetName() <<
" : [ " <<
v->getMin() <<
" , "
<<
v->getMax() <<
" ]" << std::endl;
}
}
BayesianCalculator bayesianCalc(*data, *mc);
bayesianCalc.SetConfidenceLevel(confLevel);
if (intervalType == 0)
bayesianCalc.SetShortestInterval();
if (intervalType == 1)
bayesianCalc.SetLeftSideTailFraction(0.5);
if (intervalType == 2)
bayesianCalc.SetLeftSideTailFraction(0.);
if (!integrationType.
IsNull()) {
bayesianCalc.SetIntegrationType(integrationType);
bayesianCalc.SetNumIters(nToys);
}
if (integrationType.
Contains(
"TOYMC")) {
cout << "using TOYMC integration: make nuisance pdf from the model " << std::endl;
bayesianCalc.ForceNuisancePdf(*nuisPdf);
scanPosterior = true;
}
if (scanPosterior)
bayesianCalc.SetScanOfPosterior(nScanPoints);
if (maxPOI != -999 && maxPOI > poi->
getMin())
SimpleInterval *interval = bayesianCalc.GetInterval();
cout <<
"\n>>>> RESULT : " << confLevel * 100 <<
"% interval on " << poi->
GetName() <<
" is : ["
<< interval->LowerLimit() << ", " << interval->UpperLimit() << "] " << endl;
if (!plotPosterior)
return;
cout << "\nDrawing plot of posterior function....." << endl;
bayesianCalc.SetScanOfPosterior(nScanPoints);
RooPlot *plot = bayesianCalc.GetPosteriorPlot();
}
R__EXTERN TSystem * gSystem
static int DefaultPrintLevel()
static const std::string & DefaultMinimizerType()
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
RooAbsData is the common abstract base class for binned and unbinned datasets.
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Fit PDF to given dataset.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
A RooPlot is a plot frame and a container for graphics objects within that frame.
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
RooRealVar represents a fundamental (non-derived) real valued object.
void setMax(const char *name, Double_t value)
Set maximum of name range to given value.
The RooWorkspace is a persistable container for RooFit projects.
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.
Bool_t import(const RooAbsArg &arg, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg())
Import a RooAbsArg object, e.g.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
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.
virtual const char * GetName() const
Returns name of object.
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Template specialisation used in RooAbsArg:
RooCmdArg Hesse(Bool_t flag=kTRUE)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Minimizer(const char *type, const char *alg=0)
Namespace for the RooStats classes.
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
Short_t Max(Short_t a, Short_t b)
Short_t Min(Short_t a, Short_t b)