44 using namespace RooFit;
45 using namespace RooStats;
47 struct BayesianNumericalOptions {
49 double confLevel = 0.95;
50 TString integrationType =
"";
57 bool plotPosterior =
false;
62 double nSigmaNuisance = -1;
66 BayesianNumericalOptions optBayes;
68 void StandardBayesianNumericalDemo(
const char *infile =
"",
const char *workspaceName =
"combined",
69 const char *modelConfigName =
"ModelConfig",
const char *dataName =
"obsData")
73 double confLevel = optBayes.confLevel;
74 TString integrationType = optBayes.integrationType;
75 int nToys = optBayes.nToys;
76 bool scanPosterior = optBayes.scanPosterior;
77 bool plotPosterior = optBayes.plotPosterior;
78 int nScanPoints = optBayes.nScanPoints;
79 int intervalType = optBayes.intervalType;
80 int maxPOI = optBayes.maxPOI;
81 double nSigmaNuisance = optBayes.nSigmaNuisance;
87 const char *filename =
"";
88 if (!strcmp(infile,
"")) {
89 filename =
"results/example_combined_GaussExample_model.root";
90 bool fileExist = !gSystem->AccessPathName(filename);
94 cout <<
"HistFactory file cannot be generated on Windows - exit" << endl;
98 cout <<
"will run standard hist2workspace example" << endl;
99 gROOT->ProcessLine(
".! prepareHistFactory .");
100 gROOT->ProcessLine(
".! hist2workspace config/example.xml");
101 cout <<
"\n\n---------------------" << endl;
102 cout <<
"Done creating example input" << endl;
103 cout <<
"---------------------\n\n" << endl;
110 TFile *file = TFile::Open(filename);
114 cout <<
"StandardRooStatsDemoMacro: Input file " << filename <<
" is not found" << endl;
123 RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
125 cout <<
"workspace not found" << endl;
130 ModelConfig *mc = (ModelConfig *)w->obj(modelConfigName);
133 RooAbsData *data = w->data(dataName);
138 cout <<
"data or ModelConfig was not found" << endl;
151 RooUniform prior(
"prior",
"", *mc->GetParametersOfInterest());
153 mc->SetPriorPdf(*w->pdf(
"prior"));
157 if (nSigmaNuisance > 0) {
158 RooAbsPdf *pdf = mc->GetPdf();
161 pdf->fitTo(*data, Save(
true), Minimizer(ROOT::Math::MinimizerOptions::DefaultMinimizerType().c_str()),
162 Hesse(
true), PrintLevel(ROOT::Math::MinimizerOptions::DefaultPrintLevel() - 1));
165 RooArgList nuisPar(*mc->GetNuisanceParameters());
166 for (
int i = 0; i < nuisPar.getSize(); ++i) {
167 RooRealVar *v =
dynamic_cast<RooRealVar *
>(&nuisPar[i]);
169 v->setMin(TMath::Max(v->getMin(), v->getVal() - nSigmaNuisance * v->getError()));
170 v->setMax(TMath::Min(v->getMax(), v->getVal() + nSigmaNuisance * v->getError()));
171 std::cout <<
"setting interval for nuisance " << v->GetName() <<
" : [ " << v->getMin() <<
" , "
172 << v->getMax() <<
" ]" << std::endl;
176 BayesianCalculator bayesianCalc(*data, *mc);
177 bayesianCalc.SetConfidenceLevel(confLevel);
181 if (intervalType == 0)
182 bayesianCalc.SetShortestInterval();
183 if (intervalType == 1)
184 bayesianCalc.SetLeftSideTailFraction(0.5);
185 if (intervalType == 2)
186 bayesianCalc.SetLeftSideTailFraction(0.);
188 if (!integrationType.IsNull()) {
189 bayesianCalc.SetIntegrationType(integrationType);
190 bayesianCalc.SetNumIters(nToys);
194 if (integrationType.Contains(
"TOYMC")) {
195 RooAbsPdf *nuisPdf = RooStats::MakeNuisancePdf(*mc,
"nuisance_pdf");
196 cout <<
"using TOYMC integration: make nuisance pdf from the model " << std::endl;
198 bayesianCalc.ForceNuisancePdf(*nuisPdf);
199 scanPosterior =
true;
204 bayesianCalc.SetScanOfPosterior(nScanPoints);
206 RooRealVar *poi = (RooRealVar *)mc->GetParametersOfInterest()->first();
207 if (maxPOI != -999 && maxPOI > poi->getMin())
210 SimpleInterval *interval = bayesianCalc.GetInterval();
213 cout <<
"\n>>>> RESULT : " << confLevel * 100 <<
"% interval on " << poi->GetName() <<
" is : ["
214 << interval->LowerLimit() <<
", " << interval->UpperLimit() <<
"] " << endl;
226 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::Ignore);
228 cout <<
"\nDrawing plot of posterior function....." << endl;
231 bayesianCalc.SetScanOfPosterior(nScanPoints);
233 RooPlot *plot = bayesianCalc.GetPosteriorPlot();