Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RNTupleInspector.hxx
Go to the documentation of this file.
1/// \file ROOT/RNTupleInspector.hxx
2/// \ingroup NTuple ROOT7
3/// \author Florine de Geus <florine.de.geus@cern.ch>
4/// \date 2023-01-09
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-2023, 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#ifndef ROOT7_RNTupleInspector
17#define ROOT7_RNTupleInspector
18
19#include <ROOT/RError.hxx>
21
22#include <TFile.h>
23#include <TH1D.h>
24#include <THStack.h>
25
26#include <cstdlib>
27#include <iostream>
28#include <memory>
29#include <numeric>
30#include <optional>
31#include <regex>
32#include <vector>
33
34namespace ROOT {
35class RNTuple;
36
37namespace Internal {
38class RPageSource;
39} // namespace Internal
40
41namespace Experimental {
42
45
46// clang-format off
47/**
48\class ROOT::Experimental::RNTupleInspector
49\ingroup NTuple
50\brief Inspect on-disk and storage-related information of an RNTuple.
51
52The RNTupleInspector can be used for studying an RNTuple in terms of its storage efficiency. It provides information on
53the level of the RNTuple itself, on the (sub)field level and on the column level.
54
55Example usage:
56
57~~~ {.cpp}
58#include <ROOT/RNTuple.hxx>
59#include <ROOT/RNTupleInspector.hxx>
60
61#include <iostream>
62
63using ROOT::Experimental::RNTuple;
64using ROOT::Experimental::RNTupleInspector;
65
66auto file = TFile::Open("data.rntuple");
67auto rntuple = std::unique_ptr<RNTuple>(file->Get<RNTuple>("NTupleName"));
68auto inspector = RNTupleInspector::Create(*rntuple);
69
70std::cout << "The compression factor is " << inspector->GetCompressionFactor()
71 << " using compression settings " << inspector->GetCompressionSettingsAsString()
72 << std::endl;
73~~~
74*/
75// clang-format on
77public:
78 /////////////////////////////////////////////////////////////////////////////
79 /// \brief Provides column-level storage information.
80 ///
81 /// The RColumnInspector class provides storage information for an individual column. This information is partly
82 /// collected during the construction of the RNTupleInspector object, and can partly be accessed using the
83 /// RColumnInspector that belongs to this field.
85 private:
87 const std::vector<std::uint64_t> fCompressedPageSizes = {};
88 std::uint32_t fElementSize = 0;
89 std::uint64_t fNElements = 0;
90
91 public:
98 ~RColumnInspector() = default;
99
101 const std::vector<std::uint64_t> &GetCompressedPageSizes() const { return fCompressedPageSizes; }
102 std::uint64_t GetNPages() const { return fCompressedPageSizes.size(); }
103 std::uint64_t GetCompressedSize() const
104 {
105 return std::accumulate(fCompressedPageSizes.begin(), fCompressedPageSizes.end(), static_cast<std::uint64_t>(0));
106 }
107 std::uint64_t GetUncompressedSize() const { return fElementSize * fNElements; }
108 std::uint64_t GetElementSize() const { return fElementSize; }
109 std::uint64_t GetNElements() const { return fNElements; }
111 };
112
113 /////////////////////////////////////////////////////////////////////////////
114 /// \brief Provides field-level storage information.
115 ///
116 /// The RFieldTreeInspector class provides storage information for a field **and** its subfields. This information is
117 /// partly collected during the construction of the RNTupleInspector object, and can partly be accessed using
118 /// the RFieldDescriptor that belongs to this field.
134
135private:
136 std::unique_ptr<ROOT::Internal::RPageSource> fPageSource;
138 std::optional<std::uint32_t> fCompressionSettings; ///< The compression settings are unknown for an empty ntuple
139 std::uint64_t fCompressedSize = 0;
140 std::uint64_t fUncompressedSize = 0;
141
142 std::unordered_map<int, RColumnInspector> fColumnInfo;
143 std::unordered_map<int, RFieldTreeInspector> fFieldTreeInfo;
144
145 RNTupleInspector(std::unique_ptr<ROOT::Internal::RPageSource> pageSource);
146
147 /////////////////////////////////////////////////////////////////////////////
148 /// \brief Gather column-level and RNTuple-level information.
149 ///
150 /// \note This method is called when the RNTupleInspector is initially created. This means that anything unexpected
151 /// about the RNTuple itself (e.g. inconsistent compression settings across clusters) will be detected here.
152 /// Therefore, any related exceptions will be thrown on creation of the inspector.
153 void CollectColumnInfo();
154
155 /////////////////////////////////////////////////////////////////////////////
156 /// \brief Recursively gather field-level information.
157 ///
158 /// \param[in] fieldId The ID of the field from which to start the recursive traversal. Typically this is the "zero
159 /// ID", i.e. the logical parent of all top-level fields.
160 ///
161 /// \return The RFieldTreeInspector for the provided field ID.
162 ///
163 /// This method is called when the RNTupleInspector is initially created.
165
166 /////////////////////////////////////////////////////////////////////////////
167 /// \brief Get the columns that make up the given field, including its subfields.
168 ///
169 /// \param [in] fieldId The ID of the field for which to collect the columns.
170 ///
171 /// \return A vector containing the IDs of all columns for the provided field ID.
172 std::vector<ROOT::DescriptorId_t> GetColumnsByFieldId(ROOT::DescriptorId_t fieldId) const;
173
174public:
180
181 /////////////////////////////////////////////////////////////////////////////
182 /// \brief Create a new RNTupleInspector.
183 ///
184 /// \param[in] sourceNTuple A pointer to the RNTuple to be inspected.
185 ///
186 /// \return A pointer to the newly created RNTupleInspector.
187 ///
188 /// \note When this factory method is called, all required static information is collected from the RNTuple's fields
189 /// and underlying columns are collected at ones. This means that when any inconsistencies are encountered (e.g.
190 /// inconsistent compression across clusters), it will throw an error here.
191 static std::unique_ptr<RNTupleInspector> Create(const RNTuple &sourceNTuple);
192
193 /////////////////////////////////////////////////////////////////////////////
194 /// \brief Create a new RNTupleInspector.
195 ///
196 /// \param[in] ntupleName The name of the RNTuple to be inspected.
197 /// \param[in] storage The path or URI to the RNTuple to be inspected.
198 ///
199 /// \see Create(RNTuple *sourceNTuple)
200 static std::unique_ptr<RNTupleInspector> Create(std::string_view ntupleName, std::string_view storage);
201
202 /////////////////////////////////////////////////////////////////////////////
203 /// \brief Get the descriptor for the RNTuple being inspected.
204 ///
205 /// \return A static copy of the ROOT::RNTupleDescriptor belonging to the inspected RNTuple.
207
208 /////////////////////////////////////////////////////////////////////////////
209 /// \brief Get the compression settings of the RNTuple being inspected.
210 ///
211 /// \return The integer representation (\f$algorithm * 10 + level\f$, where \f$algorithm\f$ follows
212 /// ROOT::RCompressionSetting::ELevel::EValues) of the compression settings used for the inspected RNTuple.
213 /// Empty for an empty ntuple.
214 ///
215 /// \note Here, we assume that the compression settings are consistent across all clusters and columns. If this is
216 /// not the case, an exception will be thrown when RNTupleInspector::Create is called.
217 std::optional<std::uint32_t> GetCompressionSettings() const { return fCompressionSettings; }
218
219 /////////////////////////////////////////////////////////////////////////////
220 /// \brief Get a string describing compression settings of the RNTuple being inspected.
221 ///
222 /// \return A string describing the compression used for the inspected RNTuple. The format of the string is
223 /// `"A (level L)"`, where `A` is the name of the compression algorithm and `L` the compression level.
224 ///
225 /// \note Here, we assume that the compression settings are consistent across all clusters and columns. If this is
226 /// not the case, an exception will be thrown when RNTupleInspector::Create is called.
227 std::string GetCompressionSettingsAsString() const;
228
229 /////////////////////////////////////////////////////////////////////////////
230 /// \brief Get the compressed, on-disk size of the RNTuple being inspected.
231 ///
232 /// \return The compressed size of the inspected RNTuple, in bytes, excluding the size of the header and footer.
233 std::uint64_t GetCompressedSize() const { return fCompressedSize; }
234
235 /////////////////////////////////////////////////////////////////////////////
236 /// \brief Get the uncompressed total size of the RNTuple being inspected.
237 ///
238 /// \return The uncompressed size of the inspected RNTuple, in bytes, excluding the size of the header and footer.
239 std::uint64_t GetUncompressedSize() const { return fUncompressedSize; }
240
241 /////////////////////////////////////////////////////////////////////////////
242 /// \brief Get the compression factor of the RNTuple being inspected.
243 ///
244 /// \return The compression factor of the inspected RNTuple.
245 ///
246 /// The compression factor shows how well the data present in the RNTuple is compressed by the compression settings
247 /// that were used. The compression factor is calculated as \f$size_{uncompressed} / size_{compressed}\f$.
248 float GetCompressionFactor() const { return (float)fUncompressedSize / (float)fCompressedSize; }
249
250 /////////////////////////////////////////////////////////////////////////////
251 /// \brief Get storage information for a given column.
252 ///
253 /// \param[in] physicalColumnId The physical ID of the column for which to get the information.
254 ///
255 /// \return The storage information for the provided column.
256 const RColumnInspector &GetColumnInspector(ROOT::DescriptorId_t physicalColumnId) const;
257
258 /////////////////////////////////////////////////////////////////////////////
259 /// \brief Get the number of columns of a given type present in the RNTuple.
260 ///
261 /// \param[in] colType The column type to count, as defined by ROOT::ENTupleColumnType.
262 ///
263 /// \return The number of columns present in the inspected RNTuple of the provided type.
265
266 /////////////////////////////////////////////////////////////////////////////
267 /// \brief Get the IDs of all columns with the given type.
268 ///
269 /// \param[in] colType The column type to collect, as defined by ROOT::ENTupleColumnType.
270 ///
271 /// \return A vector containing the physical IDs of columns of the provided type.
272 const std::vector<ROOT::DescriptorId_t> GetColumnsByType(ROOT::ENTupleColumnType colType);
273
274 /////////////////////////////////////////////////////////////////////////////
275 /// \brief Get all column types present in the RNTuple being inspected.
276 ///
277 /// \return A vector containing all column types present in the RNTuple.
278 const std::vector<ROOT::ENTupleColumnType> GetColumnTypes();
279
280 /////////////////////////////////////////////////////////////////////////////
281 /// \brief Print storage information per column type.
282 ///
283 /// \param[in] format Whether to print the information as a (markdown-parseable) table or in CSV format.
284 /// \param[in] output Where to write the output to. Default is `stdout`.
285 ///
286 /// The output includes for each column type its count, the total number of elements, the compressed size and the
287 /// uncompressed size.
288 ///
289 /// **Example: printing the column type information of an RNTuple as a table**
290 /// ~~~ {.cpp}
291 /// #include <ROOT/RNTupleInspector.hxx>
292 /// using ROOT::Experimental::RNTupleInspector;
293 /// using ROOT::Experimental::ENTupleInspectorPrintFormat;
294 ///
295 /// auto inspector = RNTupleInspector::Create("myNTuple", "some/file.root");
296 /// inspector->PrintColumnTypeInfo();
297 /// ~~~
298 /// Output:
299 /// ~~~
300 /// column type | count | # elements | compressed bytes | uncompressed bytes
301 /// ----------------|---------|-----------------|-------------------|--------------------
302 /// SplitIndex64 | 2 | 150 | 72 | 1200
303 /// SplitReal32 | 4 | 300 | 189 | 1200
304 /// SplitUInt32 | 3 | 225 | 123 | 900
305 /// ~~~
306 ///
307 /// **Example: printing the column type information of an RNTuple in CSV format**
308 /// ~~~ {.cpp}
309 /// #include <ROOT/RNTupleInspector.hxx>
310 /// using ROOT::Experimental::RNTupleInspector;
311 /// using ROOT::Experimental::ENTupleInspectorPrintFormat;
312 ///
313 /// auto inspector = RNTupleInspector::Create("myNTuple", "some/file.root");
314 /// inspector->PrintColumnTypeInfo();
315 /// ~~~
316 /// Output:
317 /// ~~~
318 /// columnType,count,nElements,compressedSize,uncompressedSize
319 /// SplitIndex64,2,150,72,1200
320 /// SplitReal32,4,300,189,1200
321 /// SplitUInt32,3,225,123,900
322 /// ~~~
324 std::ostream &output = std::cout);
325
326 /////////////////////////////////////////////////////////////////////////////
327 /// \brief Get a histogram showing information for each column type present,
328 ///
329 /// \param[in] histKind Which type of information should be returned.
330 /// \param[in] histName The name of the histogram. An empty string means a default name will be used.
331 /// \param[in] histTitle The title of the histogram. An empty string means a default title will be used.
332 ///
333 /// \return A pointer to a `TH1D` containing the specified kind of information.
334 ///
335 /// Get a histogram showing the count, number of elements, size on disk, or size in memory for each column
336 /// type present in the inspected RNTuple.
337 std::unique_ptr<TH1D> GetColumnTypeInfoAsHist(ENTupleInspectorHist histKind, std::string_view histName = "",
338 std::string_view histTitle = "");
339
340 /////////////////////////////////////////////////////////////////////////////
341 /// \brief Get a histogram containing the size distribution of the compressed pages for an individual column.
342 ///
343 /// \param[in] physicalColumnId The physical ID of the column for which to get the page size distribution.
344 /// \param[in] histName The name of the histogram. An empty string means a default name will be used.
345 /// \param[in] histTitle The title of the histogram. An empty string means a default title will be used.
346 /// \param[in] nBins The desired number of histogram bins.
347 ///
348 /// \return A pointer to a `TH1D` containing the page size distribution.
349 ///
350 /// The x-axis will range from the smallest page size, to the largest (inclusive).
351 std::unique_ptr<TH1D> GetPageSizeDistribution(ROOT::DescriptorId_t physicalColumnId, std::string histName = "",
352 std::string histTitle = "", size_t nBins = 64);
353
354 /////////////////////////////////////////////////////////////////////////////
355 /// \brief Get a histogram containing the size distribution of the compressed pages for all columns of a given type.
356 ///
357 /// \param[in] colType The column type for which to get the size distribution, as defined by ROOT::ENTupleColumnType.
358 /// \param[in] histName The name of the histogram. An empty string means a default name will be used.
359 /// \param[in] histTitle The title of the histogram. An empty string means a default title will be used.
360 /// \param[in] nBins The desired number of histogram bins.
361 ///
362 /// \return A pointer to a `TH1D` containing the page size distribution.
363 ///
364 /// The x-axis will range from the smallest page size, to the largest (inclusive).
365 std::unique_ptr<TH1D> GetPageSizeDistribution(ROOT::ENTupleColumnType colType, std::string histName = "",
366 std::string histTitle = "", size_t nBins = 64);
367
368 /////////////////////////////////////////////////////////////////////////////
369 /// \brief Get a histogram containing the size distribution of the compressed pages for a collection columns.
370 ///
371 /// \param[in] colIds The physical IDs of the columns for which to get the page size distribution.
372 /// \param[in] histName The name of the histogram. An empty string means a default name will be used.
373 /// \param[in] histTitle The title of the histogram. An empty string means a default title will be used.
374 /// \param[in] nBins The desired number of histogram bins.
375 ///
376 /// \return A pointer to a `TH1D` containing the (cumulative) page size distribution.
377 ///
378 /// The x-axis will range from the smallest page size, to the largest (inclusive).
379 std::unique_ptr<TH1D> GetPageSizeDistribution(std::initializer_list<ROOT::DescriptorId_t> colIds,
380 std::string histName = "", std::string histTitle = "",
381 size_t nBins = 64);
382
383 /////////////////////////////////////////////////////////////////////////////
384 /// \brief Get a histogram containing the size distribution of the compressed pages for all columns of a given list
385 /// of types.
386 ///
387 /// \param[in] colTypes The column types for which to get the size distribution, as defined by
388 /// ROOT::ENTupleColumnType. The default is an empty vector, which indicates that the distribution
389 /// for *all* physical columns will be returned.
390 /// \param[in] histName The name of the histogram. An empty string means a default name will be used. The name of
391 /// each histogram inside the `THStack` will be `histName + colType`.
392 /// \param[in] histTitle The title of the histogram. An empty string means a default title will be used.
393 /// \param[in] nBins The desired number of histogram bins.
394 ///
395 /// \return A pointer to a `THStack` with one histogram for each column type.
396 ///
397 /// The x-axis will range from the smallest page size, to the largest (inclusive).
398 ///
399 /// **Example: Drawing a non-stacked page size distribution with a legend**
400 /// ~~~ {.cpp}
401 /// auto canvas = std::make_unique<TCanvas>();
402 /// auto inspector = RNTupleInspector::Create("myNTuple", "ntuple.root");
403 ///
404 /// // We want to show the page size distributions of columns with type `kSplitReal32` and `kSplitReal64`.
405 /// auto hist = inspector->GetPageSizeDistribution(
406 /// {ROOT::ENTupleColumnType::kSplitReal32, ROOT::ENTupleColumnType::kSplitReal64});
407 /// // The "PLC" option automatically sets the line color for each histogram in the `THStack`.
408 /// // The "NOSTACK" option will draw the histograms on top of each other instead of stacked.
409 /// hist->DrawClone("PLC NOSTACK");
410 /// canvas->BuildLegend(0.7, 0.8, 0.89, 0.89);
411 /// canvas->DrawClone();
412 /// ~~~
413 std::unique_ptr<THStack> GetPageSizeDistribution(std::initializer_list<ROOT::ENTupleColumnType> colTypes = {},
414 std::string histName = "", std::string histTitle = "",
415 size_t nBins = 64);
416
417 /////////////////////////////////////////////////////////////////////////////
418 /// \brief Get storage information for a given (sub)field by ID.
419 ///
420 /// \param[in] fieldId The ID of the (sub)field for which to get the information.
421 ///
422 /// \return The storage information inspector for the provided (sub)field tree.
423 const RFieldTreeInspector &GetFieldTreeInspector(ROOT::DescriptorId_t fieldId) const;
424
425 /////////////////////////////////////////////////////////////////////////////
426 /// \brief Get a storage information inspector for a given (sub)field by name, including its subfields.
427 ///
428 /// \param[in] fieldName The name of the (sub)field for which to get the information.
429 ///
430 /// \return The storage information inspector for the provided (sub)field tree.
431 const RFieldTreeInspector &GetFieldTreeInspector(std::string_view fieldName) const;
432
433 /////////////////////////////////////////////////////////////////////////////
434 /// \brief Get the number of fields of a given type or class present in the RNTuple.
435 ///
436 /// \param[in] typeNamePattern The type or class name to count. May contain regular expression patterns for grouping
437 /// multiple kinds of types or classes.
438 /// \param[in] searchInSubfields If set to `false`, only top-level fields will be considered.
439 ///
440 /// \return The number of fields that matches the provided type.
441 size_t GetFieldCountByType(const std::regex &typeNamePattern, bool searchInSubfields = true) const;
442
443 /////////////////////////////////////////////////////////////////////////////
444 /// \brief Get the number of fields of a given type or class present in the RNTuple.
445 ///
446 /// \see GetFieldCountByType(const std::regex &typeNamePattern, bool searchInSubfields) const
447 size_t GetFieldCountByType(std::string_view typeNamePattern, bool searchInSubfields = true) const
448 {
449 return GetFieldCountByType(std::regex{std::string(typeNamePattern)}, searchInSubfields);
450 }
451
452 /////////////////////////////////////////////////////////////////////////////
453 /// \brief Get the IDs of (sub-)fields whose name matches the given string.
454 ///
455 /// \param[in] fieldNamePattern The name of the field name to get. Because field names are unique by design,
456 /// providing a single field name will return a vector containing just the ID of that field. However, regular
457 /// expression patterns are supported in order to get the IDs of all fields whose name follow a certain structure.
458 /// \param[in] searchInSubfields If set to `false`, only top-level fields will be considered.
459 ///
460 /// \return A vector containing the IDs of fields that match the provided name.
461 const std::vector<ROOT::DescriptorId_t>
462 GetFieldsByName(const std::regex &fieldNamePattern, bool searchInSubfields = true) const;
463
464 /////////////////////////////////////////////////////////////////////////////
465 /// \brief Get the IDs of (sub-)fields whose name matches the given string.
466 ///
467 /// \see GetFieldsByName(const std::regex &fieldNamePattern, bool searchInSubfields) const
468 const std::vector<ROOT::DescriptorId_t>
470 {
471 return GetFieldsByName(std::regex{std::string(fieldNamePattern)}, searchInSubfields);
472 }
473};
474} // namespace Experimental
475} // namespace ROOT
476
477#endif // ROOT7_RNTupleInspector
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format
The available trivial, native content types of a column.
Provides column-level storage information.
RColumnInspector(const ROOT::RColumnDescriptor &colDesc, const std::vector< std::uint64_t > &compressedPageSizes, std::uint32_t elemSize, std::uint64_t nElems)
const ROOT::RColumnDescriptor & GetDescriptor() const
const std::vector< std::uint64_t > fCompressedPageSizes
const std::vector< std::uint64_t > & GetCompressedPageSizes() const
RFieldTreeInspector(const ROOT::RFieldDescriptor &fieldDesc, std::uint64_t onDiskSize, std::uint64_t inMemSize)
Inspect on-disk and storage-related information of an RNTuple.
std::vector< ROOT::DescriptorId_t > GetColumnsByFieldId(ROOT::DescriptorId_t fieldId) const
Get the columns that make up the given field, including its subfields.
const std::vector< ROOT::DescriptorId_t > GetFieldsByName(std::string_view fieldNamePattern, bool searchInSubfields=true)
Get the IDs of (sub-)fields whose name matches the given string.
float GetCompressionFactor() const
Get the compression factor of the RNTuple being inspected.
RNTupleInspector & operator=(RNTupleInspector &&other)=delete
const RFieldTreeInspector & GetFieldTreeInspector(ROOT::DescriptorId_t fieldId) const
Get storage information for a given (sub)field by ID.
std::unique_ptr< TH1D > GetPageSizeDistribution(ROOT::DescriptorId_t physicalColumnId, std::string histName="", std::string histTitle="", size_t nBins=64)
Get a histogram containing the size distribution of the compressed pages for an individual column.
const ROOT::RNTupleDescriptor & GetDescriptor() const
Get the descriptor for the RNTuple being inspected.
RNTupleInspector(const RNTupleInspector &other)=delete
std::uint64_t GetCompressedSize() const
Get the compressed, on-disk size of the RNTuple being inspected.
size_t GetColumnCountByType(ROOT::ENTupleColumnType colType) const
Get the number of columns of a given type present in the RNTuple.
std::uint64_t GetUncompressedSize() const
Get the uncompressed total size of the RNTuple being inspected.
RNTupleInspector(RNTupleInspector &&other)=delete
std::optional< std::uint32_t > fCompressionSettings
The compression settings are unknown for an empty ntuple.
const std::vector< ROOT::ENTupleColumnType > GetColumnTypes()
Get all column types present in the RNTuple being inspected.
size_t GetFieldCountByType(const std::regex &typeNamePattern, bool searchInSubfields=true) const
Get the number of fields of a given type or class present in the RNTuple.
const std::vector< ROOT::DescriptorId_t > GetFieldsByName(const std::regex &fieldNamePattern, bool searchInSubfields=true) const
Get the IDs of (sub-)fields whose name matches the given string.
std::string GetCompressionSettingsAsString() const
Get a string describing compression settings of the RNTuple being inspected.
RFieldTreeInspector CollectFieldTreeInfo(ROOT::DescriptorId_t fieldId)
Recursively gather field-level information.
RNTupleInspector(std::unique_ptr< ROOT::Internal::RPageSource > pageSource)
size_t GetFieldCountByType(std::string_view typeNamePattern, bool searchInSubfields=true) const
Get the number of fields of a given type or class present in the RNTuple.
void PrintColumnTypeInfo(ENTupleInspectorPrintFormat format=ENTupleInspectorPrintFormat::kTable, std::ostream &output=std::cout)
Print storage information per column type.
std::optional< std::uint32_t > GetCompressionSettings() const
Get the compression settings of the RNTuple being inspected.
const RColumnInspector & GetColumnInspector(ROOT::DescriptorId_t physicalColumnId) const
Get storage information for a given column.
std::unique_ptr< ROOT::Internal::RPageSource > fPageSource
RNTupleInspector & operator=(const RNTupleInspector &other)=delete
static std::unique_ptr< RNTupleInspector > Create(const RNTuple &sourceNTuple)
Create a new RNTupleInspector.
std::unordered_map< int, RFieldTreeInspector > fFieldTreeInfo
const std::vector< ROOT::DescriptorId_t > GetColumnsByType(ROOT::ENTupleColumnType colType)
Get the IDs of all columns with the given type.
void CollectColumnInfo()
Gather column-level and RNTuple-level information.
std::unordered_map< int, RColumnInspector > fColumnInfo
std::unique_ptr< TH1D > GetColumnTypeInfoAsHist(ENTupleInspectorHist histKind, std::string_view histName="", std::string_view histTitle="")
Get a histogram showing information for each column type present,.
Metadata stored for every column of an RNTuple.
Metadata stored for every field of an RNTuple.
The on-storage metadata of an RNTuple.
Representation of an RNTuple data set in a ROOT file.
Definition RNTuple.hxx:65
const_iterator begin() const
const_iterator end() const
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
std::uint64_t DescriptorId_t
Distriniguishes elements of the same type within a descriptor, e.g. different fields.
static void output()