BayesianAnalysisToolkit  0.9.3
Private Attributes | List of all members
BCGoFTest Class Reference

The class for testing model hypotheses. More...

#include <BCGoFTest.h>

Inheritance diagram for BCGoFTest:
Inheritance graph
[legend]
Collaboration diagram for BCGoFTest:
Collaboration graph
[legend]

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)
 
- Public Member Functions inherited from BCIntegrate
 BCIntegrate ()
 
 BCIntegrate (const BCIntegrate &bcintegrate)
 
virtual ~BCIntegrate ()
 
BCIntegrateoperator= (const BCIntegrate &bcintegrate)
 
int ReadMarginalizedFromFile (const char *file)
 
BCH1DGetMarginalized (const BCParameter *parameter)
 
BCH1DGetMarginalized (const char *name)
 
BCH1DGetMarginalized (unsigned index)
 
BCH2DGetMarginalized (const BCParameter *parameter1, const BCParameter *parameter2)
 
BCH2DGetMarginalized (const char *name1, const char *name2)
 
BCH2DGetMarginalized (unsigned index1, unsigned index2)
 
int PrintAllMarginalized1D (const char *filebase)
 
int PrintAllMarginalized2D (const char *filebase)
 
int PrintAllMarginalized (const char *file, std::string options1d="BTsiB3CS1D0pdf0Lmeanmode", std::string options2d="BTfB3CS1meangmode", unsigned int hdiv=1, unsigned int ndiv=1)
 
double GetIntegral () const
 
BCIntegrate::BCOptimizationMethod GetOptimizationMethod () const
 
BCIntegrate::BCIntegrationMethod GetIntegrationMethod () const
 
BCIntegrate::BCMarginalizationMethod GetMarginalizationMethod () const
 
BCIntegrate::BCSASchedule GetSASchedule () const
 
void GetRandomVectorUnitHypercube (std::vector< double > &x) const
 
void GetRandomVectorInParameterSpace (std::vector< double > &x) const
 
double GetRandomPoint (std::vector< double > &x)
 
int GetNIterationsMin () const
 
int GetNIterationsMax () const
 
int GetNIterationsPrecisionCheck () const
 
int GetNIterationsOutput () const
 
int GetNIterations () const
 
double GetRelativePrecision () const
 
double GetAbsolutePrecision () const
 
BCCubaMethod GetCubaIntegrationMethod () const
 
const BCCubaOptions::VegasGetCubaVegasOptions () const
 
const BCCubaOptions::SuaveGetCubaSuaveOptions () const
 
const BCCubaOptions::DivonneGetCubaDivonneOptions () const
 
const BCCubaOptions::CuhreGetCubaCuhreOptions () const
 
BCH1DGetSlice (const BCParameter *parameter, const std::vector< double > parameters=std::vector< double >(0), int bins=0, bool flag_norm=true)
 
