Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hist004_RHist_concurrent.C File Reference

Detailed Description

Concurrent filling of RHist.

#include <ROOT/RHist.hxx>
#include <cstddef>
#include <iostream>
#include <memory>
#include <random>
#include <thread>
#include <vector>
// It is currently not possible to directly draw an RHist, so this function implements an output with ASCII characters.
{
// Get the axis object from the histogram.
auto &axis = hist.GetAxes()[0];
// Print (some of) the global statistics.
std::cout << "entries = " << hist.GetNEntries();
std::cout << ", mean = " << hist.ComputeMean();
std::cout << ", stddev = " << hist.ComputeStdDev();
std::cout << "\n";
// "Draw" the histogram with ASCII characters. The height is hard-coded to work for this tutorial.
for (int row = 15; row > 0; row--) {
auto print = [&](ROOT::Experimental::RBinIndex bin) {
auto value = hist.GetBinContent(bin);
static constexpr int Scale = 15;
std::cout << (value >= (row * Scale) ? '*' : ' ');
};
// First the underflow bin, separated by a vertical bar.
std::cout << '|';
// Now iterate the normal bins and print a '*' if the value is sufficiently large.
for (auto bin : axis.GetNormalRange()) {
print(bin);
}
// Finally the overflow bin after a separating vertical bar.
std::cout << '|';
std::cout << "\n";
}
}
static void FillHistogram(std::size_t t, ROOT::Experimental::RHistFillContext<int> &context)
{
// Create a normal distribution with mean 10.0 and stddev 4.0.
std::mt19937 gen(t);
std::normal_distribution normal(10.0, 4.0);
for (std::size_t i = 0; i < 1000; i++) {
context.Fill(normal(gen));
}
}
{
// Create an axis and a histogram that multiple threads will fill concurrently.
ROOT::Experimental::RRegularAxis axis(40, {0.0, 20.0});
auto hist = std::make_shared<ROOT::Experimental::RHist<int>>(axis);
{
// Fill histogram from multiple threads.
std::vector<std::thread> threads;
for (std::size_t t = 0; t < 4; t++) {
threads.emplace_back([t, &filler] {
// Create a fill context for this thread.
auto context = filler.CreateFillContext();
FillHistogram(t, *context);
});
}
for (auto &&t : threads) {
t.join();
}
}
// "Draw" the histogram with ASCII characters.
DrawHistogram(*hist);
}
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:23
static RBinIndex Overflow()
static RBinIndex Underflow()
A histogram filler to concurrently fill an RHist.
A context to concurrently fill an RHist.
void Fill(const std::tuple< A... > &args)
Fill an entry into the histogram.
A histogram for aggregation of data along multiple dimensions.
Definition RHist.hxx:64
double ComputeMean(std::size_t dim=0) const
Compute the arithmetic mean of unbinned values.
Definition RHist.hxx:159
double ComputeStdDev(std::size_t dim=0) const
Compute the standard deviation of unbinned values.
Definition RHist.hxx:161
const BinContentType & GetBinContent(const std::array< RBinIndex, N > &indices) const
Get the content of a single bin.
Definition RHist.hxx:182
const std::vector< RAxisVariant > & GetAxes() const
Definition RHist.hxx:151
std::uint64_t GetNEntries() const
Definition RHist.hxx:155
A regular axis with equidistant bins in the interval .
entries = 4000, mean = 10.0266, stddev = 4.01885
| |
| |
| **** |
| ***** ** * |
| * ******** * |
| ************* |
| ************* |
| ************** * |
| ******************* |
| ******************* |
| ********************** |
| ************************* |
| ***************************** |
| ********************************* |
*| *********************************** * |*
Date
February 2026
Author
The ROOT Team

Definition in file hist004_RHist_concurrent.C.