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 extrem 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 "RooBrentRootFinder.h"
58#include "RooFormulaVar.h"
59#include "RooGenericPdf.h"
60#include "RooPlot.h"
61#include "RooProdPdf.h"
62#include "RooDataSet.h"
63
64// include header file of this class
68
69#include "Math/IFunction.h"
71#include "Math/Integrator.h"
72#include "Math/RootFinder.h"
74#include "RooFunctor.h"
75#include "RooFunctor1DBinding.h"
76#include "RooTFnBinding.h"
77#include "RooMsgService.h"
78
79#include "TAxis.h"
80#include "TF1.h"
81#include "TH1.h"
82#include "TMath.h"
83
84#include <map>
85#include <cmath>
86
87#include "RConfigure.h"
88
90
91using namespace std;
92
93namespace RooStats {
94
95
96// first some utility classes and functions
97
98#ifdef R__HAS_MATHMORE
100#else
102#endif
103
104
105
106
108 LikelihoodFunction(RooFunctor & f, RooFunctor * prior = 0, double offset = 0) :
109 fFunc(f), fPrior(prior),
110 fOffset(offset), fMaxL(0) {
112 }
113
114 void SetPrior(RooFunctor * prior) { fPrior = prior; }
115
116 double operator() (const double *x ) const {
117 double nll = fFunc(x) - fOffset;
118 double likelihood = std::exp(-nll);
119
120 if (fPrior) likelihood *= (*fPrior)(x);
121
122 int nCalls = fFunc.binding().numCall();
123 if (nCalls > 0 && nCalls % 1000 == 0) {
124 ooccoutD(nullptr,Eval) << "Likelihood evaluation ncalls = " << nCalls
125 << " x0 " << x[0] << " nll = " << nll+fOffset;
126 if (fPrior) ooccoutD(nullptr,Eval) << " prior(x) = " << (*fPrior)(x);
127 ooccoutD(nullptr,Eval) << " likelihood " << likelihood
128 << " max Likelihood " << fMaxL << std::endl;
129 }
130
131 if (likelihood > fMaxL ) {
132 fMaxL = likelihood;
133 if ( likelihood > 1.E10) {
134 ooccoutW(nullptr,Eval) << "LikelihoodFunction::() WARNING - Huge likelihood value found for parameters ";
135 for (int i = 0; i < fFunc.nObs(); ++i)
136 ooccoutW(nullptr,Eval) << " x[" << i << " ] = " << x[i];
137 ooccoutW(nullptr,Eval) << " nll = " << nll << " L = " << likelihood << std::endl;
138 }
139 }
140
141 return likelihood;
142 }
143
144 // for the 1D case
145 double operator() (double x) const {
146 // just call the previous method
147 assert(fFunc.nObs() == 1); // check nobs = 1
148 double tmp = x;
149 return (*this)(&tmp);
150 }
151
152 RooFunctor & fFunc; // functor representing the nll function
153 RooFunctor * fPrior; // functor representing the prior function
154 double fOffset; // offset used to bring the nll in a reasonable range for computing the exponent
155 mutable double fMaxL;
156};
157
158
159// Posterior CDF Function class
160// for integral of posterior function in nuisance and POI
161// 1-Dim function as function of the poi
162
164
165public:
166
167 PosteriorCdfFunction(RooAbsReal & nll, RooArgList & bindParams, RooAbsReal * prior = 0, const char * integType = 0, double nllMinimum = 0) :
168 fFunctor(nll, bindParams, RooArgList() ), // functor
169 fPriorFunc(nullptr),
170 fLikelihood(fFunctor, 0, nllMinimum), // integral of exp(-nll) function
171 fIntegrator(ROOT::Math::IntegratorMultiDim::GetType(integType) ), // integrator
172 fXmin(bindParams.getSize() ), // vector of parameters (min values)
173 fXmax(bindParams.getSize() ), // vector of parameter (max values)
174 fNorm(1.0), fNormErr(0.0), fOffset(0), fMaxPOI(0),
175 fHasNorm(false), fUseOldValues(true), fError(false)
176 {
177
178 if (prior) {
179 fPriorFunc = std::make_shared<RooFunctor>(*prior, bindParams, RooArgList());
181 }
182
184
185 ooccoutD(nullptr,NumIntegration) << "PosteriorCdfFunction::Compute integral of posterior in nuisance and poi. "
186 << " nllMinimum is " << nllMinimum << std::endl;
187
188 std::vector<double> par(bindParams.getSize());
189 for (unsigned int i = 0; i < fXmin.size(); ++i) {
190 RooRealVar & var = (RooRealVar &) bindParams[i];
191 fXmin[i] = var.getMin();
192 fXmax[i] = var.getMax();
193 par[i] = var.getVal();
194 ooccoutD(nullptr,NumIntegration) << "PosteriorFunction::Integrate" << var.GetName()
195 << " in interval [ " << fXmin[i] << " , " << fXmax[i] << " ] " << std::endl;
196 }
197
198 fIntegrator.Options().Print(ooccoutD(nullptr,NumIntegration));
199
200 // store max POI value because it will be changed when evaluating the function
201 fMaxPOI = fXmax[0];
202
203 // compute first the normalization with the poi
204 fNorm = (*this)( fMaxPOI );
205 if (fError) ooccoutE(nullptr,NumIntegration) << "PosteriorFunction::Error computing normalization - norm = " << fNorm << std::endl;
206 fHasNorm = true;
207 fNormCdfValues.insert(std::make_pair(fXmin[0], 0) );
208 fNormCdfValues.insert(std::make_pair(fXmax[0], 1.0) );
209
210 }
211
212 // copy constructor (needed for Cloning the object)
213 // need special treatment because integrator
214 // has no copy constructor
216 ROOT::Math::IGenFunction(),
217 fFunctor(rhs.fFunctor),
218 //fPriorFunc(std::shared_ptr<RooFunctor>((RooFunctor*)0)),
221 fIntegrator(ROOT::Math::IntegratorMultiDim::GetType( rhs.fIntegrator.Name().c_str() ) ), // integrator
222 fXmin( rhs.fXmin),
223 fXmax( rhs.fXmax),
224 fNorm( rhs.fNorm),
225 fNormErr( rhs.fNormErr),
226 fOffset(rhs.fOffset),
227 fMaxPOI(rhs.fMaxPOI),
228 fHasNorm(rhs.fHasNorm),
230 fError(rhs.fError),
232 {
234 // need special treatment for the smart pointer
235 // if (rhs.fPriorFunc.get() ) {
236 // fPriorFunc = std::shared_ptr<RooFunctor>(new RooFunctor(*(rhs.fPriorFunc) ) );
237 // fLikelihood.SetPrior( fPriorFunc.get() );
238 // }
239 }
240
241
242 bool HasError() const { return fError; }
243
244
245 ROOT::Math::IGenFunction * Clone() const override {
246 ooccoutD(nullptr,NumIntegration) << " cloning function .........." << std::endl;
247 return new PosteriorCdfFunction(*this);
248 }
249
250 // offset value for computing the root
251 void SetOffset(double offset) { fOffset = offset; }
252
253private:
254
255 // make assignment operator private
257 return *this;
258 }
259
260 double DoEval (double x) const override {
261
262 // evaluate cdf at poi value x by integrating poi from [xmin,x] and all the nuisances
263 fXmax[0] = x;
264 if (x <= fXmin[0] ) return -fOffset;
265 // could also avoid a function evaluation at maximum
266 if (x >= fMaxPOI && fHasNorm) return 1. - fOffset; // cdf is bound to these values
267
268 // computes the integral using a previous cdf estimate
269 double normcdf0 = 0;
270 if (fHasNorm && fUseOldValues) {
271 // look in the map of the stored cdf values the closes one
272 std::map<double,double>::iterator itr = fNormCdfValues.upper_bound(x);
273 --itr; // upper bound returns a position 1 up of the value we want
274 if (itr != fNormCdfValues.end() ) {
275 fXmin[0] = itr->first;
276 normcdf0 = itr->second;
277 // ooccoutD(nullptr,NumIntegration) << "PosteriorCdfFunction: computing integral between in poi interval : "
278 // << fXmin[0] << " - " << fXmax[0] << std::endl;
279 }
280 }
281
282 fFunctor.binding().resetNumCall(); // reset number of calls for debug
283
284 double cdf = fIntegrator.Integral(&fXmin[0],&fXmax[0]);
285 double error = fIntegrator.Error();
286 double normcdf = cdf/fNorm; // normalize the cdf
287
288 ooccoutD(nullptr,NumIntegration) << "PosteriorCdfFunction: poi = [" << fXmin[0] << " , "
289 << fXmax[0] << "] integral = " << cdf << " +/- " << error
290 << " norm-integ = " << normcdf << " cdf(x) = " << normcdf+normcdf0
291 << " ncalls = " << fFunctor.binding().numCall() << std::endl;
292
293 if (TMath::IsNaN(cdf) || cdf > std::numeric_limits<double>::max()) {
294 ooccoutE(nullptr,NumIntegration) << "PosteriorFunction::Error computing integral - cdf = "
295 << cdf << std::endl;
296 fError = true;
297 }
298
299 if (cdf != 0 && error/cdf > 0.2 )
300 oocoutW(nullptr,NumIntegration) << "PosteriorCdfFunction: integration error is larger than 20 % x0 = " << fXmin[0]
301 << " x = " << x << " cdf(x) = " << cdf << " +/- " << error << std::endl;
302
303 if (!fHasNorm) {
304 oocoutI(nullptr,NumIntegration) << "PosteriorCdfFunction - integral of posterior = "
305 << cdf << " +/- " << error << std::endl;
306 fNormErr = error;
307 return cdf;
308 }
309
310 normcdf += normcdf0;
311
312 // store values in the map
313 if (fUseOldValues) {
314 fNormCdfValues.insert(std::make_pair(x, normcdf) );
315 }
316
317 double errnorm = sqrt( error*error + normcdf*normcdf * fNormErr * fNormErr )/fNorm;
318 if (normcdf > 1. + 3 * errnorm) {
319 oocoutW(nullptr,NumIntegration) << "PosteriorCdfFunction: normalized cdf values is larger than 1"
320 << " x = " << x << " normcdf(x) = " << normcdf << " +/- " << error/fNorm << std::endl;
321 }
322
323 return normcdf - fOffset; // apply an offset (for finding the roots)
324 }
325
326 mutable RooFunctor fFunctor; // functor binding nll
327 mutable std::shared_ptr<RooFunctor> fPriorFunc; // functor binding the prior
328 LikelihoodFunction fLikelihood; // likelihood function
329 mutable ROOT::Math::IntegratorMultiDim fIntegrator; // integrator (mutable because Integral() is not const
330 mutable std::vector<double> fXmin; // min value of parameters (poi+nuis) -
331 mutable std::vector<double> fXmax; // max value of parameters (poi+nuis) - max poi changes so it is mutable
332 double fNorm; // normalization value (computed in ctor)
333 mutable double fNormErr; // normalization error value (computed in ctor)
334 double fOffset; // offset for computing the root
335 double fMaxPOI; // maximum value of POI
336 bool fHasNorm; // flag to control first call to the function
337 bool fUseOldValues; // use old cdf values
338 mutable bool fError; // flag to indicate if a numerical evaluation error occurred
339 mutable std::map<double,double> fNormCdfValues;
340};
341
342//__________________________________________________________________
343// Posterior Function class
344// 1-Dim function as function of the poi
345// and it integrated all the nuisance parameters
346
348
349public:
350
351
352 PosteriorFunction(RooAbsReal & nll, RooRealVar & poi, RooArgList & nuisParams, RooAbsReal * prior = 0, const char * integType = 0, double
353 norm = 1.0, double nllOffset = 0, int niter = 0) :
354 fFunctor(nll, nuisParams, RooArgList() ),
355 fPriorFunc(nullptr),
356 fLikelihood(fFunctor, 0, nllOffset),
357 fPoi(&poi),
358 fXmin(nuisParams.getSize() ),
359 fXmax(nuisParams.getSize() ),
360 fNorm(norm),
361 fError(0)
362 {
363
364 if (prior) {
365 fPriorFunc = std::make_shared<RooFunctor>(*prior, nuisParams, RooArgList());
367 }
368
369 ooccoutD(nullptr,NumIntegration) << "PosteriorFunction::Evaluate the posterior function by integrating the nuisances: " << std::endl;
370 for (unsigned int i = 0; i < fXmin.size(); ++i) {
371 RooRealVar & var = (RooRealVar &) nuisParams[i];
372 fXmin[i] = var.getMin();
373 fXmax[i] = var.getMax();
374 ooccoutD(nullptr,NumIntegration) << "PosteriorFunction::Integrate " << var.GetName()
375 << " in interval [" << fXmin[i] << " , " << fXmax[i] << " ] " << std::endl;
376 }
377 if (fXmin.size() == 1) { // 1D case
379
380 fIntegratorOneDim->SetFunction(fLikelihood);
381 // interested only in relative tolerance
382 //fIntegratorOneDim->SetAbsTolerance(1.E-300);
383 fIntegratorOneDim->Options().Print(ooccoutD(nullptr,NumIntegration) );
384 }
385 else if (fXmin.size() > 1) { // multiDim case
386 fIntegratorMultiDim = std::make_unique<ROOT::Math::IntegratorMultiDim>(ROOT::Math::IntegratorMultiDim::GetType(integType));
387 fIntegratorMultiDim->SetFunction(fLikelihood, fXmin.size());
389 if (niter > 0) {
390 opt.SetNCalls(niter);
391 fIntegratorMultiDim->SetOptions(opt);
392 }
393 //fIntegratorMultiDim->SetAbsTolerance(1.E-300);
394 // print the options
395 opt.Print(ooccoutD(nullptr,NumIntegration) );
396 }
397 }
398
399
400 ROOT::Math::IGenFunction * Clone() const override {
401 assert(1);
402 return 0; // cannot clone this function for integrator
403 }
404
405 double Error() const { return fError;}
406
407
408private:
409 double DoEval (double x) const override {
410
411 // evaluate posterior function at a poi value x by integrating all nuisance parameters
412
413 fPoi->setVal(x);
414 fFunctor.binding().resetNumCall(); // reset number of calls for debug
415
416 double f = 0;
417 double error = 0;
418 if (fXmin.size() == 1) { // 1D case
419 f = fIntegratorOneDim->Integral(fXmin[0],fXmax[0]);
420 error = fIntegratorOneDim->Error();
421 }
422 else if (fXmin.size() > 1) { // multi-dim case
423 f = fIntegratorMultiDim->Integral(&fXmin[0],&fXmax[0]);
424 error = fIntegratorMultiDim->Error();
425 } else {
426 // no integration to be done
427 f = fLikelihood(x);
428 }
429
430 // debug
431 ooccoutD(nullptr,NumIntegration) << "PosteriorFunction: POI value = "
432 << x << "\tf(x) = " << f << " +/- " << error
433 << " norm-f(x) = " << f/fNorm
434 << " ncalls = " << fFunctor.binding().numCall() << std::endl;
435
436
437
438
439 if (f != 0 && error/f > 0.2 )
440 ooccoutW(nullptr,NumIntegration) << "PosteriorFunction::DoEval - Error from integration in "
441 << fXmin.size() << " Dim is larger than 20 % "
442 << "x = " << x << " p(x) = " << f << " +/- " << error << std::endl;
443
444 fError = error / fNorm;
445 return f / fNorm;
446 }
447
449 mutable std::shared_ptr<RooFunctor> fPriorFunc; // functor binding the prior
452 std::unique_ptr<ROOT::Math::Integrator> fIntegratorOneDim;
453 std::unique_ptr<ROOT::Math::IntegratorMultiDim> fIntegratorMultiDim;
454 std::vector<double> fXmin;
455 std::vector<double> fXmax;
456 double fNorm;
457 mutable double fError;
458};
459
460////////////////////////////////////////////////////////////////////////////////
461/// Posterior function obtaining sampling toy MC for the nuisance according to their pdf
462
464
465public:
466
467
468 PosteriorFunctionFromToyMC(RooAbsReal & nll, RooAbsPdf & pdf, RooRealVar & poi, RooArgList & nuisParams, RooAbsReal * prior = 0, double
469 nllOffset = 0, int niter = 0, bool redoToys = true ) :
470 fFunctor(nll, nuisParams, RooArgList() ),
471 fPriorFunc(nullptr),
472 fLikelihood(fFunctor, 0, nllOffset),
473 fPdf(&pdf),
474 fPoi(&poi),
475 fNuisParams(nuisParams),
476 fNumIterations(niter),
477 fError(-1),
478 fRedoToys(redoToys)
479 {
480 if (niter == 0) fNumIterations = 100; // default value
481
482 if (prior) {
483 fPriorFunc = std::make_shared<RooFunctor>(*prior, nuisParams, RooArgList());
485 }
486
487 ooccoutI(nullptr,InputArguments) << "PosteriorFunctionFromToyMC::Evaluate the posterior function by randomizing the nuisances: niter " << fNumIterations << std::endl;
488
489 ooccoutI(nullptr,InputArguments) << "PosteriorFunctionFromToyMC::Pdf used for randomizing the nuisance is " << fPdf->GetName() << std::endl;
490 // check that pdf contains the nuisance
491 std::unique_ptr<RooArgSet> vars{fPdf->getVariables()};
492 for (int i = 0; i < fNuisParams.getSize(); ++i) {
493 if (!vars->find( fNuisParams[i].GetName() ) ) {
494 ooccoutW(nullptr,InputArguments) << "Nuisance parameter " << fNuisParams[i].GetName()
495 << " is not part of sampling pdf. "
496 << "they will be treated as constant " << std::endl;
497 }
498 }
499
500 if (!fRedoToys) {
501 ooccoutI(nullptr,InputArguments) << "PosteriorFunctionFromToyMC::Generate nuisance toys only one time (for all POI points)" << std::endl;
502 GenerateToys();
503 }
504 }
505
506 // generate first n-samples of the nuisance parameters
507 void GenerateToys() const {
508 fGenParams = std::unique_ptr<RooDataSet>{fPdf->generate(fNuisParams, fNumIterations)};
509 if(fGenParams==0) {
510 ooccoutE(nullptr,InputArguments) << "PosteriorFunctionFromToyMC - failed to generate nuisance parameters" << std::endl;
511 }
512 }
513
514 double Error() const { return fError;}
515
516 ROOT::Math::IGenFunction * Clone() const override {
517 // use default copy constructor
518 //return new PosteriorFunctionFromToyMC(*this);
519 // clone not implemented
520 assert(1);
521 return 0;
522 }
523
524private:
525 // evaluate the posterior at the poi value x
526 double DoEval( double x) const override {
527
528 int npar = fNuisParams.getSize();
529 assert (npar > 0);
530
531
532 // generate the toys
533 if (fRedoToys) GenerateToys();
534 if (!fGenParams) return 0;
535
536 // evaluate posterior function at a poi value x by integrating all nuisance parameters
537
538 fPoi->setVal(x);
539
540 // loop over all of the generate data
541 double sum = 0;
542 double sum2 = 0;
543
544 for(int iter=0; iter<fNumIterations; ++iter) {
545
546 // get the set of generated parameters and set the nuisance parameters to the generated values
547 std::vector<double> p(npar);
548 for (int i = 0; i < npar; ++i) {
549 const RooArgSet* genset=fGenParams->get(iter);
550 RooAbsArg * arg = genset->find( fNuisParams[i].GetName() );
551 RooRealVar * var = dynamic_cast<RooRealVar*>(arg);
552 assert(var != 0);
553 p[i] = var->getVal();
554 ((RooRealVar &) fNuisParams[i]).setVal(p[i]);
555 }
556
557 // evaluate now the likelihood function
558 double fval = fLikelihood( &p.front() );
559
560 // likelihood already must contained the pdf we have sampled
561 // so we must divided by it. The value must be normalized on all
562 // other parameters
564 double nuisPdfVal = fPdf->getVal(&arg);
565 fval /= nuisPdfVal;
566
567
568 if( fval > std::numeric_limits<double>::max() ) {
569 ooccoutE(nullptr,Eval) << "BayesianCalculator::EvalPosteriorFunctionFromToy : "
570 << "Likelihood evaluates to infinity " << std::endl;
571 ooccoutE(nullptr,Eval) << "poi value = " << x << std::endl;
572 ooccoutE(nullptr,Eval) << "Nuisance parameter values : ";
573 for (int i = 0; i < npar; ++i)
574 ooccoutE(nullptr,Eval) << fNuisParams[i].GetName() << " = " << p[i] << " ";
575 ooccoutE(nullptr,Eval) << " - return 0 " << std::endl;
576
577 fError = 1.E30;
578 return 0;
579 }
580 if( TMath::IsNaN(fval) ) {
581 ooccoutE(nullptr,Eval) << "BayesianCalculator::EvalPosteriorFunctionFromToy : "
582 << "Likelihood is a NaN " << std::endl;
583 ooccoutE(nullptr,Eval) << "poi value = " << x << std::endl;
584 ooccoutE(nullptr,Eval) << "Nuisance parameter values : ";
585 for (int i = 0; i < npar; ++i)
586 ooccoutE(nullptr,Eval) << fNuisParams[i].GetName() << " = " << p[i] << " ";
587 ooccoutE(nullptr,Eval) << " - return 0 " << std::endl;
588 fError = 1.E30;
589 return 0;
590 }
591
592
593
594 sum += fval;
595 sum2 += fval*fval;
596 }
597
598 // compute the average and variance
599 double val = sum/double(fNumIterations);
600 double dval2 = std::max( sum2/double(fNumIterations) - val*val, 0.0);
601 fError = std::sqrt( dval2 / fNumIterations);
602
603 // debug
604 ooccoutD(nullptr,NumIntegration) << "PosteriorFunctionFromToyMC: POI value = "
605 << x << "\tp(x) = " << val << " +/- " << fError << std::endl;
606
607
608 if (val != 0 && fError/val > 0.2 ) {
609 ooccoutW(nullptr,NumIntegration) << "PosteriorFunctionFromToyMC::DoEval"
610 << " - Error in estimating posterior is larger than 20% ! "
611 << "x = " << x << " p(x) = " << val << " +/- " << fError << std::endl;
612 }
613
614
615 return val;
616 }
617
619 mutable std::shared_ptr<RooFunctor> fPriorFunc; // functor binding the prior
621 mutable RooAbsPdf * fPdf;
624 mutable std::unique_ptr<RooDataSet> fGenParams;
626 mutable double fError;
627 bool fRedoToys; // do toys every iteration
628
629};
630
631////////////////////////////////////////////////////////////////////////////////
632// Implementation of BayesianCalculator
633////////////////////////////////////////////////////////////////////////////////
634
635////////////////////////////////////////////////////////////////////////////////
636/// default constructor
637
639 fData(0),
640 fPdf(0),
641 fPriorPdf(0),
642 fNuisancePdf(0),
643 fProductPdf (0), fLikelihood (0), fIntegratedLikelihood (0), fPosteriorPdf(0),
644 fPosteriorFunction(0), fApproxPosterior(0),
645 fLower(0), fUpper(0),
646 fNLLMin(0),
647 fSize(0.05), fLeftSideFraction(0.5),
648 fBrfPrecision(0.00005),
649 fNScanBins(-1),
650 fNumIterations(0),
651 fValidInterval(false)
652{
653}
654
655////////////////////////////////////////////////////////////////////////////////
656/// Constructor from data set, model pdf, parameter of interests and prior pdf
657/// If nuisance parameters are given they will be integrated according either to the prior or
658/// their constraint term included in the model
659
660BayesianCalculator::BayesianCalculator( /* const char* name, const char* title, */
662 RooAbsPdf& pdf,
663 const RooArgSet& POI,
664 RooAbsPdf& priorPdf,
665 const RooArgSet* nuisanceParameters ) :
666 fData(&data),
667 fPdf(&pdf),
668 fPOI(POI),
669 fPriorPdf(&priorPdf),
670 fNuisancePdf(0),
671 fProductPdf (0), fLikelihood (0), fIntegratedLikelihood (0), fPosteriorPdf(0),
672 fPosteriorFunction(0), fApproxPosterior(0),
673 fLower(0), fUpper(0),
674 fNLLMin(0),
675 fSize(0.05), fLeftSideFraction(0.5),
676 fBrfPrecision(0.00005),
677 fNScanBins(-1),
678 fNumIterations(0),
679 fValidInterval(false)
680{
681
682 if (nuisanceParameters) fNuisanceParameters.add(*nuisanceParameters);
683 // remove constant nuisance parameters
685}
686
687////////////////////////////////////////////////////////////////////////////////
688/// Constructor from a data set and a ModelConfig
689/// model pdf, poi and nuisances will be taken from the ModelConfig
690
692 ModelConfig & model) :
693 fData(&data),
694 fPdf(model.GetPdf()),
695 fPriorPdf( model.GetPriorPdf()),
696 fNuisancePdf(0),
697 fProductPdf (0), fLikelihood (0), fIntegratedLikelihood (0), fPosteriorPdf(0),
698 fPosteriorFunction(0), fApproxPosterior(0),
699 fLower(0), fUpper(0),
700 fNLLMin(0),
701 fSize(0.05), fLeftSideFraction(0.5),
702 fBrfPrecision(0.00005),
703 fNScanBins(-1),
704 fNumIterations(0),
705 fValidInterval(false)
706{
707 SetModel(model);
708}
709
710
712{
713 // destructor
714 ClearAll();
715}
716
717////////////////////////////////////////////////////////////////////////////////
718/// clear all cached pdf objects
719
721 if (fProductPdf) delete fProductPdf;
722 fLogLike.reset();
723 if (fLikelihood) delete fLikelihood;
725 if (fPosteriorPdf) delete fPosteriorPdf;
728 fPosteriorPdf = 0;
730 fProductPdf = 0;
731 fLikelihood = 0;
733 fLower = 0;
734 fUpper = 0;
735 fNLLMin = 0;
736 fValidInterval = false;
737}
738
739////////////////////////////////////////////////////////////////////////////////
740/// set the model to use
741/// The model pdf, prior pdf, parameter of interest and nuisances
742/// will be taken according to the model
743
745
746 fPdf = model.GetPdf();
747 fPriorPdf = model.GetPriorPdf();
748 // assignment operator = does not do a real copy the sets (must use add method)
749 fPOI.removeAll();
753 if (model.GetParametersOfInterest()) fPOI.add( *(model.GetParametersOfInterest()) );
756 if (model.GetGlobalObservables()) fGlobalObs.add( *(model.GetGlobalObservables() ) );
757 // remove constant nuisance parameters
759
760 // invalidate the cached pointers
761 ClearAll();
762}
763
764////////////////////////////////////////////////////////////////////////////////
765/// Build and return the posterior function (not normalized) as a RooAbsReal
766/// the posterior is obtained from the product of the likelihood function and the
767/// prior pdf which is then integrated in the nuisance parameters (if existing).
768/// A prior function for the nuisance can be specified either in the prior pdf object
769/// or in the model itself. If no prior nuisance is specified, but prior parameters are then
770/// the integration is performed assuming a flat prior for the nuisance parameters.
771///
772/// NOTE: the return object is managed by the BayesianCalculator class, users do not need to delete it,
773/// but the object will be deleted when the BayesiabCalculator object is deleted
774
776{
777
779 if (fLikelihood) return fLikelihood;
780
781 // run some sanity checks
782 if (!fPdf ) {
783 coutE(InputArguments) << "BayesianCalculator::GetPosteriorPdf - missing pdf model" << std::endl;
784 return 0;
785 }
786 if (fPOI.empty()) {
787 coutE(InputArguments) << "BayesianCalculator::GetPosteriorPdf - missing parameter of interest" << std::endl;
788 return 0;
789 }
790 if (fPOI.getSize() > 1) {
791 coutE(InputArguments) << "BayesianCalculator::GetPosteriorPdf - current implementation works only on 1D intervals" << std::endl;
792 return 0;
793 }
794
795
796 std::unique_ptr<RooArgSet> constrainedParams{fPdf->getParameters(*fData)};
797 // remove the constant parameters
798 RemoveConstantParameters(&*constrainedParams);
799
800 //constrainedParams->Print("V");
801
802 // use RooFit::Constrain() to be sure constraints terms are taken into account
804
805
806
807 ccoutD(Eval) << "BayesianCalculator::GetPosteriorFunction : "
808 << " pdf value " << fPdf->getVal()
809 << " neglogLikelihood = " << fLogLike->getVal() << std::endl;
810
811 if (fPriorPdf)
812 ccoutD(Eval) << "\t\t\t priorPOI value " << fPriorPdf->getVal() << std::endl;
813
814 // check that likelihood evaluation is not infinity
815 double nllVal = fLogLike->getVal();
816 if ( nllVal > std::numeric_limits<double>::max() ) {
817 coutE(Eval) << "BayesianCalculator::GetPosteriorFunction : "
818 << " Negative log likelihood evaluates to infinity " << std::endl
819 << " Non-const Parameter values : ";
820 RooArgList p(*constrainedParams);
821 for (int i = 0; i < p.getSize(); ++i) {
822 RooRealVar * v = dynamic_cast<RooRealVar *>(&p[i] );
823 if (v!=0) ccoutE(Eval) << v->GetName() << " = " << v->getVal() << " ";
824 }
825 ccoutE(Eval) << std::endl;
826 ccoutE(Eval) << "-- Perform a full likelihood fit of the model before or set more reasonable parameter values"
827 << std::endl;
828 coutE(Eval) << "BayesianCalculator::GetPosteriorFunction : " << " cannot compute posterior function " << std::endl;
829 return 0;
830 }
831
832
833
834 // need do find minimum of log-likelihood in the range to shift function
835 // to avoid numerical errors when we compute the likelihood (overflows in the exponent)
836 // N.B.: this works for only 1 parameter of interest otherwise Minuit should be used for finding the minimum
837 RooFunctor * nllFunc = fLogLike->functor(fPOI);
838 assert(nllFunc);
839 ROOT::Math::Functor1D wnllFunc(*nllFunc);
840 RooRealVar* poi = dynamic_cast<RooRealVar*>( fPOI.first() );
841 assert(poi);
842
843 // try to reduce some error messages
844 //bool silentMode = (RooMsgService::instance().globalKillBelow() >= RooFit::ERROR || RooMsgService::instance().silentMode()) ;
846
847
848 coutI(Eval) << "BayesianCalculator::GetPosteriorFunction : "
849 << " nll value " << nllVal << " poi value = " << poi->getVal() << std::endl;
850
851
853 minim.SetFunction(wnllFunc,poi->getMin(),poi->getMax() );
854 bool ret = minim.Minimize(100,1.E-3,1.E-3);
855 fNLLMin = 0;
856 if (ret) fNLLMin = minim.FValMinimum();
857
858 coutI(Eval) << "BayesianCalculator::GetPosteriorFunction : minimum of NLL vs POI for POI = "
859 << poi->getVal() << " min NLL = " << fNLLMin << std::endl;
860
861 delete nllFunc;
862
863
864 if ( fNuisanceParameters.empty() || fIntegrationType.Contains("ROOFIT") ) {
865
866 ccoutD(Eval) << "BayesianCalculator::GetPosteriorFunction : use ROOFIT integration "
867 << std::endl;
868
869#ifdef DOLATER // (not clear why this does not work)
870 // need to make in this case a likelihood from the nll and make the product with the prior
871 TString likeName = TString("likelihood_times_prior_") + TString(fPriorPdf->GetName());
872 TString formula;
873 formula.Form("exp(-@0+%f+log(@1))",fNLLMin);
874 fLikelihood = new RooFormulaVar(likeName,formula,RooArgList(*fLogLike,*fPriorPdf));
875#else
876 // here use RooProdPdf (not very nice) but working
877
878 if (fLogLike) fLogLike.reset();
879 if (fProductPdf) {
880 delete fProductPdf;
881 fProductPdf = 0;
882 }
883
884 // // create a unique name for the product pdf
885 RooAbsPdf * pdfAndPrior = fPdf;
886 if (fPriorPdf) {
887 TString prodName = TString("product_") + TString(fPdf->GetName()) + TString("_") + TString(fPriorPdf->GetName() );
888 // save this as data member since it needs to be deleted afterwards
889 fProductPdf = new RooProdPdf(prodName,"",RooArgList(*fPdf,*fPriorPdf));
890 pdfAndPrior = fProductPdf;
891 }
892
893 std::unique_ptr<RooArgSet> constrParams{fPdf->getParameters(*fData)};
894 // remove the constant parameters
895 RemoveConstantParameters(&*constrParams);
897
898 TString likeName = TString("likelihood_times_prior_") + TString(pdfAndPrior->GetName());
899 TString formula;
900 formula.Form("exp(-@0+%f)",fNLLMin);
901 fLikelihood = new RooFormulaVar(likeName,formula,RooArgList(*fLogLike));
902#endif
903
904
905 // if no nuisance parameter we can just return the likelihood function
908 fLikelihood = 0;
909 }
910 else
911 // case of using RooFit for the integration
912 fIntegratedLikelihood = std::unique_ptr<RooAbsReal>{fLikelihood->createIntegral(fNuisanceParameters)}.release();
913
914
915 }
916
917 else if ( fIntegrationType.Contains("TOYMC") ) {
918 // compute the posterior as expectation values of the likelihood function
919 // sampling on the nuisance parameters
920
922
923 bool doToysEveryIteration = true;
924 // if type is 1-TOYMC or TOYMC-1
925 if ( fIntegrationType.Contains("1") || fIntegrationType.Contains("ONE") ) doToysEveryIteration = false;
926
927 RooAbsPdf * samplingPdf = (fNuisancePdf) ? fNuisancePdf : fPdf;
928 if (!fNuisancePdf) {
929 ccoutI(Eval) << "BayesianCalculator::GetPosteriorFunction : no nuisance pdf is provided, try using global pdf (this will be slower)"
930 << std::endl;
931 }
932 fPosteriorFunction = new PosteriorFunctionFromToyMC(*fLogLike, *samplingPdf, *poi, nuisParams, fPriorPdf, fNLLMin,
933 fNumIterations, doToysEveryIteration );
934
935 TString name = "toyposteriorfunction_from_";
936 name += fLogLike->GetName();
938
939 // need to scan likelihood in this case
940 if (fNScanBins <= 0) fNScanBins = 100;
941
942 }
943
944 else {
945
946 // use ROOT integration method if there are nuisance parameters
947
950
951 TString name = "posteriorfunction_from_";
952 name += fLogLike->GetName();
954
955 }
956
957
959 coutW(Eval) << "BayesianCalculator::GetPosteriorFunction : " << RooAbsReal::numEvalErrors() << " errors reported in evaluating log-likelihood function "
960 << std::endl;
963
965
966}
967
968////////////////////////////////////////////////////////////////////////////////
969/// Build and return the posterior pdf (i.e posterior function normalized to all range of poi)
970/// Note that an extra integration in the POI is required for the normalization
971/// NOTE: user must delete the returned object
972
974{
975
977 if (!plike) return 0;
978
979
980 // create a unique name on the posterior from the names of the components
981 TString posteriorName = this->GetName() + TString("_posteriorPdf_") + plike->GetName();
982
983 RooAbsPdf * posteriorPdf = new RooGenericPdf(posteriorName,"@0",*plike);
984
985 return posteriorPdf;
986}
987
988////////////////////////////////////////////////////////////////////////////////
989/// When am approximate posterior is computed binninig the parameter of interest (poi) range
990/// (see SetScanOfPosteriors) an histogram is created and can be returned to the user
991/// A nullptr is instead returned when the posterior is computed without binning the poi.
992///
993/// NOTE: the returned object is managed by the BayesianCalculator class,
994/// if the user wants to take ownership of the returned histogram, he needs to clone
995/// or copy the return object.
996
998{
999 return (fApproxPosterior) ? fApproxPosterior->GetHistogram() : nullptr;
1000}
1001
1002
1003////////////////////////////////////////////////////////////////////////////////
1004/// return a RooPlot with the posterior and the credibility region
1005/// NOTE: User takes ownership of the returned object
1006
1007RooPlot* BayesianCalculator::GetPosteriorPlot(bool norm, double precision ) const
1008{
1009
1011
1012 // if a scan is requested approximate the posterior
1013 if (fNScanBins > 0)
1015
1016 RooAbsReal * posterior = fIntegratedLikelihood;
1017 if (norm) {
1018 // delete and re-do always posterior pdf (could be invalid after approximating it)
1019 if (fPosteriorPdf) delete fPosteriorPdf;
1021 posterior = fPosteriorPdf;
1022 }
1023 if (!posterior) return 0;
1024
1026
1027 RooAbsRealLValue* poi = dynamic_cast<RooAbsRealLValue*>( fPOI.first() );
1028 assert(poi);
1029
1030
1031 RooPlot* plot = poi->frame();
1032 if (!plot) return 0;
1033
1034 // try to reduce some error messages
1036
1037 plot->SetTitle(TString("Posterior probability of parameter \"")+TString(poi->GetName())+TString("\""));
1039 posterior->plotOn(plot);
1040 plot->GetYaxis()->SetTitle("posterior function");
1041
1042 // reset the counts and default mode
1045
1046 return plot;
1047}
1048
1049////////////////////////////////////////////////////////////////////////////////
1050/// set the integration type (possible type are) :
1051///
1052/// - 1D integration ( used when only one nuisance and when the posterior is scanned):
1053/// adaptive , gauss, nonadaptive
1054/// - multidim:
1055/// - ADAPTIVE, adaptive numerical integration
1056/// The parameter numIters (settable with SetNumIters) is the max number of function calls.
1057/// It can be reduced to make the integration faster but it will be difficult to reach the required tolerance
1058/// - VEGAS MC integration method based on importance sampling - numIters is number of function calls
1059/// Extra Vegas parameter can be set using IntegratorMultiDimOptions class
1060/// - MISER MC integration method based on stratified sampling
1061/// See also http://en.wikipedia.org/wiki/Monte_Carlo_integration for VEGAS and MISER description
1062/// - PLAIN simple MC integration method, where the max number of calls can be specified using SetNumIters(numIters)
1063///
1064/// Extra integration types are:
1065///
1066/// - TOYMC:
1067/// evaluate posterior by generating toy MC for the nuisance parameters. It is a MC
1068/// integration, where the function is sampled according to the nuisance. It is convenient to use when all
1069/// the nuisance are uncorrelated and it is efficient to generate them
1070/// The toy are generated by default for each poi values
1071/// (this method has been proposed and provided by J.P Chou)
1072/// - 1-TOYMC : same method as before but in this case the toys are generated only one time and then used for
1073/// each poi value. It can be convenient when the generation time is much larger than the evaluation time,
1074/// otherwise it is recommended to re-generate the toy for each poi scanned point of the posterior function
1075/// - ROOFIT:
1076/// use roofit default integration methods which will produce a nested integral (not recommended for more
1077/// than 1 nuisance parameters)
1078
1080 // if type = 0 use default specified via class IntegratorMultiDimOptions::SetDefaultIntegrator
1083}
1084
1085////////////////////////////////////////////////////////////////////////////////
1086/// Compute the interval. By Default a central interval is computed
1087/// and the result is a SimpleInterval object.
1088///
1089/// Using the method (to be called before SetInterval) SetLeftSideTailFraction the user can choose the type of interval.
1090/// By default the returned interval is a central interval with the confidence level specified
1091/// previously in the constructor ( LeftSideTailFraction = 0.5).
1092/// - For lower limit use SetLeftSideTailFraction = 1
1093/// - For upper limit use SetLeftSideTailFraction = 0
1094/// - for shortest intervals use SetLeftSideTailFraction = -1 or call the method SetShortestInterval()
1095///
1096/// NOTE: The BayesianCalculator covers only the case with one
1097/// single parameter of interest
1098///
1099/// NOTE: User takes ownership of the returned object
1100
1102{
1103
1104 if (fValidInterval)
1105 coutW(Eval) << "BayesianCalculator::GetInterval - recomputing interval for the same CL and same model" << std::endl;
1106
1107 RooRealVar* poi = dynamic_cast<RooRealVar*>( fPOI.first() );
1108 if (!poi) {
1109 coutE(Eval) << "BayesianCalculator::GetInterval - no parameter of interest is set " << std::endl;
1110 return 0;
1111 }
1112
1113
1114
1115 // get integrated likelihood (posterior function)
1117
1118 //bool silentMode = (RooMsgService::instance().globalKillBelow() >= RooFit::ERROR || RooMsgService::instance().silentMode()) ;
1120
1121 if (fLeftSideFraction < 0 ) {
1122 // compute short intervals
1124 }
1125 else {
1126 // compute the other intervals
1127
1128 double lowerCutOff = fLeftSideFraction * fSize;
1129 double upperCutOff = 1. - (1.- fLeftSideFraction) * fSize;
1130
1131
1132 if (fNScanBins > 0) {
1133 ComputeIntervalFromApproxPosterior(lowerCutOff, upperCutOff);
1134 }
1135
1136 else {
1137 // use integration method if there are nuisance parameters
1138 if (fNuisanceParameters.getSize() > 0) {
1139 ComputeIntervalFromCdf(lowerCutOff, upperCutOff);
1140 }
1141 else {
1142 // case of no nuisance - just use createCdf from roofit
1143 ComputeIntervalUsingRooFit(lowerCutOff, upperCutOff);
1144 }
1145 // case cdf failed (scan then the posterior)
1146 if (!fValidInterval) {
1147 fNScanBins = 100;
1148 coutW(Eval) << "BayesianCalculator::GetInterval - computing integral from cdf failed - do a scan in "
1149 << fNScanBins << " nbins " << std::endl;
1150 ComputeIntervalFromApproxPosterior(lowerCutOff, upperCutOff);
1151 }
1152 }
1153 }
1154
1155
1156 // reset the counts and default mode
1157 if (RooAbsReal::numEvalErrors() > 0)
1158 coutW(Eval) << "BayesianCalculator::GetInterval : " << RooAbsReal::numEvalErrors() << " errors reported in evaluating log-likelihood function "
1159 << std::endl;
1160
1163
1164 if (!fValidInterval) {
1165 fLower = 1; fUpper = 0;
1166 coutE(Eval) << "BayesianCalculator::GetInterval - cannot compute a valid interval - return a dummy [1,0] interval"
1167 << std::endl;
1168 }
1169 else {
1170 coutI(Eval) << "BayesianCalculator::GetInterval - found a valid interval : [" << fLower << " , "
1171 << fUpper << " ]" << std::endl;
1172 }
1173
1174 TString interval_name = TString("BayesianInterval_a") + TString(this->GetName());
1175 SimpleInterval * interval = new SimpleInterval(interval_name,*poi,fLower,fUpper,ConfidenceLevel());
1176 interval->SetTitle("SimpleInterval from BayesianCalculator");
1177
1178 return interval;
1179}
1180
1181////////////////////////////////////////////////////////////////////////////////
1182/// Returns the value of the parameter for the point in
1183/// parameter-space that is the most likely.
1184/// How do we do if there are points that are equi-probable?
1185/// use approximate posterior
1186/// t.b.d use real function to find the mode
1187
1189
1192 return h->GetBinCenter(h->GetMaximumBin() );
1193 //return fApproxPosterior->GetMaximumX();
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// internal function compute the interval using RooFit
1198
1199void BayesianCalculator::ComputeIntervalUsingRooFit(double lowerCutOff, double upperCutOff ) const {
1200
1201 coutI(Eval) << "BayesianCalculator: Compute interval using RooFit: posteriorPdf + createCdf + RooBrentRootFinder " << std::endl;
1202
1203 RooRealVar* poi = dynamic_cast<RooRealVar*>( fPOI.first() );
1204 assert(poi);
1205
1206 fValidInterval = false;
1208 if (!fPosteriorPdf) return;
1209
1210 std::unique_ptr<RooAbsReal> cdf{fPosteriorPdf->createCdf(fPOI,RooFit::ScanNoCdf())};
1211 if (!cdf) return;
1212
1213 std::unique_ptr<RooAbsFunc> cdf_bind{cdf->bindVars(fPOI,&fPOI)};
1214 if (!cdf_bind) return;
1215
1216 RooBrentRootFinder brf(*cdf_bind);
1217 brf.setTol(fBrfPrecision); // set the brf precision
1218
1219 double tmpVal = poi->getVal(); // patch used because findRoot changes the value of poi
1220
1221 bool ret = true;
1222 if (lowerCutOff > 0) {
1223 double y = lowerCutOff;
1224 ret &= brf.findRoot(fLower,poi->getMin(),poi->getMax(),y);
1225 }
1226 else
1227 fLower = poi->getMin();
1228
1229 if (upperCutOff < 1.0) {
1230 double y=upperCutOff;
1231 ret &= brf.findRoot(fUpper,poi->getMin(),poi->getMax(),y);
1232 }
1233 else
1234 fUpper = poi->getMax();
1235 if (!ret) coutE(Eval) << "BayesianCalculator::GetInterval "
1236 << "Error returned from Root finder, estimated interval is not fully correct"
1237 << std::endl;
1238 else
1239 fValidInterval = true;
1240
1241
1242 poi->setVal(tmpVal); // patch: restore the original value of poi
1243}
1244
1245////////////////////////////////////////////////////////////////////////////////
1246/// internal function compute the interval using Cdf integration
1247
1248void BayesianCalculator::ComputeIntervalFromCdf(double lowerCutOff, double upperCutOff ) const {
1249
1250 fValidInterval = false;
1251
1252 coutI(InputArguments) << "BayesianCalculator:GetInterval Compute the interval from the posterior cdf " << std::endl;
1253
1254 RooRealVar* poi = dynamic_cast<RooRealVar*>( fPOI.first() );
1255 assert(poi);
1256 if (GetPosteriorFunction() == 0) {
1257 coutE(InputArguments) << "BayesianCalculator::GetInterval() cannot make posterior Function " << std::endl;
1258 return;
1259 }
1260
1261 // need to remove the constant parameters
1262 RooArgList bindParams;
1263 bindParams.add(fPOI);
1264 bindParams.add(fNuisanceParameters);
1265
1266 // this code could be put inside the PosteriorCdfFunction
1267
1268 //bindParams.Print("V");
1269
1271 if( cdf.HasError() ) {
1272 coutE(Eval) << "BayesianCalculator: Numerical error computing CDF integral - try a different method " << std::endl;
1273 return;
1274 }
1275
1276 //find the roots
1277
1279
1280 ccoutD(Eval) << "BayesianCalculator::GetInterval - finding roots of posterior using RF " << rf.Name()
1281 << " with precision = " << fBrfPrecision;
1282
1283 if (lowerCutOff > 0) {
1284 cdf.SetOffset(lowerCutOff);
1285 ccoutD(NumIntegration) << "Integrating posterior to get cdf and search lower limit at p =" << lowerCutOff << std::endl;
1286 bool ok = rf.Solve(cdf, poi->getMin(),poi->getMax() , 200,fBrfPrecision, fBrfPrecision);
1287 if( cdf.HasError() )
1288 coutW(Eval) << "BayesianCalculator: Numerical error integrating the CDF " << std::endl;
1289 if (!ok) {
1290 coutE(NumIntegration) << "BayesianCalculator::GetInterval - Error from root finder when searching lower limit !" << std::endl;
1291 return;
1292 }
1293 fLower = rf.Root();
1294 }
1295 else {
1296 fLower = poi->getMin();
1297 }
1298 if (upperCutOff < 1.0) {
1299 cdf.SetOffset(upperCutOff);
1300 ccoutD(NumIntegration) << "Integrating posterior to get cdf and search upper interval limit at p =" << upperCutOff << std::endl;
1301 bool ok = rf.Solve(cdf, fLower,poi->getMax() , 200, fBrfPrecision, fBrfPrecision);
1302 if( cdf.HasError() )
1303 coutW(Eval) << "BayesianCalculator: Numerical error integrating the CDF " << std::endl;
1304 if (!ok) {
1305 coutE(NumIntegration) << "BayesianCalculator::GetInterval - Error from root finder when searching upper limit !" << std::endl;
1306 return;
1307 }
1308 fUpper = rf.Root();
1309 }
1310 else {
1311 fUpper = poi->getMax();
1312 }
1313
1314 fValidInterval = true;
1315}
1316
1317////////////////////////////////////////////////////////////////////////////////
1318/// approximate posterior in nbins using a TF1
1319/// scan the poi values and evaluate the posterior at each point
1320/// and save the result in a cloned TF1
1321/// For each point the posterior is evaluated by integrating the nuisance
1322/// parameters
1323
1325
1326 if (fApproxPosterior) {
1327 // if number of bins of existing function is >= requested one - no need to redo the scan
1328 if (fApproxPosterior->GetNpx() >= fNScanBins) return;
1329 // otherwise redo the scan
1330 delete fApproxPosterior;
1331 fApproxPosterior = 0;
1332 }
1333
1334
1335 RooAbsReal * posterior = GetPosteriorFunction();
1336 if (!posterior) return;
1337
1338
1339 TF1 * tmp = posterior->asTF(fPOI);
1340 assert(tmp != 0);
1341 // binned the function in nbins and evaluate at those points
1342 if (fNScanBins > 0) tmp->SetNpx(fNScanBins); // if not use default of TF1 (which is 100)
1343
1344 coutI(Eval) << "BayesianCalculator - scan posterior function in nbins = " << tmp->GetNpx() << std::endl;
1345
1346 fApproxPosterior = (TF1*) tmp->Clone();
1347 // save this function for future reuse
1348 // I can delete now original posterior and use this approximated copy
1349 delete tmp;
1350 TString name = posterior->GetName() + TString("_approx");
1351 TString title = posterior->GetTitle() + TString("_approx");
1352 RooAbsReal * posterior2 = new RooTFnBinding(name,title,fApproxPosterior,fPOI);
1353 if (posterior == fIntegratedLikelihood) {
1354 delete fIntegratedLikelihood;
1355 fIntegratedLikelihood = posterior2;
1356 }
1357 else if (posterior == fLikelihood) {
1358 delete fLikelihood;
1359 fLikelihood = posterior2;
1360 }
1361 else {
1362 assert(1); // should never happen this case
1363 }
1364}
1365
1366////////////////////////////////////////////////////////////////////////////////
1367/// compute the interval using the approximate posterior function
1368
1369void BayesianCalculator::ComputeIntervalFromApproxPosterior(double lowerCutOff, double upperCutOff ) const {
1370
1371 ccoutD(Eval) << "BayesianCalculator: Compute interval from the approximate posterior " << std::endl;
1372
1374 if (!fApproxPosterior) return;
1375
1376 double prob[2];
1377 double limits[2] = {0,0};
1378 prob[0] = lowerCutOff;
1379 prob[1] = upperCutOff;
1380 fApproxPosterior->GetQuantiles(2,limits,prob);
1381 fLower = limits[0];
1382 fUpper = limits[1];
1383 fValidInterval = true;
1384}
1385
1386////////////////////////////////////////////////////////////////////////////////
1387/// compute the shortest interval from the histogram representing the posterior
1388
1389
1391 coutI(Eval) << "BayesianCalculator - computing shortest interval with CL = " << 1.-fSize << std::endl;
1392
1393 // compute via the approx posterior function
1395 if (!fApproxPosterior) return;
1396
1397 TH1D * h1 = dynamic_cast<TH1D*>(fApproxPosterior->GetHistogram() );
1398 assert(h1 != 0);
1400 // get bins and sort them
1401 double * bins = h1->GetArray();
1402 // exclude under/overflow bins
1403 int n = h1->GetSize()-2;
1404 std::vector<int> index(n);
1405 // exclude bins[0] (the underflow bin content)
1406 TMath::Sort(n, bins+1, &index[0]);
1407 // find cut off as test size
1408 double sum = 0;
1409 double actualCL = 0;
1410 double upper = h1->GetXaxis()->GetXmin();
1411 double lower = h1->GetXaxis()->GetXmax();
1412 double norm = h1->GetSumOfWeights();
1413
1414 for (int i = 0; i < n; ++i) {
1415 int idx = index[i];
1416 double p = bins[ idx] / norm;
1417 sum += p;
1418 if (sum > 1.-fSize ) {
1419 actualCL = sum - p;
1420 break;
1421 }
1422
1423 // histogram bin content starts from 1
1424 if ( h1->GetBinLowEdge(idx+1) < lower)
1425 lower = h1->GetBinLowEdge(idx+1);
1426 if ( h1->GetXaxis()->GetBinUpEdge(idx+1) > upper)
1427 upper = h1->GetXaxis()->GetBinUpEdge(idx+1);
1428 }
1429
1430 ccoutD(Eval) << "BayesianCalculator::ComputeShortestInterval - actual interval CL = "
1431 << actualCL << " difference from requested is " << (actualCL-(1.-fSize))/fSize*100. << "% "
1432 << " limits are [ " << lower << " , " << " upper ] " << std::endl;
1433
1434
1435 if (lower < upper) {
1436 fLower = lower;
1437 fUpper = upper;
1438
1439
1440
1441 if ( std::abs(actualCL-(1.-fSize)) > 0.1*(1.-fSize) )
1442 coutW(Eval) << "BayesianCalculator::ComputeShortestInterval - actual interval CL = "
1443 << actualCL << " differs more than 10% from desired CL value - must increase nbins "
1444 << n << " to an higher value " << std::endl;
1445 }
1446 else
1447 coutE(Eval) << "BayesianCalculator::ComputeShortestInterval " << n << " bins are not sufficient " << std::endl;
1448
1449 fValidInterval = true;
1450
1451}
1452
1453
1454
1455
1456} // end namespace RooStats
size_t fSize
#define f(i)
Definition RSha256.hxx:104
#define h(i)
Definition RSha256.hxx:106
#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)
#define ClassImp(name)
Definition Rtypes.h:377
@ kGray
Definition Rtypes.h:65
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:110
User class for performing function minimization.
void SetFunction(const ROOT::Math::IGenFunction &f, double xlow, double xup)
Sets function to be minimized.
bool Minimize(int maxIter, double absTol=1.E-8, double relTol=1.E-10) override
Find minimum position iterating until convergence specified by the absolute and relative tolerance or...
double FValMinimum() const override
Return function value at current estimate of the minimum.
Functor1D class for one-dimensional functions.
Definition Functor.h:578
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition IFunction.h:135
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
User class for performing multidimensional integration.
double Integral(const double *xmin, const double *xmax)
evaluate the integral with the previously given function between xmin[] and xmax[]
void SetFunction(Function &f, unsigned int dim)
set integration function using a generic function implementing the operator()(double *x) The dimensio...
static IntegrationMultiDim::Type GetType(const char *name)
static function to get the enumeration from a string
ROOT::Math::IntegratorMultiDimOptions Options() const
retrieve the options
double Error() const
return integration error
User Class for performing numerical integration of a function in one dimension.
Definition Integrator.h:98
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:84
const char * Name() const
Return the current and latest estimate of the lower value of the Root-finding interval (for bracketin...
Definition RootFinder.h:201
bool Solve(Function &f, Derivative &d, double start, int maxIter=100, double absTol=1E-8, double relTol=1E-10)
Solve f(x) = 0, given a derivative d.
Definition RootFinder.h:235
double Root() const
Return the current and latest estimate of the Root.
Definition RootFinder.h:171
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
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 expresssion tree)
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Storage_t const & get() const
Const access to the underlying stl container.
Int_t getSize() const
Return the number of elements in the collection.
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * first() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:59
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
RooFit::OwningPtr< RooAbsReal > createCdf(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Create a cumulative distribution function of this p.d.f in terms of the observables listed in iset.
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:60
virtual RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, const RooLinkedList &cmdList={})
Construct representation of -log(L) of PDF with given dataset.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
RooPlot * frame(const RooCmdArg &arg1, 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()) const
Create a new RooPlot on the heap with a drawing frame initialized for this object,...
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, 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()) const
Create an object that represents the integral of the function over one or more observables std::liste...
TF1 * asTF(const RooArgList &obs, const RooArgList &pars=RooArgList(), const RooArgSet &nset=RooArgSet()) const
Return a ROOT TF1,2,3 object bound to this RooAbsReal with given definition of observables and parame...
virtual RooPlot * plotOn(RooPlot *frame, 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(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
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.
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:55
Implement the abstract 1-dimensional root finding interface using the Brent-Decker method.
void setTol(double tol)
Set convergence tolerance parameter.
bool findRoot(double &result, double xlo, double xhi, double value=0) const override
Do the root finding using the Brent-Decker method.
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
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:33
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setVal(double value) override
Set value of variable to 'value'.
BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation of a...
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 ComputeIntervalUsingRooFit(double c1, double c2) const
internal function compute the interval using RooFit
void ComputeIntervalFromCdf(double c1, double c2) const
internal function compute the interval using Cdf integration
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
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)
PosteriorCdfFunction(RooAbsReal &nll, RooArgList &bindParams, RooAbsReal *prior=0, const char *integType=0, double nllMinimum=0)
std::shared_ptr< RooFunctor > fPriorFunc
PosteriorCdfFunction & operator=(const PosteriorCdfFunction &)
PosteriorCdfFunction(const PosteriorCdfFunction &rhs)
ROOT::Math::IntegratorMultiDim fIntegrator
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.
PosteriorFunctionFromToyMC(RooAbsReal &nll, RooAbsPdf &pdf, RooRealVar &poi, RooArgList &nuisParams, RooAbsReal *prior=0, double nllOffset=0, int niter=0, bool redoToys=true)
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
ROOT::Math::IGenFunction * Clone() const override
Clone a function.
std::shared_ptr< RooFunctor > fPriorFunc
std::unique_ptr< ROOT::Math::IntegratorMultiDim > fIntegratorMultiDim
PosteriorFunction(RooAbsReal &nll, RooRealVar &poi, RooArgList &nuisParams, RooAbsReal *prior=0, const char *integType=0, double norm=1.0, double nllOffset=0, int niter=0)
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:135
Double_t GetXmin() const
Definition TAxis.h:134
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:528
1-Dim function class
Definition TF1.h:213
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function Note that this histogram is managed ...
Definition TF1.cxx:1585
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition TF1.cxx:3431
virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum)
Compute Quantiles for density distribution of this function.
Definition TF1.cxx:1993
TObject * Clone(const char *newname=nullptr) const override
Make a complete copy of the underlying object.
Definition TF1.cxx:1065
virtual Int_t GetNpx() const
Definition TF1.h:495
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:620
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition TH1.cxx:9007
TAxis * GetXaxis()
Definition TH1.h:322
virtual Double_t GetBinLowEdge(Int_t bin) const
Return bin lower edge for 1D histogram.
Definition TH1.cxx:9018
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:8826
virtual Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition TH1.cxx:7827
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Basic string class.
Definition TString.h:139
void ToUpper()
Change string to upper case.
Definition TString.cxx:1183
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2334
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
RooCmdArg ScanNoCdf()
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg GlobalObservables(Args_t &&... argsOrArgSet)
RooCmdArg ConditionalObservables(Args_t &&... argsOrArgSet)
Create a RooCmdArg to declare conditional observables.
RooCmdArg FillColor(Color_t color)
RooCmdArg Precision(double prec)
RooCmdArg DrawOption(const char *opt)
RooCmdArg Range(const char *rangeName, bool adjustNorm=true)
RooCmdArg MoveToBack()
RooCmdArg VLines()
Double_t y[n]
Definition legend1.C:17
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.
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Namespace for the RooStats classes.
Definition Asimov.h:19
void RemoveConstantParameters(RooArgSet *set)
const ROOT::Math::RootFinder::EType kRootFinderType
Bool_t IsNaN(Double_t x)
Definition TMath.h:890
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:431
void SetPrior(RooFunctor *prior)
LikelihoodFunction(RooFunctor &f, RooFunctor *prior=0, double offset=0)
double operator()(const double *x) const
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345