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

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

#include <BCHistogramFitter.h>

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

Public Member Functions

Constructors and destructors
 BCHistogramFitter ()
 
 BCHistogramFitter (const char *name)
 
 BCHistogramFitter (TH1D *hist, TF1 *func)
 
 BCHistogramFitter (const char *name, TH1D *hist, TF1 *func)
 
 ~BCHistogramFitter ()
 
Member functions (get)
TH1D * GetHistogram ()
 
TH1D * GetHistogramExpected ()
 
TF1 * GetFitFunction ()
 
TGraph * GetErrorBand ()
 
TGraph * GetGraphFitFunction ()
 
Member functions (set)
int SetHistogram (TH1D *hist)
 
int SetHistogramExpected (const std::vector< double > &parameters)
 
int SetFitFunction (TF1 *func)
 
void SetFlagIntegration (bool flag)
 
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 *hist, TF1 *func)
 
void DrawFit (const char *options="HIST", bool flaglegend=false)
 
int CalculatePValueFast (const std::vector< double > &par, unsigned nIterations=100000)
 
int CalculatePValueLikelihood (const std::vector< double > &par)
 
int CalculatePValueLeastSquares (const std::vector< double > &par, bool weightExpect=true)
 
int CalculatePValueKolmogorov (const std::vector< double > &par)
 
double CDF (const std::vector< double > &parameters, int index, bool lower=false)
 
- 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)
 
- 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)
 
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 MCMCUserIterationInterface ()
 
virtual void MCMCCurrentPointInterface (std::vector< double > &, int, bool)
 

Private Attributes

TH1D * fHistogram
 
TF1 * fFitFunction
 
bool fFlagIntegration
 
TGraph * fErrorBand
 
TGraph * fGraphFitFunction
 
TH1D * fHistogramExpected
 

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 BCFitter
TH2D * fErrorBandXY
 

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 a TH1D histogram using a TF1 function.

Definition at line 35 of file BCHistogramFitter.h.

Constructor & Destructor Documentation

BCHistogramFitter::BCHistogramFitter ( )

The default constructor.

Definition at line 30 of file BCHistogramFitter.cxx.

31  : BCFitter()
32  , fHistogram(0)
33  , fFitFunction(0)
35 {
36  // set default options and values
39  fFlagIntegration = true;
40  flag_discrete = true;
41 
42  // set MCMC for marginalization
44 }
BCHistogramFitter::BCHistogramFitter ( const char *  name)

Constructor

Parameters
namename of the model

Definition at line 48 of file BCHistogramFitter.cxx.

49  : BCFitter(name)
50  , fHistogram(0)
51  , fFitFunction(0)
53 {
54  // set default options and values
57  fFlagIntegration = true;
58  flag_discrete = true;
59 
60  // set MCMC for marginalization
62 }
BCHistogramFitter::BCHistogramFitter ( TH1D *  hist,
TF1 *  func 
)

Constructor.

Parameters
histhistogram to fit
funcfit function

Definition at line 66 of file BCHistogramFitter.cxx.

67  : BCFitter()
68  , fHistogram(0)
69  , fFitFunction(0)
71 {
72  SetHistogram(hist);
73  SetFitFunction(func);
74 
77  fFlagIntegration = true;
78  flag_discrete = true;
79 
80  // set MCMC for marginalization
82 }
BCHistogramFitter::BCHistogramFitter ( const char *  name,
TH1D *  hist,
TF1 *  func 
)

Constructor

Parameters
namename of the model
histhistogram to fit
funcfit function

Definition at line 86 of file BCHistogramFitter.cxx.

87  : BCFitter(name)
88  , fHistogram(0)
89  , fFitFunction(0)
91 {
92  SetHistogram(hist);
93  SetFitFunction(func);
94 
97  fFlagIntegration = true;
98  flag_discrete = true;
99 
100  // set MCMC for marginalization
102 }
BCHistogramFitter::~BCHistogramFitter ( )

The default destructor.

Definition at line 259 of file BCHistogramFitter.cxx.

260 {
261  // todo memory leak, many members not removed
262  delete fHistogramExpected;
263 }

Member Function Documentation

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

Calculate the p-value using fast-MCMC and the likelihood as test statistic. The method is explained in the appendix of http://arxiv.org/abs/1011.1674

Note
Obtain the results with
See Also
BCModel::GetPValue and the value corrected for degrees of freedom with
BCModel::GetPValueNDoF
Parameters
parThe parameter values for which expectations are computed
nIterationsnumber of pseudo experiments generated by the Markov chain
Returns
An error code

