Estimate the error in the integral of a fitted function taking into account the errors in the parameters resulting from the fit.
The error is estimated also using the correlations values obtained from the fit
run the macro doing:
After having computed the integral and its error using the integral and the integral error using the generic functions TF1::Integral and TF1::IntegralError, we compute the integrals and its error analytically using the fact that the fitting function is .
Therefore we have:
- integral in [0,1] :
ic = p[1]* (1-std::cos(p[0]) )/p[0]
- derivative of integral with respect to p0:
c0c = p[1] * (std::cos(p[0]) + p[0]*std::sin(p[0]) -1.)/p[0]/p[0]
- derivative of integral with respect to p1:
c1c = (1-std::cos(p[0]) )/p[0]
and then we can compute the integral error using error propagation and the covariance matrix for the parameters p obtained from the fit.
integral error : sic = std::sqrt( c0c*c0c * covMatrix(0,0) + c1c*c1c * covMatrix(1,1) + 2.* c0c*c1c * covMatrix(0,1))
Note that, if possible, one should fit directly the function integral, which are the number of events of the different components (e.g. signal and background). In this way one obtains a better and more correct estimate of the integrals uncertainties, since they are obtained directly from the fit without using the approximation of error propagation. This is possible in ROOT. when using the TF1NormSum class, see the tutorial fitNormSum.C
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 49.5952
NDf = 48
Edm = 1.61787e-06
NCalls = 58
p0 = 3.13205 +/- 0.0312726
p1 = 29.7625 +/- 1.00876
Covariance matrix from the fit
2x2 matrix is as follows
| 0 | 1 |
-------------------------------
0 | 0.000978 0.009147
1 | 0.009147 1.018
Integral = 19.0047 +/- 0.616472
#include <cassert>
#include <iostream>
#include <cmath>
double f(
double *
x,
double *
p) {
}
double par[
NPAR] = { 3.14, 1.};
double integral =
fitFunc->Integral(0,1);
auto covMatrix = fitResult->GetCovarianceMatrix();
std::cout << "Covariance matrix from the fit ";
<< std::endl;
double ic =
p[1]* (1-std::cos(
p[0]) )/
p[0];
double c0c =
p[1] * (std::cos(
p[0]) +
p[0]*std::sin(
p[0]) -1.)/
p[0]/
p[0];
double c1c = (1-std::cos(
p[0]) )/
p[0];
std::cout << " ERROR: test failed : different analytical integral : "
<<
ic <<
" +/- " <<
sic << std::endl;
}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
1-D histogram with a double per channel (see TH1 documentation)
virtual void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr)
Fill histogram following distribution in function fname.
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
void Draw(Option_t *option="") override
Draw this histogram with options.
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
- Author
- Lorenzo Moneta
Definition in file ErrorIntegral.C.