ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ErrorIntegral.C
Go to the documentation of this file.
1 // Estimate the error in the integral of a fitted function
2 // taking into account the errors in the parameters resulting from the fit.
3 // The error is estimated also using the correlations values obtained from
4 // the fit
5 //
6 // run the macro doing:
7 // .x ErrorIntegral.C
8 // Author: Lorenzo Moneta
9 
10 #include "TF1.h"
11 #include "TH1D.h"
12 #include "TVirtualFitter.h"
13 #include "TMath.h"
14 #include <assert.h>
15 #include <iostream>
16 #include <cmath>
17 
18 //#define HAVE_OLD_ROOT_VERSION
19 
20 TF1 * fitFunc; // fit function pointer
21 
22 const int NPAR = 2; // number of function parameters;
23 
24 //____________________________________________________________________
25 double f(double * x, double * p) {
26  // function used to fit the data
27  return p[1]*TMath::Sin( p[0] * x[0] );
28 }
29 
30 // when TF1::IntegralError was not available
31 
32 #ifdef HAVE_OLD_ROOT_VERSION
33 //____________________________________________________________________
34 double df_dPar(double * x, double * p) {
35  // derivative of the function w.r..t parameters
36  // use calculated derivatives from TF1::GradientPar
37 
38  double grad[NPAR];
39  // p is used to specify for which parameter the derivative is computed
40  int ipar = int(p[0] );
41  assert (ipar >=0 && ipar < NPAR );
42 
43  assert(fitFunc);
44  fitFunc->GradientPar(x, grad);
45 
46  return grad[ipar];
47 }
48 
49 //____________________________________________________________________
50 double IntegralError(int npar, double * c, double * errPar,
51  double * covMatrix = 0) {
52 // calculate the error on the integral given the parameter error and
53 // the integrals of the gradient functions c[]
54 
55  double err2 = 0;
56  for (int i = 0; i < npar; ++i) {
57  if (covMatrix == 0) { // assume error are uncorrelated
58  err2 += c[i] * c[i] * errPar[i] * errPar[i];
59  } else {
60  double s = 0;
61  for (int j = 0; j < npar; ++j) {
62  s += covMatrix[i*npar + j] * c[j];
63  }
64  err2 += c[i] * s;
65  }
66  }
67 
68  return TMath::Sqrt(err2);
69 }
70 #endif
71 
72 //____________________________________________________________________
73 void ErrorIntegral() {
74  fitFunc = new TF1("f",f,0,1,NPAR);
75  TH1D * h1 = new TH1D("h1","h1",50,0,1);
76 
77  double par[NPAR] = { 3.14, 1.};
78  fitFunc->SetParameters(par);
79 
80  h1->FillRandom("f",1000); // fill histogram sampling fitFunc
81  fitFunc->SetParameter(0,3.); // vary a little the parameters
82  h1->Fit(fitFunc); // fit the histogram
83 
84  h1->Draw();
85 
86  // calculate the integral
87  double integral = fitFunc->Integral(0,1);
88 
90  assert(fitter != 0);
91  double * covMatrix = fitter->GetCovarianceMatrix();
92 
93 #ifdef HAVE_OLD_ROOT_VERSION
94 
95  // calculate now the error (needs the derivatives of the function
96  // w..r.t the parameters)
97  TF1 * deriv_par0 = new TF1("dfdp0",df_dPar,0,1,1);
98  deriv_par0->SetParameter(0,0);
99 
100  TF1 * deriv_par1 = new TF1("dfdp1",df_dPar,0,1,1);
101  deriv_par1->SetParameter(0,1.);
102 
103  double c[2];
104 
105  c[0] = deriv_par0->Integral(0,1);
106  c[1] = deriv_par1->Integral(0,1);
107 
108  double * epar = fitFunc->GetParErrors();
109 
110  // without correlations
111  double sigma_integral_0 = IntegralError(2,c,epar);
112 
113 
114 
115  // with correlations
116  double sigma_integral = IntegralError(2,c,epar,covMatrix);
117 
118 #else
119 
120  // using new function in TF1 (from 12/6/2007)
121  double sigma_integral = fitFunc->IntegralError(0,1);
122 
123 #endif
124 
125  std::cout << "Integral = " << integral << " +/- " << sigma_integral
126  << std::endl;
127 
128  // estimated integral and error analytically
129 
130  double * p = fitFunc->GetParameters();
131  double ic = p[1]* (1-std::cos(p[0]) )/p[0];
132  double c0c = p[1] * (std::cos(p[0]) + p[0]*std::sin(p[0]) -1.)/p[0]/p[0];
133  double c1c = (1-std::cos(p[0]) )/p[0];
134 
135  // estimated error with correlations
136  double sic = std::sqrt( c0c*c0c * covMatrix[0] + c1c*c1c * covMatrix[3]
137  + 2.* c0c*c1c * covMatrix[1]);
138 
139  if ( std::fabs(sigma_integral-sic) > 1.E-6*sic )
140  std::cout << " ERROR: test failed : different analytical integral : "
141  << ic << " +/- " << sic << std::endl;
142 }
double par[1]
Definition: unuranDistr.cxx:38
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
double f(double *x, double *p)
Definition: ErrorIntegral.C:25
virtual Double_t * GetCovarianceMatrix() const =0
#define assert(cond)
Definition: unittest.h:542
TF1 * fitFunc
Definition: ErrorIntegral.C:20
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
double cos(double)
double sqrt(double)
TH1F * h1
Definition: legend1.C:5
double sin(double)
virtual const Double_t * GetParErrors() const
Definition: TF1.h:370
virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps=0.01)
Compute the gradient (derivative) wrt a parameter ipar.
Definition: TF1.cxx:2117
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH1.cxx:3330
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
double IntegralError(TF1 *func, Int_t ndim, const double *a, const double *b, const double *params, const double *covmat, double epsilon)
Definition: TF1Helper.cxx:38
static TVirtualFitter * GetFitter()
static: return the current Fitter
Double_t E()
Definition: TMath.h:54
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
void ErrorIntegral()
Definition: ErrorIntegral.C:73
Abstract Base Class for Fitting.
virtual Double_t * GetParameters() const
Definition: TF1.h:358
1-Dim function class
Definition: TF1.h:149
Double_t Sin(Double_t)
Definition: TMath.h:421
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 parameteric function between a and b due to the parameter uncertainties...
Definition: TF1.cxx:2422
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:424
const int NPAR
Definition: ErrorIntegral.C:22
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
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.
Definition: TH1.cxx:3607