ROOT logo
ROOT » ROOFIT » ROOSTATS » RooStats::BayesianCalculator

class RooStats::BayesianCalculator: public RooStats::IntervalCalculator, public TNamed


BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation of a credible interval using a Bayesian method. The class works only for one single parameter of interest and it integrates the likelihood function with the given prior probability density function to compute the posterior probability. The result of the class is a one dimensional interval (class SimpleInterval ), which is obtained from inverting the cumulative posterior distribution.

The interface allows one to construct the class by passing the data set, probability density function for the model, the prior functions and then the parameter of interest to scan. The nuisance parameters can also be passed to be marginalized when computing the posterior. Alternatively, the class can be constructed by passing the data and the ModelConfig containing all the needed information (model pdf, prior pdf, parameter of interest, nuisance parameters, etc..)

After configuring the calculator, one only needs to ask GetInterval(), which will return an SimpleInterval object. By default the extrem of the integral are obtained by inverting directly the cumulative posterior distribution. By using the method SetScanOfPosterior(nbins) the interval is then obtained by scanning the posterior function in the given number of points. The firts method is in general faster but it requires an integration one extra dimension ( in the poi in addition to the nuisance parameters), therefore in some case it can be less robust.

The class can also return the posterior function (method GetPosteriorFunction) or if needed the normalized posterior function (the posterior pdf) (method GetPosteriorPdf). A posterior plot is also obtained using the GetPosteriorPlot method.

The class allows to use different integration methods for integrating in the nuisances and in the poi. All the numerical integration methods of ROOT can be used via the method SetIntegrationType (see more in the documentation of this method).



Function Members (Methods)

public:
virtual~BayesianCalculator()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
RooStats::BayesianCalculatorBayesianCalculator()
RooStats::BayesianCalculatorBayesianCalculator(const RooStats::BayesianCalculator&)
RooStats::BayesianCalculatorBayesianCalculator(RooAbsData& data, RooStats::ModelConfig& model)
RooStats::BayesianCalculatorBayesianCalculator(RooAbsData& data, RooAbsPdf& pdf, const RooArgSet& POI, RooAbsPdf& priorPOI, const RooArgSet* nuisanceParameters = 0)
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
virtual TObject*TNamed::Clone(const char* newname = "") const
virtual Int_tTNamed::Compare(const TObject* obj) const
virtual Double_tConfidenceLevel() const
virtual voidTNamed::Copy(TObject& named) const
virtual voidTObject::Delete(Option_t* option = "")MENU
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual voidTNamed::FillBuffer(char*& buffer)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
voidForceNuisancePdf(RooAbsPdf& pdf)
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
virtual RooStats::SimpleInterval*GetInterval() const
doubleGetMode() const
virtual const char*TNamed::GetName() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
RooAbsReal*GetPosteriorFunction() const
RooAbsPdf*GetPosteriorPdf() const
RooPlot*GetPosteriorPlot(bool norm = false, double precision = 0.01) const
virtual const char*TNamed::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() constMENU
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTNamed::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTNamed::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
virtual Bool_tTObject::Notify()
voidTObject::Obsolete(const char* method, const char* asOfVers, const char* removedFromVers) const
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
RooStats::BayesianCalculator&operator=(const RooStats::BayesianCalculator&)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidTNamed::Print(Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
voidSetBrfPrecision(double precision)
virtual voidSetConfidenceLevel(Double_t cl)
virtual voidSetData(RooAbsData& data)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
voidSetIntegrationType(const char* type)
voidSetLeftSideTailFraction(Double_t leftSideFraction)
virtual voidSetModel(const RooStats::ModelConfig& model)
virtual voidTNamed::SetName(const char* name)MENU
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
virtual voidSetNumIters(Int_t numIters)
static voidTObject::SetObjectStat(Bool_t stat)
voidSetScanOfPosterior(int nbin = 100)
voidSetShortestInterval()
virtual voidSetTestSize(Double_t size)
virtual voidTNamed::SetTitle(const char* title = "")MENU
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector& insp)
virtual Double_tSize() const
virtual Int_tTNamed::Sizeof() const
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
voidApproximatePosterior() const
voidClearAll() const
voidComputeIntervalFromApproxPosterior(double c1, double c2) const
voidComputeIntervalFromCdf(double c1, double c2) const
voidComputeIntervalUsingRooFit(double c1, double c2) const
voidComputeShortestInterval() const
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
voidTObject::MakeZombie()

