#include <ROOT/RNTuple.hxx>
#include <atomic>
#include <memory>
#include <mutex>
#include <thread>
#include <vector>
#include <utility>
constexpr char const *kNTupleFileName = "ntpl007_mtFill.root";
constexpr int kNWriterThreads = 4;
constexpr int kNEventsPerThread = 25000;
void FillData(std::unique_ptr<REntry> entry, RNTupleWriter *ntuple) {
static std::mutex gLock;
static std::atomic<std::uint32_t> gThreadId;
const auto threadId = ++gThreadId;
auto prng = std::make_unique<TRandom3>();
prng->SetSeed();
auto id = entry->Get<std::uint32_t>("id");
auto vpx = entry->Get<std::vector<float>>("vpx");
auto vpy = entry->Get<std::vector<float>>("vpy");
auto vpz = entry->Get<std::vector<float>>("vpz");
for (int i = 0; i < kNEventsPerThread; i++) {
vpx->clear();
vpy->clear();
vpz->clear();
*id = threadId;
int npx = static_cast<int>(prng->Rndm(1) * 15);
for (int j = 0; j < npx; ++j) {
float px, py, pz;
prng->Rannor(px, py);
pz = px*px + py*py;
vpx->emplace_back(px);
vpy->emplace_back(py);
vpz->emplace_back(pz);
}
std::lock_guard<std::mutex> guard(gLock);
ntuple->Fill(*entry);
}
}
void Write()
{
auto model = RNTupleModel::Create();
model->MakeField<std::uint32_t>("id");
model->MakeField<std::vector<float>>("vpx");
model->MakeField<std::vector<float>>("vpy");
model->MakeField<std::vector<float>>("vpz");
auto ntuple = RNTupleWriter::Recreate(std::move(model), "NTuple", kNTupleFileName);
std::vector<std::unique_ptr<REntry>> entries;
std::vector<std::thread> threads;
for (int i = 0; i < kNWriterThreads; ++i)
entries.emplace_back(ntuple->CreateEntry());
for (int i = 0; i < kNWriterThreads; ++i)
threads.emplace_back(FillData, std::move(entries[i]), ntuple.get());
for (int i = 0; i < kNWriterThreads; ++i)
threads[i].join();
}
void Read()
{
auto ntuple = RNTupleReader::Open("NTuple", kNTupleFileName);
auto viewVpx = ntuple->GetView<float>("vpx._0");
TCanvas *
c1 =
new TCanvas(
"c2",
"Multi-Threaded Filling Example", 200, 10, 1500, 500);
TH1F h(
"h",
"This is the px distribution", 100, -4, 4);
for (auto i : viewVpx.GetFieldRange())
auto nEvents = ntuple->GetNEntries();
auto viewId = ntuple->GetView<std::uint32_t>("id");
TH2F hFillSequence(
"",
"Entry Id vs Thread Id;Entry Sequence Number;Filling Thread",
100, 0, nEvents, 100, 0, kNWriterThreads);
for (auto i : ntuple->GetEntryRange())
hFillSequence.Fill(i, viewId(i));
hFillSequence.DrawCopy();
}
void ntpl007_mtFill()
{
Write();
Read();
}
#define R__LOAD_LIBRARY(LIBRARY)
R__EXTERN TStyle * gStyle
The REntry is a collection of values in an ntuple corresponding to a complete row in the data set.
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
An RNTuple that gets filled with entries (data) and writes them to storage.
1-D histogram with a float per channel (see TH1 documentation)}
2-D histogram with a float per channel (see TH1 documentation)}
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...