Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hist004_RHist_concurrent.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_histv7
3///
4/// Concurrent filling of RHist.
5///
6/// \macro_code
7/// \macro_output
8///
9/// \date February 2026
10/// \author The ROOT Team
11
12#include <ROOT/RAxisVariant.hxx>
13#include <ROOT/RBinIndex.hxx>
14#include <ROOT/RHist.hxx>
17#include <ROOT/RRegularAxis.hxx>
18
19#include <cstddef>
20#include <iostream>
21#include <memory>
22#include <random>
23#include <thread>
24#include <vector>
25
26// It is currently not possible to directly draw an RHist, so this function implements an output with ASCII characters.
27static void DrawHistogram(const ROOT::Experimental::RHist<int> &hist)
28{
29 // Get the axis object from the histogram.
30 auto &axis = hist.GetAxes()[0];
31
32 // Print (some of) the global statistics.
33 std::cout << "entries = " << hist.GetNEntries();
34 std::cout << ", mean = " << hist.ComputeMean();
35 std::cout << ", stddev = " << hist.ComputeStdDev();
36 std::cout << "\n";
37
38 // "Draw" the histogram with ASCII characters. The height is hard-coded to work for this tutorial.
39 for (int row = 15; row > 0; row--) {
40 auto print = [&](ROOT::Experimental::RBinIndex bin) {
41 auto value = hist.GetBinContent(bin);
42 static constexpr int Scale = 15;
43 std::cout << (value >= (row * Scale) ? '*' : ' ');
44 };
45
46 // First the underflow bin, separated by a vertical bar.
48 std::cout << '|';
49
50 // Now iterate the normal bins and print a '*' if the value is sufficiently large.
51 for (auto bin : axis.GetNormalRange()) {
52 print(bin);
53 }
54
55 // Finally the overflow bin after a separating vertical bar.
56 std::cout << '|';
58 std::cout << "\n";
59 }
60}
61
62static void FillHistogram(std::size_t t, ROOT::Experimental::RHistFillContext<int> &context)
63{
64 // Create a normal distribution with mean 10.0 and stddev 4.0.
65 std::mt19937 gen(t);
66 std::normal_distribution normal(10.0, 4.0);
67 for (std::size_t i = 0; i < 1000; i++) {
68 context.Fill(normal(gen));
69 }
70}
71
73{
74 // Create an axis and a histogram that multiple threads will fill concurrently.
75 ROOT::Experimental::RRegularAxis axis(40, {0.0, 20.0});
76 auto hist = std::make_shared<ROOT::Experimental::RHist<int>>(axis);
77
78 {
79 // Fill histogram from multiple threads.
81 std::vector<std::thread> threads;
82 for (std::size_t t = 0; t < 4; t++) {
83 threads.emplace_back([t, &filler] {
84 // Create a fill context for this thread.
85 auto context = filler.CreateFillContext();
86 FillHistogram(t, *context);
87 });
88 }
89
90 for (auto &&t : threads) {
91 t.join();
92 }
93 }
94
95 // "Draw" the histogram with ASCII characters.
96 DrawHistogram(*hist);
97}
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 .