Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
fitConvolution.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_fit
3 /// \notebook -js
4 /// Tutorial for convolution of two functions
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \author Aurelie Flandi
11 
12 #include <stdio.h>
13 #include <TMath.h>
14 #include <TCanvas.h>
15 #include <iostream>
16 #include <TROOT.h>
17 #include <TChain.h>
18 #include <TObject.h>
19 #include <TRandom.h>
20 #include <TFile.h>
21 #include <math.h>
22 #include <TF1Convolution.h>
23 #include <TF1.h>
24 #include <TH1F.h>
25 #include <TGraph.h>
26 #include <TStopwatch.h>
27 
28 using namespace std;
29 
30 void fitConvolution()
31 {
32  //construction of histogram to fit
33  TH1F *h_ExpGauss = new TH1F("h_ExpGauss","Exponential convoluted by gaussian",100,0.,5.);
34  for (int i=0;i<1e6;i++)
35  {
36  Double_t x = gRandom->Exp(1./0.3);//gives a alpha of -0.3 in the exp
37  x += gRandom->Gaus(0.,3.);
38  h_ExpGauss->Fill(x);//probability density function of the addition of two variables is the convolution of 2 dens. functions
39  }
40 
41  TF1Convolution *f_conv = new TF1Convolution("expo","gaus",-1,6,true);
42  f_conv->SetRange(-1.,6.);
43  f_conv->SetNofPointsFFT(1000);
44  TF1 *f = new TF1("f",*f_conv, 0., 5., f_conv->GetNpar());
45  f->SetParameters(1.,-0.3,0.,1.);
46 
47  //fit
48  new TCanvas("c","c",800,1000);
49  h_ExpGauss -> Fit("f");
50  h_ExpGauss->Draw();
51 
52 }