Definition at line 442 of file BCHistogramFitter.cxx.

443 {
444  // check size of parameter vector
445  if (par.size() != GetNParameters()) {
447  "BCHistogramFitter::CalculatePValueFast : Number of parameters is inconsistent.");
448  return 0;
449  }
450 
451  // check if histogram exists
452  if (!fHistogram) {
454  "BCHistogramFitter::CalculatePValueFast : Histogram not defined.");
455  return 0;
456  }
457 
458  //update the expected number of events for all bins
460 
461  // copy observed and expected values
462  std::vector<unsigned> observed(fHistogram->GetNbinsX());
463  std::vector<double> expected(fHistogramExpected->GetNbinsX());
464 
465  for (int ibin = 0 ; ibin < fHistogram->GetNbinsX(); ++ibin){
466  observed[ibin] = (unsigned int) fHistogram->GetBinContent(ibin + 1);
467  expected[ibin] = fHistogramExpected->GetBinContent(ibin + 1);
468  }
469 
470  // create pseudo experiments
471  fPValue = BCMath::FastPValue(observed, expected, nIterations, fRandom->GetSeed());
472 
473  // correct for fitting bias
475 
476  // no error
477  return 1;
478 }
int BCHistogramFitter::CalculatePValueKolmogorov ( const std::vector< double > &  par)

Calculate the p-value using Kolmogorov-Smirnov test. Note that the reference distribution is known only asymptotically. Some explanation is given in http://root.cern.ch/root/htmldoc/TMath.html

Parameters
parThe set of parameter values used in the model, usually the best fit parameters
pvalueThe pvalue
Returns
An error code

Definition at line 554 of file BCHistogramFitter.cxx.

555 {
556  if (!fHistogramExpected || !fHistogram) {
557  BCLog::OutError("BCHistogramFitter::CalculatePValueKolmogorov: "
558  "Please define the reference distribution by calling \n"
559  "BCHistogramFitter::SetHistogramExpected() first!");
560  return 0;
561  }
562 
563  //update expected counts with current parameters
565 
566  //perform the test
567  fPValue = fHistogramExpected->KolmogorovTest(fHistogram);
569 
570  // no error
571  return 1;
572 }
int BCHistogramFitter::CalculatePValueLeastSquares ( const std::vector< double > &  par,
bool  weightExpect = true 
)

Calculate the p-value using approximate chi^2 distribution of squared difference for conventional weights. Approximation is valid for bin contents >5 and not as as good for little data as CalculatePValueLikelihood, see eq. (32.13), PDG: Statistics, Monte Carlo, Group Theory. Physics Letters B 667, 316-339(2008).

Parameters
parThe set of parameter values used in the model, usually the best fit parameters
pvalueThe pvalue
weightuse the variance from the expected #counts (true) or the measured counts (false)
Returns
An error code

Definition at line 517 of file BCHistogramFitter.cxx.

518 {
519  // initialize test statistic chi^2
520  double chi2 = 0.0;
521 
522  //Calculate expected counts to compare with
524 
525  for (int ibin = 1; ibin <= fHistogram->GetNbinsX(); ++ibin) {
526 
527  // get the number of observed events
528  double y = fHistogram->GetBinContent(ibin);
529 
530  // get the number of expected events
531  double yexp = fHistogramExpected->GetBinContent(ibin);
532 
533  //convert 1/0.0 into 1 for weighted sum
534  double weight;
535  if (weightExpect)
536  weight = (yexp > 0) ? yexp : 1.0;
537  else
538  weight = (y > 0) ? y : 1.0;
539 
540  // get the contribution from this datapoint
541  chi2 += (y - yexp) * (y - yexp) / weight;
542  }
543 
544  // p value from chi^2 distribution
545  fPValue = TMath::Prob(chi2, GetNDataPoints());
546  fPValueNDoF = TMath::Prob(chi2, GetNDataPoints() - GetNParameters());
547 
548  // no error
549  return 1;
550 }
int BCHistogramFitter::CalculatePValueLikelihood ( const std::vector< double > &  par)

Calculate the p-value using approximate chi^2 distribution of scaled likelihood. Approximation is valid for bin contents >5, see eq. (32.12), PDG: Statistics, Monte Carlo, Group Theory. Physics Letters B 667, 316-339(2008).

Parameters
parThe set of parameter values used in the model, usually the best fit parameters
pvalueThe pvalue
Returns
An error code

