Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl007_mtFill.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Example of multi-threaded writes using multiple REntry objects
5///
6/// \macro_image
7/// \macro_code
8///
9/// \date July 2021
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
16#include <ROOT/REntry.hxx>
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 = "ntpl007_mtFill.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(std::unique_ptr<REntry> entry, RNTupleWriter *writer) {
53 // Protect the ntuple->Fill() call
54 static std::mutex gLock;
55
56 static std::atomic<std::uint32_t> gThreadId;
57 const auto threadId = ++gThreadId;
58
59 auto prng = std::make_unique<TRandom3>();
60 prng->SetSeed();
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 std::lock_guard<std::mutex> guard(gLock);
86 writer->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::Create();
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 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
101 auto writer = RNTupleWriter::Recreate(std::move(model), "NTuple", kNTupleFileName);
102
103 std::vector<std::unique_ptr<REntry>> entries;
104 std::vector<std::thread> threads;
105 for (int i = 0; i < kNWriterThreads; ++i)
106 entries.emplace_back(writer->CreateEntry());
107 for (int i = 0; i < kNWriterThreads; ++i)
108 threads.emplace_back(FillData, std::move(entries[i]), writer.get());
109 for (int i = 0; i < kNWriterThreads; ++i)
110 threads[i].join();
111
112 // The writer unique pointer goes out of scope here. On destruction, the writer flushes unwritten data to disk
113 // and closes the attached ROOT file.
114}
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
147
148void ntpl007_mtFill()
149{
150 Write();
151 Read();
152}
#define h(i)
Definition RSha256.hxx:106
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
The REntry is a collection of values in an ntuple corresponding to a complete row in the data set.
Definition REntry.hxx:45
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
An RNTuple that gets filled with entries (data) and writes them to storage.
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
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:1636
return c1
Definition legend1.C:41
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.