Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hist012_TH1_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/// \date November 2022
10/// \author Victor Perevovchikov
11
13{
14 for (Int_t j = 0; j < 3; j++)
15 c1->GetPad(j + 1)->Modified();
16
17 c1->Modified();
18 c1->Update();
19
21}
22
24{
25 // Create a new canvas.
26 TCanvas *c1 = new TCanvas("c1", "Dynamic Filling Example", 200, 10, 600, 900);
27
28 // Create a normal histogram and two TH1K histograms
29 TH1 *hpx[3];
30 hpx[0] = new TH1F("hp0", "Normal histogram", 1000, -4, 4);
31 hpx[1] = new TH1K("hk1", "Nearest Neighbour of order 3", 1000, -4, 4);
32 hpx[2] = new TH1K("hk2", "Nearest Neighbour of order 16", 1000, -4, 4, 16);
33 c1->Divide(1, 3);
34 for (Int_t j = 0; j < 3; j++) {
35 c1->cd(j + 1);
36 hpx[j]->SetFillColor(48);
37 hpx[j]->Draw();
38 }
39
40 // Fill histograms randomly
41 gRandom->SetSeed(12345);
42 Float_t px, py, pz;
43 const Int_t kUPDATE = 10;
44 for (Int_t i = 0; i <= 300; i++) {
45 gRandom->Rannor(px, py);
46 for (Int_t j = 0; j < 3; j++)
47 hpx[j]->Fill(px);
48 if (i && (i % kUPDATE) == 0)
50 }
51
52 for (Int_t j = 0; j < 3; j++)
53 hpx[j]->Fit("gaus");
54
56}
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:645
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
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