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).
virtual | ~BayesianCalculator() |
void | TObject::AbstractMethod(const char* method) const |
virtual void | TObject::AppendPad(Option_t* option = "") |
RooStats::BayesianCalculator | BayesianCalculator() |
RooStats::BayesianCalculator | BayesianCalculator(const RooStats::BayesianCalculator&) |
RooStats::BayesianCalculator | BayesianCalculator(RooAbsData& data, RooStats::ModelConfig& model) |
RooStats::BayesianCalculator | BayesianCalculator(RooAbsData& data, RooAbsPdf& pdf, const RooArgSet& POI, RooAbsPdf& priorPdf, const RooArgSet* nuisanceParameters = 0) |
virtual void | TObject::Browse(TBrowser* b) |
static TClass* | Class() |
virtual const char* | TObject::ClassName() const |
virtual void | TNamed::Clear(Option_t* option = "") |
virtual TObject* | TNamed::Clone(const char* newname = "") const |
virtual Int_t | TNamed::Compare(const TObject* obj) const |
virtual Double_t | ConfidenceLevel() const |
virtual void | TNamed::Copy(TObject& named) const |
virtual void | TObject::Delete(Option_t* option = "")MENU |
virtual Int_t | TObject::DistancetoPrimitive(Int_t px, Int_t py) |
virtual void | TObject::Draw(Option_t* option = "") |
virtual void | TObject::DrawClass() constMENU |
virtual TObject* | TObject::DrawClone(Option_t* option = "") constMENU |
virtual void | TObject::Dump() constMENU |
virtual void | TObject::Error(const char* method, const char* msgfmt) const |
virtual void | TObject::Execute(const char* method, const char* params, Int_t* error = 0) |
virtual void | TObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0) |
virtual void | TObject::ExecuteEvent(Int_t event, Int_t px, Int_t py) |
virtual void | TObject::Fatal(const char* method, const char* msgfmt) const |
virtual void | TNamed::FillBuffer(char*& buffer) |
virtual TObject* | TObject::FindObject(const char* name) const |
virtual TObject* | TObject::FindObject(const TObject* obj) const |
void | ForceNuisancePdf(RooAbsPdf& pdf) |
virtual Option_t* | TObject::GetDrawOption() const |
static Long_t | TObject::GetDtorOnly() |
virtual const char* | TObject::GetIconName() const |
virtual RooStats::SimpleInterval* | GetInterval() const |
double | GetMode() const |
virtual const char* | TNamed::GetName() const |
virtual char* | TObject::GetObjectInfo(Int_t px, Int_t py) const |
static Bool_t | TObject::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_t | TObject::GetUniqueID() const |
virtual Bool_t | TObject::HandleTimer(TTimer* timer) |
virtual ULong_t | TNamed::Hash() const |
virtual void | TObject::Info(const char* method, const char* msgfmt) const |
virtual Bool_t | TObject::InheritsFrom(const char* classname) const |
virtual Bool_t | TObject::InheritsFrom(const TClass* cl) const |
virtual void | TObject::Inspect() constMENU |
void | TObject::InvertBit(UInt_t f) |
virtual TClass* | IsA() const |
virtual Bool_t | TObject::IsEqual(const TObject* obj) const |
virtual Bool_t | TObject::IsFolder() const |
Bool_t | TObject::IsOnHeap() const |
virtual Bool_t | TNamed::IsSortable() const |
Bool_t | TObject::IsZombie() const |
virtual void | TNamed::ls(Option_t* option = "") const |
void | TObject::MayNotUse(const char* method) const |
virtual Bool_t | TObject::Notify() |
void | TObject::Obsolete(const char* method, const char* asOfVers, const char* removedFromVers) const |
static void | TObject::operator delete(void* ptr) |
static void | TObject::operator delete(void* ptr, void* vp) |
static void | TObject::operator delete[](void* ptr) |
static void | TObject::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 void | TObject::Paint(Option_t* option = "") |
virtual void | TObject::Pop() |
virtual void | TNamed::Print(Option_t* option = "") const |
virtual Int_t | TObject::Read(const char* name) |
virtual void | TObject::RecursiveRemove(TObject* obj) |
void | TObject::ResetBit(UInt_t f) |
virtual void | TObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU |
virtual void | TObject::SavePrimitive(ostream& out, Option_t* option = "") |
void | TObject::SetBit(UInt_t f) |
void | TObject::SetBit(UInt_t f, Bool_t set) |
void | SetBrfPrecision(double precision) |
virtual void | SetConditionalObservables(const RooArgSet& set) |
virtual void | SetConfidenceLevel(Double_t cl) |
virtual void | SetData(RooAbsData& data) |
virtual void | TObject::SetDrawOption(Option_t* option = "")MENU |
static void | TObject::SetDtorOnly(void* obj) |
void | SetIntegrationType(const char* type) |
void | SetLeftSideTailFraction(Double_t leftSideFraction) |
virtual void | SetModel(const RooStats::ModelConfig& model) |
virtual void | TNamed::SetName(const char* name)MENU |
virtual void | TNamed::SetNameTitle(const char* name, const char* title) |
virtual void | SetNuisanceParameters(const RooArgSet& set) |
virtual void | SetNumIters(Int_t numIters) |
static void | TObject::SetObjectStat(Bool_t stat) |
virtual void | SetParameters(const RooArgSet& set) |
virtual void | SetPriorPdf(RooAbsPdf& pdf) |
void | SetScanOfPosterior(int nbin = 100) |
void | SetShortestInterval() |
virtual void | SetTestSize(Double_t size) |
virtual void | TNamed::SetTitle(const char* title = "")MENU |
virtual void | TObject::SetUniqueID(UInt_t uid) |
virtual void | ShowMembers(TMemberInspector&) |
virtual Double_t | Size() const |
virtual Int_t | TNamed::Sizeof() const |
virtual void | Streamer(TBuffer&) |
void | StreamerNVirtual(TBuffer& ClassDef_StreamerNVirtual_b) |
virtual void | TObject::SysError(const char* method, const char* msgfmt) const |
Bool_t | TObject::TestBit(UInt_t f) const |
Int_t | TObject::TestBits(UInt_t f) const |
virtual void | TObject::UseCurrentStyle() |
virtual void | TObject::Warning(const char* method, const char* msgfmt) const |
virtual Int_t | TObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) |
virtual Int_t | TObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const |
void | ApproximatePosterior() const |
void | ClearAll() const |
void | ComputeIntervalFromApproxPosterior(double c1, double c2) const |
void | ComputeIntervalFromCdf(double c1, double c2) const |
void | ComputeIntervalUsingRooFit(double c1, double c2) const |
void | ComputeShortestInterval() const |
virtual void | TObject::DoError(int level, const char* location, const char* fmt, va_list va) const |
void | TObject::MakeZombie() |
enum TObject::EStatusBits { | kCanDelete | |
kMustCleanup | ||
kObjInCanvas | ||
kIsReferenced | ||
kHasUUID | ||
kCannotPick | ||
kNoContextMenu | ||
kInvalidObject | ||
}; | ||
enum TObject::[unnamed] { | kIsOnHeap | |
kNotDeleted | ||
kZombie | ||
kBitMask | ||
kSingleKey | ||
kOverwrite | ||
kWriteDelete | ||
}; |
TString | TNamed::fName | object identifier |
TString | TNamed::fTitle | object title |
TF1* | fApproxPosterior | TF1 representing the scanned posterior function |
double | fBrfPrecision | root finder precision |
RooArgSet | fConditionalObs | conditional observables |
RooAbsData* | fData | data set |
RooAbsReal* | fIntegratedLikelihood | integrated likelihood function, i.e - unnormalized posterior function |
TString | fIntegrationType | |
double | fLeftSideFraction | fraction of probability content on left side of interval |
RooAbsReal* | fLikelihood | internal pointer to likelihood function |
RooAbsReal* | fLogLike | internal pointer to log likelihood function |
Double_t | fLower | computer lower interval bound |
Double_t | fNLLMin | minimum value of Nll |
int | fNScanBins | number of bins to scan, if = -1 no scan is done (default) |
RooArgSet | fNuisanceParameters | nuisance parameters |
RooAbsPdf* | fNuisancePdf | nuisance pdf (needed when using nuisance sampling technique) |
int | fNumIterations | number of iterations (when using ToyMC) |
RooArgSet | fPOI | POI |
RooAbsPdf* | fPdf | model pdf (could contain the nuisance pdf as constraint term) |
ROOT::Math::IBaseFunctionOneDim* | fPosteriorFunction | function representing the posterior |
RooAbsPdf* | fPosteriorPdf | normalized (on the poi) posterior pdf |
RooAbsPdf* | fPriorPdf | prior pdf (typically for the POI) |
RooAbsPdf* | fProductPdf | internal pointer to model * prior |
double | fSize | size used for getting the interval |
Double_t | fUpper | upper interval bound |
Bool_t | fValidInterval |
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
Constructor from a data set and a ModelConfig model pdf, poi and nuisances will be taken from the ModelConfig
set the model to use The model pdf, prior pdf, parameter of interest and nuisances will be taken according to the model
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
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
return a RooPlot with the posterior and the credibility region NOTE: User takes ownership of the returned object
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
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
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
internal function compute the interval using RooFit
internal function compute the interval using Cdf integration
compute the interval using the approximate posterior function
specify the nuisance parameters (eg. the rest of the parameters)
{fNuisanceParameters.removeAll(); fNuisanceParameters.add(set);}
set the conditional observables which will be used when creating the NLL so the pdf's will not be normalized on the conditional observables when computing the NLL
{fConditionalObs.removeAll(); fConditionalObs.add(set);}
set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval)
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval)
{ SetTestSize(1.-cl); }
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;}
set the Bayesian calculator to compute the shorest interval (default is central interval) to switch off SetLeftSideTailFraction to the rght value
{ fLeftSideFraction = -1; }
set the precision of the Root Finder
{ fBrfPrecision = precision; }
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; }
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; }
force the nuisance pdf when using the toy mc sampling
{ fNuisancePdf = &pdf; }