Definition at line 481 of file BCHistogramFitter.cxx.

482 {
483  // initialize test statistic -2*lambda
484  double logLambda = 0.0;
485 
486  //Calculate expected counts to compare with
488 
489  for (int ibin = 1; ibin <= fHistogram->GetNbinsX(); ++ibin) {
490 
491  // get the number of observed events
492  double y = fHistogram->GetBinContent(ibin);
493 
494  // get the number of expected events
495  double yexp = fHistogramExpected->GetBinContent(ibin);
496 
497  // get the contribution from this datapoint
498  if (y == 0)
499  logLambda += yexp;
500  else
501  logLambda += yexp - y + y * log(y / yexp);
502  }
503 
504  // rescale
505  logLambda *= 2.0;
506 
507  //p value from chi^2 distribution, returns zero if logLambda < 0
508  fPValue = TMath::Prob(logLambda, GetNDataPoints());
509  fPValueNDoF = TMath::Prob(logLambda, GetNDataPoints() - GetNParameters());
510 
511  // no error
512  return 1;
513 }
double BCHistogramFitter::CDF ( const std::vector< double > &  ,
int  ,
bool  = false 
)
virtual

1dim cumulative distribution function of the probability to get the data f(x_i|param) for a single measurement, assumed to be of identical functional form for all measurements

Parameters
parametersThe parameter values at which point to compute the cdf
indexThe data point index starting at 0,1...N-1
loweronly needed for discrete distributions! Return the CDF for the count one less than actually observed, e.g. in Poisson process, if 3 actually observed, then CDF(2) is returned

Reimplemented from BCModel.

Definition at line 576 of file BCHistogramFitter.cxx.

577 {
578 
579  if (parameters.empty())
580  BCLog::OutWarning("BCHistogramFitter::CDF: parameter vector empty!");
581  //histogram stores underflow in bin 0, so datapoint 0 is in bin 1
582  index++;
583 
584  // get the number of observed events (should be integer)
585  double yObs = fHistogram->GetBinContent(index);
586 
587  // get function value of lower bin edge
588  double edgeLow = fHistogram->GetBinLowEdge(index);
589  double edgeHigh = edgeLow + fHistogram->GetBinWidth(index);
590 
591  // expectation value of this bin
592  double yExp = 0.0;
593 
594  // use ROOT's TH1D::Integral method
595  if (fFlagIntegration)
596  yExp = fFitFunction->Integral(edgeLow, edgeHigh, &parameters[0]);
597 
598  // use linear interpolation
599  else {
600  double dx = fHistogram->GetBinWidth(index);
601  double fEdgeLow = fFitFunction->Eval(edgeLow);
602  double fEdgeHigh = fFitFunction->Eval(edgeHigh);
603  yExp = fEdgeLow * dx + (fEdgeHigh - fEdgeLow) * dx / 2.;
604  }
605 
606  // usually Poisson bins do not agree with fixed probability bins
607  // so choose where it should belong
608 
609  if (lower) {
610  if ((int) yObs >= 1)
611  return ROOT::Math::poisson_cdf((unsigned int) (yObs - 1), yExp);
612  else
613  return 0.0;
614  }
615  // what if yObs as double doesn't reprepsent a whole number? exception?
616  if ((double) (unsigned int) yObs != yObs) {
617  BCLog::OutWarning(Form(
618  "BCHistogramFitter::CDF: Histogram values should be integer!\n"
619  " Bin %d = %f", index, yObs));
620 
621  //convert randomly to integer
622  // ex. yObs = 9.785 =>
623  // yObs --> 10 with P = 78.5%
624  // yObs --> 9 with P = 21.5%
625  double U = fRandom->Rndm();
626  double yObsLower = (unsigned int) yObs;
627  if (U > (yObs - yObsLower))
628  yObs = yObsLower;
629  else
630  yObs = yObsLower + 1;
631  }
632 
633  return ROOT::Math::poisson_cdf((unsigned int) yObs, yExp);
634 }
void BCHistogramFitter::DrawFit ( const char *  options = "HIST",
bool  flaglegend = false 
)
virtual

Draw the fit in the current pad.

Implements BCFitter.

Definition at line 389 of file BCHistogramFitter.cxx.

