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