Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RNTuple.cxx
Go to the documentation of this file.
1/// \file RNTuple.cxx
2/// \ingroup NTuple ROOT7
3/// \author Jakob Blomer <jblomer@cern.ch>
4/// \date 2018-10-04
5/// \warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback
6/// is welcome!
7
8/*************************************************************************
9 * Copyright (C) 1995-2019, Rene Brun and Fons Rademakers. *
10 * All rights reserved. *
11 * *
12 * For the licensing terms see $ROOTSYS/LICENSE. *
13 * For the list of contributors see $ROOTSYS/README/CREDITS. *
14 *************************************************************************/
15
16#include <ROOT/RNTuple.hxx>
17
19#include <ROOT/RNTupleModel.hxx>
20#include <ROOT/RPageStorage.hxx>
22#ifdef R__USE_IMT
23#include <ROOT/TTaskGroup.hxx>
24#endif
25
26#include <TError.h>
27#include <TROOT.h> // for IsImplicitMTEnabled()
28
29#include <algorithm>
30#include <exception>
31#include <functional>
32#include <iomanip>
33#include <iostream>
34#include <sstream>
35#include <string>
36#include <unordered_map>
37#include <utility>
38
39
40#ifdef R__USE_IMT
42{
43 fTaskGroup = std::make_unique<TTaskGroup>();
44}
45
46
47void ROOT::Experimental::RNTupleImtTaskScheduler::AddTask(const std::function<void(void)> &taskFunc)
48{
49 fTaskGroup->Run(taskFunc);
50}
51
52
54{
55 fTaskGroup->Wait();
56}
57#endif
58
59
60//------------------------------------------------------------------------------
61
62
64 std::unordered_map<const Detail::RFieldBase *, DescriptorId_t> fieldPtr2Id;
65 fieldPtr2Id[model.GetFieldZero()] = fSource->GetDescriptor().GetFieldZeroId();
66 for (auto &field : *model.GetFieldZero()) {
67 auto parentId = fieldPtr2Id[field.GetParent()];
68 auto fieldId = fSource->GetDescriptor().FindFieldId(field.GetName(), parentId);
70 fieldPtr2Id[&field] = fieldId;
71 Detail::RFieldFuse::Connect(fieldId, *fSource, field);
72 }
73}
74
76{
77#ifdef R__USE_IMT
78 if (IsImplicitMTEnabled()) {
79 fUnzipTasks = std::make_unique<RNTupleImtTaskScheduler>();
80 fSource->SetTaskScheduler(fUnzipTasks.get());
81 }
82#endif
83 fSource->Attach();
84 fMetrics.ObserveMetrics(fSource->GetMetrics());
85}
86
88 std::unique_ptr<ROOT::Experimental::RNTupleModel> model,
89 std::unique_ptr<ROOT::Experimental::Detail::RPageSource> source)
90 : fSource(std::move(source))
91 , fModel(std::move(model))
92 , fMetrics("RNTupleReader")
93{
96}
97
98ROOT::Experimental::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::Experimental::Detail::RPageSource> source)
99 : fSource(std::move(source))
100 , fModel(nullptr)
101 , fMetrics("RNTupleReader")
102{
104}
105
107
108std::unique_ptr<ROOT::Experimental::RNTupleReader> ROOT::Experimental::RNTupleReader::Open(
109 std::unique_ptr<RNTupleModel> model,
110 std::string_view ntupleName,
111 std::string_view storage,
112 const RNTupleReadOptions &options)
113{
114 return std::make_unique<RNTupleReader>(std::move(model), Detail::RPageSource::Create(ntupleName, storage, options));
115}
116
117std::unique_ptr<ROOT::Experimental::RNTupleReader> ROOT::Experimental::RNTupleReader::Open(
118 std::string_view ntupleName,
119 std::string_view storage,
120 const RNTupleReadOptions &options)
121{
122 return std::make_unique<RNTupleReader>(Detail::RPageSource::Create(ntupleName, storage, options));
123}
124
126{
127 if (!fModel) {
128 fModel = fSource->GetDescriptor().GenerateModel();
129 ConnectModel(*fModel);
130 }
131 return fModel.get();
132}
133
135{
136 // TODO(lesimon): In a later version, these variables may be defined by the user or the ideal width may be read out from the terminal.
137 char frameSymbol = '*';
138 int width = 80;
139 /*
140 if (width < 30) {
141 output << "The width is too small! Should be at least 30." << std::endl;
142 return;
143 }
144 */
145 std::string name = fSource->GetDescriptor().GetName();
146 switch (what) {
148 for (int i = 0; i < (width/2 + width%2 - 4); ++i)
149 output << frameSymbol;
150 output << " NTUPLE ";
151 for (int i = 0; i < (width/2 - 4); ++i)
152 output << frameSymbol;
153 output << std::endl;
154 // FitString defined in RFieldVisitor.cxx
155 output << frameSymbol << " N-Tuple : " << RNTupleFormatter::FitString(name, width-13) << frameSymbol << std::endl; // prints line with name of ntuple
156 output << frameSymbol << " Entries : " << RNTupleFormatter::FitString(std::to_string(GetNEntries()), width - 13) << frameSymbol << std::endl; // prints line with number of entries
157
158 // Traverses through all fields to gather information needed for printing.
159 RPrepareVisitor prepVisitor;
160 // Traverses through all fields to do the actual printing.
161 RPrintSchemaVisitor printVisitor(output);
162
163 // Note that we do not need to connect the model, we are only looking at its tree of fields
164 auto fullModel = fSource->GetDescriptor().GenerateModel();
165 fullModel->GetFieldZero()->AcceptVisitor(prepVisitor);
166
167 printVisitor.SetFrameSymbol(frameSymbol);
168 printVisitor.SetWidth(width);
169 printVisitor.SetDeepestLevel(prepVisitor.GetDeepestLevel());
170 printVisitor.SetNumFields(prepVisitor.GetNumFields());
171
172 for (int i = 0; i < width; ++i)
173 output << frameSymbol;
174 output << std::endl;
175 fullModel->GetFieldZero()->AcceptVisitor(printVisitor);
176 for (int i = 0; i < width; ++i)
177 output << frameSymbol;
178 output << std::endl;
179 break;
180 }
182 fSource->GetDescriptor().PrintInfo(output);
183 break;
185 fMetrics.Print(output);
186 break;
187 default:
188 // Unhandled case, internal error
189 R__ASSERT(false);
190 }
191}
192
193
195{
196 if (!fDisplayReader)
197 fDisplayReader = Clone();
198 return fDisplayReader.get();
199}
200
201
203{
204 RNTupleReader *reader = this;
205 REntry *entry = nullptr;
206 // Don't accidentally trigger loading of the entire model
207 if (fModel)
208 entry = fModel->GetDefaultEntry();
209
210 switch(format) {
212 reader = GetDisplayReader();
213 entry = reader->GetModel()->GetDefaultEntry();
214 // Fall through
216 if (!entry) {
217 output << "{}" << std::endl;
218 break;
219 }
220
221 reader->LoadEntry(index);
222 output << "{";
223 for (auto iValue = entry->begin(); iValue != entry->end(); ) {
224 output << std::endl;
225 RPrintValueVisitor visitor(*iValue, output, 1 /* level */);
226 iValue->GetField()->AcceptVisitor(visitor);
227
228 if (++iValue == entry->end()) {
229 output << std::endl;
230 break;
231 } else {
232 output << ",";
233 }
234 }
235 output << "}" << std::endl;
236 break;
237 default:
238 // Unhandled case, internal error
239 R__ASSERT(false);
240 }
241}
242
243
244//------------------------------------------------------------------------------
245
246
248 std::unique_ptr<ROOT::Experimental::RNTupleModel> model,
249 std::unique_ptr<ROOT::Experimental::Detail::RPageSink> sink)
250 : fSink(std::move(sink))
251 , fModel(std::move(model))
252 , fClusterSizeEntries(kDefaultClusterSizeEntries)
253 , fLastCommitted(0)
254 , fNEntries(0)
255{
256 fSink->Create(*fModel.get());
257}
258
260{
261 CommitCluster();
262 fSink->CommitDataset();
263}
264
265std::unique_ptr<ROOT::Experimental::RNTupleWriter> ROOT::Experimental::RNTupleWriter::Recreate(
266 std::unique_ptr<RNTupleModel> model,
267 std::string_view ntupleName,
268 std::string_view storage,
269 const RNTupleWriteOptions &options)
270{
271 return std::make_unique<RNTupleWriter>(std::move(model), Detail::RPageSink::Create(ntupleName, storage, options));
272}
273
274std::unique_ptr<ROOT::Experimental::RNTupleWriter> ROOT::Experimental::RNTupleWriter::Append(
275 std::unique_ptr<RNTupleModel> model,
276 std::string_view ntupleName,
277 TFile &file,
278 const RNTupleWriteOptions &options)
279{
280 auto sink = std::make_unique<Detail::RPageSinkFile>(ntupleName, file, options);
281 return std::make_unique<RNTupleWriter>(std::move(model), std::move(sink));
282}
283
284
286{
287 if (fNEntries == fLastCommitted) return;
288 for (auto& field : *fModel->GetFieldZero()) {
289 field.Flush();
290 field.CommitCluster();
291 }
292 fSink->CommitCluster(fNEntries);
293 fLastCommitted = fNEntries;
294}
295
296
297//------------------------------------------------------------------------------
298
299
301 : fOffset(0), fDefaultEntry(std::move(defaultEntry))
302{
303}
include TDocParser_001 C image html pict1_TDocParser_001 png width
#define R__ASSERT(e)
Definition TError.h:120
char name[80]
Definition TGX11.cxx:110
static void Connect(DescriptorId_t fieldId, RPageStorage &pageStorage, RFieldBase &field)
Definition RField.cxx:118
static std::unique_ptr< RPageSink > Create(std::string_view ntupleName, std::string_view location, const RNTupleWriteOptions &options=RNTupleWriteOptions())
Guess the concrete derived page source from the file name (location)
static std::unique_ptr< RPageSource > Create(std::string_view ntupleName, std::string_view location, const RNTupleReadOptions &options=RNTupleReadOptions())
Guess the concrete derived page source from the file name (location)
RCollectionNTuple(std::unique_ptr< REntry > defaultEntry)
Definition RNTuple.cxx:300
The REntry is a collection of values in an ntuple corresponding to a complete row in the data set.
Definition REntry.hxx:42
static std::string FitString(const std::string &str, int availableSpace)
void Reset() final
Start a new set of tasks.
Definition RNTuple.cxx:41
std::unique_ptr< TTaskGroup > fTaskGroup
Definition RNTuple.hxx:69
void AddTask(const std::function< void(void)> &taskFunc) final
Take a callable that represents a task.
Definition RNTuple.cxx:47
void Wait() final
Blocks until all scheduled tasks finished.
Definition RNTuple.cxx:53
The RNTupleModel encapulates the schema of an ntuple.
RFieldZero * GetFieldZero() const
Common user-tunable settings for reading ntuples.
An RNTuple that is used to read data from storage.
Definition RNTuple.hxx:90
RNTupleReader * GetDisplayReader()
Definition RNTuple.cxx:194
static std::unique_ptr< RNTupleReader > Open(std::unique_ptr< RNTupleModel > model, std::string_view ntupleName, std::string_view storage, const RNTupleReadOptions &options=RNTupleReadOptions())
Definition RNTuple.cxx:108
void ConnectModel(const RNTupleModel &model)
Definition RNTuple.cxx:63
RNTupleReader(std::unique_ptr< RNTupleModel > model, std::unique_ptr< Detail::RPageSource > source)
The user imposes an ntuple model, which must be compatible with the model found in the data on storag...
Definition RNTuple.cxx:87
void LoadEntry(NTupleSize_t index)
Analogous to Fill(), fills the default entry of the model.
Definition RNTuple.hxx:165
void PrintInfo(const ENTupleInfo what=ENTupleInfo::kSummary, std::ostream &output=std::cout)
Prints a detailed summary of the ntuple, including a list of fields.
Definition RNTuple.cxx:134
std::unique_ptr< RNTupleModel > fModel
Needs to be destructed before fSource.
Definition RNTuple.hxx:98
void Show(NTupleSize_t index, const ENTupleShowFormat format=ENTupleShowFormat::kCurrentModelJSON, std::ostream &output=std::cout)
Shows the values of the i-th entry/row, starting with 0 for the first entry.
Definition RNTuple.cxx:202
Common user-tunable settings for storing ntuples.
void CommitCluster()
Ensure that the data from the so far seen Fill calls has been written to storage.
Definition RNTuple.cxx:285
RNTupleWriter(std::unique_ptr< RNTupleModel > model, std::unique_ptr< Detail::RPageSink > sink)
Definition RNTuple.cxx:247
std::unique_ptr< RNTupleModel > fModel
Needs to be destructed before fSink.
Definition RNTuple.hxx:218
std::unique_ptr< Detail::RPageSink > fSink
Definition RNTuple.hxx:216
static std::unique_ptr< RNTupleWriter > Recreate(std::unique_ptr< RNTupleModel > model, std::string_view ntupleName, std::string_view storage, const RNTupleWriteOptions &options=RNTupleWriteOptions())
Definition RNTuple.cxx:265
static std::unique_ptr< RNTupleWriter > Append(std::unique_ptr< RNTupleModel > model, std::string_view ntupleName, TFile &file, const RNTupleWriteOptions &options=RNTupleWriteOptions())
Definition RNTuple.cxx:274
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.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:54
ENTupleInfo
Listing of the different options that can be printed by RNTupleReader::GetInfo()
Definition RNTuple.hxx:50
std::uint64_t NTupleSize_t
Integer type long enough to hold the maximum number of entries in a column.
ENTupleShowFormat
Listing of the different entry output formats of RNTupleReader::Show()
Definition RNTuple.hxx:59
constexpr DescriptorId_t kInvalidDescriptorId
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:556
Definition file.py:1
static const char * what
Definition stlLoader.cc:6
static void output(int code)
Definition gifencode.c:226