BayesianAnalysisToolkit  0.9.3
Classes | Private Attributes | List of all members
BCEfficiencyFitter Class Reference

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

#include <BCEfficiencyFitter.h>

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

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)
 
- 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 * fHistogram1
 
TH1D * fHistogram2
 
TGraphAsymmErrors * fHistogramRatio
 
TF1 * fFitFunction
 
bool fFlagIntegration
 
TGraph * fErrorBand
 
TGraph * fGraphFitFunction
 
TH1D * fHistogramBinomial
 
unsigned int fDataPointType
 

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 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.

33  : BCFitter()
34  , fHistogram1(0)
35  , fHistogram2(0)
36  , fFitFunction(0)
38  , fDataPointType(1)
39 {
40  // set default options and values
41  fFlagIntegration = false;
42 
43  // set MCMC for marginalization
45 }
BCEfficiencyFitter::BCEfficiencyFitter ( const char *  name)

Constructor

Parameters
namename fo the model

Definition at line 49 of file BCEfficiencyFitter.cxx.

50  : BCFitter(name)
51  , fHistogram1(0)
52  , fHistogram2(0)
53  , fFitFunction(0)
55  , fDataPointType(1)
56 {
57  // set default options and values
58  fFlagIntegration = false;
59 
60  // set MCMC for marginalization
62 }
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.

67  : BCFitter()
68  , fHistogram1(0)
69  , fHistogram2(0)
70  , fFitFunction(0)
72  , fDataPointType(1)
73 {
74  SetHistograms(hist1, hist2);
75  SetFitFunction(func);
76 
77  fFlagIntegration = false;
78 
79  // set MCMC for marginalization
81 }
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.

86  : BCFitter(name)
87  , fHistogram1(0)
88  , fHistogram2(0)
89  , fFitFunction(0)
91  , fDataPointType(1)
92 {
93  SetHistograms(hist1, hist2);
94  SetFitFunction(func);
95 
97  fFlagIntegration = false;
98 
99  // set MCMC for marginalization
101 }
BCEfficiencyFitter::~BCEfficiencyFitter ( )

The default destructor.

Definition at line 211 of file BCEfficiencyFitter.cxx.

212 {
213  delete fHistogram1;
214  delete fHistogram2;
215  delete fHistogramBinomial;
216 }

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.

454 {
455  // check size of parameter vector
456  if (par.size() != GetNParameters()) {
457  BCLog::OutError("BCEfficiencyFitter::CalculatePValueFast : Number of parameters is inconsistent.");
458  return 0;
459  }
460 
461  // check if histogram exists
462  if (!fHistogram1 || !fHistogram2) {
463  BCLog::OutError("BCEfficiencyFitter::CalculatePValueFast : Histogram not defined.");
464  return 0;
465  }
466 
467  // define temporary variables
468  int nbins = fHistogram1->GetNbinsX();
469 
470  std::vector<unsigned> histogram(nbins, 0);
471  std::vector<double> expectation(nbins, 0);
472 
473  double logp = 0;
474  double logp_start = 0;
475  int counter_pvalue = 0;
476 
477  // define starting distribution
478  for (int ibin = 0; ibin < nbins; ++ibin) {
479  // get bin boundaries
480  double xmin = fHistogram1->GetBinLowEdge(ibin+1);
481  double xmax = fHistogram1->GetBinLowEdge(ibin+2);
482 
483  // get the number of expected events
484  double yexp = fFitFunction->Integral(xmin, xmax);
485 
486  // get n
487  int n = int(fHistogram1->GetBinContent(ibin));
488 
489  // get k
490  int k = int(fHistogram2->GetBinContent(ibin));
491 
492  histogram[ibin] = k;
493  expectation[ibin] = n * yexp;
494 
495  // calculate p;
496  logp += BCMath::LogApproxBinomial(n, k, yexp);
497  }
498  logp_start = logp;
499 
500  // loop over iterations
501  for (unsigned iiter = 0; iiter < nIterations; ++iiter)
502  {
503  // loop over bins
504  for (int ibin = 0; ibin < nbins; ++ibin)
505  {
506  // get n
507  int n = int(fHistogram1->GetBinContent(ibin));
508 
509  // get k
510  int k = histogram.at(ibin);
511 
512  // get expectation
513  double yexp = 0;
514  if (n > 0)
515  yexp = expectation.at(ibin)/n;
516 
517  // random step up or down in statistics for this bin
518  double ptest = fRandom->Rndm() - 0.5;
519 
520  // continue, if efficiency is at limit
521  if (!(yexp > 0. || yexp < 1.0))
522  continue;
523 
524  // increase statistics by 1
525  if (ptest > 0 && (k < n))
526  {
527  // calculate factor of probability
528  double r = (double(n)-double(k))/(double(k)+1.) * yexp / (1. - yexp);
529 
530  // walk, or don't (this is the Metropolis part)
531  if (fRandom->Rndm() < r)
532  {
533  histogram[ibin] = k + 1;
534  logp += log(r);
535  }
536  }
537 
538  // decrease statistics by 1
539  else if (ptest <= 0 && (k > 0))
540  {
541  // calculate factor of probability
542  double r = double(k) / (double(n)-(double(k)-1)) * (1-yexp)/yexp;
543 
544  // walk, or don't (this is the Metropolis part)
545  if (fRandom->Rndm() < r)
546  {
547  histogram[ibin] = k - 1;
548  logp += log(r);
549  }
550  }
551 
552  } // end of looping over bins
553 
554  //one new toy data set is created
555  //call user interface to calculate arbitrary statistic's distribution
556  if (callback)
557  (*callback)(expectation, histogram);
558 
559  // increase counter
560  if (logp < logp_start)
561  counter_pvalue++;
562 
563  } // end of looping over iterations
564 
565  // calculate p-value
566  pvalue = double(counter_pvalue) / double(nIterations);
567 
568  // correct for fit bias
569  pvalueCorrected = BCMath::CorrectPValue(pvalue, GetNParameters(), nbins);
570 
571  // no error
572  return 1;
573 }
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.

