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 multuple 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// Until C++ runtime modules are universally used, we explicitly load the ntuple library. Otherwise
17// triggering autoloading from the use of templated types would require an exhaustive enumeration
18// of "all" template instances in the LinkDef file.
19R__LOAD_LIBRARY(ROOTNTuple)
20
21#include <ROOT/RNTuple.hxx>
22#include <ROOT/RNTupleModel.hxx>
23
24#include <TCanvas.h>
25#include <TH1F.h>
26#include <TRandom.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
37using REntry = ROOT::Experimental::REntry;
38using RNTupleModel = ROOT::Experimental::RNTupleModel;
39using RNTupleReader = ROOT::Experimental::RNTupleReader;
40using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
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 *ntuple) {
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->Get<std::uint32_t>("id");
63 auto vpx = entry->Get<std::vector<float>>("vpx");
64 auto vpy = entry->Get<std::vector<float>>("vpy");
65 auto vpz = entry->Get<std::vector<float>>("vpz");
66
67 for (int i = 0; i < kNEventsPerThread; i++) {
68 vpx->clear();
69 vpy->clear();
70 vpz->clear();
71 *id = threadId;
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 ntuple->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 ntuple = 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(ntuple->CreateEntry());
107 for (int i = 0; i < kNWriterThreads; ++i)
108 threads.emplace_back(FillData, std::move(entries[i]), ntuple.get());
109 for (int i = 0; i < kNWriterThreads; ++i)
110 threads[i].join();
111
112 // The ntuple unique pointer goes out of scope here. On destruction, the ntuple 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 ntuple = RNTupleReader::Open("NTuple", kNTupleFileName);
121 // TODO(jblomer): the "inner name" of the vector should become "vpx._0"
122 auto viewVpx = ntuple->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 = ntuple->GetNEntries();
140 auto viewId = ntuple->GetView<std::uint32_t>("id");
141 TH2F hFillSequence("","Entry Id vs Thread Id;Entry Sequence Number;Filling Thread",
142 100, 0, nEvents, 100, 0, kNWriterThreads);
143 for (auto i : ntuple->GetEntryRange())
144 hFillSequence.Fill(i, viewId(i));
145 hFillSequence.DrawCopy();
146}
147
148
149void ntpl007_mtFill()
150{
151 Write();
152 Read();
153}
#define h(i)
Definition RSha256.hxx:106
#define R__LOAD_LIBRARY(LIBRARY)
Definition Rtypes.h:491
R__EXTERN TStyle * gStyle
Definition TStyle.h:414
The REntry is a collection of values in an ntuple corresponding to a complete row in the data set.
Definition REntry.hxx:43
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
Definition RNTuple.hxx:110
An RNTuple that gets filled with entries (data) and writes them to storage.
Definition RNTuple.hxx:368
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:257
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:1589
return c1
Definition legend1.C:41
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.