Logo ROOT  
Reference Guide
 
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, interface, and data format is still subject to changes.
14// Do not use for real data!
15
17#include <ROOT/RNTupleModel.hxx>
20
21#include <TCanvas.h>
22#include <TH1F.h>
23#include <TH2F.h>
24#include <TRandom.h>
25#include <TRandom3.h>
26#include <TStyle.h>
27#include <TSystem.h>
28
29#include <atomic>
30#include <memory>
31#include <mutex>
32#include <thread>
33#include <vector>
34#include <utility>
35
36// Import classes from experimental namespace for the time being
41
42// Where to store the ntuple of this example
43constexpr char const *kNTupleFileName = "ntpl009_parallelWriter.root";
44
45// Number of parallel threads to fill the ntuple
46constexpr int kNWriterThreads = 4;
47
48// Number of events to generate is kNEventsPerThread * kNWriterThreads
49constexpr int kNEventsPerThread = 25000;
50
51// Thread function to generate and write events
52void FillData(RNTupleParallelWriter *writer)
53{
54 static std::atomic<std::uint32_t> gThreadId;
55 const auto threadId = ++gThreadId;
56
57 auto prng = std::make_unique<TRandom3>();
58 prng->SetSeed();
59
60 auto fillContext = writer->CreateFillContext();
61 auto entry = fillContext->CreateEntry();
62
63 auto id = entry->GetPtr<std::uint32_t>("id");
64 *id = threadId;
65 auto vpx = entry->GetPtr<std::vector<float>>("vpx");
66 auto vpy = entry->GetPtr<std::vector<float>>("vpy");
67 auto vpz = entry->GetPtr<std::vector<float>>("vpz");
68
69 for (int i = 0; i < kNEventsPerThread; i++) {
70 vpx->clear();
71 vpy->clear();
72 vpz->clear();
73
74 int npx = static_cast<int>(prng->Rndm(1) * 15);
75 // Set the field data for the current event
76 for (int j = 0; j < npx; ++j) {
77 float px, py, pz;
78 prng->Rannor(px, py);
79 pz = px * px + py * py;
80
81 vpx->emplace_back(px);
82 vpy->emplace_back(py);
83 vpz->emplace_back(pz);
84 }
85
86 fillContext->Fill(*entry);
87 }
88}
89
90// Generate kNEvents with multiple threads in kNTupleFileName
91void Write()
92{
93 // Create the data model
94 auto model = RNTupleModel::CreateBare();
95 model->MakeField<std::uint32_t>("id");
96 model->MakeField<std::vector<float>>("vpx");
97 model->MakeField<std::vector<float>>("vpy");
98 model->MakeField<std::vector<float>>("vpz");
99
100 // Create RNTupleWriteOptions to make the writing commit multiple clusters (so that "Entry Id vs Thread Id" shows the
101 // interleaved clusters).
102 RNTupleWriteOptions options;
103 options.SetApproxZippedClusterSize(1024 * 1024);
104
105 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
106 auto writer = RNTupleParallelWriter::Recreate(std::move(model), "NTuple", kNTupleFileName, options);
107
108 std::vector<std::thread> threads;
109 for (int i = 0; i < kNWriterThreads; ++i)
110 threads.emplace_back(FillData, writer.get());
111 for (int i = 0; i < kNWriterThreads; ++i)
112 threads[i].join();
113
114 // The writer unique pointer goes out of scope here. On destruction, the writer flushes unwritten data to disk
115 // and closes the attached ROOT file.
116}
117
118// For all of the events, histogram only one of the written vectors
119void Read()
120{
121 auto reader = RNTupleReader::Open("NTuple", kNTupleFileName);
122 auto viewVpx = reader->GetView<float>("vpx._0");
123
124 gStyle->SetOptStat(0);
125
126 TCanvas *c1 = new TCanvas("c2", "Multi-Threaded Filling Example", 200, 10, 1500, 500);
127 c1->Divide(2, 1);
128
129 c1->cd(1);
130 TH1F h("h", "This is the px distribution", 100, -4, 4);
131 h.SetFillColor(48);
132 // Iterate through all values of vpx in all events
133 for (auto i : viewVpx.GetFieldRange())
134 h.Fill(viewVpx(i));
135 // Prevent the histogram from disappearing
136 h.DrawCopy();
137
138 c1->cd(2);
139 auto nEvents = reader->GetNEntries();
140 auto viewId = reader->GetView<std::uint32_t>("id");
141 TH2F hFillSequence("", "Entry Id vs Thread Id;Entry Sequence Number;Filling Thread", 100, 0, nEvents, 100, 0,
142 kNWriterThreads + 1);
143 for (auto i : reader->GetEntryRange())
144 hFillSequence.Fill(i, viewId(i));
145 hFillSequence.DrawCopy();
146}
147
148void ntpl009_parallelWriter()
149{
150 Write();
151 Read();
152}
#define h(i)
Definition RSha256.hxx:106
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
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.