BCFitter.h
1 #ifndef __BCFITTER__H
2 #define __BCFITTER__H
3 
10 /*
11  * Copyright (C) 2007-2018, the BAT core developer team
12  * All rights reserved.
13  *
14  * For the licensing terms see doc/COPYING.
15  * For documentation see http://mpp.mpg.de/bat
16  */
17 
18 // ---------------------------------------------------------
19 
20 #include <BAT/BCAux.h>
21 #include <BAT/BCDataSet.h>
22 #include <BAT/BCModel.h>
23 
24 #include <TF1.h>
25 #include <TH2D.h>
26 
27 #include <string>
28 
29 // ---------------------------------------------------------
30 
31 class BCFitter : public BCModel
32 {
33 public:
34 
36  /* @{ */
37 
42  BCFitter(const TF1& f, const std::string& name = "fitter_model");
43 
46  virtual ~BCFitter() = 0;
47 
48  /* @} */
49 
51  /* @{ */
52 
56  {
57  return fFitFunction.at(GetCurrentChain());
58  }
59 
65  TH2* GetGraphicalErrorBandXY(double level = .68, int nsmooth = 0, bool overcoverage = true) const;
66 
69  const TH2& GetErrorBandXY()const
70  { return fErrorBandXY; };
71 
76  std::vector<double> GetErrorBand(double level) const;
77 
80  TGraph* GetErrorBandGraph(double level1, double level2) const;
81 
85  TGraph* GetFitFunctionGraph(const std::vector<double>& parameters);
86 
93 
98  TGraph* GetFitFunctionGraph(const std::vector<double>& parameters, double xmin, double xmax, int n = 1000);
99 
102  void FixDataAxis(unsigned int index, bool fixed);
103 
106  bool GetFixedDataAxis(unsigned int index) const;
107 
110  double GetPValue() const
111  {
112  return fPValue;
113  }
114 
115  /* @} */
116 
118  /* @{ */
119 
122  void SetErrorBandContinuous(bool flag);
123 
126  void SetErrorBandExtensionLowEdgeX(double extension)
127  {
128  fErrorBandExtensionLowEdgeX = extension;
129  }
130 
133  void SetErrorBandExtensionUpEdgeX(double extension)
134  {
135  fErrorBandExtensionUpEdgeX = extension;
136  }
137 
140  void SetErrorBandExtensionLowEdgeY(double extension)
141  {
142  fErrorBandExtensionLowEdgeY = extension;
143  }
144 
147  void SetErrorBandExtensionUpEdgeY(double extension)
148  {
149  fErrorBandExtensionUpEdgeY = extension;
150  }
154  void SetFillErrorBand(bool flag = true)
155  {
156  fFlagFillErrorBand = flag;
157  }
158 
163  {
164  fFlagFillErrorBand = false;
165  }
166 
170  void SetFitFunctionIndexX(int index)
171  {
172  fFitFunctionIndexX = index;
173  }
174 
178  void SetFitFunctionIndexY(int index)
179  {
180  fFitFunctionIndexY = index;
181  }
182 
187  void SetFitFunctionIndices(int indexx, int indexy)
188  {
189  SetFitFunctionIndexX(indexx);
190  SetFitFunctionIndexY(indexy);
191  }
192 
197  void SetFlagIntegration(bool flag)
198  {
199  fFlagIntegration = flag;
200  };
201 
207  virtual double FitFunction(const std::vector<double>& x, const std::vector<double>& parameters);
208 
215  double Integral(const std::vector<double>& parameters, double xmin, double xmax);
216 
217 
218  /* @} */
220  /* @{ */
223  virtual void Fit() = 0;
224 
227  virtual void DrawFit(const std::string& options, bool flaglegend = false) = 0;
228 
231  virtual void MCMCUserIterationInterface();
232 
235  virtual void MarginalizePreprocess();
236 
239  void FillErrorBand();
240 
243  void PrintShortFitSummary();
244 
245  /* @} */
246 
249  virtual void MCMCUserInitialize();
250 
251 private:
253  std::vector<TF1> fFitFunction;
254 
255 protected:
259 
266 
270 
273  std::vector<double> fErrorBandX;
274 
278 
282 
285  double fPValue;
286 
290 
293 
298 
305 
309 
313 
316  using BCModel::SetDataSet;
317  using BCModel::GetDataSet;
318 
320  virtual double GraphCorrection(unsigned /* ibin */) const
321  {
322  return 1.0;
323  }
324 };
325 
326 // ---------------------------------------------------------
327 
328 #endif
TH2D fErrorBandXY
The error band histogram.
Definition: BCFitter.h:289
void PrintShortFitSummary()
Prints a short summary of the fit results on the screen.
Definition: BCFitter.cxx:212
unsigned fErrorBandNbinsX
Number of X bins of the error band histogram.
Definition: BCFitter.h:277
double fErrorBandExtensionUpEdgeX
extends the upper edge of x range by the given value
Definition: BCFitter.h:304
double GetPValue() const
Definition: BCFitter.h:110
double fErrorBandExtensionUpEdgeY
extends the upper edge of y range by the given value
Definition: BCFitter.h:312
std::vector< double > GetErrorBand(double level) const
Returns a vector of y-values at a certain probability level.
Definition: BCFitter.cxx:223
void FixDataAxis(unsigned int index, bool fixed)
Toggle the data axis defined by index to be fixed.
Definition: BCFitter.cxx:340
TF1 & GetFitFunction()
Get fit function.
Definition: BCFitter.h:55
virtual void MarginalizePreprocess()
Overloaded from BCIntegrate.
Definition: BCFitter.cxx:89
TGraph * GetErrorBandGraph(double level1, double level2) const
Definition: BCFitter.cxx:249
TH2 * GetGraphicalErrorBandXY(double level=.68, int nsmooth=0, bool overcoverage=true) const
Definition: BCFitter.cxx:274
bool GetFixedDataAxis(unsigned int index) const
Definition: BCFitter.cxx:346
void UnsetFillErrorBand()
Turn off filling of the error band during the MCMC run.
Definition: BCFitter.h:162
void SetDataSet(BCDataSet *dataset)
Sets the data set.
Definition: BCModel.h:145
The base class for all user-defined models.
Definition: BCModel.h:39
const TH2 & GetErrorBandXY() const
Definition: BCFitter.h:69
int fFitFunctionIndexY
The index for function fits in y direction.
Definition: BCFitter.h:265
A class representing a set of data points.
Definition: BCDataSet.h:39
void SetFitFunctionIndices(int indexx, int indexy)
Sets indices of the x and y values in function fits.
Definition: BCFitter.h:187
virtual void DrawFit(const std::string &options, bool flaglegend=false)=0
Draw the fit in the current pad.
bool fFlagFillErrorBand
Flag whether or not to fill the error band.
Definition: BCFitter.h:258
void SetFlagIntegration(bool flag)
Sets the flag for integration.
Definition: BCFitter.h:197
virtual void MCMCUserIterationInterface()
Overloaded from BCEngineMCMC.
Definition: BCFitter.cxx:81
TGraph * GetFitFunctionGraph()
Definition: BCFitter.h:91
void FillErrorBand()
Fill error band histogram for current iteration.
Definition: BCFitter.cxx:131
double fErrorBandExtensionLowEdgeX
extends the lower edge of x range by the given value
Definition: BCFitter.h:301
void SetErrorBandExtensionLowEdgeY(double extension)
Extends the lower y Edge of th errorband by -extension.
Definition: BCFitter.h:140
bool fErrorBandContinuous
A flag for single point evaluation of the error "band".
Definition: BCFitter.h:269
void SetErrorBandContinuous(bool flag)
Sets the error band flag to continuous function.
Definition: BCFitter.cxx:352
void SetFitFunctionIndexY(int index)
Sets index of the y values in function fits.
Definition: BCFitter.h:178
double fErrorBandExtensionLowEdgeY
extends the upper edge of y range by the given value
Definition: BCFitter.h:308
int fFitFunctionIndexX
The index for function fits in x direction.
Definition: BCFitter.h:262
unsigned fErrorBandNbinsY
Number of Y bins of the error band histogram.
Definition: BCFitter.h:281
virtual double GraphCorrection(unsigned) const
Take care of bin width when creating a graph from the fit function.
Definition: BCFitter.h:320
void SetFitFunctionIndexX(int index)
Sets index of the x values in function fits.
Definition: BCFitter.h:170
virtual void Fit()=0
Fit the function&#39;s parameters to the data.
BCDataSet fFitterDataSet
Needed for uncertainty propagation.
Definition: BCFitter.h:292
void SetErrorBandExtensionUpEdgeY(double extension)
Extends the lower y Edge of th errorband by +extension.
Definition: BCFitter.h:147
double Integral(const std::vector< double > &parameters, double xmin, double xmax)
Compute the integral of the fit function between xmin and xmax.
Definition: BCFitter.cxx:196
void SetErrorBandExtensionUpEdgeX(double extension)
Extends the lower x Edge of th errorband by +extension.
Definition: BCFitter.h:133
virtual void MCMCUserInitialize()
Create enough TF1 copies for thread safety.
Definition: BCFitter.cxx:74
double fPValue
p value for goodness of fit
Definition: BCFitter.h:285
void SetFillErrorBand(bool flag=true)
Turn on or off the filling of the error band during the MCMC run.
Definition: BCFitter.h:154
void SetErrorBandExtensionLowEdgeX(double extension)
Extends the lower x Edge of th errorband by -extension.
Definition: BCFitter.h:126
BCFitter(const TF1 &f, const std::string &name="fitter_model")
Constructor.
Definition: BCFitter.cxx:28
std::vector< double > fErrorBandX
The x positions where the error is calculated.
Definition: BCFitter.h:273
virtual ~BCFitter()=0
The default destructor.
Definition: BCFitter.cxx:69
virtual double FitFunction(const std::vector< double > &x, const std::vector< double > &parameters)
Defines a fit function.
Definition: BCFitter.cxx:188
unsigned GetCurrentChain() const
A base class for all fitting classes.
Definition: BCFitter.h:31
BCDataSet * GetDataSet()
Definition: BCModel.h:78
virtual const std::vector< double > & GetBestFitParameters() const
bool fFlagIntegration
Flag for using the ROOT TH1::Integral method (true), or linear interpolation (false) ...
Definition: BCFitter.h:297