445 {
446  //use NULL pointer for no callback
447  return CalculatePValueFast(par, NULL, pvalue, pvalueCorrected, nIterations);
448 }
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.

350 {
351  if (!fHistogram1 || !fHistogram2) {
352  BCLog::OutError("BCEfficiencyFitter::DrawFit : Histogram(s) not defined.");
353  return;
354  }
355 
356  if (!fFitFunction) {
357  BCLog::OutError("BCEfficiencyFitter::DrawFit : Fit function not defined.");
358  return;
359  }
360 
361  // create efficiency graph
362  TGraphAsymmErrors * histRatio = new TGraphAsymmErrors();
363  histRatio->SetMarkerStyle(20);
364  histRatio->SetMarkerSize(1.5);
365 
366  int npoints = 0;
367 
368  // set points
369  for (int i = 1; i <= fHistogram1->GetNbinsX(); ++i)
370  {
371  if(int(fHistogram1->GetBinContent(i)) == 0) {
372  ++npoints;
373  continue;
374  }
375 
376  // calculate central value and uncertainties
377  double xexp, xmin, xmax;
378  GetUncertainties( int(fHistogram1->GetBinContent(i)), int(fHistogram2->GetBinContent(i)), 0.68, xexp, xmin, xmax);
379 
380  histRatio->SetPoint( npoints, fHistogram1->GetBinCenter(i), xexp);
381 
382  // no X uncertainties
383  histRatio->SetPointEXhigh(npoints, 0.);
384  histRatio->SetPointEXlow(npoints, 0.);
385 
386  // set Y uncertainties
387  histRatio->SetPointEYhigh(npoints, xmax - xexp);
388  histRatio->SetPointEYlow(npoints, xexp - xmin);
389 
390  ++npoints;
391  }
392 
393 
394  // check wheather options contain "same"
395  TString opt = options;
396  opt.ToLower();
397 
398  // if not same, draw the histogram first to get the axes
399  if(!opt.Contains("same"))
400  {
401  // create new histogram
402  TH2D * hist_axes = new TH2D("hist_axes",
403  Form(";%s;ratio", fHistogram1->GetXaxis()->GetTitle()),
404  fHistogram1->GetNbinsX(),
405  fHistogram1->GetXaxis()->GetBinLowEdge(1),
406  fHistogram1->GetXaxis()->GetBinLowEdge(fHistogram1->GetNbinsX()+1),
407  1, 0., 1.);
408  hist_axes->SetStats(kFALSE);
409  hist_axes->Draw();
410 
411  histRatio->Draw(TString::Format("%sp",opt.Data()));
412  }
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  histRatio->SetMarkerSize(.7);
420  histRatio->Draw(TString::Format("%ssamep",opt.Data()));
421 
422  // draw the fit function on top
424  fGraphFitFunction->SetLineColor(kRed);
425  fGraphFitFunction->SetLineWidth(2);
426  fGraphFitFunction->Draw("l same");
427 
428  // draw legend
429  if (flaglegend)
430  {
431  TLegend * legend = new TLegend(0.25, 0.75, 0.55, 0.95);
432  legend->SetBorderSize(0);
433  legend->SetFillColor(kWhite);
434  legend->AddEntry(histRatio, "Data", "L");
435  legend->AddEntry(fGraphFitFunction, "Best fit", "L");
436  legend->AddEntry(fErrorBand, "Error band", "F");
437  legend->Draw();
438  }
439  gPad->RedrawAxis();
440 }
int BCEfficiencyFitter::Fit ( )
virtual

