From $ROOTSYS/tutorials/fit/fitConvolution.C

//
//  Convolution.c
//  
//
//  Created by Aurélie Flandi on 09.09.14.
//
//

#include <stdio.h>
#include <TMath.h>
#include <TCanvas.h>
#include <iostream>
#include <TROOT.h>
#include <TChain.h>
#include <TObject.h>
#include <TRandom.h>
#include <TFile.h>
#include <math.h>
#include <TF1Convolution.h>
#include <TF1.h>
#include <TH1F.h>
#include <TGraph.h>
#include <TStopwatch.h>


using namespace std;


void fitConvolution()
{

   //tutorial for convolution of two functions
   
   
   //construction of histogram to fit
   TH1F *h_ExpGauss = new TH1F("h_ExpGauss","Exponential convoluted by gaussian",100,0.,5.);
   for (int i=0;i<1e6;i++)
   {
      Double_t x = gRandom->Exp(1./0.3);//gives a alpha of -0.3 in the exp
      x += gRandom->Gaus(0.,3.);
      h_ExpGauss->Fill(x);//probability density function of the addition of two variables is the convolution of 2 dens. functions
   }

   TF1Convolution *f_conv = new TF1Convolution("expo","gaus",-1,6,true);
   f_conv->SetRange(-1.,6.);
   f_conv->SetNofPointsFFT(1000);
   TF1   *f = new TF1("f",*f_conv, 0., 5., f_conv->GetNpar());
   f->SetParameters(1.,-0.3,0.,1.);
   
   //fit
   new TCanvas("c","c",800,1000);
   h_ExpGauss -> Fit("f");
   h_ExpGauss->Draw();
  
}
 fitConvolution.C:1
 fitConvolution.C:2
 fitConvolution.C:3
 fitConvolution.C:4
 fitConvolution.C:5
 fitConvolution.C:6
 fitConvolution.C:7
 fitConvolution.C:8
 fitConvolution.C:9
 fitConvolution.C:10
 fitConvolution.C:11
 fitConvolution.C:12
 fitConvolution.C:13
 fitConvolution.C:14
 fitConvolution.C:15
 fitConvolution.C:16
 fitConvolution.C:17
 fitConvolution.C:18
 fitConvolution.C:19
 fitConvolution.C:20
 fitConvolution.C:21
 fitConvolution.C:22
 fitConvolution.C:23
 fitConvolution.C:24
 fitConvolution.C:25
 fitConvolution.C:26
 fitConvolution.C:27
 fitConvolution.C:28
 fitConvolution.C:29
 fitConvolution.C:30
 fitConvolution.C:31
 fitConvolution.C:32
 fitConvolution.C:33
 fitConvolution.C:34
 fitConvolution.C:35
 fitConvolution.C:36
 fitConvolution.C:37
 fitConvolution.C:38
 fitConvolution.C:39
 fitConvolution.C:40
 fitConvolution.C:41
 fitConvolution.C:42
 fitConvolution.C:43
 fitConvolution.C:44
 fitConvolution.C:45
 fitConvolution.C:46
 fitConvolution.C:47
 fitConvolution.C:48
 fitConvolution.C:49
 fitConvolution.C:50
 fitConvolution.C:51
 fitConvolution.C:52
 fitConvolution.C:53
 fitConvolution.C:54
 fitConvolution.C:55
 fitConvolution.C:56