Logo ROOT   6.16/01
Reference Guide
BayesianCalculator.h
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#ifndef ROOSTATS_BayesianCalculator
12#define ROOSTATS_BayesianCalculator
13
14#include "TNamed.h"
15
16#include "Math/IFunctionfwd.h"
17
18#include "RooArgSet.h"
19
21
23
24class RooAbsData;
25class RooAbsPdf;
26class RooPlot;
27class RooAbsReal;
28class TF1;
29
30
31namespace RooStats {
32
33 class ModelConfig;
34 class SimpleInterval;
35
37
38 public:
39
40 // constructor
42
44 RooAbsPdf& pdf,
45 const RooArgSet& POI,
46 RooAbsPdf& priorPdf,
47 const RooArgSet* nuisanceParameters = 0 );
48
51
52 // destructor
53 virtual ~BayesianCalculator();
54
55 // get the plot with option to get it normalized
56 RooPlot* GetPosteriorPlot(bool norm = false, double precision = 0.01) const;
57
58 // return posterior pdf (object is managed by the BayesianCalculator class)
60 // return posterior function (object is managed by the BayesianCalculator class)
62
63 // compute the interval. By Default a central interval is computed
64 // By using SetLeftTileFraction can control if central/ upper/lower interval
65 // For shortest interval use SetShortestInterval(true)
66 virtual SimpleInterval* GetInterval() const ;
67
68 virtual void SetData( RooAbsData & data ) {
69 fData = &data;
70 ClearAll();
71 }
72
73
74 // set the model via the ModelConfig
75 virtual void SetModel( const ModelConfig& model );
76
77 // specify the parameters of interest in the interval
78 virtual void SetParameters(const RooArgSet& set) { fPOI.removeAll(); fPOI.add(set); }
79
80 // specify the nuisance parameters (eg. the rest of the parameters)
82
83 // Set only the Prior Pdf
84 virtual void SetPriorPdf(RooAbsPdf& pdf) { fPriorPdf = &pdf; }
85
86 // set the conditional observables which will be used when creating the NLL
87 // so the pdf's will not be normalized on the conditional observables when computing the NLL
89
90 // set the global observables which will be used when creating the NLL
91 // so the constraint pdf's will be normalized correctly on the global observables when computing the NLL
93
94 // set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval)
95 virtual void SetTestSize( Double_t size ) {
96 fSize = size;
97 fValidInterval = false;
98 }
99 // set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval)
100 virtual void SetConfidenceLevel( Double_t cl ) { SetTestSize(1.-cl); }
101 // Get the size of the test (eg. rate of Type I error)
102 virtual Double_t Size() const { return fSize; }
103 // Get the Confidence level for the test
104 virtual Double_t ConfidenceLevel() const { return 1.-fSize; }
105
106 // set the fraction of probability content on the left tail
107 // Central limits use 0.5 (default case)
108 // for upper limits it is 0 and 1 for lower limit
109 // For shortest intervals a negative value (i.e. -1) must be given
110 void SetLeftSideTailFraction(Double_t leftSideFraction ) {fLeftSideFraction = leftSideFraction;}
111
112 // set the Bayesian calculator to compute the shortest interval (default is central interval)
113 // to switch off SetLeftSideTailFraction to the right value
115
116 // set the precision of the Root Finder
117 void SetBrfPrecision( double precision ) { fBrfPrecision = precision; }
118
119 // use directly the approximate posterior function obtained by binning it in nbins
120 // by default the cdf is used by integrating the posterior
121 // if a value of nbin <= 0 the cdf function will be used
122 void SetScanOfPosterior(int nbin = 100) { fNScanBins = nbin; }
123
124 // set the number of iterations when running a MC integration algorithm
125 // If not set use default algorithmic values
126 // In case of ToyMC sampling of the nuisance the value is 100
127 // In case of using the GSL MCintegrations types the default value is
128 // defined in ROOT::Math::IntegratorMultiDimOptions::DefaultNCalls()
129 virtual void SetNumIters(Int_t numIters) { fNumIterations = numIters; }
130
131 // set the integration type (possible type are) :
132 void SetIntegrationType(const char * type);
133
134 // return the mode (most probable value of the posterior function)
135 double GetMode() const;
136
137 // force the nuisance pdf when using the toy mc sampling
139
140 protected:
141
142 void ClearAll() const;
143
144 void ApproximatePosterior() const;
145
146 void ComputeIntervalFromApproxPosterior(double c1, double c2) const;
147
148 void ComputeIntervalFromCdf(double c1, double c2) const;
149
150 void ComputeIntervalUsingRooFit(double c1, double c2) const;
151
152 void ComputeShortestInterval() const;
153
154 private:
155
156 // plan to replace the above: return a SimpleInterval integrating
157 // over all other parameters except the one specified as argument
158 // virtual SimpleInterval* GetInterval( RooRealVar* parameter ) const { return 0; }
159
160 RooAbsData* fData; // data set
161 RooAbsPdf* fPdf; // model pdf (could contain the nuisance pdf as constraint term)
163 RooAbsPdf* fPriorPdf; // prior pdf (typically for the POI)
164 RooAbsPdf* fNuisancePdf; // nuisance pdf (needed when using nuisance sampling technique)
165 RooArgSet fNuisanceParameters; // nuisance parameters
166 RooArgSet fConditionalObs ; // conditional observables
167 RooArgSet fGlobalObs; // global observables
168
169 mutable RooAbsPdf* fProductPdf; // internal pointer to model * prior
170 mutable RooAbsReal* fLogLike; // internal pointer to log likelihood function
171 mutable RooAbsReal* fLikelihood; // internal pointer to likelihood function
172 mutable RooAbsReal* fIntegratedLikelihood; // integrated likelihood function, i.e - unnormalized posterior function
173 mutable RooAbsPdf* fPosteriorPdf; // normalized (on the poi) posterior pdf
174 mutable ROOT::Math::IGenFunction * fPosteriorFunction; // function representing the posterior
175 mutable TF1 * fApproxPosterior; // TF1 representing the scanned posterior function
176 mutable Double_t fLower; // computer lower interval bound
177 mutable Double_t fUpper; // upper interval bound
178 mutable Double_t fNLLMin; // minimum value of Nll
179 double fSize; // size used for getting the interval
180 double fLeftSideFraction; // fraction of probability content on left side of interval
181 double fBrfPrecision; // root finder precision
182 mutable int fNScanBins; // number of bins to scan, if = -1 no scan is done (default)
183 int fNumIterations; // number of iterations (when using ToyMC)
185
187
188 protected:
189
190 ClassDef(BayesianCalculator,2) // BayesianCalculator class
191
192 };
193}
194
195#endif
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
#define ClassDef(name, id)
Definition: Rtypes.h:324
int type
Definition: TGX11.cxx:120
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition: IFunction.h:135
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation of a...
void SetBrfPrecision(double precision)
void ForceNuisancePdf(RooAbsPdf &pdf)
ROOT::Math::IGenFunction * fPosteriorFunction
virtual void SetConditionalObservables(const RooArgSet &set)
virtual void SetData(RooAbsData &data)
Set the DataSet ( add to the the workspace if not already there ?)
RooAbsPdf * GetPosteriorPdf() const
Build and return the posterior pdf (i.e posterior function normalized to all range of poi) Note that ...
virtual SimpleInterval * GetInterval() const
Compute the interval.
RooPlot * GetPosteriorPlot(bool norm=false, double precision=0.01) const
return a RooPlot with the posterior and the credibility region NOTE: User takes ownership of the retu...
void ClearAll() const
clear all cached pdf objects
virtual Double_t ConfidenceLevel() const
Get the Confidence level for the test.
void SetScanOfPosterior(int nbin=100)
void ComputeShortestInterval() const
compute the shortest interval
virtual Double_t Size() const
Get the size of the test (eg. rate of Type I error)
virtual void SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( e.g. 0.05 for a 95% Confidence Interval)
void ApproximatePosterior() const
approximate posterior in nbins using a TF1 scan the poi values and evaluate the posterior at each poi...
virtual void SetModel(const ModelConfig &model)
set the model to use The model pdf, prior pdf, parameter of interest and nuisances will be taken acco...
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (e.g. 0.95 for a 95% Confidence Interval)
virtual void SetNuisanceParameters(const RooArgSet &set)
double GetMode() const
Returns the value of the parameter for the point in parameter-space that is the most likely.
virtual void SetGlobalObservables(const RooArgSet &set)
virtual void SetNumIters(Int_t numIters)
virtual void SetPriorPdf(RooAbsPdf &pdf)
RooAbsReal * GetPosteriorFunction() const
Build and return the posterior function (not normalized) as a RooAbsReal the posterior is obtained fr...
virtual void SetParameters(const RooArgSet &set)
void SetIntegrationType(const char *type)
set the integration type (possible type are) :
void SetLeftSideTailFraction(Double_t leftSideFraction)
void ComputeIntervalFromApproxPosterior(double c1, double c2) const
compute the interval using the approximate posterior function
BayesianCalculator()
default constructor
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
IntervalCalculator is an interface class for a tools which produce RooStats ConfIntervals.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
SimpleInterval is a concrete implementation of the ConfInterval interface.
1-Dim function class
Definition: TF1.h:211
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
Basic string class.
Definition: TString.h:131
return c1
Definition: legend1.C:41
return c2
Definition: legend2.C:14
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20