BayesianAnalysisToolkit  0.9.3
Protected Attributes | Private Attributes | List of all members
BCFitter Class Referenceabstract

A base class for all fitting classes. More...

#include <BCFitter.h>

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

Public Member Functions

Constructors and destructors
 BCFitter ()
 
 BCFitter (const char *name)
 
 ~BCFitter ()
 
Member functions (get)
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
 
Member functions (set)
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)
 
virtual double FitFunction (const std::vector< double > &, const std::vector< double > &)
 
Member functions (miscellaneous methods)
int ReadErrorBandFromFile (const char *file)
 
virtual int Fit ()=0
 
virtual void DrawFit (const char *options, bool flaglegend=false)=0
 
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)
 
virtual double LogLikelihood (const std::vector< double > &params)=0
 
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)
 

Protected Attributes

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
 
- Protected Attributes inherited from BCIntegrate
int fID
 
TMinuit * fMinuit
 
double fMinuitArglist [2]
 
int fMinuitErrorFlag
 
bool fFlagIgnorePrevOptimization
 
double fSAT0
 
double fSATmin
 
TTree * fTreeSA
 
bool fFlagWriteSAToFile
 
int fSANIterations
 
double fSATemperature
 
double fSALogProb
 
std::vector< double > fSAx
 
std::vector< BCH1D * > fMarginalized1D
 
std::vector< BCH2D * > fMarginalized2D
 
bool fFlagMarginalized
 
- Protected Attributes inherited from BCEngineMCMC
BCParameterSet fParameters
 
unsigned fMCMCNChains
 
unsigned fMCMCNLag
 
std::vector< unsigned > fMCMCNIterations
 
int fMCMCCurrentIteration
 
int fMCMCCurrentChain
 
unsigned fMCMCNIterationsUpdate
 
unsigned fMCMCNIterationsUpdateMax
 
int fMCMCNIterationsConvergenceGlobal
 
bool fMCMCFlagConvergenceGlobal
 
unsigned fMCMCNIterationsMax
 
unsigned fMCMCNIterationsRun
 
unsigned fMCMCNIterationsPreRunMin
 
std::vector< int > fMCMCNTrialsTrue
 
unsigned fMCMCNTrials
 
bool fMCMCFlagWriteChainToFile
 
bool fMCMCFlagWritePreRunToFile
 
std::vector< double > fMCMCTrialFunctionScaleFactor
 
std::vector< double > fMCMCTrialFunctionScaleFactorStart
 
bool fMCMCFlagPreRun
 
bool fMCMCFlagRun
 
std::vector< double > fMCMCInitialPosition
 
double fMCMCEfficiencyMin
 
double fMCMCEfficiencyMax
 
int fMCMCFlagInitialPosition
 
bool fMCMCFlagOrderParameters
 
int fMCMCPhase
 
std::vector< double > fMCMCx
 
std::vector< double > fMCMCxMax
 
std::vector< double > fMCMCxMean
 
std::vector< double > fMCMCxVar
 
std::vector< double > fMCMCprob
 
std::vector< double > fMCMCprobMax
 
std::vector< double > fMCMCprobMean
 
std::vector< double > fMCMCprobVar
 
bool fMCMCRValueUseStrict
 
double fMCMCRValueCriterion
 
double fMCMCRValueParametersCriterion
 
double fMCMCRValue
 
std::vector< double > fMCMCRValueParameters
 
TRandom3 * fRandom
 
std::vector< TH1D * > fMCMCH1Marginalized
 
std::vector< TH2D * > fMCMCH2Marginalized
 
std::vector< TTree * > fMCMCTrees
 
std::vector< double > fMCMCBestFitParameters
 
double fMCMCLogMaximum
 
std::vector< double > fMarginalModes
 

Private Attributes

TGraph * fErrorBand
 
bool fFlagFillErrorBand
 
int fFitFunctionIndexX
 
int fFitFunctionIndexY
 
bool fErrorBandContinuous
 
std::vector< double > fErrorBandX
 
unsigned fErrorBandNbinsX
 
