Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hist001_RHist_basics.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_histv7
3///
4/// Basics of RHist, including filling and adding them.
5///
6/// \macro_code
7/// \macro_output
8///
9/// \date September 2025
10/// \author The ROOT Team
11
12#include <ROOT/RBinIndex.hxx>
13#include <ROOT/RHist.hxx>
14#include <ROOT/RRegularAxis.hxx>
15
16#include <cstddef>
17#include <iostream>
18#include <random>
19#include <variant>
20
21// It is currently not possible to directly draw an RHist, so this function implements an output with ASCII characters.
22static void DrawHistogram(const ROOT::Experimental::RHist<int> &hist)
23{
24 // Get the axis object from the histogram.
25 auto &axis = std::get<ROOT::Experimental::RRegularAxis>(hist.GetAxes()[0]);
26
27 // Print (some of) the global statistics.
28 std::cout << "entries = " << hist.GetNEntries();
29 std::cout << ", mean = " << hist.ComputeMean();
30 std::cout << ", stddev = " << hist.ComputeStdDev();
31 std::cout << "\n";
32
33 // "Draw" the histogram with ASCII characters. The height is hard-coded to work for this tutorial.
34 for (int row = 15; row > 0; row--) {
35 auto print = [&](ROOT::Experimental::RBinIndex bin) {
36 auto value = hist.GetBinContent(bin);
37 static constexpr int Scale = 10;
38 std::cout << (value >= (row * Scale) ? '*' : ' ');
39 };
40
41 // First the underflow bin, separated by a vertical bar.
43 std::cout << '|';
44
45 // Now iterate the normal bins and print a '*' if the value is sufficiently large.
46 for (auto bin : axis.GetNormalRange()) {
47 print(bin);
48 }
49
50 // Finally the overflow bin after a separating vertical bar.
51 std::cout << '|';
53 std::cout << "\n";
54 }
55}
56
58{
59 // Create an axis that can be used for multiple histograms.
60 ROOT::Experimental::RRegularAxis axis(40, {0.0, 20.0});
61
62 // Create a first histogram and fill with random values.
64
65 // Create a normal distribution with mean 5.0 and stddev 2.0.
66 std::mt19937 gen;
67 std::normal_distribution normal1(5.0, 2.0);
68 for (std::size_t i = 0; i < 1000; i++) {
69 hist1.Fill(normal1(gen));
70 }
71
72 // "Draw" the histogram with ASCII characters.
73 std::cout << "hist1 with expected mean = " << normal1.mean() << "\n";
75 std::cout << "\n";
76
77 // Create a second histogram and fill with random values of a different distribution.
79 std::normal_distribution normal2(13.0, 4.0);
80 for (std::size_t i = 0; i < 1500; i++) {
81 hist2.Fill(normal2(gen));
82 }
83 std::cout << "hist2 with expected mean = " << normal2.mean() << "\n";
85 std::cout << "\n";
86
87 // Create a third, merged histogram from the two previous two.
89 hist3.Add(hist1);
90 hist3.Add(hist2);
91 std::cout << "hist3 with expected entries = " << (hist1.GetNEntries() + hist2.GetNEntries()) << "\n";
93}
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 value
A bin index with special values for underflow and overflow bins.
Definition RBinIndex.hxx:22
static RBinIndex Overflow()
static RBinIndex Underflow()
A histogram for aggregation of data along multiple dimensions.
Definition RHist.hxx:55
double ComputeMean(std::size_t dim=0) const
Compute the arithmetic mean of unbinned values.
Definition RHist.hxx:113
double ComputeStdDev(std::size_t dim=0) const
Compute the standard deviation of unbinned values.
Definition RHist.hxx:115
const BinContentType & GetBinContent(const std::array< RBinIndex, N > &indices) const
Get the content of a single bin.
Definition RHist.hxx:136
const std::vector< RAxisVariant > & GetAxes() const
Definition RHist.hxx:105
std::uint64_t GetNEntries() const
Definition RHist.hxx:109
A regular axis with equidistant bins in the interval .