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
[#1] INFO:ObjectHandling -- RooWorkspace::import(combined) importing RooUniform::prior
[#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: (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 CPU computation library compiled with -mavx2
[#1] INFO:Eval -- BayesianCalculator::GetPosteriorFunction : nll value 15.6103 poi value = 1
[#1] INFO:Eval -- BayesianCalculator::GetPosteriorFunction : minimum of NLL vs POI for POI = 1.12121 min NLL = 15.5792
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter alpha_syst2. Note that passing copies of the parameters is not supported.
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter alpha_syst3. Note that passing copies of the parameters is not supported.
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter gamma_stat_channel1_bin_0. Note that passing copies of the parameters is not supported.
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter gamma_stat_channel1_bin_1. Note that passing copies of the parameters is not supported.
[#1] INFO:InputArguments -- BayesianCalculator:GetInterval Compute the interval from the posterior cdf
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter alpha_syst2. Note that passing copies of the parameters is not supported.
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter alpha_syst3. Note that passing copies of the parameters is not supported.
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter gamma_stat_channel1_bin_0. Note that passing copies of the parameters is not supported.
[#0] WARNING:InputArguments -- RooRealBinding: The function prior does not depend on the parameter gamma_stat_channel1_bin_1. Note that passing copies of the parameters is not supported.
[#1] INFO:NumericIntegration -- PosteriorCdfFunction - integral of posterior = 0.0560849 +/- 0.000423427
[#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;
if (!strcmp(infile, "")) {
filename =
"results/example_combined_GaussExample_model.root";
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
if (!file) {
cout <<
"StandardRooStatsDemoMacro: Input file " <<
filename <<
" is not found" << endl;
return;
}
cout << "workspace not found" << endl;
return;
}
cout << "data or ModelConfig was not found" << endl;
return;
}
if (nSigmaNuisance > 0) {
assert(pdf);
std::unique_ptr<RooFitResult> res{
res->Print();
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;
}
}
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())
cout <<
"\n>>>> RESULT : " << confLevel * 100 <<
"% interval on " << poi->
GetName() <<
" is : ["
if (!plotPosterior)
return;
cout << "\nDrawing plot of posterior function....." << endl;
bayesianCalc.SetScanOfPosterior(nScanPoints);
}
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
R__EXTERN TSystem * gSystem
static int DefaultPrintLevel()
static const std::string & DefaultMinimizerType()
void Print(Option_t *options=nullptr) const override
Print the object to the defaultPrintStream().
RooAbsArg * first() const
Abstract base class for binned and unbinned datasets.
Abstract interface for all probability density functions.
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
virtual double getMin(const char *name=nullptr) const
Get minimum 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.
Plot frame and a container for graphics objects within that frame.
Variable that can be changed from the outside.
void setMax(const char *name, double value)
Set maximum of name range to given value.
BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation of a...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
virtual void SetPriorPdf(const RooAbsPdf &pdf)
Set the Prior Pdf, add to the workspace if not already there.
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)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
SimpleInterval is a concrete implementation of the ConfInterval interface.
virtual double UpperLimit()
return the interval upper limit
virtual double LowerLimit()
return the interval lower limit
Persistable container for RooFit projects.
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
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.
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
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.
RooCmdArg Minimizer(const char *type, const char *alg=nullptr)
RooCmdArg Hesse(bool flag=true)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.