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

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

#include <BCGraphFitter.h>

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

Public Member Functions

Constructors and destructors
 BCGraphFitter ()
 
 BCGraphFitter (const char *name)
 
 BCGraphFitter (TGraphErrors *graph, TF1 *func)
 
 BCGraphFitter (const char *name, TGraphErrors *graph, TF1 *func)
 
 ~BCGraphFitter ()
 
Member functions (get)
TGraphErrors * GetGraph ()
 
TF1 * GetFitFunction ()
 
TGraph * GetErrorBand ()
 
TGraph * GetGraphFitFunction ()
 
Member functions (set)
int SetGraph (TGraphErrors *graph)
 
int SetFitFunction (TF1 *func)
 
Member functions (miscellaneous methods)
double LogLikelihood (const std::vector< double > &parameters)
 
double FitFunction (const std::vector< double > &x, const std::vector< double > &parameters)
 
int Fit ()
 
int Fit (TGraphErrors *graph, TF1 *func)
 
void DrawFit (const char *options="", bool flaglegend=false)
 
virtual 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

TGraphErrors * fGraph
 
TF1 * fFitFunction
 
TGraph * fErrorBand
 
TGraph * fGraphFitFunction
 

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 graphs with functions.

Author
Daniel Kollar
Kevin Kröninger
Version
1.0
Date
2008 This class allows fitting of a TGraphErrors using a TF1 function. It doeasn't take the x uncertainties into account. For that look at BCGraphXFitter (not yet implemented).

Definition at line 35 of file BCGraphFitter.h.

Constructor & Destructor Documentation

BCGraphFitter::BCGraphFitter ( )

Default constructor

Definition at line 27 of file BCGraphFitter.cxx.

28  : BCFitter()
29  , fGraph(0)
30  , fFitFunction(0)
31  , fErrorBand(0)
33 {
34  // set MCMC for marginalization
36 }
BCGraphFitter::BCGraphFitter ( const char *  name)

Constructor

Parameters
namename of the model

Definition at line 40 of file BCGraphFitter.cxx.

41  : BCFitter(name)
42  , fGraph(0)
43  , fFitFunction(0)
44  , fErrorBand(0)
46 {
47  // set MCMC for marginalization
49 }
BCGraphFitter::BCGraphFitter ( TGraphErrors *  graph,
TF1 *  func 
)

Constructor

Parameters
graphpointer to TGraphErrors
funcpointer to TF1

Definition at line 53 of file BCGraphFitter.cxx.

54  : BCFitter()
55  , fGraph(0)
56  , fFitFunction(0)
57  , fErrorBand(0)
59 {
60  SetGraph(graph);
61  SetFitFunction(func);
62 
63  // set MCMC for marginalization
65 }
BCGraphFitter::BCGraphFitter ( const char *  name,
TGraphErrors *  graph,
TF1 *  func 
)

Constructor

Parameters
namename of the model
graphpointer to TGraphErrors
funcpointer to TF1

Definition at line 69 of file BCGraphFitter.cxx.

70  : BCFitter(name)
71  , fGraph(0)
72  , fFitFunction(0)
73  , fErrorBand(0)
75 {
76  SetGraph(graph);
77  SetFitFunction(func);
78 
79  // set MCMC for marginalization
81 }
BCGraphFitter::~BCGraphFitter ( )

The default destructor.

Definition at line 206 of file BCGraphFitter.cxx.

207 {
208  // todo memory leak
209 }

Member Function Documentation

double BCGraphFitter::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 383 of file BCGraphFitter.cxx.

383  {
384 
385  //format: x y error_x error_y
386  std::vector<double> values = fDataSet->GetDataPoint(index)->GetValues();
387 
388  if (values.at(2))
389  BCLog::OutWarning("BCGraphFitter::CDF: Non-zero errors in x-direction are ignored!");
390 
391  // get the observed value
392  double yObs = values.at(1);
393 
394  // expectation value
395  double yExp = FitFunction(values, parameters);
396 
397  return ROOT::Math::normal_cdf(yObs, values.at(3), yExp);
398 }
void BCGraphFitter::DrawFit ( const char *  options = "",
bool  flaglegend = false 
)
virtual

Draw the fit in the current pad.

Implements BCFitter.

Definition at line 330 of file BCGraphFitter.cxx.