Performs the fit.

Returns
An error code.

Implements BCFitter.

Definition at line 309 of file BCEfficiencyFitter.cxx.

310 {
311  // set histogram
312  if (!fHistogram1 || !fHistogram2) {
313  BCLog::OutError("BCEfficiencyFitter::Fit : Histogram(s) not defined.");
314  return 0;
315  }
316 
317  // set function
318  if (!fFitFunction) {
319  BCLog::OutError("BCEfficiencyFitter::Fit : Fit function not defined.");
320  return 0;
321  }
322 
323  // perform marginalization
324  MarginalizeAll();
325 
326  // maximize posterior probability, using the best-fit values close
327  // to the global maximum from the MCMC
331  SetOptimizationMethod(method_temp);
332 
333  // calculate the p-value using the fast MCMC algorithm
334  double pvalue, pvalueCorrected;
335  if ( CalculatePValueFast(GetBestFitParameters(), pvalue, pvalueCorrected) )
336  fPValue = pvalue;
337  else
338  BCLog::OutError("BCEfficiencyFitter::Fit : Could not use the fast p-value evaluation.");
339 
340  // print summary to screen
342 
343  // no error
344  return 1;
345 }
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.

287 {
288  // set histogram
289  if (hist1 && hist2)
290  SetHistograms(hist1, hist2);
291  else {
292  BCLog::OutError("BCEfficiencyFitter::Fit : Histogram(s) not defined.");
293  return 0;
294  }
295 
296  // set function
297  if (func)
298  SetFitFunction(func);
299  else {
300  BCLog::OutError("BCEfficiencyFitter::Fit : Fit function not defined.");
301  return 0;
302  }
303 
304  return Fit();
305 }
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.

277 {
278  // set the parameters of the function
279  fFitFunction->SetParameters(&params[0]);
280 
281  return fFitFunction->Eval(x[0]);
282 }
TGraph* BCEfficiencyFitter::GetErrorBand ( )
inline
Returns
pointer to the error band

Definition at line 120 of file BCEfficiencyFitter.h.

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

Definition at line 115 of file BCEfficiencyFitter.h.

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

Definition at line 125 of file BCEfficiencyFitter.h.

126  { return fGraphFitFunction; };
TH1D* BCEfficiencyFitter::GetHistogram1 ( )
inline
Returns
The histogram 1

Definition at line 100 of file BCEfficiencyFitter.h.

101  { return fHistogram1; };
TH1D* BCEfficiencyFitter::GetHistogram2 ( )
inline
Returns
The histogram 2

Definition at line 105 of file BCEfficiencyFitter.h.

106  { return fHistogram2; };
TGraphAsymmErrors* BCEfficiencyFitter::GetHistogramRatio ( )
inline
Returns
The histogram ratio

Definition at line 110 of file BCEfficiencyFitter.h.

111  { return fHistogramRatio; };
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.

