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 \( f(x) = p[1]* sin(p[0]*x) \).
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
FCN=49.5952 FROM MIGRAD STATUS=CONVERGED 52 CALLS 53 TOTAL
EDM=1.22682e-09 STRATEGY= 1 ERROR MATRIX UNCERTAINTY 2.5 per cent
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 3.13201e+00 3.12699e-02 -3.64656e-05 2.15221e-03
2 p1 2.97626e+01 1.00773e+00 6.67621e-05 -4.02033e-06
Covariance matrix from the fit
2x2 matrix is as follows
| 0 | 1 |
-------------------------------
0 | 0.0009778 0.009128
1 | 0.009128 1.016
Integral = 19.005 +/- 0.6159
#include <assert.h>
#include <iostream>
#include <cmath>
const int NPAR = 2;
double f(
double *
x,
double * p) {
}
void ErrorIntegral() {
fitFunc =
new TF1(
"f",
f,0,1,NPAR);
double par[NPAR] = { 3.14, 1.};
auto fitResult =
h1->
Fit(fitFunc,
"S");
double integral = fitFunc->
Integral(0,1);
auto covMatrix = fitResult->GetCovarianceMatrix();
std::cout << "Covariance matrix from the fit ";
covMatrix.Print();
double sigma_integral = fitFunc->
IntegralError(0,1, fitResult->GetParams() , covMatrix.GetMatrixArray());
std::cout << "Integral = " << integral << " +/- " << sigma_integral
<< 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];
double sic = std::sqrt( c0c*c0c * covMatrix(0,0) + c1c*c1c * covMatrix(1,1)
+ 2.* c0c*c1c * covMatrix(0,1));
if ( std::fabs(sigma_integral-sic) > 1.E-6*sic )
std::cout << " ERROR: test failed : different analytical integral : "
<< ic << " +/- " << sic << std::endl;
}
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
virtual Double_t * GetParameters() const
virtual void SetParameters(const Double_t *params)
virtual void SetParameter(Int_t param, Double_t value)
virtual Double_t IntegralError(Double_t a, Double_t b, const Double_t *params=0, const Double_t *covmat=0, Double_t epsilon=1.E-2)
Return Error on Integral of a parametric function between a and b due to the parameter uncertainties ...
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.
virtual void Draw(Option_t *option="")
Draw this histogram with options.
- Author
- Lorenzo Moneta
Definition in file ErrorIntegral.C.