Logo ROOT  
Reference Guide
hksimple.C File Reference

Detailed Description

View in nbviewer Open in SWAN Illustrates the advantages of a TH1K histogram

void padRefresh(TPad *pad,int flag=0)
{
if (!pad) return;
pad->Modified();
pad->Update();
TList *tl = pad->GetListOfPrimitives();
if (!tl) return;
TListIter next(tl);
TObject *to;
while ((to=next())) {
if (to->InheritsFrom(TPad::Class())) padRefresh((TPad*)to,1);}
if (flag) return;
}
void hksimple()
{
// Create a new canvas.
TCanvas* c1 = new TCanvas("c1","Dynamic Filling Example",200,10,600,900);
// Create a normal histogram and two TH1K histograms
TH1 *hpx[3];
hpx[0] = new TH1F("hp0","Normal histogram",1000,-4,4);
hpx[1] = new TH1K("hk1","Nearest Neighbour of order 3",1000,-4,4);
hpx[2] = new TH1K("hk2","Nearest Neighbour of order 16",1000,-4,4,16);
c1->Divide(1,3);
Int_t j;
for (j=0;j<3;j++) {
c1->cd(j+1);
hpx[j]->SetFillColor(48);
hpx[j]->Draw();
}
// Fill histograms randomly
gRandom->SetSeed(12345);
Float_t px, py, pz;
const Int_t kUPDATE = 10;
for (Int_t i = 0; i <= 300; i++) {
gRandom->Rannor(px,py);
for (j=0;j<3;j++) {hpx[j]->Fill(px);}
if (i && (i%kUPDATE) == 0) {
padRefresh(c1);
}
}
for (j=0;j<3;j++) hpx[j]->Fit("gaus");
padRefresh(c1);
}
void Class()
Definition: Class.C:29
int Int_t
Definition: RtypesCore.h:41
float Float_t
Definition: RtypesCore.h:53
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
R__EXTERN TSystem * gSystem
Definition: TSystem.h:560
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
The Canvas class.
Definition: TCanvas.h:31
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
TH1K class supports the nearest K Neighbours method, widely used in cluster analysis.
Definition: TH1K.h:27
The TH1 histogram class.
Definition: TH1.h:56
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3275
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
Iterator of linked list.
Definition: TList.h:200
A doubly linked list.
Definition: TList.h:44
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
The most important graphics class in the ROOT system.
Definition: TPad.h:29
void Modified(Bool_t flag=1)
Definition: TPad.h:417
virtual void Update()
Update pad.
Definition: TPad.cxx:2833
TList * GetListOfPrimitives() const
Definition: TPad.h:242
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:597
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
Definition: TRandom.cxx:489
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition: TSystem.cxx:426
return c1
Definition: legend1.C:41
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:134
Author
Victor Perevovchikov

Definition in file hksimple.C.