Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
ntpl009_parallelWriter.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Example of multi-threaded writes using RNTupleParallelWriter. Adapted from the ntpl007_mtFill tutorial.
5///
6/// \macro_image
7/// \macro_code
8///
9/// \date Feburary 2024
10/// \author The ROOT Team
11
12// NOTE: The RNTuple classes are experimental at this point.
13// Functionality and interface are still subject to changes.
14
16#include <ROOT/RNTupleModel.hxx>
19
20#include <TCanvas.h>
21#include <TH1F.h>
22#include <TH2F.h>
23#include <TRandom.h>
24#include <TRandom3.h>
25#include <TStyle.h>
26#include <TSystem.h>
27
28#include <atomic>
29#include <memory>
30#include <mutex>
31#include <thread>
32#include <vector>
33#include <utility>
34
35// Import classes from experimental namespace for the time being
40
41// Where to store the ntuple of this example
42constexpr char const *kNTupleFileName = "ntpl009_parallelWriter.root";
43
44// Number of parallel threads to fill the ntuple
45constexpr int kNWriterThreads = 4;
46
47// Number of events to generate is kNEventsPerThread * kNWriterThreads
48constexpr int kNEventsPerThread = 25000;
49
50// Thread function to generate and write events
51void FillData(RNTupleParallelWriter *writer)
52{
53 static std::atomic<std::uint32_t> gThreadId;
54 const auto threadId = ++gThreadId;
55
56 auto prng = std::make_unique<TRandom3>();
57 prng->SetSeed();
58
59 auto fillContext = writer->CreateFillContext();
60 auto entry = fillContext->CreateEntry();
61
62 auto id = entry->GetPtr<std::uint32_t>("id");
63 *id = threadId;
64 auto vpx = entry->GetPtr<std::vector<float>>("vpx");
65 auto vpy = entry->GetPtr<std::vector<float>>("vpy");
66 auto vpz = entry->GetPtr<std::vector<float>>("vpz");
67
68 for (int i = 0; i < kNEventsPerThread; i++) {
69 vpx->clear();
70 vpy->clear();
71 vpz->clear();
72
73 int npx = static_cast<int>(prng->Rndm(1) * 15);
74 // Set the field data for the current event
75 for (int j = 0; j < npx; ++j) {
76 float px, py, pz;
77 prng->Rannor(px, py);
78 pz = px * px + py * py;
79
80 vpx->emplace_back(px);
81 vpy->emplace_back(py);
82 vpz->emplace_back(pz);
83 }
84
85 fillContext->Fill(*entry);
86 }
87}
88
89// Generate kNEvents with multiple threads in kNTupleFileName
90void Write()
91{
92 // Create the data model
93 auto model = RNTupleModel::CreateBare();
94 model->MakeField<std::uint32_t>("id");
95 model->MakeField<std::vector<float>>("vpx");
96 model->MakeField<std::vector<float>>("vpy");
97 model->MakeField<std::vector<float>>("vpz");
98
99 // Create RNTupleWriteOptions to make the writing commit multiple clusters (so that "Entry Id vs Thread Id" shows the
100 // interleaved clusters).
101 RNTupleWriteOptions options;
102 options.SetApproxZippedClusterSize(1024 * 1024);
103
104 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
105 auto writer = RNTupleParallelWriter::Recreate(std::move(model), "NTuple", kNTupleFileName, options);
106
107 std::vector<std::thread> threads;
108 for (int i = 0; i < kNWriterThreads; ++i)
109 threads.emplace_back(FillData, writer.get());
110 for (int i = 0; i < kNWriterThreads; ++i)
111 threads[i].join();
112
113 // The writer unique pointer goes out of scope here. On destruction, the writer flushes unwritten data to disk
114 // and closes the attached ROOT file.
115}
116
117// For all of the events, histogram only one of the written vectors
118void Read()
119{
120 auto reader = RNTupleReader::Open("NTuple", kNTupleFileName);
121 auto viewVpx = reader->GetView<float>("vpx._0");
122
123 gStyle->SetOptStat(0);
124
125 TCanvas *c1 = new TCanvas("c2", "Multi-Threaded Filling Example", 200, 10, 1500, 500);
126 c1->Divide(2, 1);
127
128 c1->cd(1);
129 TH1F h("h", "This is the px distribution", 100, -4, 4);
130 h.SetFillColor(48);
131 // Iterate through all values of vpx in all events
132 for (auto i : viewVpx.GetFieldRange())
133 h.Fill(viewVpx(i));
134 // Prevent the histogram from disappearing
135 h.DrawCopy();
136
137 c1->cd(2);
138 auto nEvents = reader->GetNEntries();
139 auto viewId = reader->GetView<std::uint32_t>("id");
140 TH2F hFillSequence("", "Entry Id vs Thread Id;Entry Sequence Number;Filling Thread", 100, 0, nEvents, 100, 0,
141 kNWriterThreads + 1);
142 for (auto i : reader->GetEntryRange())
143 hFillSequence.Fill(i, viewId(i));
144 hFillSequence.DrawCopy();
145}
146
148{
149 Write();
150 Read();
151}
#define h(i)
Definition RSha256.hxx:106
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
The RNTupleModel encapulates the schema of an ntuple.
A writer to fill an RNTuple from multiple contexts.
An RNTuple that is used to read data from storage.
Common user-tunable settings for storing ntuples.
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:622
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1640
return c1
Definition legend1.C:41
RNTupleGlobalRange GetFieldRange(const RFieldBase &field, const RPageSource &pageSource)
Helper to get the iteration space of the given field that needs to be connected to the given page sou...
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.