unsigned fErrorBandNbinsY
 

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)
 

Detailed Description

A base class for all fitting classes.

Author
Kevin Kröninger
Version
1.0
Date
08.2013 A base class for all fitting classes

Definition at line 27 of file BCFitter.h.

Constructor & Destructor Documentation

BCFitter::BCFitter ( )

The default constructor.

Definition at line 25 of file BCFitter.cxx.

25  : BCModel()
26  , fErrorBand(0)
27  , fFlagFillErrorBand(true)
28  , fFitFunctionIndexX(-1)
29  , fFitFunctionIndexY(-1)
30  , fErrorBandContinuous(true)
31  , fErrorBandX(std::vector<double>(0))
32  , fErrorBandNbinsX(100)
33  , fErrorBandNbinsY(500)
34  , fErrorBandXY(0)
35 {
36 }
BCFitter::BCFitter ( const char *  name)

Constructor

Parameters
namename of the model

Definition at line 39 of file BCFitter.cxx.

39  : BCModel(name)
40  , fErrorBand(0)
41  , fFlagFillErrorBand(true)
42  , fFitFunctionIndexX(-1)
43  , fFitFunctionIndexY(-1)
44  , fErrorBandContinuous(true)
45  , fErrorBandX(std::vector<double>(0))
46  , fErrorBandNbinsX(100)
47  , fErrorBandNbinsY(500)
48  , fErrorBandXY(0)
49 {
50 }
BCFitter::~BCFitter ( )

The default destructor.

Definition at line 53 of file BCFitter.cxx.

54 {
55 }

Member Function Documentation

virtual void BCFitter::DrawFit ( const char *  options,
bool  flaglegend = false 
)
pure virtual

Draw the fit in the current pad.

Implemented in BCEfficiencyFitter, BCHistogramFitter, and BCGraphFitter.

void BCFitter::FillErrorBand ( )

Fill error band histogram for curreent iteration. This method is called from MCMCIterationInterface()

Definition at line 102 of file BCFitter.cxx.

103 {
104  // function fitting
105  if (fFitFunctionIndexX < 0)
106  return;
107 
108  // loop over all possible x values ...
109  if (fErrorBandContinuous) {
110  double x = 0;
111  for (unsigned ix = 0; ix < fErrorBandNbinsX; ix++) {
112  // calculate x
113  x = fErrorBandXY->GetXaxis()->GetBinCenter(ix + 1);
114 
115  // calculate y
116  std::vector<double> xvec;
117  xvec.push_back(x);
118 
119  // loop over all chains
120  for (unsigned ichain = 0; ichain < MCMCGetNChains(); ++ichain) {
121  // calculate y
122  double y = FitFunction(xvec, MCMCGetx(ichain));
123 
124  // fill histogram
125  fErrorBandXY->Fill(x, y);
126  }
127 
128  xvec.clear();
129  }
130  }
131  // ... or evaluate at the data point x-values
132  else {
133  unsigned ndatapoints = fErrorBandX.size();
134  double x = 0;
135 
136  for (unsigned ix = 0; ix < ndatapoints; ++ix) {
137  // calculate x
138  x = fErrorBandX.at(ix);
139 
140  // calculate y
141  std::vector<double> xvec;
142  xvec.push_back(x);
143 
144  // loop over all chains
145  for (unsigned ichain = 0; ichain < MCMCGetNChains(); ++ichain) {
146  // calculate y
147  double y = FitFunction(xvec, MCMCGetx(ichain));
148 
149  // fill histogram
150  fErrorBandXY->Fill(x, y);
151  }
152 
153  xvec.clear();
154  }
155  }
156 }
virtual int BCFitter::Fit ( )
pure virtual

Performs the fit.

Returns
An error code.

Implemented in BCEfficiencyFitter, BCHistogramFitter, and BCGraphFitter.

virtual double BCFitter::FitFunction ( const std::vector< double > &  ,
const std::vector< double > &   
)
inlinevirtual

Defines a fit function.

Parameters
parametersA set of parameter values
xA vector of x-values
Returns
The value of the fit function at the x-values given a set of parameters

