Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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#include <ROOT/REntry.hxx>
13#include <ROOT/RNTupleModel.hxx>
16
17#include <TCanvas.h>
18#include <TH1F.h>
19#include <TH2F.h>
20#include <TRandom.h>
21#include <TRandom3.h>
22#include <TStyle.h>
23#include <TSystem.h>
24
25#include <atomic>
26#include <memory>
27#include <mutex>
28#include <thread>
29#include <vector>
30#include <utility>
31
32// Where to store the ntuple of this example
33constexpr char const *kNTupleFileName = "ntpl007_mtFill.root";
34
35// Number of parallel threads to fill the ntuple
36constexpr int kNWriterThreads = 4;
37
38// Number of events to generate is kNEventsPerThread * kNWriterThreads
39constexpr int kNEventsPerThread = 25000;
40
41// Thread function to generate and write events
42void FillData(std::unique_ptr<ROOT::REntry> entry, ROOT::RNTupleWriter *writer)
43{
44 // Protect the writer->Fill() call
45 static std::mutex gLock;
46
47 static std::atomic<std::uint32_t> gThreadId;
48 const auto threadId = ++gThreadId;
49
50 auto prng = std::make_unique<TRandom3>();
51 prng->SetSeed();
52
53 auto id = entry->GetPtr<std::uint32_t>("id");
54 *id = threadId;
55 auto vpx = entry->GetPtr<std::vector<float>>("vpx");
56 auto vpy = entry->GetPtr<std::vector<float>>("vpy");
57 auto vpz = entry->GetPtr<std::vector<float>>("vpz");
58
59 for (int i = 0; i < kNEventsPerThread; i++) {
60 vpx->clear();
61 vpy->clear();
62 vpz->clear();
63
64 int npx = static_cast<int>(prng->Rndm(1) * 15);
65 // Set the field data for the current event
66 for (int j = 0; j < npx; ++j) {
67 float px, py, pz;
68 prng->Rannor(px, py);
69 pz = px*px + py*py;
70
71 vpx->emplace_back(px);
72 vpy->emplace_back(py);
73 vpz->emplace_back(pz);
74 }
75
76 std::lock_guard<std::mutex> guard(gLock);
77 writer->Fill(*entry);
78 }
79}
80
81// Generate kNEvents with multiple threads in kNTupleFileName
82void Write()
83{
84 // Create the data model
85 auto model = ROOT::RNTupleModel::Create();
86 model->MakeField<std::uint32_t>("id");
87 model->MakeField<std::vector<float>>("vpx");
88 model->MakeField<std::vector<float>>("vpy");
89 model->MakeField<std::vector<float>>("vpz");
90
91 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
92 auto writer = ROOT::RNTupleWriter::Recreate(std::move(model), "NTuple", kNTupleFileName);
93
94 std::vector<std::unique_ptr<ROOT::REntry>> entries;
95 std::vector<std::thread> threads;
96 for (int i = 0; i < kNWriterThreads; ++i)
97 entries.emplace_back(writer->CreateEntry());
98 for (int i = 0; i < kNWriterThreads; ++i)
99 threads.emplace_back(FillData, std::move(entries[i]), writer.get());
100 for (int i = 0; i < kNWriterThreads; ++i)
101 threads[i].join();
102
103 // The writer unique pointer goes out of scope here. On destruction, the writer flushes unwritten data to disk
104 // and closes the attached ROOT file.
105}
106
107
108// For all of the events, histogram only one of the written vectors
109void Read()
110{
112 auto viewVpx = reader->GetView<float>("vpx._0");
113
114 gStyle->SetOptStat(0);
115
116 TCanvas *c1 = new TCanvas("c2", "Multi-Threaded Filling Example", 200, 10, 1500, 500);
117 c1->Divide(2, 1);
118
119 c1->cd(1);
120 TH1F h("h", "This is the px distribution", 100, -4, 4);
121 h.SetFillColor(48);
122 // Iterate through all values of vpx in all events
123 for (auto i : viewVpx.GetFieldRange())
124 h.Fill(viewVpx(i));
125 // Prevent the histogram from disappearing
126 h.DrawCopy();
127
128 c1->cd(2);
129 auto nEvents = reader->GetNEntries();
130 auto viewId = reader->GetView<std::uint32_t>("id");
131 TH2F hFillSequence("", "Entry Id vs Thread Id;Entry Sequence Number;Filling Thread", 100, 0, nEvents, 100, 0,
132 kNWriterThreads + 1);
133 for (auto i : reader->GetEntryRange())
134 hFillSequence.Fill(i, viewId(i));
135 hFillSequence.DrawCopy();
136}
137
138
139void ntpl007_mtFill()
140{
141 Write();
142 Read();
143}
#define h(i)
Definition RSha256.hxx:106
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
static std::unique_ptr< RNTupleModel > Create()
static std::unique_ptr< RNTupleReader > Open(std::string_view ntupleName, std::string_view storage, const ROOT::RNTupleReadOptions &options=ROOT::RNTupleReadOptions())
Open an RNTuple for reading.
An RNTuple that gets filled with entries (data) and writes them to storage.
static std::unique_ptr< RNTupleWriter > Recreate(std::unique_ptr< ROOT::RNTupleModel > model, std::string_view ntupleName, std::string_view storage, const ROOT::RNTupleWriteOptions &options=ROOT::RNTupleWriteOptions())
Throws an exception if the model is null.
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:650
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:1642
return c1
Definition legend1.C:41
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.
ROOT::RNTupleGlobalRange GetFieldRange(const ROOT::RFieldBase &field, const ROOT::Internal::RPageSource &pageSource)
Helper to get the iteration space of the given field that needs to be connected to the given page sou...