Data Members

protected:
TStringTNamed::fNameobject identifier
TStringTNamed::fTitleobject title
private:
TF1*fApproxPosteriorTF1 representing the scanned posterior function
doublefBrfPrecisionroot finder precision
RooAbsData*fDatadata set
RooAbsReal*fIntegratedLikelihoodintegrated likelihood function, i.e - unnormalized posterior function
TStringfIntegrationType
doublefLeftSideFractionfraction of probability content on left side of interval
RooAbsReal*fLikelihoodinternal pointer to likelihood function
RooAbsReal*fLogLikeinternal pointer to log likelihood function
Double_tfLowercomputer lower interval bound
Double_tfNLLMinminimum value of Nll
intfNScanBinsnumber of bins to scan, if = -1 no scan is done (default)
RooArgSetfNuisanceParameters
RooAbsPdf*fNuisancePdfnuisance pdf (needed when using nuisance sampling technique)
intfNumIterationsnumber of iterations (when using ToyMC)
RooArgSetfPOIPOI
RooAbsPdf*fPdfmodel pdf (could contain the nuisance pdf as constraint term)
ROOT::Math::IBaseFunctionOneDim*fPosteriorFunctionfunction representing the posterior
RooAbsPdf*fPosteriorPdfnormalized (on the poi) posterior pdf
RooAbsPdf*fPriorPOIprior pdf for POI
RooAbsPdf*fProductPdfinternal pointer to model * prior
doublefSizesize used for getting the interval
Double_tfUpperupper interval bound
Bool_tfValidInterval

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

BayesianCalculator()
 default constructor
BayesianCalculator( /* const char* name, const char* title, */ RooAbsData& data, RooAbsPdf& pdf, const RooArgSet& POI, RooAbsPdf& priorPOI, const RooArgSet* nuisanceParameters )
 Constructor from data set, model pdf, parameter of interests and prior pdf
 If nuisance parameters are given they will be integrated according either to the prior or
 their constraint term included in the model
BayesianCalculator(RooAbsData& data, RooStats::ModelConfig& model)
 Constructor from a data set and a ModelConfig
 model pdf, poi and nuisances will be taken from the ModelConfig
~BayesianCalculator()
 destructor
void ClearAll() const
 clear all cached pdf objects
void SetModel(const RooStats::ModelConfig& model)
 set the model to use
 The model pdf, prior pdf, parameter of interest and nuisances
 will be taken according to the model
RooAbsReal* GetPosteriorFunction() const
 Build and return the posterior function (not normalized) as a RooAbsReal
 the posterior is obtained from the product of the likelihood function and the
 prior pdf which is then intergated in the nuisance parameters (if existing).
 A prior function for the nuisance can be specified either in the prior pdf object
 or in the model itself. If no prior nuisance is specified, but prior parameters are then
 the integration is performed assuming a flat prior for the nuisance parameters.
 NOTE: the return object is managed by the class, users do not need to delete it