390 {
391  if (!fHistogram) {
392  BCLog::OutError("BCHistogramFitter::DrawFit : Histogram not defined.");
393  return;
394  }
395 
396  if (!fFitFunction) {
397  BCLog::OutError("BCHistogramFitter::DrawFit : Fit function not defined.");
398  return;
399  }
400 
401  if (!fErrorBandXY || GetBestFitParameters().empty()) {
402  BCLog::OutError("BCHistogramFitter::DrawFit : Fit not performed yet.");
403  return;
404  }
405 
406  // check wheather options contain "same"
407  TString opt = options;
408  opt.ToLower();
409 
410  // if not same, draw the histogram first to get the axes
411  if (!opt.Contains("same"))
412  fHistogram->Draw(opt.Data());
413 
414  // draw the error band as central 68% probability interval
415  fErrorBand = GetErrorBandGraph(0.16, 0.84);
416  fErrorBand->Draw("f same");
417 
418  // now draw the histogram again since it was covered by the band
419  fHistogram->Draw(TString::Format("%ssame", opt.Data()));
420 
421  // draw the fit function on top
423  fGraphFitFunction->SetLineColor(kRed);
424  fGraphFitFunction->SetLineWidth(2);
425  fGraphFitFunction->Draw("l same");
426 
427  // draw legend
428  if (flaglegend) {
429  TLegend * legend = new TLegend(0.25, 0.75, 0.55, 0.95);
430  legend->SetBorderSize(0);
431  legend->SetFillColor(kWhite);
432  legend->AddEntry(fHistogram, "Data", "L");
433  legend->AddEntry(fGraphFitFunction, "Best fit", "L");
434  legend->AddEntry(fErrorBand, "Error band", "F");
435  legend->Draw();
436  }
437 
438  gPad->RedrawAxis();
439 }
int BCHistogramFitter::Fit ( )
virtual

Performs the fit.

Returns
An error code.

Implements BCFitter.

Definition at line 351 of file BCHistogramFitter.cxx.

352 {
353  // set histogram
354  if (!fHistogram) {
355  BCLog::OutError("BCHistogramFitter::Fit : Histogram not defined.");
356  return 0;
357  }
358 
359  // set function
360  if (!fFitFunction) {
361  BCLog::OutError("BCHistogramFitter::Fit : Fit function not defined.");
362  return 0;
363  }
364 
365  // perform marginalization
366  MarginalizeAll();
367 
368  // maximize posterior probability, using the best-fit values close
369  // to the global maximum from the MCMC
373  SetOptimizationMethod(method_temp);
374 
375  // calculate the p-value using the fast MCMC algorithm
378  "BCHistogramFitter::Fit : Could not use the fast p-value evaluation.");
379 
380  // print summary to screen
382 
383  // no error
384  return 1;
385 }
int BCHistogramFitter::Fit ( TH1D *  hist,
TF1 *  func 
)

Performs the fit.

Parameters
histThe histogram (TH1D).
funcThe fit function.
Returns
An error code.

Definition at line 328 of file BCHistogramFitter.cxx.

329 {
330  // set histogram
331  if (hist)
332  SetHistogram(hist);
333  else {
334  BCLog::OutError("BCHistogramFitter::Fit : Histogram not defined.");
335  return 0;
336  }
337 
338  // set function
339  if (func)
340  SetFitFunction(func);
341  else {
342  BCLog::OutError("BCHistogramFitter::Fit : Fit function not defined.");
343  return 0;
344  }
345 
346  return Fit();
347 }
double BCHistogramFitter::FitFunction ( const std::vector< double > &  x,
const std::vector< double > &  parameters 
)
virtual

Plots the histogram

Parameters
optionsOptions for plotting.
filenameName of the file which the histogram is printed into. The following options are available:
F : plots the fit function on top of the data E0 : plots the fit function and the 68% prob. uncertainty band of the fit function on top of the data E1 : plots the expectation from the fit function and the uncertainty bin-by-bin as error bars. Returns the y-value of the 1-dimensional fit function at an x and for a set of parameters.
xA vector with the x-value.
parametersA set of parameters.

Reimplemented from BCFitter.

Definition at line 318 of file BCHistogramFitter.cxx.

319 {
320  // set the parameters of the function
321  fFitFunction->SetParameters(&params[0]);
322 
323  return fFitFunction->Eval(x[0]) * fHistogram->GetBinWidth(fHistogram->FindBin(x[0]));
324 }
TGraph* BCHistogramFitter::GetErrorBand ( )
inline
Returns
pointer to the error band

Definition at line 91 of file BCHistogramFitter.h.

92  { return fErrorBand; };
TF1* BCHistogramFitter::GetFitFunction ( )
inline
Returns
The fit function