577 {
578  if (n == 0)
579  {
580  xexp = 0.;
581  xmin = 0.;
582  xmax = 0.;
583  return 0;
584  }
585 
586  BCLog::OutDebug(Form("Calculating efficiency data-point of type %d for (n,k) = (%d,%d)",fDataPointType,n,k));
587 
588  // create a histogram with binomial distribution
589  if (fHistogramBinomial)
590  fHistogramBinomial->Reset();
591  else
592  fHistogramBinomial = new TH1D("hist_binomial", "", 1000, 0., 1.);
593 
594  // loop over bins and fill histogram
595  for (int i = 1; i <= 1000; ++i) {
596  double x = fHistogramBinomial->GetBinCenter(i);
597  double val = BCMath::ApproxBinomial(n, k, x);
598  fHistogramBinomial->SetBinContent(i, val);
599  }
600 
601  // normalize
602  fHistogramBinomial->Scale(1. / fHistogramBinomial->Integral());
603 
604  // calculate central value and uncertainties
605  if (fDataPointType == 0) {
606  xexp = fHistogramBinomial->GetMean();
607  double rms = fHistogramBinomial->GetRMS();
608  xmin = xexp-rms;
609  xmax = xexp+rms;
610  BCLog::OutDebug(Form(" - mean = %f , rms = %f)",xexp,rms));
611  }
612  else if (fDataPointType == 1) {
613  xexp = (double)k/(double)n;
614  BCH1D * fbh = new BCH1D((TH1D*)fHistogramBinomial->Clone("hcp"));
615  std::vector<double> intervals = fbh->GetSmallestIntervals(p);
616  int ninter = intervals.size();
617  if ( ninter<2 ) {
618  xmin = xmax = xexp = 0.;
619  }
620  else {
621  xmin = intervals[0];
622  xmax = intervals[1];
623  }
624  delete fbh;
625  }
626  else {
627  // calculate quantiles
628  int nprobSum = 3;
629  double q[3];
630  double probSum[3];
631  probSum[0] = (1.-p)/2.;
632  probSum[1] = 1.-(1.-p)/2.;
633  probSum[2] = .5;
634 
635  fHistogramBinomial->GetQuantiles(nprobSum, q, probSum);
636 
637  xexp = q[2];
638  xmin = q[0];
639  xmax = q[1];
640 
641  }
642 
643  BCLog::OutDebug(Form(" - efficiency = %f , range (%f - %f)",xexp,xmin,xmax));
644 
645  return 1;
646 }
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.

231 {
232 
233  // initialize probability
234  double loglikelihood = 0;
235 
236  // set the parameters of the function
237  fFitFunction->SetParameters(&params[0]);
238 
239  // get the number of bins
240  int nbins = fHistogram1->GetNbinsX();
241 
242  // get bin width
243  double dx = fHistogram1->GetXaxis()->GetBinWidth(0);
244 
245  // loop over all bins
246  for (int ibin = 1; ibin <= nbins; ++ibin)
247  {
248  // get n
249  int n = int(fHistogram1->GetBinContent(ibin));
250 
251  // get k
252  int k = int(fHistogram2->GetBinContent(ibin));
253 
254  // get x
255  double x = fHistogram1->GetBinCenter(ibin);
256 
257  double eff = 0;
258 
259  // use ROOT's TH1D::Integral method
260  if (fFlagIntegration)
261  eff = fFitFunction->Integral(x - dx/2., x + dx/2.) / dx;
262 
263  // use linear interpolation
264  else
265  eff = (fFitFunction->Eval(x + dx/2.) + fFitFunction->Eval(x - dx/2.)) / 2.;
266 
267  // get the value of the Poisson distribution
268  loglikelihood += BCMath::LogApproxBinomial(n, k, eff);
269  }
270 
271  return loglikelihood;
272 }
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.

221 {
222  if(type < 0 || type > 2)
223  BCLog::OutError(Form("BCEfficiencyFitter::SetDataPointType : Unknown data point type %d (should be between 0 and 2).",type));
224  else
225  fDataPointType = type;
226 }
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.

172 {
173  // check if function exists
174  if(!func) {
175  BCLog::OutError("BCEfficiencyFitter::SetFitFunction : TF1 not created.");
176  return 0;
177  }
178 
179  // set the function
180  fFitFunction = func;
181 
182  // update the model name to contain the function name
183  if(fName=="model")
184  SetName(TString::Format("BCEfficiencyFitter with %s",fFitFunction->GetName()));
185 
186  // reset parameters
187  ClearParameters(true);
188 
189  // get the new number of parameters
190  int n = func->GetNpar();
191 
192  // add parameters
193  for (int i = 0; i < n; ++i)
194  {
195  double xmin;
196  double xmax;
197 
198  func->GetParLimits(i, xmin, xmax);
199 
200  AddParameter(func->GetParName(i), xmin, xmax);
201  }
202 
203  // set flat prior for all parameters by default
205 
206  return GetNParameters();
207 }
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.