331 {
332  if (!fGraph)
333  {
334  BCLog::Out(BCLog::error, BCLog::error,"BCGraphFitter::DrawFit() : TGraphErrors not defined.");
335  return;
336  }
337 
338  if (!fFitFunction)
339  {
340  BCLog::Out(BCLog::error, BCLog::error,"BCGraphFitter::DrawFit() : Fit function not defined.");
341  return;
342  }
343 
344  // check wheather options contain "same"
345  TString opt = options;
346  opt.ToLower();
347 
348  // if not same, draw the histogram first to get the axes
349  if(!opt.Contains("same"))
350  fGraph->Draw("ap");
351 
352  // draw the error band as central 68% probability interval
353  fErrorBand = GetErrorBandGraph(0.16, 0.84);
354  fErrorBand->Draw("f same");
355 
356  // draw the fit function on top
358  fGraphFitFunction->SetLineColor(kRed);
359  fGraphFitFunction->SetLineWidth(2);
360  fGraphFitFunction->Draw("l same");
361 
362  // now draw the histogram again since it was covered by the band and
363  // the best fit
364  fGraph->Draw("p same");
365 
366  // draw legend
367  if (flaglegend)
368  {
369  TLegend * legend = new TLegend(0.25, 0.75, 0.55, 0.95);
370  legend->SetBorderSize(0);
371  legend->SetFillColor(kWhite);
372  legend->AddEntry(fGraph, "Data", "P");
373  legend->AddEntry(fGraphFitFunction, "Best fit", "L");
374  legend->AddEntry(fErrorBand, "Error band", "F");
375  legend->Draw();
376  }
377 
378  gPad->RedrawAxis();
379 }
int BCGraphFitter::Fit ( )
virtual

Performs the fit. The graph and the function has to be set beforehand.

Returns
An error code.

Implements BCFitter.

Definition at line 281 of file BCGraphFitter.cxx.

282 {
283  // set graph
284  if (!fGraph) {
285  BCLog::OutError("BCEfficiencyFitter::Fit : Graph not defined.");
286  return 0;
287  }
288 
289  // set function
290  if (!fFitFunction) {
291  BCLog::OutError("BCEfficiencyFitter::Fit : Fit function not defined.");
292  return 0;
293  }
294 
295  // check setup
297  Form("Fitting %d data points with function of %d parameters",GetNDataPoints(),GetNParameters()));
298  if(GetNDataPoints() <= GetNParameters())
299  {
301  Form("Number of parameters (%d) lower than or equal to number of points (%d)."
303  BCLog::Out(BCLog::warning,BCLog::warning,"Fit doesn't have much meaning.");
304  }
305 
306  // perform marginalization
307  MarginalizeAll();
308 
309  // maximize posterior probability, using the best-fit values close
310  // to the global maximum from the MCMC
314  SetOptimizationMethod(method_temp);
315 
316  // calculate p-value from the chi2 probability
317  // this is only valid for a product of gaussiang which is the case for
318  // the BCGraphFitter
321 
322  // print summary to screen
324 
325  return 1;
326 }
int BCGraphFitter::Fit ( TGraphErrors *  graph,
TF1 *  func 
)

Performs the fit of the graph with the function.

Parameters
graphpointer to TGraphErrors object
funcpointer to TF1 object
Returns
An error code.

Definition at line 262 of file BCGraphFitter.cxx.

263 {
264  // set graph
265  if (!SetGraph(graph)) {
266  BCLog::OutError("BCEfficiencyFitter::Fit : Graph not defined.");
267  return 0;
268  }
269 
270  // set function
271  if (!SetFitFunction(func)) {
272  return 0;
273  BCLog::OutError("BCEfficiencyFitter::Fit : Fit function not defined.");
274  }
275 
276  return Fit();
277 }
double BCGraphFitter::FitFunction ( const std::vector< double > &  x,
const std::vector< double > &  parameters 
)
virtual

Returns the value of the 1D fit function for a given set of parameters at a given x.

Parameters
xpoint to calculate the function value at
parametersparameters of the function

Reimplemented from BCFitter.

Definition at line 247 of file BCGraphFitter.cxx.

248 {
249  // set the parameters of the function
250  // passing the pointer to first element of the vector is
251  // not completely safe as there might be an implementation where
252  // the vector elements are not stored consecutively in memory.
253  // however it is much faster than copying the contents, especially
254  // for large number of parameters
255  fFitFunction->SetParameters(&params[0]);
256 
257  return fFitFunction->Eval(x[0]);
258 }
TGraph* BCGraphFitter::GetErrorBand ( )
inline
Returns
pointer to the error band

Definition at line 85 of file BCGraphFitter.h.

86  { return fErrorBand; };
TF1* BCGraphFitter::GetFitFunction ( )
inline
Returns
pointer to TF1

Definition at line 80 of file BCGraphFitter.h.

81  { return fFitFunction; };
TGraphErrors* BCGraphFitter::GetGraph ( )
inline
Returns
pointer to TGraphErrors

Definition at line 75 of file BCGraphFitter.h.

76  { return fGraph; };
TGraph* BCGraphFitter::GetGraphFitFunction ( )
inline
Returns
pointer to a graph for the fit function

Definition at line 90 of file BCGraphFitter.h.

91  { return fGraphFitFunction; };
double BCGraphFitter::LogLikelihood ( const std::vector< double > &  parameters)
virtual

The log of the prior probability. It is set to be flat in all parameters.

Parameters
parametersvector containing the parameter values The log of the conditional probability.
parametersvector containing the parameter values

Implements BCModel.

Definition at line 213 of file BCGraphFitter.cxx.

