BAT  0.9.4
The Bayesian analysis toolkit
 All Classes Namespaces Functions Variables Enumerations
BCGoFTest Class Reference

The class for testing model hypotheses. More...

#include <BCGoFTest.h>

Inheritance diagram for BCGoFTest:
Collaboration diagram for BCGoFTest:

Public Member Functions

Constructors and destructors
 BCGoFTest (const char *name)
 
 ~BCGoFTest ()
 
Member functions (get)
double GetCalculatedPValue (bool flag_histogram=false)
 
TH1D * GetHistogramLogProb ()
 
BCModelGetTestModel ()
 
Member functions (set)
void SetTestModel (BCModel *testmodel)
 
int SetTestPoint (std::vector< double > parameters)
 
Member functions (miscellaneous methods)
double LogLikelihood (const std::vector< double > &parameters)
 
double LogAPrioriProbability (const std::vector< double > &)
 
void MCMCUserIterationInterface ()
 
- Public Member Functions inherited from BCModel
 BCModel (const char *name="model")
 
 BCModel (const BCModel &bcmodel)
 
virtual ~BCModel ()
 
BCModeloperator= (const BCModel &bcmodel)
 
const std::string & GetName () const
 
double GetModelAPrioriProbability () const
 
double GetModelAPosterioriProbability () const
 
BCDataSetGetDataSet () const
 
BCDataPointGetDataPointLowerBoundaries () const
 
BCDataPointGetDataPointUpperBoundaries () const
 
double GetDataPointLowerBoundary (unsigned int index) const
 
double GetDataPointUpperBoundary (unsigned int index) const
 
bool GetFlagBoundaries () const
 
unsigned GetNDataPoints () const
 
BCDataPointGetDataPoint (unsigned int index) const
 
void SetName (const char *name)
 
void SetModelAPrioriProbability (double probability)
 
void SetModelAPosterioriProbability (double probability)
 
virtual int AddParameter (BCParameter *parameter)
 
void SetDataSet (BCDataSet *dataset)
 
void SetSingleDataPoint (BCDataPoint *datapoint)
 
void SetSingleDataPoint (BCDataSet *dataset, unsigned int index)
 
void SetDataBoundaries (unsigned int index, double lowerboundary, double upperboundary, bool fixed=false)
 
void SetDataPointLowerBoundaries (BCDataPoint *datasetlowerboundaries)
 
void SetDataPointUpperBoundaries (BCDataPoint *datasetupperboundaries)
 
void SetDataPointLowerBoundary (int index, double lowerboundary)
 
void SetDataPointUpperBoundary (int index, double upperboundary)
 
int SetPrior (int index, TF1 *f)
 
int SetPrior (const char *name, TF1 *f)
 
int SetPriorDelta (int index, double value)
 
int SetPriorDelta (const char *name, double value)
 
int SetPriorGauss (int index, double mean, double sigma)
 
int SetPriorGauss (const char *name, double mean, double sigma)
 
int SetPriorGauss (int index, double mean, double sigmadown, double sigmaup)
 
int SetPriorGauss (const char *name, double mean, double sigmadown, double sigmaup)
 
int SetPrior (int index, TH1 *h, bool flag=false)
 
int SetPrior (const char *name, TH1 *h, bool flag=false)
 
int SetPriorConstant (int index)
 
int SetPriorConstant (const char *name)
 
int SetPriorConstantAll ()
 
void Copy (const BCModel &bcmodel)
 
double APrioriProbability (const std::vector< double > &parameters)
 
virtual double Likelihood (const std::vector< double > &params)
 
double ProbabilityNN (const std::vector< double > &params)
 
double LogProbabilityNN (const std::vector< double > &parameters)
 
double Probability (const std::vector< double > &parameter)
 
double LogProbability (const std::vector< double > &parameter)
 
virtual double SamplingFunction (const std::vector< double > &parameters)
 
double Eval (const std::vector< double > &parameters)
 
virtual double LogEval (const std::vector< double > &parameters)
 
virtual void CorrelateDataPointValues (std::vector< double > &x)
 
double GetPvalueFromChi2 (const std::vector< double > &par, int sigma_index)
 
