Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fitConvolution.C File Reference

Detailed Description

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

****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 298.12
NDf = 96
Edm = 1.67196e-06
NCalls = 448
p0 = 7.32861 +/- 0.0370492
p1 = 0.0733018 +/- 0.00243973
p2 = -2.26418 +/- 0.0491372
p3 = 1.12808 +/- 0.0628185
#include <TCanvas.h>
#include <TRandom.h>
#include <TF1Convolution.h>
#include <TF1.h>
#include <TH1F.h>
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++) {
// Gives a alpha of -0.3 in the exp.
double x = gRandom->Exp(1. / 0.3);
x += gRandom->Gaus(0., 3.);
// Probability density function of the addition of two variables is the
// convolution of two density functions.
h_ExpGauss->Fill(x);
}
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.
h_ExpGauss->Fit("f");
}
#define f(i)
Definition RSha256.hxx:104
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
Class wrapping convolution of two functions.
1-Dim function class
Definition TF1.h:234
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:878
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition TRandom.cxx:275
virtual Double_t Exp(Double_t tau)
Returns an exponential deviate.
Definition TRandom.cxx:252
Double_t x[n]
Definition legend1.C:17
Author
Aurelie Flandi

Definition in file fitConvolution.C.