Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hksimple.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_hist
3/// \notebook
4/// Illustrates the advantages of a TH1K histogram
5///
6/// \macro_image
7/// \macro_code
8///
9/// \author Victor Perevovchikov
10
11void canvasRefresh(TCanvas *c1)
12{
13 for (Int_t j = 0; j < 3; j++)
14 c1->GetPad(j+1)->Modified();
15
16 c1->Modified();
17 c1->Update();
18
20}
21
22void hksimple()
23{
24// Create a new canvas.
25 TCanvas* c1 = new TCanvas("c1","Dynamic Filling Example",200,10,600,900);
26
27// Create a normal histogram and two TH1K histograms
28 TH1 *hpx[3];
29 hpx[0] = new TH1F("hp0","Normal histogram",1000,-4,4);
30 hpx[1] = new TH1K("hk1","Nearest Neighbour of order 3",1000,-4,4);
31 hpx[2] = new TH1K("hk2","Nearest Neighbour of order 16",1000,-4,4,16);
32 c1->Divide(1,3);
33 for (Int_t j = 0; j < 3; j++) {
34 c1->cd(j + 1);
35 hpx[j]->SetFillColor(48);
36 hpx[j]->Draw();
37 }
38
39// Fill histograms randomly
40 gRandom->SetSeed(12345);
41 Float_t px, py, pz;
42 const Int_t kUPDATE = 10;
43 for (Int_t i = 0; i <= 300; i++) {
44 gRandom->Rannor(px,py);
45 for (Int_t j = 0; j < 3; j++)
46 hpx[j]->Fill(px);
47 if (i && (i % kUPDATE) == 0)
48 canvasRefresh(c1);
49 }
50
51 for (Int_t j = 0; j < 3; j++)
52 hpx[j]->Fit("gaus");
53
54 canvasRefresh(c1);
55}
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
R__EXTERN TSystem * gSystem
Definition TSystem.h:555
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
TH1K class supports the nearest K Neighbours method, widely used in cluster analysis.
Definition TH1K.h:26
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:615
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:507
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:416
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:133