Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl009_parallelWriter.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Example of multi-threaded writes using RNTupleParallelWriter. Adapted from the ntpl007_mtFill tutorial.
5///
6/// \macro_image
7/// \macro_code
8///
9/// \date Feburary 2024
10/// \author The ROOT Team
11
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 = "ntpl009_parallelWriter.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
43{
44 static std::atomic<std::uint32_t> gThreadId;
45 const auto threadId = ++gThreadId;
46
47 auto prng = std::make_unique<TRandom3>();
48 prng->SetSeed();
49
50 auto fillContext = writer->CreateFillContext();
51 auto entry = fillContext->CreateEntry();
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 fillContext->Fill(*entry);
77 }
78}
79
80// Generate kNEvents with multiple threads in kNTupleFileName
81void Write()
82{
83 // Create the data model
84 auto model = ROOT::RNTupleModel::CreateBare();
85 model->MakeField<std::uint32_t>("id");
86 model->MakeField<std::vector<float>>("vpx");
87 model->MakeField<std::vector<float>>("vpy");
88 model->MakeField<std::vector<float>>("vpz");
89
90 // Create RNTupleWriteOptions to make the writing commit multiple clusters (so that "Entry Id vs Thread Id" shows the
91 // interleaved clusters).
93 options.SetApproxZippedClusterSize(1024 * 1024);
94
95 // We hand-over the data model to a newly created ntuple of name "NTuple", stored in kNTupleFileName
96 auto writer = ROOT::RNTupleParallelWriter::Recreate(std::move(model), "NTuple", kNTupleFileName, options);
97
98 std::vector<std::thread> threads;
99 for (int i = 0; i < kNWriterThreads; ++i)
100 threads.emplace_back(FillData, writer.get());
101 for (int i = 0; i < kNWriterThreads; ++i)
102 threads[i].join();
103
104 // The writer unique pointer goes out of scope here. On destruction, the writer flushes unwritten data to disk
105 // and closes the attached ROOT file.
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())
135 hFillSequence.DrawCopy();
136}
137
139{
140 Write();
141 Read();
142}
#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.
Int_t Fill(Double_t x) override
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
static std::unique_ptr< RNTupleModel > CreateBare()
Creates a "bare model", i.e. an RNTupleModel with no default entry.
A writer to fill an RNTuple from multiple contexts.
static std::unique_ptr< RNTupleParallelWriter > Recreate(std::unique_ptr< ROOT::RNTupleModel > model, std::string_view ntupleName, std::string_view storage, const ROOT::RNTupleWriteOptions &options=ROOT::RNTupleWriteOptions())
Recreate a new file and return a writer to write an RNTuple.
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.
Common user-tunable settings for storing RNTuples.
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:879
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:345
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:1641
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...