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 and interface are still subject to changes.
14
15#include <ROOT/REntry.hxx>
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 = "ntpl007_mtFill.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(std::unique_ptr<REntry> entry, RNTupleWriter *writer) {
52 // Protect the ntuple->Fill() call
53 static std::mutex gLock;
54
55 static std::atomic<std::uint32_t> gThreadId;
56 const auto threadId = ++gThreadId;
57
58 auto prng = std::make_unique<TRandom3>();
59 prng->SetSeed();
60
61 auto id = entry->GetPtr<std::uint32_t>("id");
62 *id = threadId;
63 auto vpx = entry->GetPtr<std::vector<float>>("vpx");
64 auto vpy = entry->GetPtr<std::vector<float>>("vpy");
65 auto vpz = entry->GetPtr<std::vector<float>>("vpz");
66
67 for (int i = 0; i < kNEventsPerThread; i++) {
68 vpx->clear();
69 vpy->clear();
70 vpz->clear();
71
72 int npx = static_cast<int>(prng->Rndm(1) * 15);
73 // Set the field data for the current event
74 for (int j = 0; j < npx; ++j) {
75 float px, py, pz;
76 prng->Rannor(px, py);
77 pz = px*px + py*py;
78
79 vpx->emplace_back(px);
80 vpy->emplace_back(py);
81 vpz->emplace_back(pz);
82 }
83
84 std::lock_guard<std::mutex> guard(gLock);
85 writer->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::Create();
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 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
100 auto writer = RNTupleWriter::Recreate(std::move(model), "NTuple", kNTupleFileName);
101
102 std::vector<std::unique_ptr<REntry>> entries;
103 std::vector<std::thread> threads;
104 for (int i = 0; i < kNWriterThreads; ++i)
105 entries.emplace_back(writer->CreateEntry());
106 for (int i = 0; i < kNWriterThreads; ++i)
107 threads.emplace_back(FillData, std::move(entries[i]), writer.get());
108 for (int i = 0; i < kNWriterThreads; ++i)
109 threads[i].join();
110
111 // The writer unique pointer goes out of scope here. On destruction, the writer flushes unwritten data to disk
112 // and closes the attached ROOT file.
113}
114
115
116// For all of the events, histogram only one of the written vectors
117void Read()
118{
119 auto reader = RNTupleReader::Open("NTuple", kNTupleFileName);
120 auto viewVpx = reader->GetView<float>("vpx._0");
121
122 gStyle->SetOptStat(0);
123
124 TCanvas *c1 = new TCanvas("c2", "Multi-Threaded Filling Example", 200, 10, 1500, 500);
125 c1->Divide(2, 1);
126
127 c1->cd(1);
128 TH1F h("h", "This is the px distribution", 100, -4, 4);
129 h.SetFillColor(48);
130 // Iterate through all values of vpx in all events
131 for (auto i : viewVpx.GetFieldRange())
132 h.Fill(viewVpx(i));
133 // Prevent the histogram from disappearing
134 h.DrawCopy();
135
136 c1->cd(2);
137 auto nEvents = reader->GetNEntries();
138 auto viewId = reader->GetView<std::uint32_t>("id");
139 TH2F hFillSequence("", "Entry Id vs Thread Id;Entry Sequence Number;Filling Thread", 100, 0, nEvents, 100, 0,
140 kNWriterThreads + 1);
141 for (auto i : reader->GetEntryRange())
142 hFillSequence.Fill(i, viewId(i));
143 hFillSequence.DrawCopy();
144}
145
146
147void ntpl007_mtFill()
148{
149 Write();
150 Read();
151}
#define h(i)
Definition RSha256.hxx:106
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
The REntry is a collection of values in an ntuple corresponding to a complete row in the data set.
Definition REntry.hxx:51
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:623
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:308
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.