Logo ROOT   6.10/02
Reference Guide
fitConvolution.C File Reference

Detailed Description

Tutorial for convolution of two functions.

pict1_fitConvolution.C.png
Processing /builddir/build/BUILD/root-6.10.02/tutorials/fit/fitConvolution.C...
FCN=298.12 FROM MIGRAD STATUS=CONVERGED 446 CALLS 447 TOTAL
EDM=7.46173e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 7.32859e+00 3.70657e-02 1.23436e-05 -5.01424e-01
2 p1 7.33043e-02 2.44038e-03 3.62160e-06 -1.15101e+00
3 p2 -2.26420e+00 4.91619e-02 5.24449e-05 -1.41239e-01
4 p3 1.12810e+00 6.28532e-02 1.94710e-05 -3.53890e-01
#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()
{
//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();
}
Author
Aurelie Flandi

Definition in file fitConvolution.C.