159  { fFlagIntegration = flag; };
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.

106 {
107  // check if histogram exists
108  if (!hist1 || !hist2) {
109  BCLog::OutError("BCEfficiencyFitter::SetHistograms : TH1D not created.");
110  return 0;
111  }
112 
113  // check compatibility of both histograms : number of bins
114  if (hist1->GetNbinsX() != hist2->GetNbinsX()) {
115  BCLog::OutError("BCEfficiencyFitter::SetHistograms : Histograms do not have the same number of bins.");
116  return 0;
117  }
118 
119  // check compatibility of both histograms : bin content
120  for (int i = 1; i <= hist1->GetNbinsX(); ++i)
121  if (hist1->GetBinContent(i) < hist2->GetBinContent(i)) {
122  BCLog::OutError("BCEfficiencyFitter::SetHistograms : Histogram 1 has fewer entries than histogram 2.");
123  return 0;
124  }
125 
126  // set pointer to histograms
127  fHistogram1 = hist1;
128  fHistogram2 = hist2;
129 
130  // create a data set. this is necessary in order to calculate the
131  // error band. the data set contains as many data points as there
132  // are bins. for now, the data points are empty.
133  BCDataSet * ds = new BCDataSet();
134 
135  // create data points and add them to the data set.
136  int nbins = fHistogram1->GetNbinsX();
137  for (int i = 0; i < nbins; ++i) {
138  BCDataPoint* dp = new BCDataPoint(2);
139  ds->AddDataPoint(dp);
140  }
141 
142  // set the new data set.
143  SetDataSet(ds);
144 
145  // calculate the lower and upper edge in x.
146  double xmin = hist1->GetBinLowEdge(1);
147  double xmax = hist1->GetBinLowEdge(nbins+1);
148 
149 // // calculate the minimum and maximum range in y.
150 // double histymin = hist2->GetMinimum();
151 // double histymax = hist1->GetMaximum();
152 
153 // // calculate the minimum and maximum of the function value based on
154 // // the minimum and maximum value in y.
155 // double ymin = TMath::Max(0., histymin - 5.*sqrt(histymin));
156 // double ymax = histymax + 5.*sqrt(histymax);
157 
158  // set the data boundaries for x and y values.
159  SetDataBoundaries(0, xmin, xmax);
160  SetDataBoundaries(1, 0.0, 1.0);
161 
162  // set the indeces for fitting.
163  SetFitFunctionIndices(0, 1);
164 
165  // no error
166  return 1;
167 }

Member Data Documentation

unsigned int BCEfficiencyFitter::fDataPointType
private

Type of point to plot for efficiency data 0 - mean + rms 1 - mode + smallest 68% 2 - median + central 68%

Definition at line 274 of file BCEfficiencyFitter.h.

TGraph* BCEfficiencyFitter::fErrorBand
private

Pointer to the error band (for legend)

Definition at line 260 of file BCEfficiencyFitter.h.

TF1* BCEfficiencyFitter::fFitFunction
private

The fit function

Definition at line 251 of file BCEfficiencyFitter.h.

bool BCEfficiencyFitter::fFlagIntegration
private

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

Definition at line 256 of file BCEfficiencyFitter.h.

TGraph* BCEfficiencyFitter::fGraphFitFunction
private

Pointer to a graph for displaying the fit function

Definition at line 264 of file BCEfficiencyFitter.h.

TH1D* BCEfficiencyFitter::fHistogram1
private

The histogram containing the larger numbers.

Definition at line 239 of file BCEfficiencyFitter.h.

TH1D* BCEfficiencyFitter::fHistogram2
private

The histogram containing the smaller numbers.

Definition at line 243 of file BCEfficiencyFitter.h.

TH1D* BCEfficiencyFitter::fHistogramBinomial
private

Temporary histogram for calculating the binomial qunatiles

Definition at line 268 of file BCEfficiencyFitter.h.

TGraphAsymmErrors* BCEfficiencyFitter::fHistogramRatio
private

The efficiency histogram.

Definition at line 247 of file BCEfficiencyFitter.h.


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