Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
fitConvolution.C File Reference

Detailed Description

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

FCN=298.12 FROM MIGRAD STATUS=CONVERGED 457 CALLS 458 TOTAL
EDM=1.08093e-08 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 7.32859e+00 3.70795e-02 1.23437e-05 -3.46193e-02
2 p1 7.33040e-02 2.44083e-03 3.62176e-06 -7.16223e-02
3 p2 -2.26420e+00 4.91803e-02 5.24021e-05 -1.27917e-02
4 p3 1.12811e+00 6.28810e-02 1.94847e-05 -2.72591e-02
#include <TCanvas.h>
#include <TRandom.h>
#include <TF1Convolution.h>
#include <TF1.h>
#include <TH1F.h>
{
// 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.
new TCanvas("c", "c", 800, 1000);
h_ExpGauss->Fit("f");
h_ExpGauss->Draw();
}
#define f(i)
Definition RSha256.hxx:104
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
The Canvas class.
Definition TCanvas.h:23
Class wrapping convolution of two functions.
Int_t GetNpar() const
void SetRange(Double_t a, Double_t b) override
Set the actual range used for the convolution.
void SetNofPointsFFT(Int_t n)
Set the number of points used for the FFT convolution.
1-Dim function class
Definition TF1.h:213
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition TH1.cxx:3894
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3338
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3060
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:274
virtual Double_t Exp(Double_t tau)
Returns an exponential deviate.
Definition TRandom.cxx:251
Double_t x[n]
Definition legend1.C:17
Author
Aurelie Flandi

Definition in file fitConvolution.C.