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>
20
21#include <TBrowser.h>
22#include <TObject.h>
23#include <TPad.h>
24#include <TText.h>
25
26#include <memory>
27#include <string>
28
29namespace {
30
31class RFieldBrowsable final : public TObject {
32private:
33 std::shared_ptr<ROOT::RNTupleReader> fReader;
36 bool fIsLeaf = false;
37 std::unique_ptr<TH1> fHistogram;
38 std::string fFieldName;
39 std::string fTypeName;
40
41public:
42 RFieldBrowsable(std::shared_ptr<ROOT::RNTupleReader> reader, ROOT::DescriptorId_t fieldId)
43 : fReader(reader), fFieldId(fieldId)
44 {
45 const auto &desc = fReader->GetDescriptor();
46 fBrowsableFieldId = ROOT::Internal::GetNextBrowsableField(fFieldId, desc);
47 fIsLeaf = desc.GetFieldDescriptor(fBrowsableFieldId).GetLinkIds().empty();
48 fFieldName = desc.GetFieldDescriptor(fFieldId).GetFieldName();
49 fTypeName = desc.GetFieldDescriptor(fFieldId).GetTypeName();
50 }
51
52 void Browse(TBrowser *b) final
53 {
54 if (!b)
55 return;
56
57 const auto &desc = fReader->GetDescriptor();
58
59 if (fIsLeaf) {
60 if (!gPad)
61 return;
62
63 auto view = fReader->GetView<void>(desc.GetQualifiedFieldName(fBrowsableFieldId));
64
65 ROOT::Internal::RNTupleDrawVisitor drawVisitor(fReader, desc.GetFieldDescriptor(fFieldId).GetFieldName());
66 view.GetField().AcceptVisitor(drawVisitor);
67 fHistogram = std::unique_ptr<TH1>(drawVisitor.MoveHist());
68 if (fHistogram->GetEntries() == 0) {
69 gPad->DrawFrame(-1., -1., 1., 1.);
70 TText *textEmpty = new TText(0., 0., "Empty");
71 textEmpty->SetTextAlign(22);
72 textEmpty->SetTextFont(42);
73 textEmpty->SetTextSize(0.1);
74 textEmpty->SetTextColor(1);
75 textEmpty->Draw();
76 } else {
77 fHistogram->Draw();
78 }
79 gPad->Update();
80 } else {
81 for (const auto &f : desc.GetFieldIterable(fBrowsableFieldId)) {
82 b->Add(new RFieldBrowsable(fReader, f.GetId()), f.GetFieldName().c_str());
83 }
84 }
85 }
86
87 bool IsFolder() const final { return !fIsLeaf; }
88 const char *GetIconName() const final { return IsFolder() ? "RNTuple-folder" : "RNTuple-leaf"; }
89
90 const char *GetName() const final { return fFieldName.c_str(); }
91 const char *GetTitle() const final { return fTypeName.c_str(); }
92};
93
94} // anonymous namespace
95
97{
98 if (!b)
99 return;
100
101 std::shared_ptr<ROOT::RNTupleReader> reader = RNTupleReader::Open(*static_cast<const ROOT::RNTuple *>(ntuple));
102 const auto &desc = reader->GetDescriptor();
103 for (const auto &f : desc.GetTopLevelFields()) {
104 b->Add(new RFieldBrowsable(reader, f.GetId()), f.GetFieldName().c_str());
105 }
106}
#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
struct void * fTypeName
Definition cppyy.h:9
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