Logo ROOT   6.08/07
Reference Guide
fitConvolution.C File Reference

Detailed Description

View in nbviewer Open in SWAN Tutorial for convolution of two functions

pict1_fitConvolution.C.png
Processing /mnt/build/workspace/root-makedoc-v608/rootspi/rdoc/src/v6-08-00-patches/tutorials/fit/fitConvolution.C...
FCN=298.12 FROM MIGRAD STATUS=CONVERGED 457 CALLS 458 TOTAL
EDM=1.08095e-08 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 7.32859e+00 3.70797e-02 1.23437e-05 -3.46156e-02
2 p1 7.33040e-02 2.44084e-03 3.62176e-06 -7.16205e-02
3 p2 -2.26420e+00 4.91805e-02 5.24021e-05 -1.27912e-02
4 p3 1.12811e+00 6.28812e-02 1.94847e-05 -2.72567e-02
#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.