43using namespace std::chrono_literals;
46constexpr char const *kNTupleFileName =
"ntpl017_shared_reader.root";
61std::ostream &
operator<<(std::ostream &os,
const EReadMode &
e)
64 case EReadMode::kNaive: os <<
"naive";
break;
65 case EReadMode::kInformed: os <<
"informed";
break;
72static std::atomic<int> gNEntriesDone;
73static std::atomic<int> gThreadId;
76void ProcessEntries(
ROOT::RNTupleReader &reader, EReadMode readMode,
const std::chrono::microseconds &usPerEvent,
77 std::vector<int> &sumLoadedClusters, std::vector<int> &sumUnsealedPages,
78 std::vector<int> &nClusters, std::vector<int> &nPages)
80 static std::mutex gLock;
82 const int threadId = gThreadId++;
85 for (
unsigned int i = threadId; i < reader.
GetNEntries(); i += 2) {
87 std::lock_guard<std::mutex> guard(gLock);
92 case EReadMode::kInformed: token.SetEntryNumber(i);
break;
93 case EReadMode::kNaive:
break;
94 default: std::terminate();
100 std::this_thread::sleep_for(usPerEvent);
102 const int entriesDone = gNEntriesDone++;
103 sumLoadedClusters.at(entriesDone) =
105 sumUnsealedPages[entriesDone] =
107 nClusters[entriesDone] =
109 nPages[entriesDone] =
114void Read(EReadMode readMode, std::vector<int> &sumLoadedClusters, std::vector<int> &sumUnsealedPages,
115 std::vector<int> &nClusters, std::vector<int> &nPages)
124 sumLoadedClusters.resize(
N);
125 sumUnsealedPages.resize(
N);
129 std::array<std::thread, 2> threads;
130 threads[0] = std::thread(ProcessEntries, std::ref(*reader), readMode, 100us, std::ref(sumLoadedClusters),
131 std::ref(sumUnsealedPages), std::ref(nClusters), std::ref(nPages));
132 threads[1] = std::thread(ProcessEntries, std::ref(*reader), readMode, 200us, std::ref(sumLoadedClusters),
133 std::ref(sumUnsealedPages), std::ref(nClusters), std::ref(nPages));
134 for (
auto &t : threads) {
138 std::cout <<
"Reading in mode '" << readMode <<
"':" << std::endl;
139 std::cout <<
"===========================" << std::endl;
141 std::cout <<
"===========================" << std::endl << std::endl;
147 auto ptrPoint = model->MakeField<Point2D>(
"point");
151 for (
int i = 0; i < 1000; ++i) {
155 auto prng = std::make_unique<TRandom3>();
158 ptrPoint->fX = prng->Rndm(1);
159 ptrPoint->fY = prng->Rndm(1);
164TGraph *GetGraph(
const std::vector<int> &counts)
166 auto graph =
new TGraph();
167 for (
unsigned int i = 0; i < counts.size(); ++i) {
168 graph->SetPoint(i, i, counts[i]);
170 graph->GetXaxis()->SetTitle(
"Number of processed entries");
174void ntpl017_shared_reader()
181 std::vector<int> sumLoadedClusters;
182 std::vector<int> sumUnsealedPages;
183 std::vector<int> nClusters;
184 std::vector<int> nPages;
189 gStyle->SetMarkerStyle(8);
190 TCanvas *canvas =
new TCanvas(
"",
"Shared Reader Example", 200, 10, 1500, 1000);
194 Read(EReadMode::kNaive, sumLoadedClusters, sumUnsealedPages, nClusters, nPages);
197 auto graph1 = GetGraph(sumUnsealedPages);
198 graph1->SetLineColor(
kRed);
200 auto graph2 = GetGraph(sumLoadedClusters);
201 graph2->SetLineColor(
kBlue);
202 graph2->Draw(
"SAME L");
204 auto legend =
new TLegend(0.125, 0.725, 0.625, 0.875);
205 legend->AddEntry(graph1,
"Number of decompressed pages",
"l");
206 legend->AddEntry(graph2,
"Number of loaded clusters",
"l");
214 auto graph3 = GetGraph(nPages);
215 graph3->SetMarkerColor(
kRed);
216 graph3->GetYaxis()->SetNdivisions(8);
217 graph3->GetYaxis()->SetRangeUser(-0.5, 14);
220 auto graph4 = GetGraph(nClusters);
221 graph4->SetMarkerColor(
kBlue);
222 graph4->Draw(
"SAME P");
224 legend =
new TLegend(0.35, 0.725, 0.85, 0.875);
225 legend->AddEntry(graph3,
"Number of currently cached pages",
"p");
226 legend->AddEntry(graph4,
"Number of currently cached clusters",
"p");
231 Read(EReadMode::kInformed, sumLoadedClusters, sumUnsealedPages, nClusters, nPages);
234 auto graph5 = GetGraph(sumUnsealedPages);
235 graph5->SetLineColor(
kRed);
238 auto graph6 = GetGraph(sumLoadedClusters);
239 graph6->SetLineColor(
kBlue);
240 graph6->Draw(
"SAME L");
247 auto graph7 = GetGraph(nPages);
248 graph7->SetMarkerColor(
kRed);
249 graph7->GetYaxis()->SetNdivisions(8);
250 graph7->GetYaxis()->SetRangeUser(-0.5, 14);
253 auto graph8 = GetGraph(nClusters);
254 graph8->SetMarkerColor(
kBlue);
255 graph8->Draw(
"SAME P");
TBuffer & operator<<(TBuffer &buf, const Tmpl *obj)
const RNTuplePerfCounter * GetCounter(std::string_view name) const
Searches this object and all the observed sub metrics.
void Print(std::ostream &output, const std::string &prefix="") const
virtual std::int64_t GetValueAsInt() const =0
static std::unique_ptr< RNTupleModel > Create()
Reads RNTuple data from storage.
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.
const Experimental::Detail::RNTupleMetrics & GetMetrics() const
ROOT::NTupleSize_t GetNEntries() const
Returns the number of entries in this RNTuple.
RActiveEntryToken CreateActiveEntryToken()
Create a new active entry token, which will not be bound to any entry number initially.
void EnableMetrics()
Enable performance measurements (decompression time, bytes read from storage, etc....
void LoadEntry(ROOT::NTupleSize_t index)
Fills the default entry of the model.
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())
Creates an RNTupleWriter backed by storage, overwriting it if one with the same URI exists.
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
TLatex * DrawLatexNDC(Double_t x, Double_t y, const char *text)
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...