BCH1DGetSlice (const char *name, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
BCH2DGetSlice (const BCParameter *parameter1, const BCParameter *parameter2, const std::vector< double > parameters=std::vector< double >(0), int bins=0, bool flag_norm=true)
 
BCH2DGetSlice (const char *name1, const char *name2, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
BCH2DGetSlice (unsigned index1, unsigned index2, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
double GetError () const
 
TMinuit * GetMinuit ()
 
int GetMinuitErrorFlag () const
 
double GetSAT0 () const
 
double GetSATmin () const
 
double GetBestFitParameter (unsigned index) const
 
double GetBestFitParameterError (unsigned index) const
 
double GetLogMaximum ()
 
const std::vector< double > & GetBestFitParameters () const
 
const std::vector< double > & GetBestFitParameterErrors () const
 
void SetMinuitArlist (double *arglist)
 
void SetFlagIgnorePrevOptimization (bool flag)
 
void SetOptimizationMethod (BCIntegrate::BCOptimizationMethod method)
 
void SetIntegrationMethod (BCIntegrate::BCIntegrationMethod method)
 
void SetMarginalizationMethod (BCIntegrate::BCMarginalizationMethod method)
 
void SetSASchedule (BCIntegrate::BCSASchedule schedule)
 
void SetNIterationsMin (int niterations)
 
void SetNIterationsMax (int niterations)
 
void SetNIterationsPrecisionCheck (int niterations)
 
void SetNIterationsOutput (int niterations)
 
void SetRelativePrecision (double relprecision)
 
void SetAbsolutePrecision (double absprecision)
 
void SetCubaIntegrationMethod (BCCubaMethod type)
 
void SetCubaOptions (const BCCubaOptions::Vegas &options)
 
void SetCubaOptions (const BCCubaOptions::Suave &options)
 
void SetCubaOptions (const BCCubaOptions::Divonne &options)
 
void SetCubaOptions (const BCCubaOptions::Cuhre &options)
 
void SetSAT0 (double T0)
 
void SetSATmin (double Tmin)
 
void SetFlagWriteSAToFile (bool flag)
 
TTree * GetSATree ()
 
void InitializeSATree ()
 
double Normalize ()
 
double Integrate (BCIntegrationMethod intmethod)
 
double Integrate ()
 
double Integrate (BCIntegrationMethod type, tRandomizer randomizer, tEvaluator evaluator, tIntegralUpdater updater, std::vector< double > &sums)
 
double EvaluatorMC (std::vector< double > &sums, const std::vector< double > &point, bool &accepted)
 
int MarginalizeAll ()
 
int MarginalizeAll (BCMarginalizationMethod margmethod)
 
virtual void MarginalizePreprocess ()
 
virtual void MarginalizePostprocess ()
 
void SAInitialize ()
 
std::vector< double > FindMode (std::vector< double > start=std::vector< double >())
 
std::vector< double > FindMode (BCIntegrate::BCOptimizationMethod optmethod, std::vector< double > start=std::vector< double >())
 
double SATemperature (double t)
 
double SATemperatureBoltzmann (double t)
 
double SATemperatureCauchy (double t)
 
virtual double SATemperatureCustom (double t)
 
std::vector< double > GetProposalPointSA (const std::vector< double > &x, int t)
 
std::vector< double > GetProposalPointSABoltzmann (const std::vector< double > &x, int t)
 
std::vector< double > GetProposalPointSACauchy (const std::vector< double > &x, int t)
 
virtual std::vector< double > GetProposalPointSACustom (const std::vector< double > &x, int t)
 
std::vector< double > SAHelperGetRandomPointOnHypersphere ()
 
double SAHelperGetRadialCauchy ()
 
double SAHelperSinusToNIntegral (int dim, double theta)
 
virtual void ResetResults ()
 
std::string DumpIntegrationMethod (BCIntegrationMethod type)
 
std::string DumpCurrentIntegrationMethod ()
 
std::string DumpUsedIntegrationMethod ()
 
std::string DumpMarginalizationMethod (BCMarginalizationMethod type)
 
std::string DumpCurrentMarginalizationMethod ()
 
std::string DumpUsedMarginalizationMethod ()
 
std::string DumpOptimizationMethod (BCOptimizationMethod type)
 
std::string DumpCurrentOptimizationMethod ()
 
std::string DumpUsedOptimizationMethod ()
 
std::string DumpCubaIntegrationMethod (BCCubaMethod type)
 
std::string DumpCubaIntegrationMethod ()
 
void SetBestFitParameters (const std::vector< double > &x)
 
void SetBestFitParameters (const std::vector< double > &x, const double &new_value, double &old_value)
 
unsigned GetNIntegrationVariables ()
 
double CalculateIntegrationVolume ()
 
bool CheckMarginalizationAvailability (BCMarginalizationMethod type)
 
bool CheckMarginalizationIndices (TH1 *hist, const std::vector< unsigned > &index)
 
- Public Member Functions inherited from BCEngineMCMC
void WriteMarkovChain (bool flag)
 
 BCEngineMCMC ()
 
 BCEngineMCMC (const BCEngineMCMC &enginemcmc)
 
virtual ~BCEngineMCMC ()
 
BCEngineMCMCoperator= (const BCEngineMCMC &engineMCMC)
 
unsigned MCMCGetNChains () const
 
unsigned MCMCGetNLag () const
 
const std::vector< unsigned > & MCMCGetNIterations () const
 
unsigned MCMCGetCurrentIteration () const
 
unsigned MCMCGetCurrentChain () const
 
unsigned MCMCGetNIterationsConvergenceGlobal () const
 
bool MCMCGetFlagConvergenceGlobal () const
 
unsigned MCMCGetNIterationsMax () const
 
unsigned MCMCGetNIterationsRun () const
 
unsigned MCMCGetNIterationsPreRunMin () const
 
unsigned MCMCGetNIterationsUpdate () const
 
unsigned MCMCGetNIterationsUpdateMax () const
 
std::vector< int > MCMCGetNTrialsTrue () const
 
int MCMCGetNTrials () const
 
const std::vector< double > & MCMCGetprobMean () const
 
const std::vector< double > & MCMCGetVariance () const
 
const std::vector< double > & MCMCGetTrialFunctionScaleFactor () const
 
std::vector< double > MCMCGetTrialFunctionScaleFactor (unsigned ichain) const
 
double MCMCGetTrialFunctionScaleFactor (unsigned ichain, unsigned ipar)
 
const std::vector< double > & MCMCGetx () const
 
std::vector< double > MCMCGetx (unsigned ichain)
 
double MCMCGetx (unsigned ichain, unsigned ipar) const
 
const std::vector< double > & MCMCGetLogProbx () const
 
double MCMCGetLogProbx (unsigned ichain)
 
int MCMCGetPhase () const
 
const std::vector< double > & MCMCGetMaximumPoints () const
 
std::vector< double > MCMCGetMaximumPoint (unsigned i) const
 
const std::vector< double > & MCMCGetMaximumLogProb () const
 
int MCMCGetFlagInitialPosition () const
 
double MCMCGetRValueCriterion () const
 
double MCMCGetRValueParametersCriterion () const
 
double MCMCGetRValue () const
 
double MCMCGetRValueParameters (unsigned i)
 
bool MCMCGetRValueStrict () const
 
bool MCMCGetFlagRun () const
 
TTree * MCMCGetMarkovChainTree (unsigned i)
 
BCH1DMCMCGetH1Marginalized (unsigned i)
 
BCH2DMCMCGetH2Marginalized (unsigned i, unsigned j)
 
BCParameterGetParameter (int index) const
 
BCParameterGetParameter (const char *name) const
 
unsigned int GetNParameters () const
 
unsigned int GetNFixedParameters ()
 
unsigned int GetNFreeParameters ()
 
const std::vector< double > & GetBestFitParametersMarginalized () const
 
void MCMCSetTrialFunctionScaleFactor (std::vector< double > scale)
 
void MCMCSetNChains (unsigned n)
 
void MCMCSetNLag (unsigned n)
 
void MCMCSetNIterationsMax (unsigned n)
 
void MCMCSetNIterationsRun (unsigned n)
 
void MCMCSetNIterationsPreRunMin (unsigned n)
 
void MCMCSetNIterationsUpdate (unsigned n)
 
void MCMCSetNIterationsUpdateMax (unsigned n)
 
void MCMCSetMinimumEfficiency (double efficiency)
 
void MCMCSetMaximumEfficiency (double efficiency)
 
void MCMCSetRandomSeed (unsigned seed)
 
void MCMCSetWriteChainToFile (bool flag)
 
void MCMCSetWritePreRunToFile (bool flag)
 
void MCMCSetInitialPositions (const std::vector< double > &x0s)
 
void MCMCSetInitialPositions (std::vector< std::vector< double > > x0s)
 
void MCMCSetFlagInitialPosition (int flag)
 
void MCMCSetFlagOrderParameters (bool flag)
 
void MCMCSetFlagFillHistograms (bool flag)
 
void MCMCSetFlagPreRun (bool flag)
 
void MCMCSetRValueCriterion (double r)
 
void MCMCSetRValueParametersCriterion (double r)
 
void MCMCSetRValueStrict (bool strict=true)
 
void MCMCSetMarkovChainTrees (const std::vector< TTree * > &trees)
 
void MCMCInitializeMarkovChainTrees ()
 
int SetMarginalized (unsigned index, TH1D *h)
 
int SetMarginalized (unsigned index1, unsigned index2, TH2D *h)
 
void MCMCSetValuesDefault ()
 
void MCMCSetValuesQuick ()
 
void MCMCSetValuesDetail ()
 
void MCMCSetPrecision (BCEngineMCMC::Precision precision)
 
void SetNbins (unsigned int nbins)
 
void Copy (const BCEngineMCMC &enginemcmc)
 
virtual int AddParameter (const char *name, double min, double max, const char *latexname="")
 
virtual void MCMCTrialFunction (unsigned ichain, std::vector< double > &x)
 
virtual double MCMCTrialFunctionSingle (unsigned ichain, unsigned ipar)
 
bool MCMCGetProposalPointMetropolis (unsigned chain, std::vector< double > &x)
 
bool MCMCGetProposalPointMetropolis (unsigned chain, unsigned parameter, std::vector< double > &x)
 
bool MCMCGetNewPointMetropolis (unsigned chain=0)
 
bool MCMCGetNewPointMetropolis (unsigned chain, unsigned parameter)
 
void MCMCInChainCheckMaximum ()
 
void MCMCInChainUpdateStatistics ()
 
void MCMCInChainFillHistograms ()
 
void MCMCInChainTestConvergenceAllChains ()
 
void MCMCInChainWriteChains ()
 
int MCMCMetropolis ()
 
int MCMCMetropolisPreRun ()
 
void MCMCResetRunStatistics ()
 
void MCMCInitializeMarkovChains ()
 
int MCMCInitialize ()
 
void ClearParameters (bool hard=false)
 
virtual void MCMCIterationInterface ()
 
virtual void MCMCCurrentPointInterface (std::vector< double > &, int, bool)
 

Private Attributes

std::vector< int > fMapDataPoint
 
std::vector< int > fMapDataValue
 
int fPValueBelow
 
int fPValueAbove
 
BCModelfTestModel
 
BCDataSetfTemporaryDataSet
 
double fLogLikelihood
 
double fLogLikelihoodMin
 
double fLogLikelihoodMax
 
TH1D * fHistogramLogProb
 

Additional Inherited Members

- Public Types inherited from BCIntegrate
enum  BCOptimizationMethod {
  kOptEmpty, kOptSimAnn, kOptMetropolis, kOptMinuit,
  kOptDefault, NOptMethods
}
 
enum  BCIntegrationMethod {
  kIntEmpty, kIntMonteCarlo, kIntCuba, kIntGrid,
  kIntDefault, NIntMethods
}
 
enum  BCMarginalizationMethod {
  kMargEmpty, kMargMetropolis, kMargMonteCarlo, kMargGrid,
  kMargDefault, NMargMethods
}
 
enum  BCSASchedule { kSACauchy, kSABoltzmann, kSACustom, NSAMethods }
 
enum  BCCubaMethod {
  kCubaVegas, kCubaSuave, kCubaDivonne, kCubaCuhre,
  NCubaMethods
}
 
typedef void(BCIntegrate::* tRandomizer )(std::vector< double > &) const
 
typedef double(BCIntegrate::* tEvaluator )(std::vector< double > &, const std::vector< double > &, bool &)
 
typedef void(* tIntegralUpdater )(const std::vector< double > &, const int &, double &, double &)
 
- Static Public Member Functions inherited from BCIntegrate
static void IntegralUpdaterMC (const std::vector< double > &sums, const int &nIterations, double &integral, double &absprecision)
 
static int CubaIntegrand (const int *ndim, const double xx[], const int *ncomp, double ff[], void *userdata)
 
static void FCNLikelihood (int &npar, double *grad, double &fval, double *par, int flag)
 
- Protected Member Functions inherited from BCIntegrate
unsigned IntegrationOutputFrequency () const
 
void LogOutputAtStartOfIntegration (BCIntegrationMethod type, BCCubaMethod cubatype)
 
void LogOutputAtIntegrationStatusUpdate (BCIntegrationMethod type, double integral, double absprecision, int nIterations)
 
void LogOutputAtEndOfIntegration (double integral, double absprecision, double relprecision, int nIterations)
 
void Copy (const BCIntegrate &bcintegrate)
 
- 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.

23  : BCModel(name)
24 {
25  // set original data set to zero
27 
28  // set test mode to zero
29  fTestModel = 0;
30 
31  // reset pvalue and counter
32  fPValue = 0;
33  fPValueAbove = 0;
34  fPValueBelow = 0;
35 
36  // reset loglikelihood and range
37  fLogLikelihood = 0;
38  fLogLikelihoodMin = 1e99;
39  fLogLikelihoodMax = -1e99;
40 
41  // define new histogram
43 
44  // set defaults for the MCMC
45  MCMCSetNChains(5);
46  MCMCSetNIterationsMax(100000);
48 }
BCGoFTest::~BCGoFTest ( )

Default destructor.

Definition at line 52 of file BCGoFTest.cxx.

53 {
54  // restore original data set
55 
56  // get number of data points and values
57  int ndatapoints = fTemporaryDataSet->GetNDataPoints();
58  int ndatavalues = fTemporaryDataSet->GetDataPoint(0)->GetNValues();
59 
60  for (int i = 0; i < ndatapoints; ++i)
61  for (int j = 0; j < ndatavalues; ++j)
63 
64  // restore data point limits
65  for (unsigned int i = 0; i < GetNParameters(); ++i)
67  fMapDataValue[i],
68  GetParameter(i)->GetLowerLimit(),
69  GetParameter(i)->GetUpperLimit());
70 
71  // delete temporary data set
72  delete fTemporaryDataSet;
73 }

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.

228 {
229  // set histogram point to null
230  fHistogramLogProb = 0;
231 
232  if (flag_histogram)
233  {
234  // perform first run to obtain limits for the log(likelihood)
235  MarginalizeAll();
236 
237  // create histogram
239  fHistogramLogProb = new TH1D(Form("hist_%s_logprob", GetName().data()), ";ln(prob);N", 100, fLogLikelihoodMin - 0.1*D, fLogLikelihoodMax + 0.1*D);
240  fHistogramLogProb->SetStats(kFALSE);
241  }
242  else
243  {
244  }
245 
246  // run MCMC
247  MarginalizeAll();
248 
249  // check for convergence
251  {
252  BCLog::OutDetail(" --> MCMC did not converge in evaluation of the p-value.");
253  return -1;
254  }
255 
256  // calculate p-value
257  fPValue = double(fPValueBelow) / double(fPValueBelow + fPValueAbove);
258 
259  // return p-value
260  return fPValue;
261 }
TH1D* BCGoFTest::GetHistogramLogProb ( )
inline
Returns
distribution of log(likelihood)

Definition at line 62 of file BCGoFTest.h.

63  { return fHistogramLogProb; };
BCModel* BCGoFTest::GetTestModel ( )
inline
Returns
pointer to the tested model

Definition at line 67 of file BCGoFTest.h.

68  { return fTestModel; };
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.

93  { return 0; };
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.

78 {
79  // set the original data set to the new parameters
80  for (int i = 0; i < int(parameters.size()); ++i)
82 
83  // calculate likelihood at the point of the original parameters
84  double loglikelihood = fTestModel->LogLikelihood(fDataSet->GetDataPoint(0)->GetValues());
85 
86  // return likelihood
87  return loglikelihood;
88 }
void BCGoFTest::MCMCUserIterationInterface ( )
virtual

Method executed for every iteration of the MCMC. User's code should be provided via overloading in the derived class

Reimplemented from BCEngineMCMC.

Definition at line 92 of file BCGoFTest.cxx.

93 {
94  int nchains = MCMCGetNChains();
95 
96  for (int i = 0; i < nchains; ++i)
97  {
98  // get likelihood at the point of the original parameters
99  double loglikelihood = MCMCGetLogProbx(i);
100 
101  // calculate pvalue
102  if (loglikelihood < fLogLikelihood)
103  fPValueBelow++;
104  else
105  fPValueAbove++;
106 
107  // if histogram exists already, then fill it ...
108  if (fHistogramLogProb)
109  fHistogramLogProb->Fill(loglikelihood);
110  // ...otherwise find range
111  else
112  {
113  if (loglikelihood > fLogLikelihoodMax)
114  fLogLikelihoodMax = loglikelihood;
115  else if (loglikelihood < fLogLikelihoodMin)
116  fLogLikelihoodMin = loglikelihood;
117  }
118  }
119 }
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.

78  { fTestModel = testmodel; };
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.

124 {
125  // check if the boundaries of the original data set exist.
127  {
128  BCLog::OutError("BCGoFTest::SetTestDataPoint : Boundaries of the original data set are not defined.");
129  return 0;
130  }
131 
132  // reset histogram
133  if (fHistogramLogProb)
134  {
135  delete fHistogramLogProb;
136  fHistogramLogProb = 0;
137  }
138 
139  // reset variables
140  fPValue = 0;
141  fPValueAbove = 0;
142  fPValueBelow = 0;
143 
144  // create temporary data set ...
146 
147  // ... and fill with the original one
148 
149  // get number of data points and values
150  int ndatapoints = fTestModel->GetDataSet()->GetNDataPoints();
151  int ndatavalues = fTestModel->GetDataSet()->GetDataPoint(0)->GetNValues();
152 
153  for (int i = 0; i < ndatapoints; ++i)
154  {
157  }
158 
159  // clear maps
160  fMapDataPoint.clear();
161  fMapDataValue.clear();
162 
163  int counter = 0;
164 
165  // remove parameters, but doesn't clear up memory
167 
168  // loop through data points and values
169  for (int i = 0; i < ndatapoints; ++i)
170  for (int j = 0; j < ndatavalues; ++j)
171  {
172  // debugKK
173  // needs to be fixed
174  // if (fTestModel->GetFixedDataAxis(j))
175  // continue;
176 
177  // add parameter to this model
178  std::string parName = Form("parameter_%i", counter);
179  AddParameter(
180  parName.c_str(),
183 
184  // add another element to the maps
185  fMapDataPoint.push_back(i);
186  fMapDataValue.push_back(j);
187 
188  // increase counter
189  counter ++;
190  }
191 
192  // check if there are any non-fixed data values left
193  if (counter == 0)
194  {
195  BCLog::OutError("BCGoFTest::SetTestDataPoint : No non-fixed data values left.");
196  return 0;
197  }
198 
199  // create a new data set containing the vector of parameters which
200  // are to be tested
201  BCDataPoint * datapoint = new BCDataPoint(parameters);
202  BCDataSet * dataset = new BCDataSet();
203  dataset->AddDataPoint(datapoint);
204 
205  // calculate likelihood of the original data set
206  fLogLikelihood = fTestModel->LogLikelihood(parameters);
207 
208  // if data set has been set before, delete
209  if (fDataSet)
210  delete fDataSet;
211 
212  // set data set of this model
213  fDataSet = dataset;
214 
215  // put proper range to new data set
216  for (int i = 0; i < int(parameters.size()); ++i)
218  i,
221 
222  return 1;
223 }

Member Data Documentation

TH1D* BCGoFTest::fHistogramLogProb
private

The distribution of log(likelihood).

Definition at line 127 of file BCGoFTest.h.

double BCGoFTest::fLogLikelihood
private

The log(likelihood) and its range.

Definition at line 121 of file BCGoFTest.h.

double BCGoFTest::fLogLikelihoodMax
private

Definition at line 123 of file BCGoFTest.h.

double BCGoFTest::fLogLikelihoodMin
private

Definition at line 122 of file BCGoFTest.h.

std::vector<int> BCGoFTest::fMapDataPoint
private

A map of data points and data values.

Definition at line 103 of file BCGoFTest.h.

std::vector<int> BCGoFTest::fMapDataValue
private

Definition at line 104 of file BCGoFTest.h.

int BCGoFTest::fPValueAbove
private

Definition at line 109 of file BCGoFTest.h.

int BCGoFTest::fPValueBelow
private

Counter for the evaluation of the p-value.

Definition at line 108 of file BCGoFTest.h.

BCDataSet* BCGoFTest::fTemporaryDataSet
private

A data set used for temporary storage.

Definition at line 117 of file BCGoFTest.h.

BCModel* BCGoFTest::fTestModel
private

A pointer to the model which is tested.

Definition at line 113 of file BCGoFTest.h.


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