ROOT 6.13/01 Reference Guide |
#include "TFumili.h"
#include "Riostream.h"
#include "TGraphAsymmErrors.h"
#include "TF1.h"
#include "TF2.h"
#include "TF3.h"
#include "TH1.h"
#include "TMath.h"
#include "TROOT.h"
#include "TVirtualFitter.h"
Functions | |
void | GraphFitChisquareFumili (Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag) |
Minimization function for Graphs using a Chisquare method. More... | |
void | H1FitChisquareFumili (Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag) |
Minimization function for H1s using a Chisquare method. More... | |
void | H1FitLikelihoodFumili (Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag) |
Minimization function for H1s using a Likelihood method. More... | |
Variables | |
TFumili * | gFumili =0 |
static const Double_t | gMAXDOUBLE =1e300 |
static const Double_t | gMINDOUBLE =-1e300 |
void GraphFitChisquareFumili | ( | Int_t & | npar, |
Double_t * | gin, | ||
Double_t & | f, | ||
Double_t * | u, | ||
Int_t | flag | ||
) |
Minimization function for Graphs using a Chisquare method.
In case of a TGraphErrors object, ex, the error along x, is projected along the y-direction by calculating the function at the points x-exlow and x+exhigh.
The chisquare is computed as the sum of the quantity below at each point:
(y - f(x))**2 ----------------------------------- ey**2 + (0.5*(exl + exh)*f'(x))**2
where x and y are the point coordinates and f'(x) is the derivative of function f(x). This method to approximate the uncertainty in y because of the errors in x, is called "effective variance" method. The improvement, compared to the previously used method (f(x+ exhigh) - f(x-exlow))/2 is of (error of x)**2 order. NOTE: 1) By using the "effective variance" method a simple linear regression becomes a non-linear case , which takes several iterations instead of 0 as in the linear case .
2) The effective variance technique assumes that there is no correlation between the x and y coordinate .
In case the function lies below (above) the data point, ey is ey_low (ey_high).
Definition at line 2145 of file TFumili.cxx.
void H1FitChisquareFumili | ( | Int_t & | npar, |
Double_t * | gin, | ||
Double_t & | f, | ||
Double_t * | u, | ||
Int_t | flag | ||
) |
Minimization function for H1s using a Chisquare method.
Definition at line 2091 of file TFumili.cxx.
void H1FitLikelihoodFumili | ( | Int_t & | npar, |
Double_t * | gin, | ||
Double_t & | f, | ||
Double_t * | u, | ||
Int_t | flag | ||
) |
Minimization function for H1s using a Likelihood method.
Basically, it forms the likelihood by determining the Poisson probability that given a number of entries in a particular bin, the fit would predict it's value. This is then done for each bin, and the sum of the logs is taken as the likelihood. PDF: P=exp(-f(x_i))/[F_i]!*(f(x_i))^[F_i] where F_i - experimental value, f(x_i) - expected theoretical value [F_i] - integer part of F_i. drawback is that if F_i>Int_t - GetSumLog will fail for big F_i is faster to use Euler's Gamma-function
Definition at line 2109 of file TFumili.cxx.
TFumili* gFumili =0 |
Definition at line 120 of file TFumili.cxx.
|
static |
Definition at line 123 of file TFumili.cxx.
|
static |
Definition at line 124 of file TFumili.cxx.