214 {
215  // initialize probability
216  double logl = 0.;
217 
218  // set the parameters of the function
219  // passing the pointer to first element of the vector is
220  // not completely safe as there might be an implementation where
221  // the vector elements are not stored consecutively in memory.
222  // however it is much faster than copying the contents, especially
223  // for large number of parameters
224  fFitFunction->SetParameters(&params[0]);
225 
226  // loop over all data points
227  for (unsigned i = 0; i < GetNDataPoints(); i++)
228  {
229  std::vector<double> x = GetDataPoint(i)->GetValues();
230 
231  // her we ignore the errors on x even when they're available
232  // i.e. we treat them just as the region specifiers
233  double y = x[1];
234  double yerr = x[3];
235  double yexp = fFitFunction->Eval(x[0]);
236 
237  // calculate log of probability assuming
238  // a Gaussian distribution for each point
239  logl += BCMath::LogGaus(y, yexp, yerr, true);
240  }
241 
242  return logl;
243 }
int BCGraphFitter::SetFitFunction ( TF1 *  func)
Parameters
funcpointer to TF1 object

Definition at line 162 of file BCGraphFitter.cxx.

163 {
164  if(!func)
165  {
166  BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetFitFunction() : TF1 not created.");
167  return 0;
168  }
169 
170  // get the new number of parameters
171  int npar = func->GetNpar();
172  if(!npar)
173  {
174  BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetFitFunction() : TF1 has zero parameters. Not able to fit.");
175  return 0;
176  }
177 
178  // set the function
179  fFitFunction = func;
180 
181  // update the model name to contain the function name
182  if(fName=="model")
183  SetName(TString::Format("GraphFitter with %s",fFitFunction->GetName()));
184 
185  // reset parameters
186  ClearParameters(true);
187 
188  // add parameters
189  for (int i = 0; i < npar; ++i)
190  {
191  double xmin;
192  double xmax;
193  fFitFunction->GetParLimits(i, xmin, xmax);
194 
195  AddParameter(fFitFunction->GetParName(i), xmin, xmax);
196  }
197 
198  // set flat prior for all parameters by default
200 
201  return GetNParameters();
202 }
int BCGraphFitter::SetGraph ( TGraphErrors *  graph)
Parameters
graphpointer to TGraphErrors object

Definition at line 85 of file BCGraphFitter.cxx.

86 {
87  if(!graph) {
88  BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors not created.");
89  return 0;
90  }
91 
92  int npoints = graph->GetN();
93  if(!npoints)
94  {
95  BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors is empty.");
96  return 0;
97  }
98  else if(npoints==1)
99  {
100  BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors has only one point. Not able to fit.");
101  return 0;
102  }
103 
104  fGraph = graph;
105 
106  double * x = fGraph->GetX();
107  double * y = fGraph->GetY();
108  double * ex = fGraph->GetEX();
109  double * ey = fGraph->GetEY();
110 
111  if(!ey)
112  {
113  BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors has NO errors set on Y. Not able to fit.");
114  return 0;
115  }
116 
117  BCDataSet * ds = new BCDataSet();
118 
119  // fill the dataset
120  // find x and y boundaries for the error band calculation
121  double xmin=x[0];
122  double xmax=x[0];
123  double ymin=y[0];
124  double ymax=y[0];
125  for (int i = 0; i < npoints; ++i)
126  {
127  // if x errors are not set, set them to zero
128  double errx = ex ? ex[i] : 0.;
129 
130  // create the data point
131  BCDataPoint * dp = new BCDataPoint(4);
132  dp->SetValue(0, x[i]);
133  dp->SetValue(1, y[i]);
134  dp->SetValue(2, errx);
135  dp->SetValue(3, ey[i]);
136  ds->AddDataPoint(dp);
137 
138  if(x[i]-errx < xmin)
139  xmin = x[i]-errx;
140  else if(x[i]+errx > xmax)
141  xmax = x[i]+errx;
142 
143  if(y[i] - 5.*ey[i] < ymin)
144  ymin = y[i] - 5.*ey[i];
145  else if(y[i] + 5.*ey[i] > ymax)
146  ymax = y[i] + 5.*ey[i];
147  }
148 
149  SetDataSet(ds);
150 
151  // set boundaries for the error band calculation
152  SetDataBoundaries(0, xmin, xmax);
153  SetDataBoundaries(1, ymin, ymax);
154 
155  SetFitFunctionIndices(0, 1);
156 
157  return GetNDataPoints();
158 }

Member Data Documentation

TGraph* BCGraphFitter::fErrorBand
private

Pointer to the error band (for legend)

Definition at line 159 of file BCGraphFitter.h.

TF1* BCGraphFitter::fFitFunction
private

The fit function

Definition at line 155 of file BCGraphFitter.h.

TGraphErrors* BCGraphFitter::fGraph
private

The graph containing the data.

Definition at line 151 of file BCGraphFitter.h.

TGraph* BCGraphFitter::fGraphFitFunction
private

Pointer to a graph for displaying the fit function

Definition at line 163 of file BCGraphFitter.h.


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