RooAbsPdf* GetPosteriorPdf() const
 Build and return the posterior pdf (i.e posterior function normalized to all range of poi)
 Note that an extra integration in the POI is required for the normalization
 NOTE: user must delete the returned object
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 returned object
void SetIntegrationType(const char* type)
 set the integration type (possible type are) :
 1D integration ( used when only one nuisance and when the posterior is scanned):
    adaptive , gauss, nonadaptive
 multidim:
       ADAPTIVE,   adaptive numerical integration
                    The parameter numIters (settable with SetNumIters) is  the max number of function calls.
                    It can be reduced to make teh integration faster but it will be difficult to reach the required tolerance
       VEGAS   MC integration method based on importance sampling - numIters is number of function calls
                Extra Vegas parameter can be set using  IntegratorMultiDimOptions class
       MISER    MC integration method based on stratified sampling
                See also http://en.wikipedia.org/wiki/Monte_Carlo_integration for VEGAS and MISER description
       PLAIN    simple MC integration method, where the max  number of calls can be specified using SetNumIters(numIters)

 Extra integration types are:

    TOYMC:
       evaluate posterior by generating toy MC for the nuisance parameters. It is a MC
       integration, where the function is sampled according to the nuisance. It is convenient to use when all
       the nuisance are uncorrelated and it is efficient to generate them
       The toy are generated by default for each  poi values
       (this method has been proposed and provided by J.P Chou)
    1-TOYMC  : same method as before but in this case the toys are generated only one time and then used for
       each poi value. It can be convenient when the generation time is much larger than the evaluation time,
       otherwise it is recoomended to re-generate the toy for each poi scanned point of the posterior function


    ROOFIT:
       use roofit default integration methods which will produce a nested integral (not reccomended for more
       than 1 nuisance parameters)

 if type = 0 use default specified via class IntegratorMultiDimOptions::SetDefaultIntegrator
SimpleInterval* GetInterval() const
 Compute the interval. By Default a central interval is computed
 and the result is a SimpleInterval object.
 Using the method (to be called before SetInterval) SetLeftSideTailFraction the user can choose the type of interval.
 By default the returned interval is a central interval with the confidence level specified
 previously in the constructor ( LeftSideTailFraction = 0.5).
  For lower limit use SetLeftSideTailFraction = 1
  For upper limit use SetLeftSideTailFraction = 0
  for shortest intervals use SetLeftSideTailFraction = -1 or call the method SetShortestInterval()
 NOTE: The BayesianCalculator covers only the case with one
 single parameter of interest
 NOTE: User takes ownership of the returned object
double GetMode() const
 Returns the value of the parameter for the point in
 parameter-space that is the most likely.
  How do we do if there are points that are equi-probable?
 use approximate posterior
 t.b.d use real function to find the mode
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
void ApproximatePosterior() const
 approximate posterior in nbins using a TF1
 scan the poi values and evaluate the posterior at each point
 and save the result in a cloned TF1
 For each point the posterior is evaluated by integrating the nuisance
 parameters
void ComputeIntervalFromApproxPosterior(double c1, double c2) const
 compute the interval using the approximate posterior function
void ComputeShortestInterval() const
 compute the shortest interval
BayesianCalculator()
 constructor
void SetData(RooAbsData& data)
void SetTestSize(Double_t size)
 set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval)
void SetConfidenceLevel(Double_t cl)
 set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval)
{ SetTestSize(1.-cl); }
Double_t Size() const
 Get the size of the test (eg. rate of Type I error)
{ return fSize; }
Double_t ConfidenceLevel() const
 Get the Confidence level for the test
{ return 1.-fSize; }
void SetLeftSideTailFraction(Double_t leftSideFraction)
 set the fraction of probability content on the left tail
 Central limits use 0.5 (default case)
 for upper limits it is 0 and 1 for lower limit
 For shortest intervals a negative value (i.e. -1) must be given
{fLeftSideFraction = leftSideFraction;}
void SetShortestInterval()
 set the Bayesian calculator to compute the shorest interval (default is central interval)
 to switch off SetLeftSideTailFraction to the rght value
void SetBrfPrecision(double precision)
 set the precision of the Root Finder
{ fBrfPrecision = precision; }
void SetScanOfPosterior(int nbin = 100)
 use directly the approximate posterior function obtained by binning it in nbins
 by default the cdf is used by integrating the posterior
 if a value of nbin <= 0 the cdf function will be used
{ fNScanBins = nbin; }
void SetNumIters(Int_t numIters)
 set the number of iterations when running a MC integration algorithm
 If not set use default algorithmic values
 In case of ToyMC sampling of the nuisance the value is 100
 In case of using the GSL MCintegrations types the default value is
 defined in ROOT::Math::IntegratorMultiDimOptions::DefaultNCalls()
{ fNumIterations = numIters; }
void ForceNuisancePdf(RooAbsPdf& pdf)
 force the nuisance pdf when using the toy mc sampling
{ fNuisancePdf = &pdf; }