Reimplemented in BCEfficiencyFitter, BCHistogramFitter, and BCGraphFitter.

Definition at line 135 of file BCFitter.h.

136  { return 0; }
void BCFitter::FixDataAxis ( unsigned int  index,
bool  fixed 
)

Definition at line 324 of file BCFitter.cxx.

325 {
326  // check if index is within range
327  if (index > fDataSet->GetDataPoint(0)->GetNValues()) {
328  BCLog::OutError("BCFitter::FixDataAxis : Index out of range.");
329  return;
330  }
331 
332  if (fDataFixedValues.size() == 0)
334  false);
335 
336  fDataFixedValues[index] = fixed;
337 }
TGraph* BCFitter::GetErrorBand ( )
inline
Returns
pointer to the error band

Definition at line 54 of file BCFitter.h.

55  { return fErrorBand; };
std::vector< double > BCFitter::GetErrorBand ( double  level) const

Returns a vector of y-values at a certain probability level.

Parameters
levelThe level of probability
Returns
vector of y-values

Definition at line 159 of file BCFitter.cxx.

160 {
161  std::vector<double> errorband;
162 
163  if (!fErrorBandXY)
164  return errorband;
165 
166  int nx = fErrorBandXY->GetNbinsX();
167  errorband.assign(nx, 0.);
168 
169  // loop over x and y bins
170  for (int ix = 1; ix <= nx; ix++) {
171  TH1D * temphist = fErrorBandXY->ProjectionY("temphist", ix, ix);
172 
173  int nprobSum = 1;
174  double q[1];
175  double probSum[1];
176  probSum[0] = level;
177 
178  temphist->GetQuantiles(nprobSum, q, probSum);
179 
180  errorband[ix - 1] = q[0];
181  }
182 
183  return errorband;
184 }
TGraph * BCFitter::GetErrorBandGraph ( double  level1,
double  level2 
) const

Definition at line 187 of file BCFitter.cxx.

188 {
189  if (!fErrorBandXY)
190  return 0;
191 
192  // define new graph
193  int nx = fErrorBandXY->GetNbinsX();
194 
195  TGraph * graph = new TGraph(2 * nx);
196  graph->SetFillStyle(1001);
197  graph->SetFillColor(kYellow);
198 
199  // get error bands
200  std::vector<double> ymin = GetErrorBand(level1);
201  std::vector<double> ymax = GetErrorBand(level2);
202 
203  for (int i = 0; i < nx; i++) {
204  graph->SetPoint(i, fErrorBandXY->GetXaxis()->GetBinCenter(i + 1), ymin[i]);
205  graph->SetPoint(nx + i, fErrorBandXY->GetXaxis()->GetBinCenter(nx - i), ymax[nx - i - 1]);
206  }
207 
208  return graph;
209 }
TH2D* BCFitter::GetErrorBandXY ( ) const
inline

const BCParameter * GetParameter(const char * name);

Returns
The 2-d histogram of the error band.

Definition at line 60 of file BCFitter.h.

61  { return fErrorBandXY; }
TH2D * BCFitter::GetErrorBandXY_yellow ( double  level = .68,
int  nsmooth = 0 
) const

Definition at line 212 of file BCFitter.cxx.

213 {
214  if (!fErrorBandXY)
215  return 0;
216 
217  int nx = fErrorBandXY->GetNbinsX();
218  int ny = fErrorBandXY->GetNbinsY();
219 
220  // copy existing histogram
221  TH2D * hist_tempxy = (TH2D*) fErrorBandXY->Clone(
222  TString::Format("%s_sub_%f.2", fErrorBandXY->GetName(), level));
223  hist_tempxy->Reset();
224  hist_tempxy->SetFillColor(kYellow);
225 
226  // loop over x bins
227  for (int ix = 1; ix < nx; ix++) {
228  BCH1D * hist_temp = new BCH1D();
229 
230  TH1D * hproj = fErrorBandXY->ProjectionY("temphist", ix, ix);
231  if (nsmooth > 0)
232  hproj->Smooth(nsmooth);
233 
234  hist_temp->SetHistogram(hproj);
235 
236  TH1D * hist_temp_yellow = hist_temp->GetSmallestIntervalHistogram(level);
237 
238  for (int iy = 1; iy <= ny; ++iy)
239  hist_tempxy->SetBinContent(ix, iy, hist_temp_yellow->GetBinContent(iy));
240 
241  delete hist_temp_yellow;
242  delete hist_temp;
243  }
244 
245  return hist_tempxy;
246 }
TGraph * BCFitter::GetFitFunctionGraph ( const std::vector< double > &  parameters)

