Basics of RHist, including filling and adding them.
#include <cstddef>
#include <iostream>
#include <random>
#include <variant>
{
std::cout << "\n";
for (int row = 15; row > 0; row--) {
static constexpr int Scale = 10;
std::cout << (value >= (row * Scale) ? '*' : ' ');
};
std::cout << '|';
for (auto bin : axis.GetNormalRange()) {
print(bin);
}
std::cout << '|';
std::cout << "\n";
}
}
void hist001_RHist_basics()
{
std::mt19937 gen;
std::normal_distribution normal1(5.0, 2.0);
for (std::size_t i = 0; i < 1000; i++) {
hist1.Fill(normal1(gen));
}
std::cout << "hist1 with expected mean = " << normal1.mean() << "\n";
DrawHistogram(hist1);
std::cout << "\n";
std::normal_distribution normal2(13.0, 4.0);
for (std::size_t i = 0; i < 1500; i++) {
hist2.Fill(normal2(gen));
}
std::cout << "hist2 with expected mean = " << normal2.mean() << "\n";
DrawHistogram(hist2);
std::cout << "\n";
hist3.Add(hist1);
hist3.Add(hist2);
std::cout << "hist3 with expected entries = " << (hist1.GetNEntries() + hist2.GetNEntries()) << "\n";
DrawHistogram(hist3);
}
A bin index with special values for underflow and overflow bins.
static RBinIndex Overflow()
static RBinIndex Underflow()
A histogram for aggregation of data along multiple dimensions.
double ComputeMean(std::size_t dim=0) const
Compute the arithmetic mean of unbinned values.
double ComputeStdDev(std::size_t dim=0) const
Compute the standard deviation of unbinned values.
const BinContentType & GetBinContent(const std::array< RBinIndex, N > &indices) const
Get the content of a single bin.
const std::vector< RAxisVariant > & GetAxes() const
std::uint64_t GetNEntries() const
A regular axis with equidistant bins in the interval .
hist1 with expected mean = 5
entries = 1000, mean = 4.98577, stddev = 1.98663
| |
| |
| |
| |
| |
| * |
| ***** |
| ***** |
| ****** * |
| ******** |
| ********** |
| ********** |
| * *********** |
| *************** |
*| ***************** |
hist2 with expected mean = 13
entries = 1500, mean = 12.9495, stddev = 4.00479
| |
| |
| |
| |
| |
| |
| |
| * * |
| * * * |
| *********** |*
| * ************** |*
| ****************** |*
| * ******************** * |*
| ** ********************** |*
| * *****************************|*
hist3 with expected entries = 2500
entries = 2500, mean = 9.76403, stddev = 5.14032
| |
| |
| |
| |
| * |
| ***** * |
| ***** * |
| ******** * * |
| ********* * ** * * |
| ********* *** *********** |*
| ***************************** |*
| ***************************** |*
| ********************************* * |*
| *********************************** |*
*| **************************************|*
- Date
- September 2025
- Author
- The ROOT Team
Definition in file hist001_RHist_basics.C.