Definition at line 86 of file BCHistogramFitter.h.

87  { return fFitFunction; };
TGraph* BCHistogramFitter::GetGraphFitFunction ( )
inline
Returns
pointer to a graph for the fit function

Definition at line 96 of file BCHistogramFitter.h.

97  { return fGraphFitFunction; };
TH1D* BCHistogramFitter::GetHistogram ( )
inline
Returns
The data histogram

Definition at line 75 of file BCHistogramFitter.h.

76  { return fHistogram; };
TH1D* BCHistogramFitter::GetHistogramExpected ( )
inline
Returns
The histogram of expected counts

Definition at line 80 of file BCHistogramFitter.h.

81  { return fHistogramExpected; };
double BCHistogramFitter::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 267 of file BCHistogramFitter.cxx.

268 {
269  // initialize probability
270  double loglikelihood = 0;
271 
272  // set the parameters of the function
273  fFitFunction->SetParameters(&params[0]);
274 
275  // get the number of bins
276  int nbins = fHistogram->GetNbinsX();
277 
278  // get bin width
279  double dx = fHistogram->GetBinWidth(1);
280 
281  // get function value of lower bin edge
282  double fedgelow = fFitFunction->Eval(fHistogram->GetBinLowEdge(1));
283 
284  // loop over all bins
285  for (int ibin = 1; ibin <= nbins; ++ibin) {
286  // get upper bin edge
287  double xedgehi = fHistogram->GetBinLowEdge(ibin) + dx;
288 
289  // get function value at upper bin edge
290  double fedgehi = fFitFunction->Eval(xedgehi);
291 
292  // get the number of observed events
293  double y = fHistogram->GetBinContent(ibin);
294 
295  double yexp = 0.;
296 
297  // use ROOT's TH1D::Integral method
298  if (fFlagIntegration)
299  yexp = fFitFunction->Integral(xedgehi - dx, xedgehi);
300 
301  // use linear interpolation
302  else {
303  yexp = fedgelow * dx + (fedgehi - fedgelow) * dx / 2.;
304 
305  // make the upper edge the lower edge for the next iteration
306  fedgelow = fedgehi;
307  }
308 
309  // get the value of the Poisson distribution
310  loglikelihood += BCMath::LogPoisson(y, yexp);
311  }
312 
313  return loglikelihood;
314 }
int BCHistogramFitter::SetFitFunction ( TF1 *  func)
Parameters
funcThe fit function
Returns
An error code (1:pass, 0:fail).

Definition at line 220 of file BCHistogramFitter.cxx.

221 {
222  // check if function exists
223  if (!func) {
224  BCLog::OutError("BCHistogramFitter::SetFitFunction : TF1 not created.");
225  return 0;
226  }
227 
228  // set the function
229  fFitFunction = func;
230 
231  // update the model name to contain the function name
232  if(fName=="model")
233  SetName(TString::Format("HistogramFitter with %s", fFitFunction->GetName()));
234 
235  // reset parameters
236  ClearParameters(true);
237 
238  // get the new number of parameters
239  int n = func->GetNpar();
240 
241  // add parameters
242  for (int i = 0; i < n; ++i) {
243  double xmin;
244  double xmax;
245 
246  func->GetParLimits(i, xmin, xmax);
247 
248  AddParameter(func->GetParName(i), xmin, xmax);
249  }
250 
251  // set flat prior for all parameters by default
253 
254  return GetNParameters();
255 }
void BCHistogramFitter::SetFlagIntegration ( bool  flag)
inline

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

Definition at line 126 of file BCHistogramFitter.h.

127  { fFlagIntegration = flag; };
int BCHistogramFitter::SetHistogram ( TH1D *  hist)
Parameters
histThe histogram containing the data
Returns
An error code (1:pass, 0:fail).

Definition at line 106 of file BCHistogramFitter.cxx.

