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
26{
29 // We must not use the descriptor guard to prevent recursive locking in field.ConnectPageSource
30 ROOT::DescriptorId_t fieldZeroId = fSource->GetSharedDescriptorGuard()->GetFieldZeroId();
31 fieldZero.SetOnDiskId(fieldZeroId);
32 // Iterate only over fieldZero's direct subfields; their descendants are recursively handled in
33 // RFieldBase::ConnectPageSource
34 for (auto &field : fieldZero.GetMutableSubfields()) {
35 // If the model has been created from the descriptor, the on-disk IDs are already set.
36 // User-provided models instead need to find their corresponding IDs in the descriptor.
37 if (field->GetOnDiskId() == ROOT::kInvalidDescriptorId) {
38 field->SetOnDiskId(fSource->GetSharedDescriptorGuard()->FindFieldId(field->GetFieldName(), fieldZeroId));
39 }
41 }
43}
44
46{
47#ifdef R__USE_IMT
48 if (IsImplicitMTEnabled() &&
49 fSource->GetReadOptions().GetUseImplicitMT() == ROOT::RNTupleReadOptions::EImplicitMT::kDefault) {
50 fUnzipTasks = std::make_unique<Experimental::Internal::RNTupleImtTaskScheduler>();
51 fSource->SetTaskScheduler(fUnzipTasks.get());
52 }
53#endif
54 fMetrics.ObserveMetrics(fSource->GetMetrics());
55 if (enableMetrics)
56 EnableMetrics();
57 fSource->Attach();
58 fNEntries = fSource->GetNEntries();
59}
60
61ROOT::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::RNTupleModel> model,
62 std::unique_ptr<ROOT::Internal::RPageSource> source,
63 const ROOT::RNTupleReadOptions &options)
64 : fSource(std::move(source)), fModel(std::move(model)), fMetrics("RNTupleReader")
65{
66 // TODO(jblomer): properly support projected fields
68 if (!projectedFields.IsEmpty()) {
69 throw RException(R__FAIL("model has projected fields, which is incompatible with providing a read model"));
70 }
71 fModel->Freeze();
73 ConnectModel(*fModel, false /* allowFieldSubstitutions */);
74}
75
76ROOT::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::Internal::RPageSource> source,
77 const ROOT::RNTupleReadOptions &options)
78 : fSource(std::move(source)), fModel(nullptr), fMetrics("RNTupleReader")
79{
81}
82
84
85std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::unique_ptr<ROOT::RNTupleModel> model,
86 std::string_view ntupleName, std::string_view storage,
87 const ROOT::RNTupleReadOptions &options)
88{
89 return std::unique_ptr<RNTupleReader>(
90 new RNTupleReader(std::move(model), Internal::RPageSource::Create(ntupleName, storage, options), options));
91}
92
93std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::string_view ntupleName, std::string_view storage,
94 const ROOT::RNTupleReadOptions &options)
95{
96 return std::unique_ptr<RNTupleReader>(
98}
99
100std::unique_ptr<ROOT::RNTupleReader>
102{
103 return std::unique_ptr<RNTupleReader>(
105}
106
107std::unique_ptr<ROOT::RNTupleReader> ROOT::RNTupleReader::Open(std::unique_ptr<ROOT::RNTupleModel> model,
108 const ROOT::RNTuple &ntuple,
109 const ROOT::RNTupleReadOptions &options)
110{
111 return std::unique_ptr<RNTupleReader>(
112 new RNTupleReader(std::move(model), Internal::RPageSourceFile::CreateFromAnchor(ntuple, options), options));
113}
114
115std::unique_ptr<ROOT::RNTupleReader>
117 std::string_view ntupleName, std::string_view storage,
118 const ROOT::RNTupleReadOptions &options)
119{
120 auto reader = std::unique_ptr<RNTupleReader>(
122 reader->fCreateModelOptions = createModelOpts;
123 return reader;
124}
125
126std::unique_ptr<ROOT::RNTupleReader>
128 const ROOT::RNTuple &ntuple, const ROOT::RNTupleReadOptions &options)
129{
130 auto reader = std::unique_ptr<RNTupleReader>(
132 reader->fCreateModelOptions = createModelOpts;
133 return reader;
134}
135
137{
138 if (!fModel) {
139 fModel = fSource->GetSharedDescriptorGuard()->CreateModel(
140 fCreateModelOptions.value_or(ROOT::RNTupleDescriptor::RCreateModelOptions{}));
141 ConnectModel(*fModel, true /* allowFieldSubstitutions */);
142 }
143 return *fModel;
144}
145
146std::unique_ptr<ROOT::REntry> ROOT::RNTupleReader::CreateEntry()
147{
148 return GetModel().CreateEntry();
149}
150
152{
153 using namespace ROOT::Internal;
154
155 // TODO(lesimon): In a later version, these variables may be defined by the user or the ideal width may be read out
156 // from the terminal.
157 char frameSymbol = '*';
158 int width = 80;
159 /*
160 if (width < 30) {
161 output << "The width is too small! Should be at least 30." << std::endl;
162 return;
163 }
164 */
165 switch (what) {
167 std::string name;
168 std::unique_ptr<ROOT::RNTupleModel> fullModel;
169 {
170 auto descriptorGuard = fSource->GetSharedDescriptorGuard();
171 name = descriptorGuard->GetName();
173 opts.SetCreateBare(true);
174 // When printing the schema we always try to reconstruct the whole thing even when we are missing the
175 // dictionaries.
176 opts.SetEmulateUnknownTypes(true);
177 fullModel = descriptorGuard->CreateModel(opts);
178 }
179
180 for (int i = 0; i < (width / 2 + width % 2 - 4); ++i)
182 output << " NTUPLE ";
183 for (int i = 0; i < (width / 2 - 4); ++i)
185 output << "\n";
186 // FitString defined in RFieldVisitor.cxx
187 output << frameSymbol << " N-Tuple : " << RNTupleFormatter::FitString(name, width - 13) << frameSymbol
188 << "\n"; // prints line with name of ntuple
189 output << frameSymbol << " Entries : " << RNTupleFormatter::FitString(std::to_string(GetNEntries()), width - 13)
190 << frameSymbol << "\n"; // prints line with number of entries
191
192 // Traverses through all fields to gather information needed for printing.
194 // Traverses through all fields to do the actual printing.
196
197 // Note that we do not need to connect the model, we are only looking at its tree of fields
198 fullModel->GetConstFieldZero().AcceptVisitor(prepVisitor);
199
200 printVisitor.SetFrameSymbol(frameSymbol);
201 printVisitor.SetWidth(width);
202 printVisitor.SetDeepestLevel(prepVisitor.GetDeepestLevel());
203 printVisitor.SetNumFields(prepVisitor.GetNumFields());
204
205 for (int i = 0; i < width; ++i)
207 output << "\n";
208 fullModel->GetConstFieldZero().AcceptVisitor(printVisitor);
209 for (int i = 0; i < width; ++i)
211 output << std::endl;
212 break;
213 }
214 case ENTupleInfo::kStorageDetails: fSource->GetSharedDescriptorGuard()->PrintInfo(output); break;
215 case ENTupleInfo::kMetrics: fMetrics.Print(output); break;
216 default:
217 // Unhandled case, internal error
218 R__ASSERT(false);
219 }
220}
221
223{
224 if (!fDisplayReader) {
226 opts.SetEmulateUnknownTypes(true);
227 auto fullModel = fSource->GetSharedDescriptorGuard()->CreateModel(opts);
228 fDisplayReader = std::unique_ptr<RNTupleReader>(
229 new RNTupleReader(std::move(fullModel), fSource->Clone(), ROOT::RNTupleReadOptions{}));
230 }
231 return fDisplayReader.get();
232}
233
235{
236 auto reader = GetDisplayReader();
237 const auto &entry = reader->GetModel().GetDefaultEntry();
238
239 reader->LoadEntry(index);
240 output << "{";
241 for (auto iValue = entry.begin(); iValue != entry.end();) {
242 output << std::endl;
244 iValue->GetField().AcceptVisitor(visitor);
245
246 if (++iValue == entry.end()) {
247 output << std::endl;
248 break;
249 } else {
250 output << ",";
251 }
252 }
253 output << "}" << std::endl;
254}
255
257{
258 auto descriptorGuard = fSource->GetSharedDescriptorGuard();
259 if (!fCachedDescriptor || fCachedDescriptor->GetGeneration() != descriptorGuard->GetGeneration())
260 fCachedDescriptor = descriptorGuard->Clone();
261 return *fCachedDescriptor;
262}
263
265{
266 auto fieldId = fSource->GetSharedDescriptorGuard()->FindFieldId(fieldName);
268 throw RException(R__FAIL("no field named '" + std::string(fieldName) + "' in RNTuple '" +
269 fSource->GetSharedDescriptorGuard()->GetName() + "'"));
270 }
271 return fieldId;
272}
#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
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.
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< 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()
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()
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:34
ROOT::RFieldZero & GetFieldZeroOfModel(RNTupleModel &model)
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.
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()