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 <TH2F.h>
27#include <TRandom.h>
28#include <TRandom3.h>
29#include <TStyle.h>
30#include <TSystem.h>
31
32#include <atomic>
33#include <memory>
34#include <mutex>
35#include <thread>
36#include <vector>
37#include <utility>
38
39// Import classes from experimental namespace for the time being
40using REntry = ROOT::Experimental::REntry;
41using RNTupleModel = ROOT::Experimental::RNTupleModel;
42using RNTupleReader = ROOT::Experimental::RNTupleReader;
43using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
44
45// Where to store the ntuple of this example
46constexpr char const *kNTupleFileName = "ntpl007_mtFill.root";
47
48// Number of parallel threads to fill the ntuple
49constexpr int kNWriterThreads = 4;
50
51// Number of events to generate is kNEventsPerThread * kNWriterThreads
52constexpr int kNEventsPerThread = 25000;
53
54// Thread function to generate and write events
55void FillData(std::unique_ptr<REntry> entry, RNTupleWriter *ntuple) {
56 // Protect the ntuple->Fill() call
57 static std::mutex gLock;
58
59 static std::atomic<std::uint32_t> gThreadId;
60 const auto threadId = ++gThreadId;
61
62 auto prng = std::make_unique<TRandom3>();
63 prng->SetSeed();
64
65 auto id = entry->Get<std::uint32_t>("id");
66 auto vpx = entry->Get<std::vector<float>>("vpx");
67 auto vpy = entry->Get<std::vector<float>>("vpy");
68 auto vpz = entry->Get<std::vector<float>>("vpz");
69
70 for (int i = 0; i < kNEventsPerThread; i++) {
71 vpx->clear();
72 vpy->clear();
73 vpz->clear();
74 *id = threadId;
75
76 int npx = static_cast<int>(prng->Rndm(1) * 15);
77 // Set the field data for the current event
78 for (int j = 0; j < npx; ++j) {
79 float px, py, pz;
80 prng->Rannor(px, py);
81 pz = px*px + py*py;
82
83 vpx->emplace_back(px);
84 vpy->emplace_back(py);
85 vpz->emplace_back(pz);
86 }
87
88 std::lock_guard<std::mutex> guard(gLock);
89 ntuple->Fill(*entry);
90 }
91}
92
93// Generate kNEvents with multiple threads in kNTupleFileName
94void Write()
95{
96 // Create the data model
97 auto model = RNTupleModel::Create();
98 model->MakeField<std::uint32_t>("id");
99 model->MakeField<std::vector<float>>("vpx");
100 model->MakeField<std::vector<float>>("vpy");
101 model->MakeField<std::vector<float>>("vpz");
102
103 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
104 auto ntuple = RNTupleWriter::Recreate(std::move(model), "NTuple", kNTupleFileName);
105
106 std::vector<std::unique_ptr<REntry>> entries;
107 std::vector<std::thread> threads;
108 for (int i = 0; i < kNWriterThreads; ++i)
109 entries.emplace_back(ntuple->CreateEntry());
110 for (int i = 0; i < kNWriterThreads; ++i)
111 threads.emplace_back(FillData, std::move(entries[i]), ntuple.get());
112 for (int i = 0; i < kNWriterThreads; ++i)
113 threads[i].join();
114
115 // The ntuple unique pointer goes out of scope here. On destruction, the ntuple flushes unwritten data to disk
116 // and closes the attached ROOT file.
117}
118
119
120// For all of the events, histogram only one of the written vectors
121void Read()
122{
123 auto ntuple = RNTupleReader::Open("NTuple", kNTupleFileName);
124 // TODO(jblomer): the "inner name" of the vector should become "vpx._0"
125 auto viewVpx = ntuple->GetView<float>("vpx._0");
126
127 gStyle->SetOptStat(0);
128
129 TCanvas *c1 = new TCanvas("c2", "Multi-Threaded Filling Example", 200, 10, 1500, 500);
130 c1->Divide(2, 1);
131
132 c1->cd(1);
133 TH1F h("h", "This is the px distribution", 100, -4, 4);
134 h.SetFillColor(48);
135 // Iterate through all values of vpx in all events
136 for (auto i : viewVpx.GetFieldRange())
137 h.Fill(viewVpx(i));
138 // Prevent the histogram from disappearing
139 h.DrawCopy();
140
141 c1->cd(2);
142 auto nEvents = ntuple->GetNEntries();
143 auto viewId = ntuple->GetView<std::uint32_t>("id");
144 TH2F hFillSequence("","Entry Id vs Thread Id;Entry Sequence Number;Filling Thread",
145 100, 0, nEvents, 100, 0, kNWriterThreads);
146 for (auto i : ntuple->GetEntryRange())
147 hFillSequence.Fill(i, viewId(i));
148 hFillSequence.DrawCopy();
149}
150
151
152void ntpl007_mtFill()
153{
154 Write();
155 Read();
156}
#define h(i)
Definition RSha256.hxx:106
#define R__LOAD_LIBRARY(LIBRARY)
Definition Rtypes.h:491
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:42
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
Definition RNTuple.hxx:101
An RNTuple that gets filled with entries (data) and writes them to storage.
Definition RNTuple.hxx:358
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:258
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.