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

A class for fitting histograms with functions. More...

#include <BCEfficiencyFitter.h>

Inheritance diagram for BCEfficiencyFitter:
Collaboration diagram for BCEfficiencyFitter:

Classes

class  ToyDataInterface
 

Public Member Functions

Constructors and destructors
 BCEfficiencyFitter ()
 
 BCEfficiencyFitter (const char *name)
 
 BCEfficiencyFitter (TH1D *hist1, TH1D *hist2, TF1 *func)
 
 BCEfficiencyFitter (const char *name, TH1D *hist1, TH1D *hist2, TF1 *func)
 
 ~BCEfficiencyFitter ()
 
Member functions (get)
TH1D * GetHistogram1 ()
 
TH1D * GetHistogram2 ()
 
TGraphAsymmErrors * GetHistogramRatio ()
 
TF1 * GetFitFunction ()
 
TGraph * GetErrorBand ()
 
TGraph * GetGraphFitFunction ()
 
int GetUncertainties (int n, int k, double p, double &xexp, double &xmin, double &xmax)
 
Member functions (set)
int SetHistograms (TH1D *hist1, TH1D *hist2)
 
int SetFitFunction (TF1 *func)
 
void SetFlagIntegration (bool flag)
 
void SetDataPointType (int type)
 
Member functions (miscellaneous methods)
virtual double LogLikelihood (const std::vector< double > &parameters)
 
double FitFunction (const std::vector< double > &x, const std::vector< double > &parameters)
 
int Fit ()
 
int Fit (TH1D *hist1, TH1D *hist2, TF1 *func)
 
void DrawFit (const char *options="", bool flaglegend=false)
 
int CalculatePValueFast (const std::vector< double > &par, BCEfficiencyFitter::ToyDataInterface *callback, double &pvalue, double &pvalueCorrected, unsigned nIterations=100000)
 
int CalculatePValueFast (const std::vector< double > &par, double &pvalue, double &pvalueCorrected, unsigned nIterations=100000)
 
- Public Member Functions inherited from BCFitter
 BCFitter ()
 
 BCFitter (const char *name)
 
 ~BCFitter ()
 
TGraph * GetErrorBand ()
 
TH2D * GetErrorBandXY () const
 
TH2D * GetErrorBandXY_yellow (double level=.68, int nsmooth=0) const
 
std::vector< double > GetErrorBand (double level) const
 
TGraph * GetErrorBandGraph (double level1, double level2) const
 
TGraph * GetFitFunctionGraph (const std::vector< double > &parameters)
 
TGraph * GetFitFunctionGraph ()
 
TGraph * GetFitFunctionGraph (const std::vector< double > &parameters, double xmin, double xmax, int n=1000)
 
void FixDataAxis (unsigned int index, bool fixed)
 
bool GetFixedDataAxis (unsigned int index) const
 
void SetFillErrorBand (bool flag=true)
 
void SetErrorBandHisto (TH2D *h)
 
void UnsetFillErrorBand ()
 
void SetFitFunctionIndexX (int index)
 
void SetFitFunctionIndexY (int index)
 
void SetFitFunctionIndices (int indexx, int indexy)
 
void SetErrorBandContinuous (bool flag)
 
int ReadErrorBandFromFile (const char *file)
 
void MCMCIterationInterface ()
 
void MarginalizePreprocess ()
 
void MarginalizePostprocess ()
 
void FillErrorBand ()
 
