ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
exampleTKDE.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// Example of using the TKDE class (kernel density estimator)
4 ///
5 /// \macro_image
6 /// \macro_code
7 ///
8 /// \authors Lorenzo Moneta, Bartolomeu Rabacal (Dec 2010)
9 
10 #include "TH1.h"
11 #include "TF1.h"
12 #include "TKDE.h"
13 #include "TCanvas.h"
14 //#include "TStopwatch.h"
15 #include "TRandom.h"
16 #ifndef __CINT__
17 #include "Math/DistFunc.h"
18 #endif
19 #include "TLegend.h"
20 
21 // test TKDE
22 
23 void exampleTKDE(int n = 1000) {
24 
25  // generate some gaussian points
26 
27  int nbin = 100;
28  double xmin = 0;
29  double xmax = 10;
30 
31  TH1D * h1 = new TH1D("h1","h1",nbin,xmin,xmax);
32 
33  // generate some points with bi- gaussian distribution
34 
35  std::vector<double> data(n);
36  for (int i = 0; i < n; ++i) {
37  if (i < 0.4*n) {
38  data[i] = gRandom->Gaus(2,1);
39  h1->Fill(data[i]);
40  }
41  else {
42  data[i] = gRandom->Gaus(7,1.5);
43  h1->Fill(data[i]);
44  }
45  }
46 
47  // scale histogram
48  h1->Scale(1./h1->Integral(),"width" );
49  h1->SetStats(false);
50  h1->SetTitle("Bi-Gaussian");
51  h1->Draw();
52 
53  // drawn true normalized density
54  TF1 * f1 = new TF1("f1","0.4*ROOT::Math::normal_pdf(x,1,2)+0.6*ROOT::Math::normal_pdf(x,1.5,7)",xmin,xmax);
55  f1->SetLineColor(kGreen+2);
56  f1->Draw("SAME");
57 
58  // create TKDE class
59  double rho = 1.0; //default value
60  TKDE * kde = new TKDE(n, &data[0], xmin,xmax, "", rho);
61  //kde->Draw("ConfidenceInterval@0.95 Same");
62  kde->Draw("SAME");
63 
64  TLegend * legend = new TLegend(0.6,0.7,0.9,0.95);
65  legend->AddEntry(f1,"True function");
66  legend->AddEntry(kde->GetDrawnFunction(),"TKDE");
67  legend->AddEntry(kde->GetDrawnLowerFunction(),"TKDE - #sigma");
68  legend->AddEntry(kde->GetDrawnUpperFunction(),"TKDE + #sigma");
69  legend->Draw();
70 }
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
TF1 * GetDrawnLowerFunction()
Definition: TKDE.h:129
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
Kernel Density Estimation class.
Definition: TKDE.h:37
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
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1052
Definition: Rtypes.h:61
virtual void Draw(const Option_t *option="")
Definition: TKDE.cxx:777
TF1 * GetDrawnUpperFunction()
Definition: TKDE.h:128
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual Double_t Integral(Option_t *option="") const
Return integral of bin contents.
Definition: TH1.cxx:7378
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
TF1 * GetDrawnFunction()
Definition: TKDE.h:127
1-Dim function class
Definition: TF1.h:149
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:6268
const Int_t n
Definition: legend1.C:16
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320