double GetPvalueFromKolmogorov (const std::vector< double > &par, int index)
 
double GetPvalueFromChi2NDoF (std::vector< double > par, int sigma_index)
 
BCH1DCalculatePValue (std::vector< double > par, bool flag_histogram=false)
 
double GetPValue ()
 
double GetPValueNDoF ()
 
double GetChi2NDoF ()
 
std::vector< double > GetChi2Runs (int dataIndex, int sigmaIndex)
 
void SetGoFNIterationsMax (int n)
 
void SetGoFNIterationsRun (int n)
 
void SetGoFNChains (int n)
 
double HessianMatrixElement (const BCParameter *parameter1, const BCParameter *parameter2, std::vector< double > point)
 
void PrintSummary ()
 
void PrintResults (const char *file)
 
void PrintShortFitSummary (int chi2flag=0)
 
void PrintHessianMatrix (std::vector< double > parameters)
 
virtual double CDF (const std::vector< double > &, int, bool)
 

Additional Inherited Members

- Protected Attributes inherited from BCModel
std::string fName
 
double fModelAPriori
 
double fModelAPosteriori
 
BCDataSetfDataSet
 
BCDataPointfDataPointLowerBoundaries
 
BCDataPointfDataPointUpperBoundaries
 
std::vector< bool > fDataFixedValues
 
double fPValue
 
double fChi2NDoF
 
double fPValueNDoF
 
bool flag_discrete
 
int fGoFNIterationsMax
 
int fGoFNIterationsRun
 
int fGoFNChains
 
std::vector< TNamed * > fPriorContainer
 
bool fPriorConstantAll
 
std::vector< bool > fPriorContainerConstant
 
std::vector< bool > fPriorContainerInterpolate
 

Detailed Description

The class for testing model hypotheses.

Author
Daniel Kollar
Kevin Kröninger
Version
1.0
Date
08.2008 This class is used for calculating the p-value of a model.

Definition at line 34 of file BCGoFTest.h.

Constructor & Destructor Documentation

BCGoFTest::BCGoFTest ( const char *  name)

Default constructor.

Definition at line 23 of file BCGoFTest.cxx.

BCGoFTest::~BCGoFTest ( )

Default destructor.

Definition at line 52 of file BCGoFTest.cxx.

Member Function Documentation

double BCGoFTest::GetCalculatedPValue ( bool  flag_histogram = false)

Calculated the p-value.

Parameters
flag_histogramA histogram is either filled or not.
Returns
p-value

Definition at line 227 of file BCGoFTest.cxx.

TH1D* BCGoFTest::GetHistogramLogProb ( )
inline
Returns
distribution of log(likelihood)

Definition at line 62 of file BCGoFTest.h.

BCModel* BCGoFTest::GetTestModel ( )
inline
Returns
pointer to the tested model

Definition at line 67 of file BCGoFTest.h.

double BCGoFTest::LogAPrioriProbability ( const std::vector< double > &  parameters)
inlinevirtual

Returns natural logarithm of the prior probability. Method needs to be overloaded by the user.

Parameters
parametersA set of parameter values
Returns
The prior probability p(parameters)
See Also
GetPrior(std::vector<double> parameters)

Reimplemented from BCModel.

Definition at line 92 of file BCGoFTest.h.

double BCGoFTest::LogLikelihood ( const std::vector< double > &  params)
virtual

Calculates natural logarithm of the likelihood. Method needs to be overloaded by the user.

Parameters
paramsA set of parameter values
Returns
Natural logarithm of the likelihood

Implements BCModel.

Definition at line 77 of file BCGoFTest.cxx.

void BCGoFTest::SetTestModel ( BCModel testmodel)
inline

Set the model to be tested.

Parameters
testmodelpointer to the model to be tested

Definition at line 77 of file BCGoFTest.h.

int BCGoFTest::SetTestPoint ( std::vector< double >  parameters)

Sets the set of parameters which the p-values is calculated for.

Parameters
parametersparameters
Returns
error code

Definition at line 123 of file BCGoFTest.cxx.


The documentation for this class was generated from the following files: