Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RNTupleTreeMap.cxx
Go to the documentation of this file.
1/// \file ROOT/RNTupleTreeMap.cxx
2/// \ingroup NTuple
3/// \author Patryk Tymoteusz Pilichowski <patryk.tymoteusz.pilichowski@cern.ch>
4/// \date 2025-09-15
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-2025, 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
18
19#include <queue>
20
23 std::uint64_t childrenIdx, std::uint64_t nChildren, ROOT::DescriptorId_t rootId, size_t rootSize)
24{
25 uint64_t size =
26 (rootId != fldDesc.GetId()) ? insp.GetFieldTreeInspector(fldDesc.GetId()).GetCompressedSize() : rootSize;
27 return {fldDesc.GetFieldName(), "", size, childrenIdx, nChildren};
28}
29
36
37std::unique_ptr<ROOT::Experimental::RTreeMapPainter>
39{
40 auto treemap = std::make_unique<ROOT::Experimental::RTreeMapPainter>();
41 const auto &descriptor = insp.GetDescriptor();
42 const auto rootId = descriptor.GetFieldZero().GetId();
43 size_t rootSize = 0;
44 for (const auto &childId : descriptor.GetFieldDescriptor(rootId).GetLinkIds()) {
45 rootSize += insp.GetFieldTreeInspector(childId).GetCompressedSize();
46 }
47
48 std::queue<std::pair<uint64_t, bool>> queue; // (columnid/fieldid, isfield)
49 queue.emplace(rootId, true);
50 while (!queue.empty()) {
51 size_t levelSize = queue.size();
52 size_t levelChildrenStart = treemap->fNodes.size() + levelSize;
53 for (size_t i = 0; i < levelSize; ++i) {
54 const auto &current = queue.front();
55 queue.pop();
56
57 size_t nChildren = 0;
58 if (current.second) {
59 std::vector<uint64_t> children;
60 const auto &fldDesc = descriptor.GetFieldDescriptor(current.first);
61 children = fldDesc.GetLinkIds();
62 for (const auto childId : children) {
63 queue.emplace(childId, 1);
64 }
65 for (const auto &columnDesc : descriptor.GetColumnIterable(fldDesc.GetId())) {
66 const auto &columnId = columnDesc.GetPhysicalId();
67 children.push_back(columnId);
68 queue.emplace(columnId, 0);
69 }
70 nChildren = children.size();
72 treemap->fNodes.push_back(node);
73 } else {
74 const auto &colInsp = insp.GetColumnInspector(current.first);
75 const auto &node = CreateNode(colInsp, levelChildrenStart);
76 treemap->fNodes.push_back(node);
77 }
78
80 }
81 }
82 return treemap;
83}
84
85std::unique_ptr<ROOT::Experimental::RTreeMapPainter>
static ROOT::Experimental::RTreeMapBase::Node CreateNode(const ROOT::Experimental::RNTupleInspector &insp, const ROOT::RFieldDescriptor &fldDesc, std::uint64_t childrenIdx, std::uint64_t nChildren, ROOT::DescriptorId_t rootId, size_t rootSize)
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Provides column-level storage information.
Inspect on-disk and storage-related information of an RNTuple.
static std::unique_ptr< RNTupleInspector > Create(const RNTuple &sourceNTuple)
Create a new RNTupleInspector.
static const char * GetColumnTypeName(ROOT::ENTupleColumnType type)
Metadata stored for every field of an RNTuple.
std::unique_ptr< RTreeMapPainter > CreateTreeMapFromRNTuple(const RNTupleInspector &insp)
Logic for converting an RNTuple to RTreeMapPainter given RNTupleInspector.
std::uint64_t DescriptorId_t
Distriniguishes elements of the same type within a descriptor, e.g. different fields.