Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RNTupleClassicBrowse.cxx
Go to the documentation of this file.
1/// \file RNTupleClassicBrowse.cxx
2/// \ingroup NTuple
3/// \author Jakob Blomer <jblomer@cern.ch>
4/// \date 2025-06-30
5
6/*************************************************************************
7 * Copyright (C) 1995-2025, 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
14#include <ROOT/RNTuple.hxx>
21
22#include <TBrowser.h>
23#include <TObject.h>
24#include <TPad.h>
25#include <TText.h>
26
27namespace {
28
29class RFieldBrowsable final : public TObject {
30private:
31 std::shared_ptr<ROOT::RNTupleReader> fReader;
34 bool fIsLeaf = false;
35 std::unique_ptr<TH1> fHistogram;
36 std::string fFieldName;
37 std::string fTypeName;
38
39public:
40 RFieldBrowsable(std::shared_ptr<ROOT::RNTupleReader> reader, ROOT::DescriptorId_t fieldId)
41 : fReader(reader), fFieldId(fieldId)
42 {
43 const auto &desc = fReader->GetDescriptor();
44 fBrowsableFieldId = ROOT::Internal::GetNextBrowsableField(fFieldId, desc);
45 fIsLeaf = desc.GetFieldDescriptor(fBrowsableFieldId).GetLinkIds().empty();
46 fFieldName = desc.GetFieldDescriptor(fFieldId).GetFieldName();
47 fTypeName = desc.GetFieldDescriptor(fFieldId).GetTypeName();
48 }
49
50 void Browse(TBrowser *b) final
51 {
52 if (!b)
53 return;
54
55 const auto &desc = fReader->GetDescriptor();
56
57 if (fIsLeaf) {
58 if (!gPad)
59 return;
60
61 auto view = fReader->GetView<void>(desc.GetQualifiedFieldName(fBrowsableFieldId));
62
63 ROOT::Internal::RNTupleDrawVisitor drawVisitor(fReader, desc.GetFieldDescriptor(fFieldId).GetFieldName());
64 view.GetField().AcceptVisitor(drawVisitor);
65 fHistogram = std::unique_ptr<TH1>(drawVisitor.MoveHist());
66 if (fHistogram->GetEntries() == 0) {
67 gPad->DrawFrame(-1., -1., 1., 1.);
68 TText *textEmpty = new TText(0., 0., "Empty");
69 textEmpty->SetTextAlign(22);
70 textEmpty->SetTextFont(42);
71 textEmpty->SetTextSize(0.1);
72 textEmpty->SetTextColor(1);
73 textEmpty->Draw();
74 } else {
75 fHistogram->Draw();
76 }
77 gPad->Update();
78 } else {
79 for (const auto &f : desc.GetFieldIterable(fBrowsableFieldId)) {
80 b->Add(new RFieldBrowsable(fReader, f.GetId()), f.GetFieldName().c_str());
81 }
82 }
83 }
84
85 bool IsFolder() const final { return !fIsLeaf; }
86 const char *GetIconName() const final { return IsFolder() ? "RNTuple-folder" : "RNTuple-leaf"; }
87
88 const char *GetName() const final { return fFieldName.c_str(); }
89 const char *GetTitle() const final { return fTypeName.c_str(); }
90};
91
92class RVisualizationBrowsable : public TObject {
93private:
94 std::unique_ptr<ROOT::Experimental::RNTupleInspector> fInspector;
95 std::unique_ptr<ROOT::Experimental::RTreeMapPainter> fTreeMap;
96
97public:
98 RVisualizationBrowsable(const ROOT::RNTuple &ntuple)
99 : fInspector(ROOT::Experimental::RNTupleInspector::Create(ntuple))
100 {
101 }
102 void Browse(TBrowser *b) final
103 {
104 if (!b || !gPad)
105 return;
106 gPad->GetListOfPrimitives()->Clear();
107 fTreeMap = ROOT::Experimental::CreateTreeMapFromRNTuple(*fInspector);
108 fTreeMap->Paint("");
109 gPad->Update();
110 }
111
112 const char *GetIconName() const final { return "RNTuple-visualization"; }
113 bool IsFolder() const final { return false; }
114 const char *GetName() const final { return "Visualization"; }
115 const char *GetTitle() const final { return "TreeMap visualization of RNTuple structure and disk usage"; }
116};
117
118} // anonymous namespace
119
121{
122 if (!b)
123 return;
124
125 std::shared_ptr<ROOT::RNTupleReader> reader = RNTupleReader::Open(*static_cast<const ROOT::RNTuple *>(ntuple));
126 const auto &desc = reader->GetDescriptor();
127 b->Add(new RVisualizationBrowsable(*static_cast<const ROOT::RNTuple *>(ntuple)), "Visualization");
128 for (const auto &f : desc.GetTopLevelFields()) {
129 b->Add(new RFieldBrowsable(reader, f.GetId()), f.GetFieldName().c_str());
130 }
131}
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gPad
Representation of an RNTuple data set in a ROOT file.
Definition RNTuple.hxx:67
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
Mother of all ROOT objects.
Definition TObject.h:41
Base class for several text objects.
Definition TText.h:22
std::unique_ptr< RTreeMapPainter > CreateTreeMapFromRNTuple(const RNTupleInspector &insp)
Logic for converting an RNTuple to RTreeMapPainter given RNTupleInspector.
void BrowseRNTuple(const void *ntuple, TBrowser *b)
DescriptorId_t GetNextBrowsableField(DescriptorId_t fieldId, const RNTupleDescriptor &desc)
std::uint64_t DescriptorId_t
Distriniguishes elements of the same type within a descriptor, e.g. different fields.
constexpr DescriptorId_t kInvalidDescriptorId