33constexpr char const *kNTupleFileName =
"ntpl009_parallelWriter.root";
36constexpr int kNWriterThreads = 4;
39constexpr int kNEventsPerThread = 25000;
44 static std::atomic<std::uint32_t> gThreadId;
45 const auto threadId = ++gThreadId;
47 auto prng = std::make_unique<TRandom3>();
50 auto fillContext =
writer->CreateFillContext();
51 auto entry = fillContext->CreateEntry();
53 auto id = entry->GetPtr<std::uint32_t>(
"id");
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");
59 for (
int i = 0; i < kNEventsPerThread; i++) {
64 int npx =
static_cast<int>(prng->Rndm(1) * 15);
66 for (
int j = 0; j < npx; ++j) {
69 pz = px * px + py * py;
71 vpx->emplace_back(px);
72 vpy->emplace_back(py);
73 vpz->emplace_back(pz);
76 fillContext->Fill(*entry);
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");
93 options.SetApproxZippedClusterSize(1024 * 1024);
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)
112 auto viewVpx = reader->GetView<
float>(
"vpx._0");
116 TCanvas *
c1 =
new TCanvas(
"c2",
"Multi-Threaded Filling Example", 200, 10, 1500, 500);
120 TH1F h(
"h",
"This is the px distribution", 100, -4, 4);
123 for (
auto i : viewVpx.GetFieldRange())
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();
138void ntpl009_parallelWriter()
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.
1-D histogram with a float per channel (see TH1 documentation)
2-D histogram with a float per channel (see TH1 documentation)
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.