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

Detailed Description

View in nbviewer Open in SWAN Example of using the TKDE class (kernel density estimator)

#include "TH1.h"
#include "TF1.h"
#include "TKDE.h"
#include "TCanvas.h"
/*#include "TStopwatch.h"*/
#include "TRandom.h"
#include "Math/DistFunc.h"
#include "TLegend.h"
// test TKDE
void exampleTKDE(int n = 1000) {
// generate some gaussian points
int nbin = 100;
double xmin = 0;
double xmax = 10;
TH1D * h1 = new TH1D("h1","h1",nbin,xmin,xmax);
// generate some points with bi- gaussian distribution
std::vector<double> data(n);
for (int i = 0; i < n; ++i) {
if (i < 0.4*n) {
data[i] = gRandom->Gaus(2,1);
h1->Fill(data[i]);
}
else {
data[i] = gRandom->Gaus(7,1.5);
h1->Fill(data[i]);
}
}
// scale histogram
h1->Scale(1./h1->Integral(),"width" );
h1->SetStats(false);
h1->SetTitle("Bi-Gaussian");
h1->Draw();
// drawn true normalized density
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);
f1->Draw("SAME");
// create TKDE class
double rho = 1.0; //default value
TKDE * kde = new TKDE(n, &data[0], xmin,xmax, "", rho);
//kde->Draw("ConfidenceInterval@0.95 Same");
kde->Draw("SAME");
TLegend * legend = new TLegend(0.6,0.7,0.9,0.95);
legend->AddEntry(f1,"True function");
legend->AddEntry(kde->GetDrawnFunction(),"TKDE");
legend->AddEntry(kde->GetDrawnLowerFunction(),"TKDE - #sigma");
legend->AddEntry(kde->GetDrawnUpperFunction(),"TKDE + #sigma");
legend->Draw();
}
@ kGreen
Definition Rtypes.h:66
float xmin
float xmax
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
1-Dim function class
Definition TF1.h:213
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition TF1.cxx:1322
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:618
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition TH1.cxx:6678
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3350
virtual Double_t Integral(Option_t *option="") const
Return integral of bin contents.
Definition TH1.cxx:7834
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3073
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition TH1.cxx:6564
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8830
Kernel Density Estimation class.
Definition TKDE.h:34
TF1 * GetDrawnUpperFunction()
Definition TKDE.h:148
TF1 * GetDrawnLowerFunction()
Definition TKDE.h:149
TF1 * GetDrawnFunction()
Definition TKDE.h:147
virtual void Draw(const Option_t *option="")
Draws either the KDE functions or its errors.
Definition TKDE.cxx:897
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition TLegend.cxx:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition TLegend.cxx:423
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
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5
TF1 * f1
Definition legend1.C:11
Authors
Lorenzo Moneta, Bartolomeu Rabacal (Dec 2010)

Definition in file exampleTKDE.C.