#include <BCGraphFitter.h>
Definition at line 34 of file BCGraphFitter.h.
Public Member Functions | |
Constructors and destructors | |
BCGraphFitter (TGraphErrors *graph, TF1 *func) | |
BCGraphFitter () | |
~BCGraphFitter () | |
Member functions (miscellaneous methods) | |
void | DrawFit (const char *options="", bool flaglegend=false) |
int | Fit (TGraphErrors *graph, TF1 *func) |
int | Fit () |
double | FitFunction (std::vector< double > x, std::vector< double > parameters) |
double | LogAPrioriProbability (std::vector< double > parameters) |
double | LogLikelihood (std::vector< double > parameters) |
Member functions (get) | |
TGraph * | GetErrorBand () |
TF1 * | GetFitFunction () |
TGraphErrors * | GetGraph () |
TGraph * | GetGraphFitFunction () |
Member functions (set) | |
int | SetFitFunction (TF1 *func) |
int | SetGraph (TGraphErrors *graph) |
Private Attributes | |
TGraph * | fErrorBand |
TF1 * | fFitFunction |
TGraphErrors * | fGraph |
TGraph * | fGraphFitFunction |
BCGraphFitter::BCGraphFitter | ( | ) |
Default constructor
Definition at line 25 of file BCGraphFitter.cxx.
00025 : BCModel("GraphFitter") 00026 { 00027 fGraph = 0; 00028 fFitFunction = 0; 00029 fErrorBand = 0; 00030 fGraphFitFunction = 0; 00031 00032 this -> MCMCSetNIterationsRun(2000); 00033 00034 this -> SetFillErrorBand(); 00035 }
BCGraphFitter::BCGraphFitter | ( | TGraphErrors * | graph, | |
TF1 * | func | |||
) |
Constructor
graph | pointer to TGraphErrors | |
func | pointer to TF1 |
Definition at line 39 of file BCGraphFitter.cxx.
00039 : BCModel("GraphFitter") 00040 { 00041 fGraph = 0; 00042 fFitFunction = 0; 00043 fErrorBand = 0; 00044 fGraphFitFunction = 0; 00045 00046 this -> MCMCSetNIterationsRun(2000); 00047 00048 this -> SetGraph(graph); 00049 this -> SetFitFunction(func); 00050 00051 this -> SetFillErrorBand(); 00052 }
BCGraphFitter::~BCGraphFitter | ( | ) |
void BCGraphFitter::DrawFit | ( | const char * | options = "" , |
|
bool | flaglegend = false | |||
) |
Draw the fit in the current pad.
Definition at line 270 of file BCGraphFitter.cxx.
00271 { 00272 if (!fGraph) 00273 { 00274 BCLog::Out(BCLog::error, BCLog::error,"BCGraphFitter::DrawFit() : TGraphErrors not defined."); 00275 return; 00276 } 00277 00278 if (!fFitFunction) 00279 { 00280 BCLog::Out(BCLog::error, BCLog::error,"BCGraphFitter::DrawFit() : Fit function not defined."); 00281 return; 00282 } 00283 00284 // check wheather options contain "same" 00285 TString opt = options; 00286 opt.ToLower(); 00287 00288 // if not same, draw the histogram first to get the axes 00289 if(!opt.Contains("same")) 00290 fGraph -> Draw("ap"); 00291 00292 // draw the error band as central 68% probability interval 00293 fErrorBand = this -> GetErrorBandGraph(0.16, 0.84); 00294 fErrorBand -> Draw("f same"); 00295 00296 // draw the fit function on top 00297 fGraphFitFunction = this -> GetFitFunctionGraph( this->GetBestFitParameters() ); 00298 fGraphFitFunction -> SetLineColor(kRed); 00299 fGraphFitFunction -> SetLineWidth(2); 00300 fGraphFitFunction -> Draw("l same"); 00301 00302 // now draw the histogram again since it was covered by the band and 00303 // the best fit 00304 fGraph -> Draw("p same"); 00305 00306 // draw legend 00307 if (flaglegend) 00308 { 00309 TLegend * legend = new TLegend(0.25, 0.75, 0.55, 0.95); 00310 legend -> SetBorderSize(0); 00311 legend -> SetFillColor(kWhite); 00312 legend -> AddEntry(fGraph, "Data", "P"); 00313 legend -> AddEntry(fGraphFitFunction, "Best fit", "L"); 00314 legend -> AddEntry(fErrorBand, "Error band", "F"); 00315 legend -> Draw(); 00316 } 00317 00318 gPad -> RedrawAxis(); 00319 }
int BCGraphFitter::Fit | ( | TGraphErrors * | graph, | |
TF1 * | func | |||
) |
Performs the fit of the graph with the function.
graph | pointer to TGraphErrors object | |
func | pointer to TF1 object |
Definition at line 240 of file BCGraphFitter.cxx.
00241 { 00242 // set graph 00243 if (!this -> SetGraph(graph)) 00244 return 0; 00245 00246 // set function 00247 if (!this -> SetFitFunction(func)) 00248 return 0; 00249 00250 // perform marginalization 00251 this -> MarginalizeAll(); 00252 00253 // maximize posterior probability, using the best-fit values close 00254 // to the global maximum from the MCMC 00255 this -> FindModeMinuit( this -> GetBestFitParameters(), -1); 00256 00257 // calculate p-value from the chi2 probability 00258 // this is only valid for a product of gaussiang which is the case for 00259 // the BCGraphFitter 00260 this -> GetPvalueFromChi2(this -> GetBestFitParameters(), 3); 00261 00262 // print summary to screen 00263 this -> PrintShortFitSummary(); 00264 00265 return 1; 00266 }
int BCGraphFitter::Fit | ( | ) | [inline] |
Performs the fit. The graph and the function has to beset beforehand.
Definition at line 117 of file BCGraphFitter.h.
00118 { return this -> Fit(fGraph, fFitFunction); };
double BCGraphFitter::FitFunction | ( | std::vector< double > | x, | |
std::vector< double > | parameters | |||
) | [virtual] |
Returns the value of the 1D fit function for a given set of parameters at a given x.
x | point to calculate the function value at | |
parameters | parameters of the function |
Reimplemented from BCIntegrate.
Definition at line 230 of file BCGraphFitter.cxx.
00231 { 00232 // set the parameters of the function 00233 fFitFunction -> SetParameters(¶ms[0]); 00234 00235 return fFitFunction -> Eval(x[0]); 00236 }
TGraph* BCGraphFitter::GetErrorBand | ( | ) | [inline] |
Definition at line 72 of file BCGraphFitter.h.
00073 { return fErrorBand; };
TF1* BCGraphFitter::GetFitFunction | ( | ) | [inline] |
Definition at line 67 of file BCGraphFitter.h.
00068 { return fFitFunction; };
TGraphErrors* BCGraphFitter::GetGraph | ( | ) | [inline] |
Definition at line 62 of file BCGraphFitter.h.
00063 { return fGraph; };
TGraph* BCGraphFitter::GetGraphFitFunction | ( | ) | [inline] |
Definition at line 77 of file BCGraphFitter.h.
00078 { return fGraphFitFunction; };
double BCGraphFitter::LogAPrioriProbability | ( | std::vector< double > | parameters | ) | [virtual] |
The log of the prior probability. It is set to be flat in all parameters.
parameters | vector containing the parameter values |
Reimplemented from BCModel.
Definition at line 184 of file BCGraphFitter.cxx.
00185 { 00186 // using flat probability in all parameters 00187 double logprob = 0.; 00188 for(unsigned int i=0; i < this -> GetNParameters(); i++) 00189 logprob -= log(this -> GetParameter(i) -> GetRangeWidth()); 00190 00191 return logprob; 00192 }
double BCGraphFitter::LogLikelihood | ( | std::vector< double > | parameters | ) | [virtual] |
The log of the conditional probability.
parameters | vector containing the parameter values |
Reimplemented from BCModel.
Definition at line 196 of file BCGraphFitter.cxx.
00197 { 00198 // initialize probability 00199 double logl = 0.; 00200 00201 // set the parameters of the function 00202 // passing the pointer to first element of the vector is 00203 // not completely safe as there might be an implementation where 00204 // the vector elements are not stored consecutively in memory. 00205 // however it is much faster than copying the contents, especially 00206 // for large number of parameters 00207 fFitFunction -> SetParameters(¶ms[0]); 00208 00209 // loop over all data points 00210 for (int i = 0; i < this -> GetNDataPoints(); i++) 00211 { 00212 std::vector <double> x = GetDataPoint(i) -> GetValues(); 00213 00214 // her we ignore the errors on x even when they're available 00215 // i.e. we treat them just as the region specifiers 00216 double y = x[1]; 00217 double yerr = x[3]; 00218 double yexp = this -> FitFunction(x,params); 00219 00220 // calculate log of probability assuming 00221 // a Gaussian distribution for each point 00222 logl += BCMath::LogGaus(y, yexp, yerr, true); 00223 } 00224 00225 return logl; 00226 }
int BCGraphFitter::SetFitFunction | ( | TF1 * | func | ) |
func | pointer to TF1 object |
Definition at line 139 of file BCGraphFitter.cxx.
00140 { 00141 if(!func) 00142 { 00143 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetFitFunction() : TF1 not created."); 00144 return 0; 00145 } 00146 00147 // get the new number of parameters 00148 int npar = func -> GetNpar(); 00149 if(!npar) 00150 { 00151 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetFitFunction() : TF1 has zero parameters. Not able to fit."); 00152 return 0; 00153 } 00154 00155 // set the function 00156 fFitFunction = func; 00157 00158 // update the model name to contain the function name 00159 this -> SetName(TString::Format("GraphFitter with %s",fFitFunction->GetName())); 00160 00161 // reset parameters 00162 fParameterSet -> clear(); 00163 00164 // add parameters 00165 for (int i = 0; i < npar; ++i) 00166 { 00167 double xmin; 00168 double xmax; 00169 fFitFunction -> GetParLimits(i, xmin, xmax); 00170 00171 this -> AddParameter(fFitFunction->GetParName(i), xmin, xmax); 00172 } 00173 00174 return this -> GetNParameters(); 00175 }
int BCGraphFitter::SetGraph | ( | TGraphErrors * | graph | ) |
graph | pointer to TGraphErrors object |
Definition at line 56 of file BCGraphFitter.cxx.
00057 { 00058 if(!graph) 00059 { 00060 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors not created."); 00061 return 0; 00062 } 00063 00064 int npoints = graph -> GetN(); 00065 if(!npoints) 00066 { 00067 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors is empty."); 00068 return 0; 00069 } 00070 else if(npoints==1) 00071 { 00072 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors has only one point. Not able to fit."); 00073 return 0; 00074 } 00075 else if(npoints==2) 00076 { 00077 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors has only two points. Fit has no meaning."); 00078 return 0; 00079 } 00080 00081 fGraph = graph; 00082 00083 double * x = fGraph -> GetX(); 00084 double * y = fGraph -> GetY(); 00085 double * ex = fGraph -> GetEX(); 00086 double * ey = fGraph -> GetEY(); 00087 00088 if(!ey) 00089 { 00090 BCLog::Out(BCLog::error,BCLog::error,"BCGraphFitter::SetGraph() : TGraphErrors has no errors set on Y. Not able to fit."); 00091 return 0; 00092 } 00093 00094 BCDataSet * ds = new BCDataSet(); 00095 00096 // fill the dataset 00097 // find x and y boundaries for the error band calculation 00098 double xmin=x[0]; 00099 double xmax=x[0]; 00100 double ymin=y[0]; 00101 double ymax=y[0]; 00102 for (int i = 0; i < npoints; ++i) 00103 { 00104 // if x errors are not set, set them to zero 00105 double errx = ex ? ex[i] : 0.; 00106 00107 // create the data point 00108 BCDataPoint * dp = new BCDataPoint(4); 00109 dp -> SetValue(0, x[i]); 00110 dp -> SetValue(1, y[i]); 00111 dp -> SetValue(2, errx); 00112 dp -> SetValue(3, ey[i]); 00113 ds -> AddDataPoint(dp); 00114 00115 if(x[i]-errx < xmin) 00116 xmin = x[i]-errx; 00117 else if(x[i]+errx > xmax) 00118 xmax = x[i]+errx; 00119 00120 if(y[i] - 5.*ey[i] < ymin) 00121 ymin = y[i] - 5.*ey[i]; 00122 else if(y[i] + 5.*ey[i] > ymax) 00123 ymax = y[i] + 5.*ey[i]; 00124 } 00125 00126 this -> SetDataSet(ds); 00127 00128 // set boundaries for the error band calculation 00129 this -> SetDataBoundaries(0, xmin, xmax); 00130 this -> SetDataBoundaries(1, ymin, ymax); 00131 00132 this -> SetFitFunctionIndices(0, 1); 00133 00134 return this -> GetNDataPoints(); 00135 }
TGraph* BCGraphFitter::fErrorBand [private] |
Pointer to the error band (for legend)
Definition at line 145 of file BCGraphFitter.h.
TF1* BCGraphFitter::fFitFunction [private] |
The fit function
Definition at line 141 of file BCGraphFitter.h.
TGraphErrors* BCGraphFitter::fGraph [private] |
The graph containing the data.
Definition at line 137 of file BCGraphFitter.h.
TGraph* BCGraphFitter::fGraphFitFunction [private] |
Pointer to a graph for displaying the fit function
Definition at line 149 of file BCGraphFitter.h.