Logo ROOT  
Reference Guide
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"
21
22#include <algorithm>
23#include <iomanip>
24#include <iostream>
25#include <sstream>
26#include <string>
27#include <unordered_map>
28#include <utility>
29
30#include <TFile.h>
32
33
34ROOT::Experimental::Detail::RNTuple::RNTuple(std::unique_ptr<ROOT::Experimental::RNTupleModel> model)
35 : fModel(std::move(model))
36 , fNEntries(0)
37{
38}
39
41{
42}
43
44//------------------------------------------------------------------------------
45
47 std::unordered_map<const Detail::RFieldBase *, DescriptorId_t> fieldPtr2Id;
48 fieldPtr2Id[fModel->GetRootField()] = fSource->GetDescriptor().FindFieldId("", kInvalidDescriptorId);
49 for (auto &field : *fModel->GetRootField()) {
50 auto parentId = fieldPtr2Id[field.GetParent()];
51 auto fieldId = fSource->GetDescriptor().FindFieldId(field.GetName(), parentId);
53 fieldPtr2Id[&field] = fieldId;
54 Detail::RFieldFuse::Connect(fieldId, *fSource, field);
55 }
56}
57
59 std::unique_ptr<ROOT::Experimental::RNTupleModel> model,
60 std::unique_ptr<ROOT::Experimental::Detail::RPageSource> source)
61 : ROOT::Experimental::Detail::RNTuple(std::move(model))
62 , fSource(std::move(source))
63 , fMetrics("RNTupleReader")
64{
65 fSource->Attach();
67 fNEntries = fSource->GetNEntries();
68 fMetrics.ObserveMetrics(fSource->GetMetrics());
69}
70
71ROOT::Experimental::RNTupleReader::RNTupleReader(std::unique_ptr<ROOT::Experimental::Detail::RPageSource> source)
72 : ROOT::Experimental::Detail::RNTuple(nullptr)
73 , fSource(std::move(source))
74 , fMetrics("RNTupleReader")
75{
76 fSource->Attach();
77 fModel = fSource->GetDescriptor().GenerateModel();
79 fNEntries = fSource->GetNEntries();
80 fMetrics.ObserveMetrics(fSource->GetMetrics());
81}
82
84{
85 // needs to be destructed before the page source
86 fModel = nullptr;
87}
88
89std::unique_ptr<ROOT::Experimental::RNTupleReader> ROOT::Experimental::RNTupleReader::Open(
90 std::unique_ptr<RNTupleModel> model,
91 std::string_view ntupleName,
92 std::string_view storage)
93{
94 return std::make_unique<RNTupleReader>(std::move(model), Detail::RPageSource::Create(ntupleName, storage));
95}
96
97std::unique_ptr<ROOT::Experimental::RNTupleReader> ROOT::Experimental::RNTupleReader::Open(
98 std::string_view ntupleName,
99 std::string_view storage)
100{
101 return std::make_unique<RNTupleReader>(Detail::RPageSource::Create(ntupleName, storage));
102}
103
105{
106 // 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.
107 char frameSymbol = '*';
108 int width = 80;
109 /*
110 if (width < 30) {
111 output << "The width is too small! Should be at least 30." << std::endl;
112 return;
113 }
114 */
115 std::string name = fSource->GetDescriptor().GetName();
116 //prepVisitor traverses through all fields to gather information needed for printing.
117 RPrepareVisitor prepVisitor;
118 //printVisitor traverses through all fields to do the actual printing.
119 RPrintVisitor printVisitor(output);
120 switch (what) {
122 for (int i = 0; i < (width/2 + width%2 - 4); ++i)
123 output << frameSymbol;
124 output << " NTUPLE ";
125 for (int i = 0; i < (width/2 - 4); ++i)
126 output << frameSymbol;
127 output << std::endl;
128 // FitString defined in RFieldVisitor.cxx
129 output << frameSymbol << " N-Tuple : " << RNTupleFormatter::FitString(name, width-13) << frameSymbol << std::endl; // prints line with name of ntuple
130 output << frameSymbol << " Entries : " << RNTupleFormatter::FitString(std::to_string(GetNEntries()), width - 13) << frameSymbol << std::endl; // prints line with number of entries
131 GetModel()->GetRootField()->TraverseVisitor(prepVisitor);
132
133 printVisitor.SetFrameSymbol(frameSymbol);
134 printVisitor.SetWidth(width);
135 printVisitor.SetDeepestLevel(prepVisitor.GetDeepestLevel());
136 printVisitor.SetNumFields(prepVisitor.GetNumFields());
137 GetModel()->GetRootField()->TraverseVisitor(printVisitor);
138
139 for (int i = 0; i < width; ++i)
140 output << frameSymbol;
141 output << std::endl;
142 break;
144 fSource->GetDescriptor().PrintInfo(output);
145 break;
147 fMetrics.Print(output);
148 break;
149 default:
150 // Unhandled case, internal error
151 assert(false);
152 }
153}
154//------------------------------------------------------------------------------
155
157 std::unique_ptr<ROOT::Experimental::RNTupleModel> model,
158 std::unique_ptr<ROOT::Experimental::Detail::RPageSink> sink)
159 : ROOT::Experimental::Detail::RNTuple(std::move(model))
160 , fSink(std::move(sink))
161 , fClusterSizeEntries(kDefaultClusterSizeEntries)
162 , fLastCommitted(0)
163{
164 fSink->Create(*fModel.get());
165}
166
168{
169 CommitCluster();
170 fSink->CommitDataset();
171 // needs to be destructed before the page sink
172 fModel = nullptr;
173}
174
175std::unique_ptr<ROOT::Experimental::RNTupleWriter> ROOT::Experimental::RNTupleWriter::Recreate(
176 std::unique_ptr<RNTupleModel> model,
177 std::string_view ntupleName,
178 std::string_view storage,
179 const RNTupleWriteOptions &options)
180{
181 return std::make_unique<RNTupleWriter>(std::move(model), Detail::RPageSink::Create(ntupleName, storage, options));
182}
183
184
186{
187 if (fNEntries == fLastCommitted) return;
188 for (auto& field : *fModel->GetRootField()) {
189 field.Flush();
190 field.CommitCluster();
191 }
192 fSink->CommitCluster(fNEntries);
193 fLastCommitted = fNEntries;
194}
195
196
197//------------------------------------------------------------------------------
198
199
201 : fOffset(0), fDefaultEntry(std::move(defaultEntry))
202{
203}
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
static void Connect(DescriptorId_t fieldId, RPageStorage &pageStorage, RFieldBase &field)
Definition: RField.cxx:76
void ObserveMetrics(RNTupleMetrics &observee)
NTupleSize_t fNEntries
The number of entries is constant for reading and reflects the sum of Fill() operations when writing.
Definition: RNTuple.hxx:62
RNTuple(std::unique_ptr< RNTupleModel > model)
Only the derived RNTupleReader and RNTupleWriter can be instantiated.
Definition: RNTuple.cxx:34
std::unique_ptr< RNTupleModel > fModel
Definition: RNTuple.hxx:60
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:200
static std::string FitString(const std::string &str, int availableSpace)
static std::unique_ptr< RNTupleReader > Open(std::unique_ptr< RNTupleModel > model, std::string_view ntupleName, std::string_view storage)
Definition: RNTuple.cxx:89
Detail::RNTupleMetrics fMetrics
Definition: RNTuple.hxx:103
std::unique_ptr< Detail::RPageSource > fSource
Definition: RNTuple.hxx:102
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:58
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:104
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:185
RNTupleWriter(std::unique_ptr< RNTupleModel > model, std::unique_ptr< Detail::RPageSink > sink)
Definition: RNTuple.cxx:156
std::unique_ptr< Detail::RPageSink > fSink
Definition: RNTuple.hxx:191
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:175
The RNTuple represents a live dataset, whose structure is defined by an RNTupleModel.
Visitor used for a prepare run to collect information needed by another visitor class.
Contains settings for printing and prints a summary of an RField instance.
basic_string_view< char > string_view
ENTupleInfo
Listing of the different options that can be returned by RNTupleReader::GetInfo()
Definition: RNTuple.hxx:81
constexpr DescriptorId_t kInvalidDescriptorId
Definition: RNTupleUtil.hxx:80
VSD Structures.
Definition: StringConv.hxx:21
static void output(int code)
Definition: gifencode.c:226