Definition at line 249 of file BCFitter.cxx.

250 {
251  if (!fErrorBandXY)
252  return 0;
253 
254  // define new graph
255  int nx = fErrorBandXY->GetNbinsX();
256  TGraph * graph = new TGraph(nx);
257 
258  // loop over x values
259  for (int i = 0; i < nx; i++) {
260  double x = fErrorBandXY->GetXaxis()->GetBinCenter(i + 1);
261 
262  std::vector<double> xvec;
263  xvec.push_back(x);
264  double y = FitFunction(xvec, parameters);
265  xvec.clear();
266 
267  graph->SetPoint(i, x, y);
268  }
269 
270  return graph;
271 }
TGraph* BCFitter::GetFitFunctionGraph ( )
inline

Definition at line 75 of file BCFitter.h.

TGraph * BCFitter::GetFitFunctionGraph ( const std::vector< double > &  parameters,
double  xmin,
double  xmax,
int  n = 1000 
)

Definition at line 274 of file BCFitter.cxx.

275 {
276  // define new graph
277  TGraph * graph = new TGraph(n + 1);
278 
279  double dx = (xmax - xmin) / (double) n;
280 
281  // loop over x values
282  for (int i = 0; i <= n; i++) {
283  double x = (double) i * dx + xmin;
284  std::vector<double> xvec;
285  xvec.push_back(x);
286  double y = FitFunction(xvec, parameters);
287 
288  xvec.clear();
289 
290  graph->SetPoint(i, x, y);
291  }
292 
293  return graph;
294 }
bool BCFitter::GetFixedDataAxis ( unsigned int  index) const

Definition at line 340 of file BCFitter.cxx.

341 {
342  // check if index is within range
343  if (index > fDataSet->GetDataPoint(0)->GetNValues()) {
344  BCLog::OutError("BCFitter::GetFixedDataAxis : Index out of range.");
345  return false;
346  }
347 
348  return fDataFixedValues[index];
349 }
void BCFitter::MarginalizePostprocess ( )
inlinevirtual

Overloaded from BCIntegrate.

Reimplemented from BCIntegrate.

Definition at line 165 of file BCFitter.h.

166  {;};
void BCFitter::MarginalizePreprocess ( )
virtual

Overloaded from BCIntegrate.

Reimplemented from BCIntegrate.

Definition at line 69 of file BCFitter.cxx.

70 {
71  // prepare function fitting
72  double dx = 0.;
73  double dy = 0.;
74 
75  if (fFitFunctionIndexX >= 0) {
78  / (double) fErrorBandNbinsX;
79 
82  / (double) fErrorBandNbinsY;
83 
85  = new TH2D(TString::Format("errorbandxy_%d", BCLog::GetHIndex()), "",
92  fErrorBandXY->SetStats(kFALSE);
93 
94  for (unsigned ix = 1; ix <= fErrorBandNbinsX; ++ix)
95  for (unsigned iy = 1; iy <= fErrorBandNbinsX; ++iy)
96  fErrorBandXY->SetBinContent(ix, iy, 0.);
97  }
98 
99 }
void BCFitter::MCMCIterationInterface ( )
virtual

Overloaded from BCEngineMCMC

Reimplemented from BCEngineMCMC.

Definition at line 58 of file BCFitter.cxx.

59 {
60  // call base interface
62 
63  // fill error band
65  FillErrorBand();
66 }
int BCFitter::ReadErrorBandFromFile ( const char *  file)

Read

Definition at line 297 of file BCFitter.cxx.

298 {
299  TFile * froot = TFile::Open(file);
300  if (!froot->IsOpen()) {
301  BCLog::OutError(Form("BCFitter::ReadErrorBandFromFile. Couldn't open file %s.", file));
302  return 0;
303  }
304 
305  int r = 0;
306 
307  TH2D * h2 = (TH2D*) froot->Get("errorbandxy");
308  if (h2) {
309  h2->SetDirectory(0);
310  h2->SetName(TString::Format("errorbandxy_%d", BCLog::GetHIndex()));
311  SetErrorBandHisto(h2);
312  r = 1;
313  }
314  else
316  Form("BCFitter::ReadErrorBandFromFile : Couldn't read histogram \"errorbandxy\" from file %s.",file));
317 
318  froot->Close();
319 
320  return r;
321 }
void BCFitter::SetErrorBandContinuous ( bool  flag)