- 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 LogAPrioriProbability (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 BCFitter
TH2D * fErrorBandXY
 
- 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

A class for fitting histograms with functions.

Author
Daniel Kollar
Kevin Kröninger
Version
1.0
Date
11.2008 This class allows fitting of efficiencies defined as a ratio of two TH1D histograms using a TF1 function. It uses binomial probabilities calculated based on the number of entries in histograms. This is only applicable if the numerator is a subset of the denominator.

Definition at line 39 of file BCEfficiencyFitter.h.

Constructor & Destructor Documentation

BCEfficiencyFitter::BCEfficiencyFitter ( )

The default constructor.

Definition at line 32 of file BCEfficiencyFitter.cxx.

BCEfficiencyFitter::BCEfficiencyFitter ( const char *  name)

Constructor

Parameters
namename fo the model

Definition at line 49 of file BCEfficiencyFitter.cxx.

BCEfficiencyFitter::BCEfficiencyFitter ( TH1D *  hist1,
TH1D *  hist2,
TF1 *  func 
)

A constructor.

Parameters
hist1The histogram with the larger numbers
hist2The histogram with the smaller numbers
funcThe fit function.

Definition at line 66 of file BCEfficiencyFitter.cxx.

BCEfficiencyFitter::BCEfficiencyFitter ( const char *  name,
TH1D *  hist1,
TH1D *  hist2,
TF1 *  func 
)

Constructor.

Parameters
namename fo the model
hist1The histogram with the larger numbers
hist2The histogram with the smaller numbers
funcThe fit function.

Definition at line 85 of file BCEfficiencyFitter.cxx.

BCEfficiencyFitter::~BCEfficiencyFitter ( )

The default destructor.

Definition at line 211 of file BCEfficiencyFitter.cxx.

Member Function Documentation

int BCEfficiencyFitter::CalculatePValueFast ( const std::vector< double > &  par,
BCEfficiencyFitter::ToyDataInterface callback,
double &  pvalue,
double &  pvalueCorrected,
unsigned  nIterations = 100000 
)

Calculate the p-value using fast-MCMC. In every iteration, a new toy data set is created. By providing a suitable implementation of ToyDataInterface, the user can calculate the distribution of an arbitrary statistic. Each toy data set as well as the expected values for the parameter values are passed on to the interface.

Parameters
parA set of parameter values
pvalueThe p-value for the default likelihood statistic
pvalueCorrectedThe p-value for the default likelihood statistic with corrections for DoF to reduce fitting bias, in all nonpathological cases the correction lowers the p-value.
callbackrequires class with operator(...) defined.
nIterationsnumber of toy data sets generated by the Markov chain
Returns
An error code

Definition at line 451 of file BCEfficiencyFitter.cxx.

int BCEfficiencyFitter::CalculatePValueFast ( const std::vector< double > &  par,
double &  pvalue,
double &  pvalueCorrected,
unsigned  nIterations = 100000 
)

Calculate the p-value using fast-MCMC.

Parameters
parA set of parameter values
pvalueThe pvalue
pvalueCorrectedThe p-value for the default likelihood statistic with corrections for DoF to reduce fitting bias
nIterationsnumber of pseudo experiments generated by the Markov chain
Returns
An error code

Definition at line 443 of file BCEfficiencyFitter.cxx.

void BCEfficiencyFitter::DrawFit ( const char *  options = "",
bool  flaglegend = false 
)
virtual

Draw the fit in the current pad.

Implements BCFitter.

Definition at line 349 of file BCEfficiencyFitter.cxx.

int BCEfficiencyFitter::Fit ( )
virtual

Performs the fit.

Returns
An error code.

Implements BCFitter.

Definition at line 309 of file BCEfficiencyFitter.cxx.

int BCEfficiencyFitter::Fit ( TH1D *  hist1,
TH1D *  hist2,
TF1 *  func 
)

Performs the fit.

Parameters
hist1The histogram with the larger number.
hist2The histogram with the smaller number.
funcThe fit function.
Returns
An error code.

Definition at line 286 of file BCEfficiencyFitter.cxx.

double BCEfficiencyFitter::FitFunction ( const std::vector< double > &  x,
const std::vector< double > &  parameters 
)
virtual

Returns the y-value of the 1-dimensional fit function at an x and for a set of parameters.

Parameters
xA vector with the x-value.
parametersA set of parameters.

Reimplemented from BCFitter.

Definition at line 276 of file BCEfficiencyFitter.cxx.

TGraph* BCEfficiencyFitter::GetErrorBand ( )
inline
Returns
pointer to the error band

Definition at line 120 of file BCEfficiencyFitter.h.

TF1* BCEfficiencyFitter::GetFitFunction ( )
inline
Returns
The fit function

Definition at line 115 of file BCEfficiencyFitter.h.

TGraph* BCEfficiencyFitter::GetGraphFitFunction ( )
inline
Returns
pointer to a graph for the fit function

Definition at line 125 of file BCEfficiencyFitter.h.

TH1D* BCEfficiencyFitter::GetHistogram1 ( )
inline
Returns
The histogram 1

Definition at line 100 of file BCEfficiencyFitter.h.

TH1D* BCEfficiencyFitter::GetHistogram2 ( )
inline
Returns
The histogram 2

Definition at line 105 of file BCEfficiencyFitter.h.

TGraphAsymmErrors* BCEfficiencyFitter::GetHistogramRatio ( )
inline
Returns
The histogram ratio

Definition at line 110 of file BCEfficiencyFitter.h.

int BCEfficiencyFitter::GetUncertainties ( int  n,
int  k,
double  p,
double &  xexp,
double &  xmin,
double &  xmax 
)

Calculates the central value and the lower and upper limits for a given probability.

Parameters
nn for the binomial.
kk for the binomial.
pThe central probability defining the limits.
xminThe central value.
xminThe lower limit.
xmaxThe upper limit.
Returns
A flag (=1 plot point, !=1 do not plot point).

Definition at line 576 of file BCEfficiencyFitter.cxx.

double BCEfficiencyFitter::LogLikelihood ( const std::vector< double > &  parameters)
virtual

The log of the prior probability. Overloaded from BCModel.

Parameters
parametersA vector of doubles containing the parameter values. The log of the conditional probability. Overloaded from BCModel.
parametersA vector of doubles containing the parameter values.

Implements BCModel.

Definition at line 230 of file BCEfficiencyFitter.cxx.

void BCEfficiencyFitter::SetDataPointType ( int  type)

Set type of point to be used to plot the efficiency data 0 - mean + rms 1 - mode + smallest 68% 2 - median + central 68%

Definition at line 220 of file BCEfficiencyFitter.cxx.

int BCEfficiencyFitter::SetFitFunction ( TF1 *  func)
Parameters
funcThe fit function
Returns
An error code (1:pass, 0:fail).

Definition at line 171 of file BCEfficiencyFitter.cxx.

void BCEfficiencyFitter::SetFlagIntegration ( bool  flag)
inline

Sets the flag for integration.
true: use ROOT's TH1D::Integrate()
false: use linear interpolation

Definition at line 158 of file BCEfficiencyFitter.h.

int BCEfficiencyFitter::SetHistograms ( TH1D *  hist1,
TH1D *  hist2 
)
Parameters
histThe histogram 1
histThe histogram 2
Returns
An error code (1:pass, 0:fail).

Definition at line 105 of file BCEfficiencyFitter.cxx.


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