Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RNTupleReader.cxx
Go to the documentation of this file.
1/// \file RNTupleReader.cxx
2/// \ingroup NTuple
3/// \author Jakob Blomer <jblomer@cern.ch>
4/// \date 2024-02-20
5
6/*************************************************************************
7 * Copyright (C) 1995-2024, Rene Brun and Fons Rademakers. *
8 * All rights reserved. *
9 * *
10 * For the licensing terms see $ROOTSYS/LICENSE. *
11 * For the list of contributors see $ROOTSYS/README/CREDITS. *
12 *************************************************************************/
13
15
16#include <ROOT/RField.hxx>
19#include <ROOT/RNTuple.hxx>
20#include <ROOT/RNTupleModel.hxx>
22
23#include <TROOT.h>
24
25#include <cassert>
26
28{
29 const auto descGuard = fPtrControlBlock->fPageSource->GetSharedDescriptorGuard();
32 throw RException(R__FAIL(std::string("entry number ") + std::to_string(entryNumber) + " out of range"));
33
34 auto [itr, _] = fPtrControlBlock->fActiveClusters.try_emplace(clusterId, 0);
35 if (itr->second++ == 0)
36 fPtrControlBlock->fPageSource->PinCluster(clusterId);
37}
38
40{
41 const auto descGuard = fPtrControlBlock->fPageSource->GetSharedDescriptorGuard();
43
45 throw RException(R__FAIL(std::string("entry number ") + std::to_string(entryNumber) + " out of range"));
46
47 auto itr = fPtrControlBlock->fActiveClusters.find(clusterId);
48 assert(itr != fPtrControlBlock->fActiveClusters.end());
49
50 if (--(itr->second) == 0) {
51 fPtrControlBlock->fActiveClusters.erase(itr);
52 fPtrControlBlock->fPageSource->UnpinCluster(clusterId);
53 }
54}
55
57{
58 if (entryNumber == fEntryNumber || entryNumber == kInvalidNTupleIndex)
59 return;
60
61 const std::lock_guard<std::mutex> _(fPtrControlBlock->fLock);
62 if (fPtrControlBlock->fPageSource) {
63 ActivateEntry(entryNumber);
64 if (fEntryNumber != kInvalidNTupleIndex)
65 DeactivateEntry(fEntryNumber);
66 }
67 fEntryNumber = entryNumber;
68}
69
71{
72 if (fEntryNumber == kInvalidNTupleIndex)
73 return;
74
75 const std::lock_guard<std::mutex> _(fPtrControlBlock->fLock);
76 if (fPtrControlBlock->fPageSource) {
77 DeactivateEntry(fEntryNumber);
78 }
79 fEntryNumber = kInvalidNTupleIndex;
80}
81
83{
84 fPtrControlBlock = other.fPtrControlBlock;
85 SetEntryNumber(other.fEntryNumber);
86}
87
89{
90 std::swap(fEntryNumber, other.fEntryNumber);
91 std::swap(fPtrControlBlock, other.fPtrControlBlock);
92 assert(other.fEntryNumber == kInvalidNTupleIndex);
93 assert(!other.fPtrControlBlock);
94}
95
98{
99 if (this == &other)
100 return *this;
101 if (other.fPtrControlBlock.get() != fPtrControlBlock.get()) {
102 Reset();
103 fPtrControlBlock = other.fPtrControlBlock;
104 }
105 SetEntryNumber(other.fEntryNumber);
106 return *this;
107}
108
110{
111 std::swap(fEntryNumber, other.fEntryNumber);
112 std::swap(fPtrControlBlock, other.fPtrControlBlock);
113 return *this;
114}
115
117{
120 // We must not use the descriptor guard to prevent recursive locking in field.ConnectPageSource
121 ROOT::DescriptorId_t fieldZeroId = fSource->GetSharedDescriptorGuard()->GetFieldZeroId();
122 fieldZero.SetOnDiskId(fieldZeroId);
123 // Iterate only over fieldZero's direct subfields; their descendants are recursively handled in
124 // RFieldBase::ConnectPageSource
125 for (auto &field : fieldZero.GetMutableSubfields()) {
126 // If the model has been created from the descriptor, the on-disk IDs are already set.
127 // User-provided models instead need to find their corresponding IDs in the descriptor.
128 if (field->GetOnDiskId() == ROOT::kInvalidDescriptorId) {
129 field->SetOnDiskId(fSource->GetSharedDescriptorGuard()->FindFieldId(field->GetFieldName(), fieldZeroId));
130 }
132 }
134}
135
137{
138#ifdef R__USE_IMT
139 if (IsImplicitMTEnabled() &&
140 fSource->GetReadOptions().GetUseImplicitMT() == ROOT::RNTupleReadOptions::EImplicitMT::kDefault) {
141 fUnzipTasks = std::make_unique<Experimental::Internal::RNTupleImtTaskScheduler>();
142 fSource->SetTaskScheduler(fUnzipTasks.get());
143 }
144#endif
145 fMetrics.ObserveMetrics(fSource->GetMetrics());
146 if (enableMetrics)
148 fSource->Attach();
149 fNEntries = fSource->GetNEntries();
150 fActiveEntriesControlBlock = std::make_shared<RActiveEntriesControlBlock>(fSource.get());
151}
152
153ROOT::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::RNTupleModel> model,
154 std::unique_ptr<ROOT::Internal::RPageSource> source,
155 const ROOT::RNTupleReadOptions &options)
156 : fSource(std::move(source)), fModel(std::move(model)), fMetrics("RNTupleReader")
157{
158 // TODO(jblomer): properly support projected fields
160 if (!projectedFields.IsEmpty()) {
161 throw RException(R__FAIL("model has projected fields, which is incompatible with providing a read model"));
162 }
163 fModel->Freeze();
165 ConnectModel(*fModel, false /* allowFieldSubstitutions */);
166}
167
168ROOT::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::Internal::RPageSource> source,
169 const ROOT::RNTupleReadOptions &options)
170 : fSource(std::move(source)), fModel(nullptr), fMetrics("RNTupleReader")
171{
173}
174
176{
177 const std::lock_guard<std::mutex> _(fActiveEntriesControlBlock->fLock);
178 fActiveEntriesControlBlock->fPageSource = nullptr;
179}
180
181std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::unique_ptr<ROOT::RNTupleModel> model,
182 std::string_view ntupleName, std::string_view storage,
183 const ROOT::RNTupleReadOptions &options)
184{
185 return std::unique_ptr<RNTupleReader>(
186 new RNTupleReader(std::move(model), Internal::RPageSource::Create(ntupleName, storage, options), options));
187}
188
189std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::string_view ntupleName, std::string_view storage,
190 const ROOT::RNTupleReadOptions &options)
191{
192 return std::unique_ptr<RNTupleReader>(
194}
195
196std::unique_ptr<ROOT::RNTupleReader>
198{
199 return std::unique_ptr<RNTupleReader>(
201}
202
203std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::unique_ptr<ROOT::RNTupleModel> model,
204 const ROOT::RNTuple &ntuple,
205 const ROOT::RNTupleReadOptions &options)
206{
207 return std::unique_ptr<RNTupleReader>(
208 new RNTupleReader(std::move(model), Internal::RPageSourceFile::CreateFromAnchor(ntuple, options), options));
209}
210
211std::unique_ptr<ROOT::RNTupleReader>
213 std::string_view ntupleName, std::string_view storage,
214 const ROOT::RNTupleReadOptions &options)
215{
216 auto reader = std::unique_ptr<RNTupleReader>(
218 reader->fCreateModelOptions = createModelOpts;
219 return reader;
220}
221
222std::unique_ptr<ROOT::RNTupleReader>
224 const ROOT::RNTuple &ntuple, const ROOT::RNTupleReadOptions &options)
225{
226 auto reader = std::unique_ptr<RNTupleReader>(
228 reader->fCreateModelOptions = createModelOpts;
229 return reader;
230}
231
233{
234 if (!fModel) {
235 fModel = fSource->GetSharedDescriptorGuard()->CreateModel(
236 fCreateModelOptions.value_or(ROOT::RNTupleDescriptor::RCreateModelOptions{}));
237 ConnectModel(*fModel, true /* allowFieldSubstitutions */);
238 }
239 return *fModel;
240}
241
242std::unique_ptr<ROOT::REntry> ROOT::RNTupleReader::CreateEntry()
243{
244 return GetModel().CreateEntry();
245}
246
248{
249 using namespace ROOT::Internal;
250
251 // TODO(lesimon): In a later version, these variables may be defined by the user or the ideal width may be read out
252 // from the terminal.
253 char frameSymbol = '*';
254 int width = 80;
255 /*
256 if (width < 30) {
257 output << "The width is too small! Should be at least 30." << std::endl;
258 return;
259 }
260 */
261 switch (what) {
263 std::string name;
264 std::unique_ptr<ROOT::RNTupleModel> fullModel;
265 {
266 auto descriptorGuard = fSource->GetSharedDescriptorGuard();
267 name = descriptorGuard->GetName();
269 opts.SetCreateBare(true);
270 // When printing the schema we always try to reconstruct the whole thing even when we are missing the
271 // dictionaries.
272 opts.SetEmulateUnknownTypes(true);
273 fullModel = descriptorGuard->CreateModel(opts);
274 }
275
276 for (int i = 0; i < (width / 2 + width % 2 - 4); ++i)
278 output << " NTUPLE ";
279 for (int i = 0; i < (width / 2 - 4); ++i)
281 output << "\n";
282 // FitString defined in RFieldVisitor.cxx
283 output << frameSymbol << " N-Tuple : " << RNTupleFormatter::FitString(name, width - 13) << frameSymbol
284 << "\n"; // prints line with name of ntuple
285 output << frameSymbol << " Entries : " << RNTupleFormatter::FitString(std::to_string(GetNEntries()), width - 13)
286 << frameSymbol << "\n"; // prints line with number of entries
287
288 // Traverses through all fields to gather information needed for printing.
290 // Traverses through all fields to do the actual printing.
292
293 // Note that we do not need to connect the model, we are only looking at its tree of fields
294 fullModel->GetConstFieldZero().AcceptVisitor(prepVisitor);
295
296 printVisitor.SetFrameSymbol(frameSymbol);
297 printVisitor.SetWidth(width);
298 printVisitor.SetDeepestLevel(prepVisitor.GetDeepestLevel());
299 printVisitor.SetNumFields(prepVisitor.GetNumFields());
300
301 for (int i = 0; i < width; ++i)
303 output << "\n";
304 fullModel->GetConstFieldZero().AcceptVisitor(printVisitor);
305 for (int i = 0; i < width; ++i)
307 output << std::endl;
308 break;
309 }
310 case ENTupleInfo::kStorageDetails: fSource->GetSharedDescriptorGuard()->PrintInfo(output); break;
311 case ENTupleInfo::kMetrics: fMetrics.Print(output); break;
312 default:
313 // Unhandled case, internal error
314 R__ASSERT(false);
315 }
316}
317
319{
320 if (!fDisplayReader) {
322 opts.SetEmulateUnknownTypes(true);
323 auto fullModel = fSource->GetSharedDescriptorGuard()->CreateModel(opts);
324 fDisplayReader = std::unique_ptr<RNTupleReader>(
325 new RNTupleReader(std::move(fullModel), fSource->Clone(), ROOT::RNTupleReadOptions{}));
326 }
327 return fDisplayReader.get();
328}
329
331{
332 auto reader = GetDisplayReader();
333 const auto &entry = reader->GetModel().GetDefaultEntry();
334
335 reader->LoadEntry(index);
336 output << "{";
337 for (auto iValue = entry.begin(); iValue != entry.end();) {
338 output << std::endl;
340 iValue->GetField().AcceptVisitor(visitor);
341
342 if (++iValue == entry.end()) {
343 output << std::endl;
344 break;
345 } else {
346 output << ",";
347 }
348 }
349 output << "}" << std::endl;
350}
351
353{
354 auto descriptorGuard = fSource->GetSharedDescriptorGuard();
355 if (!fCachedDescriptor || fCachedDescriptor->GetGeneration() != descriptorGuard->GetGeneration())
356 fCachedDescriptor = descriptorGuard->Clone();
357 return *fCachedDescriptor;
358}
359
361{
362 auto fieldId = fSource->GetSharedDescriptorGuard()->FindFieldId(fieldName);
364 throw RException(R__FAIL("no field named '" + std::string(fieldName) + "' in RNTuple '" +
365 fSource->GetSharedDescriptorGuard()->GetName() + "'"));
366 }
367 return fieldId;
368}
#define R__FAIL(msg)
Short-hand to return an RResult<T> in an error state; the RError is implicitly converted into RResult...
Definition RError.hxx:300
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t width
char name[80]
Definition TGX11.cxx:110
#define _(A, B)
Definition cfortran.h:108
void ObserveMetrics(RNTupleMetrics &observee)
static std::string FitString(const std::string &str, int availableSpace)
static std::unique_ptr< RPageSourceFile > CreateFromAnchor(const RNTuple &anchor, const ROOT::RNTupleReadOptions &options=ROOT::RNTupleReadOptions())
Used from the RNTuple class to build a datasource if the anchor is already available.
static std::unique_ptr< RPageSource > Create(std::string_view ntupleName, std::string_view location, const ROOT::RNTupleReadOptions &options=ROOT::RNTupleReadOptions())
Guess the concrete derived page source from the file name (location)
Visitor used for a pre-processing run to collect information needed by another visitor class.
Contains settings for printing and prints a summary of an RField instance.
Renders a JSON value corresponding to the field.
Base class for all ROOT issued exceptions.
Definition RError.hxx:79
The on-storage metadata of an RNTuple.
The RNTupleModel encapulates the schema of an RNTuple.
Common user-tunable settings for reading RNTuples.
An active entry token is a pledge for the data of a certain entry number not to be evicted from the p...
std::shared_ptr< RActiveEntriesControlBlock > fPtrControlBlock
RActiveEntryToken & operator=(const RActiveEntryToken &other)
void ActivateEntry(NTupleSize_t entryNumber)
void Reset()
Release the entry number, i.e.
RActiveEntryToken(std::shared_ptr< RActiveEntriesControlBlock > ptrControlBlock)
void DeactivateEntry(NTupleSize_t entryNumber)
void SetEntryNumber(NTupleSize_t entryNumber)
Set or replace the entry number.
Reads RNTuple data from storage.
std::unique_ptr< ROOT::REntry > CreateEntry()
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.
void InitPageSource(bool enableMetrics)
std::unique_ptr< Internal::RPageStorage::RTaskScheduler > fUnzipTasks
Set as the page source's scheduler for parallel page decompression if implicit multi-threading (IMT) ...
ROOT::NTupleSize_t fNEntries
We know that the RNTupleReader is always reading a single RNTuple, so the number of entries is fixed.
std::unique_ptr< ROOT::RNTupleModel > fModel
Needs to be destructed before fSource.
const ROOT::RNTupleDescriptor & GetDescriptor()
Returns a cached copy of the page source descriptor.
void PrintInfo(const ENTupleInfo what=ENTupleInfo::kSummary, std::ostream &output=std::cout) const
Prints a detailed summary of the RNTuple, including a list of fields.
std::unique_ptr< Internal::RPageSource > fSource
RNTupleReader(std::unique_ptr< ROOT::RNTupleModel > model, std::unique_ptr< Internal::RPageSource > source, const ROOT::RNTupleReadOptions &options)
const ROOT::RNTupleModel & GetModel()
std::shared_ptr< RActiveEntriesControlBlock > fActiveEntriesControlBlock
Initialized when the page source is connected.
void ConnectModel(ROOT::RNTupleModel &model, bool allowFieldSubstitutions)
ROOT::DescriptorId_t RetrieveFieldId(std::string_view fieldName) const
void Show(ROOT::NTupleSize_t index, std::ostream &output=std::cout)
Shows the values of the i-th entry/row, starting with 0 for the first entry.
RNTupleReader * GetDisplayReader()
Experimental::Detail::RNTupleMetrics fMetrics
void EnableMetrics()
Enable performance measurements (decompression time, bytes read from storage, etc....
Representation of an RNTuple data set in a ROOT file.
Definition RNTuple.hxx:67
const_iterator begin() const
const_iterator end() const
void SetAllowFieldSubstitutions(RFieldZero &fieldZero, bool val)
Definition RField.cxx:36
ROOT::RFieldZero & GetFieldZeroOfModel(RNTupleModel &model)
ROOT::DescriptorId_t CallFindClusterIdOn(const ROOT::RNTupleDescriptor &desc, ROOT::NTupleSize_t entryIdx)
void CallConnectPageSourceOnField(RFieldBase &, ROOT::Internal::RPageSource &)
RProjectedFields & GetProjectedFieldsOfModel(RNTupleModel &model)
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:600
std::uint64_t DescriptorId_t
Distriniguishes elements of the same type within a descriptor, e.g. different fields.
constexpr NTupleSize_t kInvalidNTupleIndex
ENTupleInfo
Listing of the different options that can be printed by RNTupleReader::GetInfo()
std::uint64_t NTupleSize_t
Integer type long enough to hold the maximum number of entries in a column.
constexpr DescriptorId_t kInvalidDescriptorId
static const char * what
Definition stlLoader.cc:5
static void output()