Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hist009_TH1_normalize.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_hist
3/// \notebook -js
4/// Normalizing a Histogram
5///
6/// Image produced by `.x NormalizeHistogram.C`
7/// Two different methods of normalizing histograms
8/// are shown, each with the original histogram.
9/// next to the normalized one.
10/// \macro_image
11/// \macro_code
12///
13/// \date November 2024
14/// \author Advait Dhingra
15
17{
18 const std::array<double, 6> binsx{0, 5, 10, 20, 50, 100};
19 TH1D *orig = new TH1D("orig", "Original histogram before normalization", binsx.size() - 1, binsx.data());
20
22
23 // Filling histogram with random entries
25 for (int i = 0; i < 100000; ++i) {
26 double r = rand.Rndm() * 100;
27 orig->Fill(r);
28 }
29
30 TH1D *norm = static_cast<TH1D *>(orig->Clone("norm"));
31 norm->SetTitle("Normalized Histogram");
32
33 // Normalizing the Histogram by scaling by 1 / the integral and taking width into account
34 norm->Scale(1. / norm->Integral(), "width");
35
36 // Drawing everything
37 TCanvas *c1 = new TCanvas("c1", "Histogram Normalization", 700, 900);
38 c1->Divide(1, 2);
39
40 c1->cd(1);
41 orig->Draw();
42 c1->cd(2);
43 norm->Draw();
44}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:693
Random number generator class based on the maximally quidistributed combined Tausworthe generator by ...
Definition TRandom2.h:27
void SetTitleFontSize(Float_t size=0)
Definition TStyle.h:411
return c1
Definition legend1.C:41