Sets the error band flag to continuous function

Definition at line 352 of file BCFitter.cxx.

353 {
354  fErrorBandContinuous = flag;
355 
356  if (flag)
357  return;
358 
359  // clear x-values
360  fErrorBandX.clear();
361 
362  // copy data x-values
363  for (unsigned int i = 0; i < fDataSet->GetNDataPoints(); ++i)
365 }
void BCFitter::SetErrorBandHisto ( TH2D *  h)
inline

Sets errorband histogram

Definition at line 97 of file BCFitter.h.

98  { fErrorBandXY = h; }
void BCFitter::SetFillErrorBand ( bool  flag = true)
inline

Turn on or off the filling of the error band during the MCMC run.

Parameters
flagset to true for turning on the filling

Definition at line 92 of file BCFitter.h.

93  { fFlagFillErrorBand=flag; }
void BCFitter::SetFitFunctionIndexX ( int  index)
inline

Sets index of the x values in function fits.

Parameters
indexIndex of the x values

Definition at line 109 of file BCFitter.h.

110  { fFitFunctionIndexX = index; }
void BCFitter::SetFitFunctionIndexY ( int  index)
inline

Sets index of the y values in function fits.

Parameters
indexIndex of the y values

Definition at line 115 of file BCFitter.h.

116  { fFitFunctionIndexY = index; }
void BCFitter::SetFitFunctionIndices ( int  indexx,
int  indexy 
)
inline

Sets indices of the x and y values in function fits.

Parameters
indexxIndex of the x values
indexyIndex of the y values

Definition at line 122 of file BCFitter.h.

123  { SetFitFunctionIndexX(indexx);
124  SetFitFunctionIndexY(indexy); }
void BCFitter::UnsetFillErrorBand ( )
inline

Turn off filling of the error band during the MCMC run. This method is equivalent to SetFillErrorBand(false)

Definition at line 103 of file BCFitter.h.

104  { fFlagFillErrorBand=false; }

Member Data Documentation

TGraph* BCFitter::fErrorBand
private

Pointer to the error band (for legend)

Definition at line 178 of file BCFitter.h.

bool BCFitter::fErrorBandContinuous
private

A flag for single point evaluation of the error "band"

Definition at line 191 of file BCFitter.h.

unsigned BCFitter::fErrorBandNbinsX
private

Number of X bins of the error band histogram

Definition at line 196 of file BCFitter.h.

unsigned BCFitter::fErrorBandNbinsY
private

Number of Y bins of the error band histogram

Definition at line 200 of file BCFitter.h.

std::vector<double> BCFitter::fErrorBandX
private

Definition at line 192 of file BCFitter.h.

TH2D* BCFitter::fErrorBandXY
protected

The error band histogram

Definition at line 206 of file BCFitter.h.

int BCFitter::fFitFunctionIndexX
private

The indices for function fits

Definition at line 186 of file BCFitter.h.

int BCFitter::fFitFunctionIndexY
private

Definition at line 187 of file BCFitter.h.

bool BCFitter::fFlagFillErrorBand
private

Flag whether or not to fill the error band

Definition at line 182 of file BCFitter.h.


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