107 {
108  // check if histogram exists
109  if (!hist) {
110  BCLog::OutError("BCHistogramFitter::SetHistogram : TH1D not created.");
111  return 0;
112  }
113 
114  // set pointer to histogram
115  fHistogram = hist;
116 
117  // create a data set. this is necessary in order to calculate the
118  // error band. the data set contains as many data points as there
119  // are bins.
120  BCDataSet * ds = new BCDataSet();
121 
122  // create data points and add them to the data set.
123  // the x value is the lower edge of the bin, and
124  // the y value is the bin count
125  int nbins = fHistogram->GetNbinsX();
126  for (int i = 0; i < nbins; ++i) {
127  BCDataPoint* dp = new BCDataPoint(2);
128  dp ->SetValue(0, fHistogram->GetBinLowEdge(i + 1));
129  dp ->SetValue(1, fHistogram->GetBinContent(i + 1));
130  ds->AddDataPoint(dp);
131  }
132 
133  // set the new data set.
134  SetDataSet(ds);
135 
136  // calculate the lower and upper edge in x.
137  double xmin = hist->GetBinLowEdge(1);
138  double xmax = hist->GetBinLowEdge(nbins + 1);
139 
140  // calculate the minimum and maximum range in y.
141  double histymin = hist->GetMinimum();
142  double histymax = hist->GetMaximum();
143 
144  // calculate the minimum and maximum of the function value based on
145  // the minimum and maximum value in y.
146  double ymin = TMath::Max(0., histymin - 5. * sqrt(histymin));
147  double ymax = histymax + 5. * sqrt(histymax);
148 
149  // set the data boundaries for x and y values.
150  SetDataBoundaries(0, xmin, xmax);
151  SetDataBoundaries(1, ymin, ymax);
152 
153  // set the indeces for fitting.
154  SetFitFunctionIndices(0, 1);
155 
156  // no error
157  return 1;
158 }
int BCHistogramFitter::SetHistogramExpected ( const std::vector< double > &  parameters)
Parameters
histThe histogram with the expected counts (typically non-integer values!)
Returns
An error code (1:pass, 0:fail).

Definition at line 162 of file BCHistogramFitter.cxx.

163 {
164  //clear up previous memory
165  if (fHistogramExpected) {
166  delete fHistogramExpected;
167  }
168  //copy all properties from the data histogram
169  fHistogramExpected = new TH1D(*fHistogram);
170 
171  // get the number of bins
172  int nBins = fHistogramExpected->GetNbinsX();
173 
174  // get bin width
175  double dx = fHistogramExpected->GetBinWidth(1);
176 
177  //set the parameters of fit function
178  fFitFunction->SetParameters(&parameters[0]);
179 
180  // get function value of lower bin edge
181  double fedgelow = fFitFunction->Eval(fHistogramExpected->GetBinLowEdge(1));
182 
183  // loop over all bins, skip underflow
184  for (int ibin = 1; ibin <= nBins; ++ibin) {
185  // get upper bin edge
186  double xedgehi = fHistogramExpected->GetBinLowEdge(ibin) + dx;
187 
188  //expected count
189  double yexp = 0.;
190 
191  // use ROOT's TH1D::Integral method
192  if (fFlagIntegration)
193  yexp = fFitFunction->Integral(xedgehi - dx, xedgehi);
194 
195  // use linear interpolation
196  else {
197  // get function value at upper bin edge
198  double fedgehi = fFitFunction->Eval(xedgehi);
199  yexp = fedgelow * dx + (fedgehi - fedgelow) * dx / 2.;
200 
201  // make the upper edge the lower edge for the next iteration
202  fedgelow = fedgehi;
203  }
204 
205  //write the expectation for the bin
206  fHistogramExpected->SetBinContent(ibin, yexp);
207 
208  //avoid automatic error as sqrt(yexp), used e.g. in Kolmogorov correction factor
209  fHistogramExpected->SetBinError(ibin, 0.0);
210 
211  // but the data under this model have that sqrt(yexp) uncertainty
212  fHistogram->SetBinError(ibin, sqrt(yexp));
213 
214  }
215  return 1;
216 }

Member Data Documentation

TGraph* BCHistogramFitter::fErrorBand
private

Pointer to the error band (for legend)

Definition at line 240 of file BCHistogramFitter.h.

TF1* BCHistogramFitter::fFitFunction
private

The fit function

Definition at line 231 of file BCHistogramFitter.h.

bool BCHistogramFitter::fFlagIntegration
private

Flag for using the ROOT TH1D::Integral method (true), or linear interpolation (false)

Definition at line 236 of file BCHistogramFitter.h.

TGraph* BCHistogramFitter::fGraphFitFunction
private

Pointer to a graph for displaying the fit function

Definition at line 244 of file BCHistogramFitter.h.

TH1D* BCHistogramFitter::fHistogram
private

The histogram containing the data.

Definition at line 227 of file BCHistogramFitter.h.

TH1D* BCHistogramFitter::fHistogramExpected
private

The histogram containing the expected data.

Definition at line 249 of file BCHistogramFitter.h.


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