ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
fitConvolution.C
Go to the documentation of this file.
1 //
2 // Convolution.c
3 //
4 //
5 // Created by AurĂ©lie Flandi on 09.09.14.
6 //
7 //
8 
9 #include <stdio.h>
10 #include <TMath.h>
11 #include <TCanvas.h>
12 #include <iostream>
13 #include <TROOT.h>
14 #include <TChain.h>
15 #include <TObject.h>
16 #include <TRandom.h>
17 #include <TFile.h>
18 #include <math.h>
19 #include <TF1Convolution.h>
20 #include <TF1.h>
21 #include <TH1F.h>
22 #include <TGraph.h>
23 #include <TStopwatch.h>
24 
25 
26 using namespace std;
27 
28 
30 {
31 
32  //tutorial for convolution of two functions
33 
34 
35  //construction of histogram to fit
36  TH1F *h_ExpGauss = new TH1F("h_ExpGauss","Exponential convoluted by gaussian",100,0.,5.);
37  for (int i=0;i<1e6;i++)
38  {
39  Double_t x = gRandom->Exp(1./0.3);//gives a alpha of -0.3 in the exp
40  x += gRandom->Gaus(0.,3.);
41  h_ExpGauss->Fill(x);//probability density function of the addition of two variables is the convolution of 2 dens. functions
42  }
43 
44  TF1Convolution *f_conv = new TF1Convolution("expo","gaus",-1,6,true);
45  f_conv->SetRange(-1.,6.);
46  f_conv->SetNofPointsFFT(1000);
47  TF1 *f = new TF1("f",*f_conv, 0., 5., f_conv->GetNpar());
48  f->SetParameters(1.,-0.3,0.,1.);
49 
50  //fit
51  new TCanvas("c","c",800,1000);
52  h_ExpGauss -> Fit("f");
53  h_ExpGauss->Draw();
54 
55 }
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:235
Int_t GetNpar() const
TFile * f
Class wrapping convolution of two functions.
Double_t x[n]
Definition: legend1.C:17
void SetRange(Double_t a, Double_t b)
void SetNofPointsFFT(Int_t n)
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
void fitConvolution()
double Double_t
Definition: RtypesCore.h:55
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:133
virtual Double_t Exp(Double_t tau)
Returns an exponential deviate.
Definition: TRandom.cxx:212