Logo ROOT   6.07/09
Reference Guide
exampleTKDE.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// \notebook
4 /// Example of using the TKDE class (kernel density estimator)
5 ///
6 /// \macro_image
7 /// \macro_code
8 ///
9 /// \authors Lorenzo Moneta, Bartolomeu Rabacal (Dec 2010)
10 
11 #include "TH1.h"
12 #include "TF1.h"
13 #include "TKDE.h"
14 #include "TCanvas.h"
15 /*#include "TStopwatch.h"*/
16 #include "TRandom.h"
17 #include "Math/DistFunc.h"
18 #include "TLegend.h"
19 
20 // test TKDE
21 
22 void exampleTKDE(int n = 1000) {
23 
24  // generate some gaussian points
25 
26  int nbin = 100;
27  double xmin = 0;
28  double xmax = 10;
29 
30  TH1D * h1 = new TH1D("h1","h1",nbin,xmin,xmax);
31 
32  // generate some points with bi- gaussian distribution
33 
34  std::vector<double> data(n);
35  for (int i = 0; i < n; ++i) {
36  if (i < 0.4*n) {
37  data[i] = gRandom->Gaus(2,1);
38  h1->Fill(data[i]);
39  }
40  else {
41  data[i] = gRandom->Gaus(7,1.5);
42  h1->Fill(data[i]);
43  }
44  }
45 
46  // scale histogram
47  h1->Scale(1./h1->Integral(),"width" );
48  h1->SetStats(false);
49  h1->SetTitle("Bi-Gaussian");
50  h1->Draw();
51 
52  // drawn true normalized density
53  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);
54  f1->SetLineColor(kGreen+2);
55  f1->Draw("SAME");
56 
57  // create TKDE class
58  double rho = 1.0; //default value
59  TKDE * kde = new TKDE(n, &data[0], xmin,xmax, "", rho);
60  //kde->Draw("ConfidenceInterval@0.95 Same");
61  kde->Draw("SAME");
62 
63  TLegend * legend = new TLegend(0.6,0.7,0.9,0.95);
64  legend->AddEntry(f1,"True function");
65  legend->AddEntry(kde->GetDrawnFunction(),"TKDE");
66  legend->AddEntry(kde->GetDrawnLowerFunction(),"TKDE - #sigma");
67  legend->AddEntry(kde->GetDrawnUpperFunction(),"TKDE + #sigma");
68  legend->Draw();
69 }
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:5893
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3127
TF1 * GetDrawnLowerFunction()
Definition: TKDE.h:129
float xmin
Definition: THbookFile.cxx:93
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:27
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:1081
Definition: Rtypes.h:61
virtual void Draw(const Option_t *option="")
Definition: TKDE.cxx:777
TLegend * legend
Definition: pirndm.C:35
TF1 * GetDrawnUpperFunction()
Definition: TKDE.h:128
TH1F * h1
Definition: legend1.C:5
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
virtual Double_t Integral(Option_t *option="") const
Return integral of bin contents.
Definition: TH1.cxx:7081
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2853
float xmax
Definition: THbookFile.cxx:93
R__EXTERN TRandom * gRandom
Definition: TRandom.h:66
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:618
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
TF1 * f1
Definition: legend1.C:11
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:301
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:5984
const Int_t n
Definition: legend1.C:16
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8058