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>
20#include <ROOT/RNTuple.hxx>
21#include <ROOT/RNTupleModel.hxx>
23
24#include <TROOT.h>
25
26#include <cassert>
27
29{
30 const auto descGuard = fPtrControlBlock->fPageSource->GetSharedDescriptorGuard();
33 throw RException(R__FAIL(std::string("entry number ") + std::to_string(entryNumber) + " out of range"));
34
35 auto [itr, _] = fPtrControlBlock->fActiveClusters.try_emplace(clusterId, 0);
36 if (itr->second++ == 0)
37 fPtrControlBlock->fPageSource->PinCluster(clusterId);
38}
39
41{
42 const auto descGuard = fPtrControlBlock->fPageSource->GetSharedDescriptorGuard();
44
46 throw RException(R__FAIL(std::string("entry number ") + std::to_string(entryNumber) + " out of range"));
47
48 auto itr = fPtrControlBlock->fActiveClusters.find(clusterId);
49 assert(itr != fPtrControlBlock->fActiveClusters.end());
50
51 if (--(itr->second) == 0) {
52 fPtrControlBlock->fActiveClusters.erase(itr);
53 fPtrControlBlock->fPageSource->UnpinCluster(clusterId);
54 }
55}
56
58{
59 if (entryNumber == fEntryNumber || entryNumber == kInvalidNTupleIndex)
60 return;
61
62 const std::lock_guard<std::mutex> _(fPtrControlBlock->fLock);
63 if (fPtrControlBlock->fPageSource) {
64 ActivateEntry(entryNumber);
65 if (fEntryNumber != kInvalidNTupleIndex)
66 DeactivateEntry(fEntryNumber);
67 }
68 fEntryNumber = entryNumber;
69}
70
72{
73 if (fEntryNumber == kInvalidNTupleIndex)
74 return;
75
76 const std::lock_guard<std::mutex> _(fPtrControlBlock->fLock);
77 if (fPtrControlBlock->fPageSource) {
78 DeactivateEntry(fEntryNumber);
79 }
80 fEntryNumber = kInvalidNTupleIndex;
81}
82
84{
85 fPtrControlBlock = other.fPtrControlBlock;
86 SetEntryNumber(other.fEntryNumber);
87}
88
90{
91 std::swap(fEntryNumber, other.fEntryNumber);
92 std::swap(fPtrControlBlock, other.fPtrControlBlock);
93 assert(other.fEntryNumber == kInvalidNTupleIndex);
94 assert(!other.fPtrControlBlock);
95}
96
99{
100 if (this == &other)
101 return *this;
102 if (other.fPtrControlBlock.get() != fPtrControlBlock.get()) {
103 Reset();
104 fPtrControlBlock = other.fPtrControlBlock;
105 }
106 SetEntryNumber(other.fEntryNumber);
107 return *this;
108}
109
111{
112 std::swap(fEntryNumber, other.fEntryNumber);
113 std::swap(fPtrControlBlock, other.fPtrControlBlock);
114 return *this;
115}
116
118{
121 // We must not use the descriptor guard to prevent recursive locking in field.ConnectPageSource
122 ROOT::DescriptorId_t fieldZeroId = fSource->GetSharedDescriptorGuard()->GetFieldZeroId();
123 fieldZero.SetOnDiskId(fieldZeroId);
124 // Iterate only over fieldZero's direct subfields; their descendants are recursively handled in
125 // RFieldBase::ConnectPageSource
126 for (auto &field : fieldZero.GetMutableSubfields()) {
127 // If the model has been created from the descriptor, the on-disk IDs are already set.
128 // User-provided models instead need to find their corresponding IDs in the descriptor.
129 if (field->GetOnDiskId() == ROOT::kInvalidDescriptorId) {
130 field->SetOnDiskId(fSource->GetSharedDescriptorGuard()->FindFieldId(field->GetFieldName(), fieldZeroId));
131 }
133 }
135}
136
138{
139#ifdef R__USE_IMT
140 if (IsImplicitMTEnabled() &&
141 fSource->GetReadOptions().GetUseImplicitMT() == ROOT::RNTupleReadOptions::EImplicitMT::kDefault) {
142 fUnzipTasks = std::make_unique<Experimental::Internal::RNTupleImtTaskScheduler>();
143 fSource->SetTaskScheduler(fUnzipTasks.get());
144 }
145#endif
146 fMetrics.ObserveMetrics(fSource->GetMetrics());
147 if (enableMetrics)
149 fSource->Attach();
150 fNEntries = fSource->GetNEntries();
151 fActiveEntriesControlBlock = std::make_shared<RActiveEntriesControlBlock>(fSource.get());
152}
153
154ROOT::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::RNTupleModel> model,
155 std::unique_ptr<ROOT::Internal::RPageSource> source,
156 const ROOT::RNTupleReadOptions &options)
157 : fSource(std::move(source)), fModel(std::move(model)), fMetrics("RNTupleReader")
158{
159 // TODO(jblomer): properly support projected fields
161 if (!projectedFields.IsEmpty()) {
162 throw RException(R__FAIL("model has projected fields, which is incompatible with providing a read model"));
163 }
164 fModel->Freeze();
166 ConnectModel(*fModel, false /* allowFieldSubstitutions */);
167}
168
169ROOT::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::Internal::RPageSource> source,
170 const ROOT::RNTupleReadOptions &options)
171 : fSource(std::move(source)), fModel(nullptr), fMetrics("RNTupleReader")
172{
174}
175
177{
178 const std::lock_guard<std::mutex> _(fActiveEntriesControlBlock->fLock);
179 fActiveEntriesControlBlock->fPageSource = nullptr;
180}
181
182std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::unique_ptr<ROOT::RNTupleModel> model,
183 std::string_view ntupleName, std::string_view storage,
184 const ROOT::RNTupleReadOptions &options)
185{
186 return std::unique_ptr<RNTupleReader>(
187 new RNTupleReader(std::move(model), Internal::RPageSource::Create(ntupleName, storage, options), options));
188}
189
190std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::string_view ntupleName, std::string_view storage,
191 const ROOT::RNTupleReadOptions &options)
192{
193 return std::unique_ptr<RNTupleReader>(
195}
196
197std::unique_ptr<ROOT::RNTupleReader>
199{
200 return std::unique_ptr<RNTupleReader>(
202}
203
204std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::unique_ptr<ROOT::RNTupleModel> model,
205 const ROOT::RNTuple &ntuple,
206 const ROOT::RNTupleReadOptions &options)
207{
208 return std::unique_ptr<RNTupleReader>(
209 new RNTupleReader(std::move(model), Internal::RPageSourceFile::CreateFromAnchor(ntuple, options), options));
210}
211
212std::unique_ptr<ROOT::RNTupleReader>
214 std::string_view ntupleName, std::string_view storage,
215 const ROOT::RNTupleReadOptions &options)
216{
217 auto reader = std::unique_ptr<RNTupleReader>(
219 reader->fCreateModelOptions = createModelOpts;
220 return reader;
221}
222
223std::unique_ptr<ROOT::RNTupleReader>
225 const ROOT::RNTuple &ntuple, const ROOT::RNTupleReadOptions &options)
226{
227 auto reader = std::unique_ptr<RNTupleReader>(
229 reader->fCreateModelOptions = createModelOpts;
230 return reader;
231}
232
234{
235 if (!fModel) {
236 fModel = fSource->GetSharedDescriptorGuard()->CreateModel(
237 fCreateModelOptions.value_or(ROOT::RNTupleDescriptor::RCreateModelOptions{}));
238 ConnectModel(*fModel, true /* allowFieldSubstitutions */);
239 }
240 return *fModel;
241}
242
243std::unique_ptr<ROOT::REntry> ROOT::RNTupleReader::CreateEntry()
244{
245 return GetModel().CreateEntry();
246}
247
248void ROOT::RNTupleReader::PrintInfo(const ENTupleInfo what, std::ostream &output) const
249{
250 using namespace ROOT::Internal;
251
252 // TODO(lesimon): In a later version, these variables may be defined by the user or the ideal width may be read out
253 // from the terminal.
254 char frameSymbol = '*';
255 int width = 80;
256 /*
257 if (width < 30) {
258 output << "The width is too small! Should be at least 30." << std::endl;
259 return;
260 }
261 */
262 switch (what) {
264 std::string name;
265 std::unique_ptr<ROOT::RNTupleModel> fullModel;
266 {
267 auto descriptorGuard = fSource->GetSharedDescriptorGuard();
268 name = descriptorGuard->GetName();
270 opts.SetCreateBare(true);
271 // When printing the schema we always try to reconstruct the whole thing even when we are missing the
272 // dictionaries.
273 opts.SetEmulateUnknownTypes(true);
274 fullModel = descriptorGuard->CreateModel(opts);
275 }
276
277 for (int i = 0; i < (width / 2 + width % 2 - 4); ++i)
278 output << frameSymbol;
279 output << " NTUPLE ";
280 for (int i = 0; i < (width / 2 - 4); ++i)
281 output << frameSymbol;
282 output << "\n";
283 // FitString defined in RFieldVisitor.cxx
284 output << frameSymbol << " N-Tuple : " << RNTupleFormatter::FitString(name, width - 13) << frameSymbol
285 << "\n"; // prints line with name of ntuple
286 output << frameSymbol << " Entries : " << RNTupleFormatter::FitString(std::to_string(GetNEntries()), width - 13)
287 << frameSymbol << "\n"; // prints line with number of entries
288
289 // Traverses through all fields to gather information needed for printing.
291 // Traverses through all fields to do the actual printing.
293
294 // Note that we do not need to connect the model, we are only looking at its tree of fields
295 fullModel->GetConstFieldZero().AcceptVisitor(prepVisitor);
296
297 printVisitor.SetFrameSymbol(frameSymbol);
298 printVisitor.SetWidth(width);
299 printVisitor.SetDeepestLevel(prepVisitor.GetDeepestLevel());
300 printVisitor.SetNumFields(prepVisitor.GetNumFields());
301
302 for (int i = 0; i < width; ++i)
303 output << frameSymbol;
304 output << "\n";
305 fullModel->GetConstFieldZero().AcceptVisitor(printVisitor);
306 for (int i = 0; i < width; ++i)
307 output << frameSymbol;
308 output << std::endl;
309 break;
310 }
311 case ENTupleInfo::kStorageDetails: fSource->GetSharedDescriptorGuard()->PrintInfo(output); break;
312 case ENTupleInfo::kMetrics: fMetrics.Print(output); break;
313 default:
314 // Unhandled case, internal error
315 R__ASSERT(false);
316 }
317}
318
320{
321 if (!fDisplayReader) {
323 opts.SetEmulateUnknownTypes(true);
324 auto fullModel = fSource->GetSharedDescriptorGuard()->CreateModel(opts);
325 fDisplayReader = std::unique_ptr<RNTupleReader>(
326 new RNTupleReader(std::move(fullModel), fSource->Clone(), ROOT::RNTupleReadOptions{}));
327 }
328 return fDisplayReader.get();
329}
330
332{
333 auto reader = GetDisplayReader();
334 const auto &entry = reader->GetModel().GetDefaultEntry();
335
336 reader->LoadEntry(index);
337 output << "{";
338 for (auto iValue = entry.begin(); iValue != entry.end();) {
339 output << std::endl;
340 ROOT::Internal::RPrintValueVisitor visitor(*iValue, output, 1 /* level */);
341 iValue->GetField().AcceptVisitor(visitor);
342
343 if (++iValue == entry.end()) {
344 output << std::endl;
345 break;
346 } else {
347 output << ",";
348 }
349 }
350 output << "}" << std::endl;
351}
352
354{
355 auto descriptorGuard = fSource->GetSharedDescriptorGuard();
356 if (!fCachedDescriptor || fCachedDescriptor->GetGeneration() != descriptorGuard->GetGeneration())
357 fCachedDescriptor = descriptorGuard->Clone();
358 return *fCachedDescriptor;
359}
360
362{
363 auto fieldId = fSource->GetSharedDescriptorGuard()->FindFieldId(fieldName);
365 throw RException(R__FAIL("no field named '" + std::string(fieldName) + "' in RNTuple '" +
366 fSource->GetSharedDescriptorGuard()->GetName() + "'"));
367 }
368 return fieldId;
369}
370
371std::unique_ptr<ROOT::Experimental::RNTupleAttrSetReader>
373{
374 auto attrSets = GetDescriptor().GetAttrSetIterable();
375 const auto it =
376 std::find_if(attrSets.begin(), attrSets.end(), [&](const auto &d) { return d.GetName() == attrSetName; });
377 if (it == attrSets.end())
378 throw ROOT::RException(R__FAIL(std::string("No such attribute set: ") + std::string(attrSetName)));
379
380 auto attrSource = fSource->OpenWithDifferentAnchor({it->GetAnchorLocator(), it->GetAnchorLength()}, readOpts);
381 auto newReader = std::unique_ptr<RNTupleReader>(new RNTupleReader(std::move(attrSource), readOpts));
383 auto attrSetReader = std::unique_ptr<ROOT::Experimental::RNTupleAttrSetReader>(
384 new ROOT::Experimental::RNTupleAttrSetReader(std::move(newReader), it->GetSchemaVersionMajor()));
385 return attrSetReader;
386}
#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
#define d(i)
Definition RSha256.hxx:102
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:148
#define _(A, B)
Definition cfortran.h:108
void ObserveMetrics(RNTupleMetrics &observee)
Class used to read a RNTupleAttrSet in the context of a RNTupleReader.
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::unique_ptr< Experimental::RNTupleAttrSetReader > OpenAttributeSet(std::string_view attrSetName, const ROOT::RNTupleReadOptions &options={})
Looks for an attribute set with the given name and creates an RNTupleAttrSetReader for it,...
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:68
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:675
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