Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
BayesianCalculator.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11
12/** \class RooStats::BayesianCalculator
13 \ingroup Roostats
14
15BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation
16of a credible interval using a Bayesian method.
17The class works only for one single parameter of interest and it integrates the likelihood function with the given prior
18probability density function to compute the posterior probability. The result of the class is a one dimensional interval
19(class SimpleInterval ), which is obtained from inverting the cumulative posterior distribution.
20This calculator works then only for model with a single parameter of interest.
21The model can instead have several nuisance parameters which are integrated (marginalized) in the computation of the posterior function.
22The integration and normalization of the posterior is computed using numerical integration methods provided by ROOT.
23See the MCMCCalculator for model with multiple parameters of interest.
24
25The interface allows one to construct the class by passing the data set, probability density function for the model, the prior
26functions and then the parameter of interest to scan. The nuisance parameters can also be passed to be marginalized when
27computing the posterior. Alternatively, the class can be constructed by passing the data and the ModelConfig containing
28all the needed information (model pdf, prior pdf, parameter of interest, nuisance parameters, etc..)
29
30After configuring the calculator, one only needs to ask GetInterval(), which
31will return an SimpleInterval object. By default the extreme of the integral are obtained by inverting directly the
32cumulative posterior distribution. By using the method SetScanOfPosterior(nbins) the interval is then obtained by
33scanning the posterior function in the given number of points. The first method is in general faster but it requires an
34integration one extra dimension ( in the poi in addition to the nuisance parameters), therefore in some case it can be
35less robust.
36
37The class can also return the posterior function (method GetPosteriorFunction) or if needed the normalized
38posterior function (the posterior pdf) (method GetPosteriorPdf). A posterior plot is also obtained using
39the GetPosteriorPlot method.
40
41The class allows to use different integration methods for integrating in (marginalizing) the nuisances and in the poi. All the numerical
42integration methods of ROOT can be used via the method SetIntegrationType (see more in the documentation of
43this method).
44
45Calculator estimating a credible interval using the Bayesian procedure.
46The calculator computes given the model the posterior distribution and estimates the
47credible interval from the given function.
48*/
49
50
51// include other header files
52
53#include "RooAbsFunc.h"
54#include "RooAbsReal.h"
55#include "RooRealVar.h"
56#include "RooArgSet.h"
57#include "RooFormulaVar.h"
58#include "RooGenericPdf.h"
59#include "RooPlot.h"
60#include "RooProdPdf.h"
61#include "RooDataSet.h"
62
63// include header file of this class
67
68#include "Math/IFunction.h"
70#include "Math/Integrator.h"
71#include "Math/RootFinder.h"
73#include "RooFunctor.h"
74#include "RooFunctor1DBinding.h"
75#include "RooTFnBinding.h"
76#include "RooMsgService.h"
77
78#include "TAxis.h"
79#include "TF1.h"
80#include "TH1.h"
81#include "TMath.h"
82
83#include <map>
84#include <cmath>
85#include <memory>
86
87#include "RConfigure.h"
88
89
90namespace RooStats {
91
92
93// first some utility classes and functions
94
95#ifdef R__HAS_MATHMORE
97#else
99#endif
100
101
102
103
105 LikelihoodFunction(RooFunctor & f, RooFunctor * prior = nullptr, double offset = 0) :
106 fFunc(f), fPrior(prior),
108 {
110 }
111
113
114 double operator() (const double *x ) const {
115 double nll = fFunc(x) - fOffset;
116 double likelihood = std::exp(-nll);
117
118 if (fPrior) likelihood *= (*fPrior)(x);
119
120 int nCalls = fFunc.binding().numCall();
121 if (nCalls > 0 && nCalls % 1000 == 0) {
122 ooccoutD(nullptr,Eval) << "Likelihood evaluation ncalls = " << nCalls
123 << " x0 " << x[0] << " nll = " << nll+fOffset;
124 if (fPrior) ooccoutD(nullptr,Eval) << " prior(x) = " << (*fPrior)(x);
125 ooccoutD(nullptr,Eval) << " likelihood " << likelihood
126 << " max Likelihood " << fMaxL << std::endl;
127 }
128
129 if (likelihood > fMaxL ) {
130 fMaxL = likelihood;
131 if ( likelihood > 1.E10) {
132 ooccoutW(nullptr,Eval) << "LikelihoodFunction::() WARNING - Huge likelihood value found for parameters ";
133 for (int i = 0; i < fFunc.nObs(); ++i)
134 ooccoutW(nullptr,Eval) << " x[" << i << " ] = " << x[i];
135 ooccoutW(nullptr,Eval) << " nll = " << nll << " L = " << likelihood << std::endl;
136 }
137 }
138
139 return likelihood;
140 }
141
142 // for the 1D case
143 double operator() (double x) const {
144 // just call the previous method
145 assert(fFunc.nObs() == 1); // check nobs = 1
146 double tmp = x;
147 return (*this)(&tmp);
148 }
149
150 RooFunctor & fFunc; // functor representing the nll function
151 RooFunctor * fPrior; // functor representing the prior function
152 double fOffset; // offset used to bring the nll in a reasonable range for computing the exponent
153 mutable double fMaxL = 0;
154};
155
156
157// Posterior CDF Function class
158// for integral of posterior function in nuisance and POI
159// 1-Dim function as function of the poi
160
162
163public:
164
165 PosteriorCdfFunction(RooAbsReal & nll, RooArgList & bindParams, RooAbsReal * prior = nullptr, const char * integType = nullptr, double nllMinimum = 0) :
166 fFunctor(nll, bindParams, RooArgList() ), // functor
167 fPriorFunc(nullptr),
168 fLikelihood(fFunctor, nullptr, nllMinimum), // integral of exp(-nll) function
169 fXmin(bindParams.size()), // vector of parameters (min values)
171 {
172 if (prior) {
173 fPriorFunc = std::make_shared<RooFunctor>(*prior, bindParams, RooArgList());
175 }
176
177 // ROOT::Math::IntegratorMultiDim does not support dim = 1, so fall back to
178 // the 1D integrator when there are no nuisance parameters.
179 if (bindParams.size() == 1) {
180 fIntegratorOneDim = std::make_unique<ROOT::Math::Integrator>(ROOT::Math::IntegratorOneDim::GetType(integType));
181 fIntegratorOneDim->SetFunction(fLikelihood);
182 } else {
183 fIntegrator = std::make_unique<ROOT::Math::IntegratorMultiDim>(ROOT::Math::IntegratorMultiDim::GetType(integType));
184 fIntegrator->SetFunction(fLikelihood, bindParams.size());
185 }
186
187 ooccoutD(nullptr,NumericIntegration) << "PosteriorCdfFunction::Compute integral of posterior in nuisance and poi. "
188 << " nllMinimum is " << nllMinimum << std::endl;
189
190 std::vector<double> par(bindParams.size());
191 for (unsigned int i = 0; i < fXmin.size(); ++i) {
192 RooRealVar & var = static_cast<RooRealVar &>( bindParams[i]);
193 fXmin[i] = var.getMin();
194 fXmax[i] = var.getMax();
195 par[i] = var.getVal();
196 ooccoutD(nullptr,NumericIntegration) << "PosteriorFunction::Integrate" << var.GetName()
197 << " in interval [ " << fXmin[i] << " , " << fXmax[i] << " ] " << std::endl;
198 }
199
200 if (fIntegrator)
201 fIntegrator->Options().Print(ooccoutD(nullptr,NumericIntegration));
202 else
203 fIntegratorOneDim->Options().Print(ooccoutD(nullptr,NumericIntegration));
204
205 // store max POI value because it will be changed when evaluating the function
206 fMaxPOI = fXmax[0];
207
208 // compute first the normalization with the poi
209 fNorm = (*this)( fMaxPOI );
210 if (fError) ooccoutE(nullptr,NumericIntegration) << "PosteriorFunction::Error computing normalization - norm = " << fNorm << std::endl;
211 fHasNorm = true;
212 fNormCdfValues.insert(std::make_pair(fXmin[0], 0) );
213 fNormCdfValues.insert(std::make_pair(fXmax[0], 1.0) );
214
215 }
216
217 // copy constructor (needed for Cloning the object)
218 // need special treatment because integrator
219 // has no copy constructor
221 ROOT::Math::IGenFunction(),
223 //fPriorFunc(std::shared_ptr<RooFunctor>((RooFunctor*)0)),
226 fXmin( rhs.fXmin),
227 fXmax( rhs.fXmax),
228 fNorm( rhs.fNorm),
236 {
237 if (rhs.fIntegratorOneDim) {
238 fIntegratorOneDim = std::make_unique<ROOT::Math::Integrator>(ROOT::Math::IntegratorOneDim::GetType(rhs.fIntegratorOneDim->Name().c_str()));
239 fIntegratorOneDim->SetFunction(fLikelihood);
240 } else {
241 fIntegrator = std::make_unique<ROOT::Math::IntegratorMultiDim>(ROOT::Math::IntegratorMultiDim::GetType(rhs.fIntegrator->Name().c_str()));
242 fIntegrator->SetFunction(fLikelihood, fXmin.size());
243 }
244 // need special treatment for the smart pointer
245 // if (rhs.fPriorFunc.get() ) {
246 // fPriorFunc = std::shared_ptr<RooFunctor>(new RooFunctor(*(rhs.fPriorFunc) ) );
247 // fLikelihood.SetPrior( fPriorFunc.get() );
248 // }
249 }
250
251
252 bool HasError() const { return fError; }
253
254
255 ROOT::Math::IGenFunction * Clone() const override {
256 ooccoutD(nullptr,NumericIntegration) << " cloning function .........." << std::endl;
257 return new PosteriorCdfFunction(*this);
258 }
259
260 // offset value for computing the root
261 void SetOffset(double offset) { fOffset = offset; }
262
263private:
264
265 // make assignment operator private
267 return *this;
268 }
269
270 double DoEval (double x) const override {
271
272 // evaluate cdf at poi value x by integrating poi from [xmin,x] and all the nuisances
273 fXmax[0] = x;
274 if (x <= fXmin[0] ) return -fOffset;
275 // could also avoid a function evaluation at maximum
276 if (x >= fMaxPOI && fHasNorm) return 1. - fOffset; // cdf is bound to these values
277
278 // computes the integral using a previous cdf estimate
279 double normcdf0 = 0;
280 if (fHasNorm && fUseOldValues) {
281 // look in the map of the stored cdf values the closes one
282 std::map<double,double>::iterator itr = fNormCdfValues.upper_bound(x);
283 --itr; // upper bound returns a position 1 up of the value we want
284 if (itr != fNormCdfValues.end() ) {
285 fXmin[0] = itr->first;
286 normcdf0 = itr->second;
287 // ooccoutD(nullptr,NumericIntegration) << "PosteriorCdfFunction: computing integral between in poi interval : "
288 // << fXmin[0] << " - " << fXmax[0] << std::endl;
289 }
290 }
291
292 fFunctor.binding().resetNumCall(); // reset number of calls for debug
293
294 double cdf;
295 double error;
296 if (fIntegratorOneDim) {
297 cdf = fIntegratorOneDim->Integral(fXmin[0], fXmax[0]);
298 error = fIntegratorOneDim->Error();
299 } else {
300 cdf = fIntegrator->Integral(&fXmin[0], &fXmax[0]);
301 error = fIntegrator->Error();
302 }
303 double normcdf = cdf/fNorm; // normalize the cdf
304
305 ooccoutD(nullptr,NumericIntegration) << "PosteriorCdfFunction: poi = [" << fXmin[0] << " , "
306 << fXmax[0] << "] integral = " << cdf << " +/- " << error
307 << " norm-integ = " << normcdf << " cdf(x) = " << normcdf+normcdf0
308 << " ncalls = " << fFunctor.binding().numCall() << std::endl;
309
310 if (TMath::IsNaN(cdf) || cdf > std::numeric_limits<double>::max()) {
311 ooccoutE(nullptr,NumericIntegration) << "PosteriorFunction::Error computing integral - cdf = "
312 << cdf << std::endl;
313 fError = true;
314 }
315
316 if (cdf != 0 && error / cdf > 0.2) {
317 oocoutW(nullptr, NumericIntegration)
318 << "PosteriorCdfFunction: integration error is larger than 20 % x0 = " << fXmin[0] << " x = " << x
319 << " cdf(x) = " << cdf << " +/- " << error << std::endl;
320 }
321
322 if (!fHasNorm) {
323 oocoutI(nullptr,NumericIntegration) << "PosteriorCdfFunction - integral of posterior = "
324 << cdf << " +/- " << error << std::endl;
325 fNormErr = error;
326 return cdf;
327 }
328
329 normcdf += normcdf0;
330
331 // store values in the map
332 if (fUseOldValues) {
333 fNormCdfValues.insert(std::make_pair(x, normcdf) );
334 }
335
336 double errnorm = sqrt( error*error + normcdf*normcdf * fNormErr * fNormErr )/fNorm;
337 if (normcdf > 1. + 3 * errnorm) {
338 oocoutW(nullptr,NumericIntegration) << "PosteriorCdfFunction: normalized cdf values is larger than 1"
339 << " x = " << x << " normcdf(x) = " << normcdf << " +/- " << error/fNorm << std::endl;
340 }
341
342 return normcdf - fOffset; // apply an offset (for finding the roots)
343 }
344
345 mutable RooFunctor fFunctor; // functor binding nll
346 mutable std::shared_ptr<RooFunctor> fPriorFunc; // functor binding the prior
347 LikelihoodFunction fLikelihood; // likelihood function
348 mutable std::unique_ptr<ROOT::Math::IntegratorMultiDim> fIntegrator; // multi-dim integrator (for >=2 dims)
349 mutable std::unique_ptr<ROOT::Math::Integrator> fIntegratorOneDim; // 1D integrator (for POI-only case)
350 mutable std::vector<double> fXmin; // min value of parameters (poi+nuis) -
351 mutable std::vector<double> fXmax; // max value of parameters (poi+nuis) - max poi changes so it is mutable
352 double fNorm = 1.0; // normalization value (computed in constructor)
353 mutable double fNormErr = 0.0; // normalization error value (computed in constructor)
354 double fOffset = 0; // offset for computing the root
355 double fMaxPOI = 0; // maximum value of POI
356 bool fHasNorm = false; // flag to control first call to the function
357 bool fUseOldValues = true; // use old cdf values
358 mutable bool fError = false; // flag to indicate if a numerical evaluation error occurred
359 mutable std::map<double,double> fNormCdfValues;
360};
361
362//__________________________________________________________________
363// Posterior Function class
364// 1-Dim function as function of the poi
365// and it integrated all the nuisance parameters
366
368
369public:
370
371
372 PosteriorFunction(RooAbsReal & nll, RooRealVar & poi, RooArgList & nuisParams, RooAbsReal * prior = nullptr, const char * integType = nullptr, double
373 norm = 1.0, double nllOffset = 0, int niter = 0) :
375 fPriorFunc(nullptr),
376 fLikelihood(fFunctor, nullptr, nllOffset),
377 fPoi(&poi),
380 fNorm(norm)
381 {
382
383 if (prior) {
384 fPriorFunc = std::make_shared<RooFunctor>(*prior, nuisParams, RooArgList());
386 }
387
388 ooccoutD(nullptr,NumericIntegration) << "PosteriorFunction::Evaluate the posterior function by integrating the nuisances: " << std::endl;
389 for (unsigned int i = 0; i < fXmin.size(); ++i) {
390 RooRealVar & var = static_cast<RooRealVar &>( nuisParams[i]);
391 fXmin[i] = var.getMin();
392 fXmax[i] = var.getMax();
393 ooccoutD(nullptr,NumericIntegration) << "PosteriorFunction::Integrate " << var.GetName()
394 << " in interval [" << fXmin[i] << " , " << fXmax[i] << " ] " << std::endl;
395 }
396 if (fXmin.size() == 1) { // 1D case
397 fIntegratorOneDim = std::make_unique<ROOT::Math::Integrator>( ROOT::Math::IntegratorOneDim::GetType(integType) );
398
399 fIntegratorOneDim->SetFunction(fLikelihood);
400 // interested only in relative tolerance
401 //fIntegratorOneDim->SetAbsTolerance(1.E-300);
402 fIntegratorOneDim->Options().Print(ooccoutD(nullptr,NumericIntegration) );
403 }
404 else if (fXmin.size() > 1) { // multiDim case
405 fIntegratorMultiDim = std::make_unique<ROOT::Math::IntegratorMultiDim>(ROOT::Math::IntegratorMultiDim::GetType(integType));
406 fIntegratorMultiDim->SetFunction(fLikelihood, fXmin.size());
408 if (niter > 0) {
409 opt.SetNCalls(niter);
410 fIntegratorMultiDim->SetOptions(opt);
411 }
412 //fIntegratorMultiDim->SetAbsTolerance(1.E-300);
413 // print the options
414 opt.Print(ooccoutD(nullptr,NumericIntegration) );
415 }
416 }
417
418
419 ROOT::Math::IGenFunction * Clone() const override {
420 assert(1);
421 return nullptr; // cannot clone this function for integrator
422 }
423
424 double Error() const { return fError;}
425
426
427private:
428 double DoEval (double x) const override {
429
430 // evaluate posterior function at a poi value x by integrating all nuisance parameters
431
432 fPoi->setVal(x);
433 fFunctor.binding().resetNumCall(); // reset number of calls for debug
434
435 double f = 0;
436 double error = 0;
437 if (fXmin.size() == 1) { // 1D case
438 f = fIntegratorOneDim->Integral(fXmin[0],fXmax[0]);
439 error = fIntegratorOneDim->Error();
440 }
441 else if (fXmin.size() > 1) { // multi-dim case
442 f = fIntegratorMultiDim->Integral(&fXmin[0],&fXmax[0]);
443 error = fIntegratorMultiDim->Error();
444 } else {
445 // no integration to be done
446 f = fLikelihood(x);
447 }
448
449 // debug
450 ooccoutD(nullptr,NumericIntegration) << "PosteriorFunction: POI value = "
451 << x << "\tf(x) = " << f << " +/- " << error
452 << " norm-f(x) = " << f/fNorm
453 << " ncalls = " << fFunctor.binding().numCall() << std::endl;
454
455 if (f != 0 && error / f > 0.2) {
456 ooccoutW(nullptr, NumericIntegration)
457 << "PosteriorFunction::DoEval - Error from integration in " << fXmin.size() << " Dim is larger than 20 % "
458 << "x = " << x << " p(x) = " << f << " +/- " << error << std::endl;
459 }
460
461 fError = error / fNorm;
462 return f / fNorm;
463 }
464
466 mutable std::shared_ptr<RooFunctor> fPriorFunc; // functor binding the prior
469 std::unique_ptr<ROOT::Math::Integrator> fIntegratorOneDim;
470 std::unique_ptr<ROOT::Math::IntegratorMultiDim> fIntegratorMultiDim;
471 std::vector<double> fXmin;
472 std::vector<double> fXmax;
473 double fNorm;
474 mutable double fError = 0;
475};
476
477////////////////////////////////////////////////////////////////////////////////
478/// Posterior function obtaining sampling toy MC for the nuisance according to their pdf
479
481
482public:
484 RooAbsReal *prior = nullptr, double nllOffset = 0, int niter = 0, bool redoToys = true)
485 : fFunctor(nll, nuisParams, RooArgList()),
486 fPriorFunc(nullptr),
487 fLikelihood(fFunctor, nullptr, nllOffset),
488 fPdf(&pdf),
489 fPoi(&poi),
492
494 {
495 if (niter == 0) fNumIterations = 100; // default value
496
497 if (prior) {
498 fPriorFunc = std::make_shared<RooFunctor>(*prior, nuisParams, RooArgList());
500 }
501
502 ooccoutI(nullptr,InputArguments) << "PosteriorFunctionFromToyMC::Evaluate the posterior function by randomizing the nuisances: niter " << fNumIterations << std::endl;
503
504 ooccoutI(nullptr,InputArguments) << "PosteriorFunctionFromToyMC::Pdf used for randomizing the nuisance is " << fPdf->GetName() << std::endl;
505 // check that pdf contains the nuisance
506 std::unique_ptr<RooArgSet> vars{fPdf->getVariables()};
507 for (std::size_t i = 0; i < fNuisParams.size(); ++i) {
508 if (!vars->find( fNuisParams[i].GetName() ) ) {
509 ooccoutW(nullptr,InputArguments) << "Nuisance parameter " << fNuisParams[i].GetName()
510 << " is not part of sampling pdf. "
511 << "they will be treated as constant " << std::endl;
512 }
513 }
514
515 if (!fRedoToys) {
516 ooccoutI(nullptr,InputArguments) << "PosteriorFunctionFromToyMC::Generate nuisance toys only one time (for all POI points)" << std::endl;
517 GenerateToys();
518 }
519 }
520
521 // generate first n-samples of the nuisance parameters
522 void GenerateToys() const {
523 fGenParams = std::unique_ptr<RooDataSet>{fPdf->generate(fNuisParams, fNumIterations)};
524 if(fGenParams==nullptr) {
525 ooccoutE(nullptr,InputArguments) << "PosteriorFunctionFromToyMC - failed to generate nuisance parameters" << std::endl;
526 }
527 }
528
529 double Error() const { return fError;}
530
531 ROOT::Math::IGenFunction * Clone() const override {
532 // use default copy constructor
533 //return new PosteriorFunctionFromToyMC(*this);
534 // clone not implemented
535 assert(1);
536 return nullptr;
537 }
538
539private:
540 // evaluate the posterior at the poi value x
541 double DoEval( double x) const override {
542
543 int npar = fNuisParams.size();
544 assert (npar > 0);
545
546
547 // generate the toys
548 if (fRedoToys) GenerateToys();
549 if (!fGenParams) return 0;
550
551 // evaluate posterior function at a poi value x by integrating all nuisance parameters
552
553 fPoi->setVal(x);
554
555 // loop over all of the generate data
556 double sum = 0;
557 double sum2 = 0;
558
559 for(int iter=0; iter<fNumIterations; ++iter) {
560
561 // get the set of generated parameters and set the nuisance parameters to the generated values
562 std::vector<double> p(npar);
563 for (int i = 0; i < npar; ++i) {
564 const RooArgSet* genset=fGenParams->get(iter);
565 RooAbsArg * arg = genset->find( fNuisParams[i].GetName() );
566 RooRealVar * var = dynamic_cast<RooRealVar*>(arg);
567 assert(var != nullptr);
568 p[i] = var->getVal();
569 (static_cast<RooRealVar &>( fNuisParams[i])).setVal(p[i]);
570 }
571
572 // evaluate now the likelihood function
573 double fval = fLikelihood( &p.front() );
574
575 // likelihood already must contained the pdf we have sampled
576 // so we must divided by it. The value must be normalized on all
577 // other parameters
579 double nuisPdfVal = fPdf->getVal(&arg);
580 fval /= nuisPdfVal;
581
582
583 if( fval > std::numeric_limits<double>::max() ) {
584 ooccoutE(nullptr,Eval) << "BayesianCalculator::EvalPosteriorFunctionFromToy : "
585 << "Likelihood evaluates to infinity " << std::endl;
586 ooccoutE(nullptr,Eval) << "poi value = " << x << std::endl;
587 ooccoutE(nullptr,Eval) << "Nuisance parameter values : ";
588 for (int i = 0; i < npar; ++i)
589 ooccoutE(nullptr,Eval) << fNuisParams[i].GetName() << " = " << p[i] << " ";
590 ooccoutE(nullptr,Eval) << " - return 0 " << std::endl;
591
592 fError = 1.E30;
593 return 0;
594 }
595 if( TMath::IsNaN(fval) ) {
596 ooccoutE(nullptr,Eval) << "BayesianCalculator::EvalPosteriorFunctionFromToy : "
597 << "Likelihood is a NaN " << std::endl;
598 ooccoutE(nullptr,Eval) << "poi value = " << x << std::endl;
599 ooccoutE(nullptr,Eval) << "Nuisance parameter values : ";
600 for (int i = 0; i < npar; ++i)
601 ooccoutE(nullptr,Eval) << fNuisParams[i].GetName() << " = " << p[i] << " ";
602 ooccoutE(nullptr,Eval) << " - return 0 " << std::endl;
603 fError = 1.E30;
604 return 0;
605 }
606
607
608
609 sum += fval;
610 sum2 += fval*fval;
611 }
612
613 // compute the average and variance
614 double val = sum/double(fNumIterations);
615 double dval2 = std::max( sum2/double(fNumIterations) - val*val, 0.0);
616 fError = std::sqrt( dval2 / fNumIterations);
617
618 // debug
619 ooccoutD(nullptr,NumericIntegration) << "PosteriorFunctionFromToyMC: POI value = "
620 << x << "\tp(x) = " << val << " +/- " << fError << std::endl;
621
622
623 if (val != 0 && fError/val > 0.2 ) {
624 ooccoutW(nullptr,NumericIntegration) << "PosteriorFunctionFromToyMC::DoEval"
625 << " - Error in estimating posterior is larger than 20% ! "
626 << "x = " << x << " p(x) = " << val << " +/- " << fError << std::endl;
627 }
628
629
630 return val;
631 }
632
634 mutable std::shared_ptr<RooFunctor> fPriorFunc; // functor binding the prior
636 mutable RooAbsPdf * fPdf;
639 mutable std::unique_ptr<RooDataSet> fGenParams;
641 mutable double fError = -1;
642 bool fRedoToys; // do toys every iteration
643
644};
645
646////////////////////////////////////////////////////////////////////////////////
647// Implementation of BayesianCalculator
648////////////////////////////////////////////////////////////////////////////////
649
650////////////////////////////////////////////////////////////////////////////////
651/// default constructor
652
654 fData(nullptr),
655 fPdf(nullptr),
656 fPriorPdf(nullptr),
657 fNuisancePdf(nullptr),
658 fProductPdf (nullptr), fLikelihood (nullptr), fIntegratedLikelihood (nullptr), fPosteriorPdf(nullptr),
659 fPosteriorFunction(nullptr), fApproxPosterior(nullptr),
660 fLower(0), fUpper(0),
661 fNLLMin(0),
662 fSize(0.05), fLeftSideFraction(0.5),
663 fBrfPrecision(0.00005),
664 fNScanBins(-1),
665 fNumIterations(0),
666 fValidInterval(false)
667{
668}
669
670////////////////////////////////////////////////////////////////////////////////
671/// Constructor from data set, model pdf, parameter of interests and prior pdf
672/// If nuisance parameters are given they will be integrated according either to the prior or
673/// their constraint term included in the model
674
675BayesianCalculator::BayesianCalculator( /* const char* name, const char* title, */
677 RooAbsPdf& pdf,
678 const RooArgSet& POI,
681 fData(&data),
682 fPdf(&pdf),
683 fPOI(POI),
684 fPriorPdf(&priorPdf),
685 fNuisancePdf(nullptr),
686 fProductPdf (nullptr), fLikelihood (nullptr), fIntegratedLikelihood (nullptr), fPosteriorPdf(nullptr),
687 fPosteriorFunction(nullptr), fApproxPosterior(nullptr),
688 fLower(0), fUpper(0),
689 fNLLMin(0),
690 fSize(0.05), fLeftSideFraction(0.5),
691 fBrfPrecision(0.00005),
692 fNScanBins(-1),
693 fNumIterations(0),
694 fValidInterval(false)
695{
696
698 // remove constant nuisance parameters
700}
701
702////////////////////////////////////////////////////////////////////////////////
703/// Constructor from a data set and a ModelConfig
704/// model pdf, poi and nuisances will be taken from the ModelConfig
705
707 ModelConfig & model) :
708 fData(&data),
709 fPdf(model.GetPdf()),
710 fPriorPdf( model.GetPriorPdf()),
711 fNuisancePdf(nullptr),
712 fProductPdf (nullptr), fLikelihood (nullptr), fIntegratedLikelihood (nullptr), fPosteriorPdf(nullptr),
713 fPosteriorFunction(nullptr), fApproxPosterior(nullptr),
714 fLower(0), fUpper(0),
715 fNLLMin(0),
716 fSize(0.05), fLeftSideFraction(0.5),
717 fBrfPrecision(0.00005),
718 fNScanBins(-1),
719 fNumIterations(0),
720 fValidInterval(false)
721{
722 SetModel(model);
723}
724
725
727{
728 // destructor
729 ClearAll();
730}
731
732////////////////////////////////////////////////////////////////////////////////
733/// clear all cached pdf objects
734
736 if (fProductPdf) delete fProductPdf;
737 fLogLike.reset();
738 if (fLikelihood) delete fLikelihood;
740 if (fPosteriorPdf) delete fPosteriorPdf;
743 fPosteriorPdf = nullptr;
744 fPosteriorFunction = nullptr;
745 fProductPdf = nullptr;
746 fLikelihood = nullptr;
747 fIntegratedLikelihood = nullptr;
748 fLower = 0;
749 fUpper = 0;
750 fNLLMin = 0;
751 fValidInterval = false;
752}
753
754////////////////////////////////////////////////////////////////////////////////
755/// set the model to use
756/// The model pdf, prior pdf, parameter of interest and nuisances
757/// will be taken according to the model
758
760
761 fPdf = model.GetPdf();
762 fPriorPdf = model.GetPriorPdf();
763 // assignment operator = does not do a real copy the sets (must use add method)
764 fPOI.removeAll();
768 if (model.GetParametersOfInterest()) fPOI.add( *(model.GetParametersOfInterest()) );
771 if (model.GetGlobalObservables()) fGlobalObs.add( *(model.GetGlobalObservables() ) );
772 // remove constant nuisance parameters
774
775 // invalidate the cached pointers
776 ClearAll();
777}
778
779////////////////////////////////////////////////////////////////////////////////
780/// Build and return the posterior function (not normalized) as a RooAbsReal
781/// the posterior is obtained from the product of the likelihood function and the
782/// prior pdf which is then integrated in the nuisance parameters (if existing).
783/// A prior function for the nuisance can be specified either in the prior pdf object
784/// or in the model itself. If no prior nuisance is specified, but prior parameters are then
785/// the integration is performed assuming a flat prior for the nuisance parameters.
786///
787/// NOTE: the return object is managed by the BayesianCalculator class, users do not need to delete it,
788/// but the object will be deleted when the BayesiabCalculator object is deleted
789
791{
792
794 if (fLikelihood) return fLikelihood;
795
796 // run some sanity checks
797 if (!fPdf ) {
798 coutE(InputArguments) << "BayesianCalculator::GetPosteriorPdf - missing pdf model" << std::endl;
799 return nullptr;
800 }
801 if (fPOI.empty()) {
802 coutE(InputArguments) << "BayesianCalculator::GetPosteriorPdf - missing parameter of interest" << std::endl;
803 return nullptr;
804 }
805 if (fPOI.size() > 1) {
806 coutE(InputArguments) << "BayesianCalculator::GetPosteriorPdf - current implementation works only on 1D intervals" << std::endl;
807 return nullptr;
808 }
809
810
811 std::unique_ptr<RooArgSet> constrainedParams{fPdf->getParameters(*fData)};
812 // remove the constant parameters
814
815 //constrainedParams->Print("V");
816
817 // use RooFit::Constrain() to be sure constraints terms are taken into account
819
820
821
822 ccoutD(Eval) << "BayesianCalculator::GetPosteriorFunction : "
823 << " pdf value " << fPdf->getVal()
824 << " neglogLikelihood = " << fLogLike->getVal() << std::endl;
825
826 if (fPriorPdf)
827 ccoutD(Eval) << "\t\t\t priorPOI value " << fPriorPdf->getVal() << std::endl;
828
829 // check that likelihood evaluation is not infinity
830 double nllVal = fLogLike->getVal();
831 if ( nllVal > std::numeric_limits<double>::max() ) {
832 coutE(Eval) << "BayesianCalculator::GetPosteriorFunction : "
833 << " Negative log likelihood evaluates to infinity " << std::endl
834 << " Non-const Parameter values : ";
836 for (std::size_t i = 0; i < p.size(); ++i) {
837 RooRealVar * v = dynamic_cast<RooRealVar *>(&p[i] );
838 if (v!=nullptr) ccoutE(Eval) << v->GetName() << " = " << v->getVal() << " ";
839 }
840 ccoutE(Eval) << std::endl;
841 ccoutE(Eval) << "-- Perform a full likelihood fit of the model before or set more reasonable parameter values"
842 << std::endl;
843 coutE(Eval) << "BayesianCalculator::GetPosteriorFunction : " << " cannot compute posterior function " << std::endl;
844 return nullptr;
845 }
846
847
848
849 // need do find minimum of log-likelihood in the range to shift function
850 // to avoid numerical errors when we compute the likelihood (overflows in the exponent)
851 // N.B.: this works for only 1 parameter of interest otherwise Minuit should be used for finding the minimum
852 RooFunctor * nllFunc = fLogLike->functor(fPOI);
855 RooRealVar* poi = dynamic_cast<RooRealVar*>(fPOI.first() );
856 assert(poi);
857
858 // try to reduce some error messages
859 //bool silentMode = (RooMsgService::instance().globalKillBelow() >= RooFit::ERROR || RooMsgService::instance().silentMode()) ;
861
862
863 coutI(Eval) << "BayesianCalculator::GetPosteriorFunction : "
864 << " nll value " << nllVal << " poi value = " << poi->getVal() << std::endl;
865
866
868 minim.SetFunction(wnllFunc,poi->getMin(),poi->getMax() );
869 bool ret = minim.Minimize(100,1.E-3,1.E-3);
870 fNLLMin = 0;
871 if (ret) fNLLMin = minim.FValMinimum();
872
873 coutI(Eval) << "BayesianCalculator::GetPosteriorFunction : minimum of NLL vs POI for POI = "
874 << poi->getVal() << " min NLL = " << fNLLMin << std::endl;
875
876 delete nllFunc;
877
878
879 if ( fNuisanceParameters.empty() || fIntegrationType.Contains("ROOFIT") ) {
880
881 ccoutD(Eval) << "BayesianCalculator::GetPosteriorFunction : use ROOFIT integration "
882 << std::endl;
883
884#ifdef DOLATER // (not clear why this does not work)
885 // need to make in this case a likelihood from the nll and make the product with the prior
886 std::string likeName = std::string{"likelihood_times_prior_"} + fPriorPdf->GetName();
887 std::stringstream formula;
888 formula << "std::exp(-@0+" << fNllMin << "+log(@1))";
889 fLikelihood = new RooFormulaVar(likeName.c_str(), formula, RooArgList(*fLogLike, *fPriorPdf));
890#else
891 // here use RooProdPdf (not very nice) but working
892
893 if (fLogLike) fLogLike.reset();
894 if (fProductPdf) {
895 delete fProductPdf;
896 fProductPdf = nullptr;
897 }
898
899 // // create a unique name for the product pdf
901 if (fPriorPdf) {
902 std::string prodName = std::string{"product_"} + fPdf->GetName() + "_" + fPriorPdf->GetName();
903 // save this as data member since it needs to be deleted afterwards
904 fProductPdf = new RooProdPdf(prodName.c_str(), "", RooArgList(*fPdf, *fPriorPdf));
906 }
907
908 std::unique_ptr<RooArgSet> constrParams{fPdf->getParameters(*fData)};
909 // remove the constant parameters
912
913 std::string likeName = std::string{"likelihood_times_prior_"} + pdfAndPrior->GetName();
914 std::stringstream formula;
915 formula << "exp(-@0+" << fNLLMin << ")";
916 fLikelihood = new RooFormulaVar(likeName.c_str(), formula.str().c_str(), RooArgList(*fLogLike));
917#endif
918
919
920 // if no nuisance parameter we can just return the likelihood function
923 fLikelihood = nullptr;
924 } else {
925 // case of using RooFit for the integration
927 std::unique_ptr<RooAbsReal>{fLikelihood->createIntegral(fNuisanceParameters)}.release();
928 }
929
930 }
931
932 else if ( fIntegrationType.Contains("TOYMC") ) {
933 // compute the posterior as expectation values of the likelihood function
934 // sampling on the nuisance parameters
935
937
938 bool doToysEveryIteration = true;
939 // if type is 1-TOYMC or TOYMC-1
941
943 if (!fNuisancePdf) {
944 ccoutI(Eval) << "BayesianCalculator::GetPosteriorFunction : no nuisance pdf is provided, try using global pdf (this will be slower)"
945 << std::endl;
946 }
949
950 TString name = "toyposteriorfunction_from_";
951 name += fLogLike->GetName();
953
954 // need to scan likelihood in this case
955 if (fNScanBins <= 0) fNScanBins = 100;
956
957 }
958
959 else {
960
961 // use ROOT integration method if there are nuisance parameters
962
965
966 TString name = "posteriorfunction_from_";
967 name += fLogLike->GetName();
969
970 }
971
972 if (RooAbsReal::numEvalErrors() > 0) {
973 coutW(Eval) << "BayesianCalculator::GetPosteriorFunction : " << RooAbsReal::numEvalErrors()
974 << " errors reported in evaluating log-likelihood function " << std::endl;
975 }
978
980
981}
982
983////////////////////////////////////////////////////////////////////////////////
984/// Build and return the posterior pdf (i.e posterior function normalized to all range of poi)
985/// Note that an extra integration in the POI is required for the normalization
986/// NOTE: user must delete the returned object
987
989{
990
992 if (!plike) return nullptr;
993
994
995 // create a unique name on the posterior from the names of the components
996 TString posteriorName = this->GetName() + TString("_posteriorPdf_") + plike->GetName();
997
999
1000 return posteriorPdf;
1001}
1002
1003////////////////////////////////////////////////////////////////////////////////
1004/// When am approximate posterior is computed binninig the parameter of interest (poi) range
1005/// (see SetScanOfPosteriors) an histogram is created and can be returned to the user
1006/// A nullptr is instead returned when the posterior is computed without binning the poi.
1007///
1008/// NOTE: the returned object is managed by the BayesianCalculator class,
1009/// if the user wants to take ownership of the returned histogram, he needs to clone
1010/// or copy the return object.
1011
1016
1017
1018////////////////////////////////////////////////////////////////////////////////
1019/// return a RooPlot with the posterior and the credibility region
1020/// NOTE: User takes ownership of the returned object
1021
1023{
1024
1026
1027 // if a scan is requested approximate the posterior
1028 if (fNScanBins > 0)
1030
1032 if (norm) {
1033 // delete and re-do always posterior pdf (could be invalid after approximating it)
1034 if (fPosteriorPdf) delete fPosteriorPdf;
1037 }
1038 if (!posterior) return nullptr;
1039
1041
1042 RooAbsRealLValue* poi = dynamic_cast<RooAbsRealLValue*>(fPOI.first() );
1043 assert(poi);
1044
1045
1046 RooPlot* plot = poi->frame();
1047 if (!plot) return nullptr;
1048
1049 // try to reduce some error messages
1051
1052 plot->SetTitle(TString("Posterior probability of parameter \"")+TString(poi->GetName())+TString("\""));
1054 posterior->plotOn(plot);
1055 plot->GetYaxis()->SetTitle("posterior function");
1056
1057 // reset the counts and default mode
1060
1061 return plot;
1062}
1063
1064////////////////////////////////////////////////////////////////////////////////
1065/// set the integration type (possible type are) :
1066///
1067/// - 1D integration ( used when only one nuisance and when the posterior is scanned):
1068/// adaptive , gauss, nonadaptive
1069/// - multidim:
1070/// - ADAPTIVE, adaptive numerical integration
1071/// The parameter numIters (settable with SetNumIters) is the max number of function calls.
1072/// It can be reduced to make the integration faster but it will be difficult to reach the required tolerance
1073/// - VEGAS MC integration method based on importance sampling - numIters is number of function calls
1074/// Extra Vegas parameter can be set using IntegratorMultiDimOptions class
1075/// - MISER MC integration method based on stratified sampling
1076/// See also http://en.wikipedia.org/wiki/Monte_Carlo_integration for VEGAS and MISER description
1077/// - PLAIN simple MC integration method, where the max number of calls can be specified using SetNumIters(numIters)
1078///
1079/// Extra integration types are:
1080///
1081/// - TOYMC:
1082/// evaluate posterior by generating toy MC for the nuisance parameters. It is a MC
1083/// integration, where the function is sampled according to the nuisance. It is convenient to use when all
1084/// the nuisance are uncorrelated and it is efficient to generate them
1085/// The toy are generated by default for each poi values
1086/// (this method has been proposed and provided by J.P Chou)
1087/// - 1-TOYMC : same method as before but in this case the toys are generated only one time and then used for
1088/// each poi value. It can be convenient when the generation time is much larger than the evaluation time,
1089/// otherwise it is recommended to re-generate the toy for each poi scanned point of the posterior function
1090/// - ROOFIT:
1091/// use roofit default integration methods which will produce a nested integral (not recommended for more
1092/// than 1 nuisance parameters)
1093
1095 // if type = 0 use default specified via class IntegratorMultiDimOptions::SetDefaultIntegrator
1098}
1099
1100////////////////////////////////////////////////////////////////////////////////
1101/// Compute the interval. By Default a central interval is computed
1102/// and the result is a SimpleInterval object.
1103///
1104/// Using the method (to be called before SetInterval) SetLeftSideTailFraction the user can choose the type of interval.
1105/// By default the returned interval is a central interval with the confidence level specified
1106/// previously in the constructor ( LeftSideTailFraction = 0.5).
1107/// - For lower limit use SetLeftSideTailFraction = 1
1108/// - For upper limit use SetLeftSideTailFraction = 0
1109/// - for shortest intervals use SetLeftSideTailFraction = -1 or call the method SetShortestInterval()
1110///
1111/// NOTE: The BayesianCalculator covers only the case with one
1112/// single parameter of interest
1113///
1114/// NOTE: User takes ownership of the returned object
1115
1117{
1118
1119 if (fValidInterval)
1120 coutW(Eval) << "BayesianCalculator::GetInterval - recomputing interval for the same CL and same model" << std::endl;
1121
1122 RooRealVar* poi = dynamic_cast<RooRealVar*>(fPOI.first() );
1123 if (!poi) {
1124 coutE(Eval) << "BayesianCalculator::GetInterval - no parameter of interest is set " << std::endl;
1125 return nullptr;
1126 }
1127
1128
1129
1130 // get integrated likelihood (posterior function)
1132
1133 //bool silentMode = (RooMsgService::instance().globalKillBelow() >= RooFit::ERROR || RooMsgService::instance().silentMode()) ;
1135
1136 if (fLeftSideFraction < 0 ) {
1137 // compute short intervals
1139 }
1140 else {
1141 // compute the other intervals
1142
1144 double upperCutOff = 1. - (1.- fLeftSideFraction) * fSize;
1145
1146
1147 if (fNScanBins > 0) {
1149 }
1150
1151 else {
1153 // case cdf failed (scan then the posterior)
1154 if (!fValidInterval) {
1155 fNScanBins = 100;
1156 coutW(Eval) << "BayesianCalculator::GetInterval - computing integral from cdf failed - do a scan in "
1157 << fNScanBins << " nbins " << std::endl;
1159 }
1160 }
1161 }
1162
1163
1164 // reset the counts and default mode
1165 if (RooAbsReal::numEvalErrors() > 0) {
1166 coutW(Eval) << "BayesianCalculator::GetInterval : " << RooAbsReal::numEvalErrors()
1167 << " errors reported in evaluating log-likelihood function " << std::endl;
1168 }
1169
1172
1173 if (!fValidInterval) {
1174 fLower = 1; fUpper = 0;
1175 coutE(Eval) << "BayesianCalculator::GetInterval - cannot compute a valid interval - return a dummy [1,0] interval"
1176 << std::endl;
1177 }
1178 else {
1179 coutI(Eval) << "BayesianCalculator::GetInterval - found a valid interval : [" << fLower << " , "
1180 << fUpper << " ]" << std::endl;
1181 }
1182
1183 TString interval_name = TString("BayesianInterval_a") + TString(this->GetName());
1185 interval->SetTitle("SimpleInterval from BayesianCalculator");
1186
1187 return interval;
1188}
1189
1190////////////////////////////////////////////////////////////////////////////////
1191/// Returns the value of the parameter for the point in
1192/// parameter-space that is the most likely.
1193/// How do we do if there are points that are equi-probable?
1194/// use approximate posterior
1195/// t.b.d use real function to find the mode
1196
1198
1201 return h->GetBinCenter(h->GetMaximumBin() );
1202 //return fApproxPosterior->GetMaximumX();
1203}
1204
1205////////////////////////////////////////////////////////////////////////////////
1206/// internal function compute the interval using Cdf integration
1207
1209
1210 fValidInterval = false;
1211
1212 coutI(InputArguments) << "BayesianCalculator:GetInterval Compute the interval from the posterior cdf " << std::endl;
1213
1214 RooRealVar* poi = dynamic_cast<RooRealVar*>(fPOI.first() );
1215 assert(poi);
1216 if (GetPosteriorFunction() == nullptr) {
1217 coutE(InputArguments) << "BayesianCalculator::GetInterval() cannot make posterior Function " << std::endl;
1218 return;
1219 }
1220
1221 // need to remove the constant parameters
1223 bindParams.add(fPOI);
1225
1226 // this code could be put inside the PosteriorCdfFunction
1227
1228 //bindParams.Print("V");
1229
1231 if( cdf.HasError() ) {
1232 coutE(Eval) << "BayesianCalculator: Numerical error computing CDF integral - try a different method " << std::endl;
1233 return;
1234 }
1235
1236 //find the roots
1237
1239
1240 ccoutD(Eval) << "BayesianCalculator::GetInterval - finding roots of posterior using RF " << rf.Name()
1241 << " with precision = " << fBrfPrecision;
1242
1243 if (lowerCutOff > 0) {
1244 cdf.SetOffset(lowerCutOff);
1245 ccoutD(NumericIntegration) << "Integrating posterior to get cdf and search lower limit at p =" << lowerCutOff << std::endl;
1246 bool ok = rf.Solve(cdf, poi->getMin(),poi->getMax() , 200,fBrfPrecision, fBrfPrecision);
1247 if( cdf.HasError() )
1248 coutW(Eval) << "BayesianCalculator: Numerical error integrating the CDF " << std::endl;
1249 if (!ok) {
1250 coutE(NumericIntegration) << "BayesianCalculator::GetInterval - Error from root finder when searching lower limit !" << std::endl;
1251 return;
1252 }
1253 fLower = rf.Root();
1254 }
1255 else {
1256 fLower = poi->getMin();
1257 }
1258 if (upperCutOff < 1.0) {
1259 cdf.SetOffset(upperCutOff);
1260 ccoutD(NumericIntegration) << "Integrating posterior to get cdf and search upper interval limit at p =" << upperCutOff << std::endl;
1261 bool ok = rf.Solve(cdf, fLower,poi->getMax() , 200, fBrfPrecision, fBrfPrecision);
1262 if( cdf.HasError() )
1263 coutW(Eval) << "BayesianCalculator: Numerical error integrating the CDF " << std::endl;
1264 if (!ok) {
1265 coutE(NumericIntegration) << "BayesianCalculator::GetInterval - Error from root finder when searching upper limit !" << std::endl;
1266 return;
1267 }
1268 fUpper = rf.Root();
1269 }
1270 else {
1271 fUpper = poi->getMax();
1272 }
1273
1274 fValidInterval = true;
1275}
1276
1277////////////////////////////////////////////////////////////////////////////////
1278/// approximate posterior in nbins using a TF1
1279/// scan the poi values and evaluate the posterior at each point
1280/// and save the result in a cloned TF1
1281/// For each point the posterior is evaluated by integrating the nuisance
1282/// parameters
1283
1285
1286 if (fApproxPosterior) {
1287 // if number of bins of existing function is >= requested one - no need to redo the scan
1288 if (fApproxPosterior->GetNpx() >= fNScanBins) return;
1289 // otherwise redo the scan
1290 delete fApproxPosterior;
1291 fApproxPosterior = nullptr;
1292 }
1293
1294
1296 if (!posterior) return;
1297
1298
1299 TF1 * tmp = posterior->asTF(fPOI);
1300 assert(tmp != nullptr);
1301 // binned the function in nbins and evaluate at those points
1302 if (fNScanBins > 0) tmp->SetNpx(fNScanBins); // if not use default of TF1 (which is 100)
1303
1304 coutI(Eval) << "BayesianCalculator - scan posterior function in nbins = " << tmp->GetNpx() << std::endl;
1305
1306 fApproxPosterior = static_cast<TF1*>(tmp->Clone());
1307 // save this function for future reuse
1308 // I can delete now original posterior and use this approximated copy
1309 delete tmp;
1310 TString name = posterior->GetName() + TString("_approx");
1311 TString title = posterior->GetTitle() + TString("_approx");
1314 delete fIntegratedLikelihood;
1316 }
1317 else if (posterior == fLikelihood) {
1318 delete fLikelihood;
1320 }
1321 else {
1322 assert(1); // should never happen this case
1323 }
1324}
1325
1326////////////////////////////////////////////////////////////////////////////////
1327/// compute the interval using the approximate posterior function
1328
1330
1331 ccoutD(Eval) << "BayesianCalculator: Compute interval from the approximate posterior " << std::endl;
1332
1334 if (!fApproxPosterior) return;
1335
1336 double prob[2];
1337 double limits[2] = {0,0};
1338 prob[0] = lowerCutOff;
1339 prob[1] = upperCutOff;
1341 fLower = limits[0];
1342 fUpper = limits[1];
1343 fValidInterval = true;
1344}
1345
1346////////////////////////////////////////////////////////////////////////////////
1347/// compute the shortest interval from the histogram representing the posterior
1348
1349
1351 coutI(Eval) << "BayesianCalculator - computing shortest interval with CL = " << 1.-fSize << std::endl;
1352
1353 // compute via the approx posterior function
1355 if (!fApproxPosterior) return;
1356
1357 TH1D * h1 = dynamic_cast<TH1D*>(fApproxPosterior->GetHistogram() );
1358 assert(h1 != nullptr);
1360 // get bins and sort them
1361 double * bins = h1->GetArray();
1362 // exclude under/overflow bins
1363 int n = h1->GetSize()-2;
1364 std::vector<int> index(n);
1365 // exclude bins[0] (the underflow bin content)
1366 TMath::Sort(n, bins+1, &index[0]);
1367 // find cut off as test size
1368 double sum = 0;
1369 double actualCL = 0;
1370 double upper = h1->GetXaxis()->GetXmin();
1371 double lower = h1->GetXaxis()->GetXmax();
1372 double norm = h1->GetSumOfWeights();
1373
1374 for (int i = 0; i < n; ++i) {
1375 int idx = index[i];
1376 double p = bins[ idx] / norm;
1377 sum += p;
1378 if (sum > 1.-fSize ) {
1379 actualCL = sum - p;
1380 break;
1381 }
1382
1383 // histogram bin content starts from 1
1384 if ( h1->GetBinLowEdge(idx+1) < lower)
1385 lower = h1->GetBinLowEdge(idx+1);
1386 if ( h1->GetXaxis()->GetBinUpEdge(idx+1) > upper)
1387 upper = h1->GetXaxis()->GetBinUpEdge(idx+1);
1388 }
1389
1390 ccoutD(Eval) << "BayesianCalculator::ComputeShortestInterval - actual interval CL = "
1391 << actualCL << " difference from requested is " << (actualCL-(1.-fSize))/fSize*100. << "% "
1392 << " limits are [ " << lower << " , " << " upper ] " << std::endl;
1393
1394
1395 if (lower < upper) {
1396 fLower = lower;
1397 fUpper = upper;
1398
1399 if (std::abs(actualCL - (1. - fSize)) > 0.1 * (1. - fSize)) {
1400 coutW(Eval) << "BayesianCalculator::ComputeShortestInterval - actual interval CL = " << actualCL
1401 << " differs more than 10% from desired CL value - must increase nbins " << n
1402 << " to an higher value " << std::endl;
1403 }
1404 }
1405 else
1406 coutE(Eval) << "BayesianCalculator::ComputeShortestInterval " << n << " bins are not sufficient " << std::endl;
1407
1408 fValidInterval = true;
1409
1410}
1411
1412
1413
1414
1415} // end namespace RooStats
dim_t fSize
#define f(i)
Definition RSha256.hxx:104
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutI(a)
#define ccoutE(a)
#define oocoutW(o, a)
#define coutW(a)
#define oocoutI(o, a)
#define coutE(a)
#define ooccoutW(o, a)
#define ooccoutI(o, a)
#define ooccoutD(o, a)
#define ccoutI(a)
#define ccoutD(a)
#define ooccoutE(o, a)
@ kGray
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:148
User class for performing function minimization.
Functor1D class for one-dimensional functions.
Definition Functor.h:97
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition IFunction.h:157
Numerical multi dimensional integration options.
void Print(std::ostream &os=std::cout) const
print all the options
void SetNCalls(unsigned int calls)
set maximum number of function calls
static IntegrationMultiDim::Type GetType(const char *name)
static function to get the enumeration from a string
static IntegrationOneDim::Type GetType(const char *name)
static function to get the enumeration from a string
User Class to find the Root of one dimensional functions.
Definition RootFinder.h:74
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * first() const
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:56
void resetNumCall() const
Reset function call counter.
Definition RooAbsFunc.h:52
Int_t numCall() const
Return number of function calls since last reset.
Definition RooAbsFunc.h:47
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Construct representation of -log(L) of PDF with given dataset.
Definition RooAbsPdf.h:155
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:49
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
RooPlot * frame(const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create a new RooPlot on the heap with a drawing frame initialized for this object,...
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
static Int_t numEvalErrors()
Return the number of logged evaluation errors since the last clearing.
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
static void clearEvalErrorLog()
Clear the stack of evaluation error messages.
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooCFunction1Binding is a templated implementation of class RooAbsReal that binds generic C(++) funct...
Lightweight interface adaptor that exports a RooAbsPdf as a functor.
Definition RooFunctor.h:25
Int_t nObs() const
Definition RooFunctor.h:34
RooAbsFunc & binding()
Definition RooFunctor.h:51
Implementation of a probability density function that takes a RooArgList of servers and a C++ express...
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:36
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
ROOT::Math::IGenFunction * fPosteriorFunction
function representing the posterior
RooAbsReal * fLikelihood
internal pointer to likelihood function
double fNLLMin
minimum value of Nll
int fNumIterations
number of iterations (when using ToyMC)
RooAbsPdf * GetPosteriorPdf() const
return posterior pdf (object is managed by the user)
RooPlot * GetPosteriorPlot(bool norm=false, double precision=0.01) const
get the plot with option to get it normalized
int fNScanBins
number of bins to scan, if = -1 no scan is done (default)
void ClearAll() const
clear all cached pdf objects
void ComputeShortestInterval() const
compute the shortest interval from the histogram representing the posterior
RooAbsPdf * fNuisancePdf
nuisance pdf (needed when using nuisance sampling technique)
RooArgSet fConditionalObs
conditional observables
double fBrfPrecision
root finder precision
RooAbsReal * fIntegratedLikelihood
integrated likelihood function, i.e - unnormalized posterior function
void ApproximatePosterior() const
approximate posterior in nbins using a TF1 scan the poi values and evaluate the posterior at each poi...
double fSize
size used for getting the interval
RooArgSet fNuisanceParameters
nuisance parameters
double fLeftSideFraction
fraction of probability content on left side of interval
RooArgSet fGlobalObs
global observables
double GetMode() const
return the mode (most probable value of the posterior function)
RooAbsPdf * fPdf
model pdf (could contain the nuisance pdf as constraint term)
SimpleInterval * GetInterval() const override
compute the interval.
RooAbsPdf * fProductPdf
internal pointer to model * prior
TF1 * fApproxPosterior
TF1 representing the scanned posterior function.
RooAbsReal * GetPosteriorFunction() const
return posterior function (object is managed by the BayesianCalculator class)
void SetIntegrationType(const char *type)
set the integration type (possible type are) :
RooAbsPdf * fPosteriorPdf
normalized (on the poi) posterior pdf
double fUpper
upper interval bound
double fLower
computer lower interval bound
TH1 * GetPosteriorHistogram() const
return the approximate posterior as histogram (TH1 object). Note the object is managed by the Bayesia...
void ComputeIntervalFromApproxPosterior(double c1, double c2) const
compute the interval using the approximate posterior function
void SetModel(const ModelConfig &model) override
set the model via the ModelConfig
RooAbsPdf * fPriorPdf
prior pdf (typically for the POI)
std::unique_ptr< RooAbsReal > fLogLike
internal pointer to log likelihood function
double ConfidenceLevel() const override
Get the Confidence level for the test.
~BayesianCalculator() override
destructor
void ComputeIntervalFromCdf(double c1, double c2) const
internal function compute the interval using Cdf integration
< A class that holds configuration information for a model using a workspace as a store
Definition ModelConfig.h:34
const RooArgSet * GetConditionalObservables() const
get RooArgSet for conditional observables (return nullptr if not existing)
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
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)
RooAbsPdf * GetPriorPdf() const
get parameters prior pdf (return nullptr if not existing)
std::shared_ptr< RooFunctor > fPriorFunc
PosteriorCdfFunction & operator=(const PosteriorCdfFunction &)
std::unique_ptr< ROOT::Math::IntegratorMultiDim > fIntegrator
PosteriorCdfFunction(const PosteriorCdfFunction &rhs)
PosteriorCdfFunction(RooAbsReal &nll, RooArgList &bindParams, RooAbsReal *prior=nullptr, const char *integType=nullptr, double nllMinimum=0)
std::unique_ptr< ROOT::Math::Integrator > fIntegratorOneDim
double DoEval(double x) const override
implementation of the evaluation function. Must be implemented by derived classes
std::map< double, double > fNormCdfValues
ROOT::Math::IGenFunction * Clone() const override
Clone a function.
Posterior function obtaining sampling toy MC for the nuisance according to their pdf.
std::shared_ptr< RooFunctor > fPriorFunc
std::unique_ptr< RooDataSet > fGenParams
double DoEval(double x) const override
implementation of the evaluation function. Must be implemented by derived classes
PosteriorFunctionFromToyMC(RooAbsReal &nll, RooAbsPdf &pdf, RooRealVar &poi, RooArgList &nuisParams, RooAbsReal *prior=nullptr, double nllOffset=0, int niter=0, bool redoToys=true)
ROOT::Math::IGenFunction * Clone() const override
Clone a function.
PosteriorFunction(RooAbsReal &nll, RooRealVar &poi, RooArgList &nuisParams, RooAbsReal *prior=nullptr, const char *integType=nullptr, double norm=1.0, double nllOffset=0, int niter=0)
std::shared_ptr< RooFunctor > fPriorFunc
std::unique_ptr< ROOT::Math::IntegratorMultiDim > fIntegratorMultiDim
ROOT::Math::IGenFunction * Clone() const override
Clone a function.
double DoEval(double x) const override
implementation of the evaluation function. Must be implemented by derived classes
std::unique_ptr< ROOT::Math::Integrator > fIntegratorOneDim
SimpleInterval is a concrete implementation of the ConfInterval interface.
Use TF1, TF2, TF3 functions as RooFit objects.
const Float_t * GetArray() const
Definition TArrayF.h:43
Int_t GetSize() const
Definition TArray.h:47
Double_t GetXmax() const
Definition TAxis.h:142
Double_t GetXmin() const
Definition TAxis.h:141
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:532
1-Dim function class
Definition TF1.h:182
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function Note that this histogram is managed ...
Definition TF1.cxx:1634
virtual Int_t GetQuantiles(Int_t n, Double_t *xp, const Double_t *p)
Compute Quantiles for density distribution of this function.
Definition TF1.cxx:2042
virtual Int_t GetNpx() const
Definition TF1.h:455
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:926
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
TAxis * GetXaxis()
Definition TH1.h:571
virtual Double_t GetSumOfWeights() const
Return the sum of weights across all bins excluding under/overflows.
Definition TH1.h:559
virtual Double_t GetBinLowEdge(Int_t bin) const
Return bin lower edge for 1D histogram.
Definition TH1.cxx:9286
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:9097
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Basic string class.
Definition TString.h:138
void ToUpper()
Change string to upper case.
Definition TString.cxx:1202
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:641
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg GlobalObservables(Args_t &&... argsOrArgSet)
RooCmdArg ConditionalObservables(Args_t &&... argsOrArgSet)
Create a RooCmdArg to declare conditional observables.
RooCmdArg FillColor(TColorNumber color)
RooCmdArg Precision(double prec)
RooCmdArg DrawOption(const char *opt)
RooCmdArg Range(const char *rangeName, bool adjustNorm=true)
RooCmdArg MoveToBack()
RooCmdArg VLines()
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5
Namespace for new Math classes and functions.
Namespace for the RooStats classes.
Definition CodegenImpl.h:66
void RemoveConstantParameters(RooArgSet *set)
const ROOT::Math::RootFinder::EType kRootFinderType
Bool_t IsNaN(Double_t x)
Definition TMath.h:903
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
Definition TMathBase.h:413
void SetPrior(RooFunctor *prior)
LikelihoodFunction(RooFunctor &f, RooFunctor *prior=nullptr, double offset=0)
double operator()(const double *x) const
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2338