Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RNTupleBrowseProvider.cxx
Go to the documentation of this file.
1/*************************************************************************
2 * Copyright (C) 1995-2021, Rene Brun and Fons Rademakers. *
3 * All rights reserved. *
4 * *
5 * For the licensing terms see $ROOTSYS/LICENSE. *
6 * For the list of contributors see $ROOTSYS/README/CREDITS. *
7 *************************************************************************/
8
13
16
17#include "TClass.h"
18#include "RFieldHolder.hxx"
19
20
21using namespace std::string_literals;
22
23using namespace ROOT::Browsable;
24
25
26// ==============================================================================================
27
28/** \class RFieldElement
29\ingroup rbrowser
30\brief Browsing element representing of RField
31\author Sergey Linev <S.Linev@gsi.de>
32\date 2021-03-08
33\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
34*/
35
36class RFieldElement : public RElement {
37protected:
38 std::shared_ptr<ROOT::Experimental::RNTupleReader> fNtplReader;
39
40 std::string fParentName;
41
43
44public:
45 RFieldElement(std::shared_ptr<ROOT::Experimental::RNTupleReader> ntplReader, const std::string &parent_name,
46 const ROOT::DescriptorId_t id)
48 {
49 }
50
51 virtual ~RFieldElement() = default;
52
53 /** Name of RField */
54 std::string GetName() const override
55 {
56 return fNtplReader->GetDescriptor().GetFieldDescriptor(fFieldId).GetFieldName();
57 }
58
59 /** Title of RField */
60 std::string GetTitle() const override
61 {
62 auto &fld = fNtplReader->GetDescriptor().GetFieldDescriptor(fFieldId);
63 return "RField name "s + fld.GetFieldName() + " type "s + fld.GetTypeName();
64 }
65
66 std::unique_ptr<RLevelIter> GetChildsIter() override;
67
68 /** Return class of field - for a moment using RNTuple class as dummy */
69 const TClass *GetClass() const { return TClass::GetClass<ROOT::RNTuple>(); }
70
71 std::unique_ptr<RHolder> GetObject() override
72 {
73 return std::make_unique<RFieldHolder>(fNtplReader, fParentName, fFieldId);
74 }
75
77 {
78 auto range = fNtplReader->GetDescriptor().GetFieldIterable(fFieldId);
79 if (range.begin() != range.end()) return kActNone;
80 return kActDraw7;
81 }
82
83 bool IsCapable(EActionKind kind) const override
84 {
85 if ((kind == kActDraw6) || (kind == kActDraw7))
86 return GetDefaultAction() == kActDraw7;
87
88 return false;
89 }
90
91};
92
93// ==============================================================================================
94
95/** \class RNTupleElement
96\ingroup rbrowser
97\brief Browsing element representing of RNTuple
98\author Sergey Linev <S.Linev@gsi.de>
99\date 2021-03-08
100\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
101*/
102
103class RNTupleElement : public RElement {
104protected:
105 std::shared_ptr<ROOT::Experimental::RNTupleReader> fNtplReader;
106
107public:
108 RNTupleElement(const std::string &ntplName, const std::string &filename)
109 {
111 }
112
113 virtual ~RNTupleElement() = default;
114
115 /** Returns true if no ntuple found */
116 bool IsNull() const { return !fNtplReader; }
117
118 /** Name of NTuple */
119 std::string GetName() const override { return fNtplReader->GetDescriptor().GetName(); }
120
121 /** Title of NTuple */
122 std::string GetTitle() const override { return "RNTuple title"s; }
123
124 /** Create iterator for childs elements if any */
125 std::unique_ptr<RLevelIter> GetChildsIter() override;
126
127 const TClass *GetClass() const { return TClass::GetClass<ROOT::RNTuple>(); }
128
129 std::unique_ptr<RItem> CreateItem() const override
130 {
131 auto item = std::make_unique<RItem>(GetName(), -1, "sap-icon://table-chart");
132 item->SetTitle(GetTitle());
133 return item;
134 }
135
136 //EActionKind GetDefaultAction() const override;
137
138 //bool IsCapable(EActionKind) const override;
139};
140
141
142// ==============================================================================================
143
144/** \class RFieldsIterator
145\ingroup rbrowser
146\brief Iterator over RNTuple fields
147\author Sergey Linev <S.Linev@gsi.de>
148\date 2021-03-08
149\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
150*/
151
152
154
155 std::shared_ptr<ROOT::Experimental::RNTupleReader> fNtplReader;
156 std::vector<ROOT::DescriptorId_t> fFieldIds;
157 std::string fParentName;
158 int fCounter{-1};
159
160public:
161 RFieldsIterator(std::shared_ptr<ROOT::Experimental::RNTupleReader> ntplReader,
162 std::vector<ROOT::DescriptorId_t> &&ids, const std::string &parent_name = ""s)
164 {
165 }
166
167 virtual ~RFieldsIterator() = default;
168
169 bool Next() override
170 {
171 return ++fCounter < (int) fFieldIds.size();
172 }
173
174 std::string GetItemName() const override
175 {
176 return fNtplReader->GetDescriptor().GetFieldDescriptor(fFieldIds[fCounter]).GetFieldName();
177 }
178
179 bool CanItemHaveChilds() const override
180 {
181 auto subrange = fNtplReader->GetDescriptor().GetFieldIterable(fFieldIds[fCounter]);
182 return subrange.begin() != subrange.end();
183 }
184
185 /** Create element for the browser */
186 std::unique_ptr<RItem> CreateItem() override
187 {
188 int nchilds = 0;
189 for (auto &sub : fNtplReader->GetDescriptor().GetFieldIterable(fFieldIds[fCounter])) {
190 (void)sub;
191 nchilds++;
192 }
193
194 const auto &field = fNtplReader->GetDescriptor().GetFieldDescriptor(fFieldIds[fCounter]);
195
196 auto item =
197 std::make_unique<RItem>(field.GetFieldName(), nchilds, nchilds > 0 ? "sap-icon://split" : "sap-icon://e-care");
198
199 item->SetTitle("RField name "s + field.GetFieldName() + " type "s + field.GetTypeName());
200
201 return item;
202 }
203
204 std::shared_ptr<RElement> GetElement() override
205 {
206 return std::make_shared<RFieldElement>(fNtplReader, fParentName, fFieldIds[fCounter]);
207 }
208};
209
210
211std::unique_ptr<RLevelIter> RFieldElement::GetChildsIter()
212{
213 std::vector<ROOT::DescriptorId_t> ids;
214 std::string prefix;
215
216 for (auto &f : fNtplReader->GetDescriptor().GetFieldIterable(fFieldId))
217 ids.emplace_back(f.GetId());
218
219 if (ids.size() == 0)
220 return nullptr;
221
222 prefix = fParentName;
223 const auto &fld = fNtplReader->GetDescriptor().GetFieldDescriptor(fFieldId);
224 prefix.append(fld.GetFieldName());
225 prefix.append(".");
226
227 return std::make_unique<RFieldsIterator>(fNtplReader, std::move(ids), prefix);
228}
229
230std::unique_ptr<RLevelIter> RNTupleElement::GetChildsIter()
231{
232 std::vector<ROOT::DescriptorId_t> ids;
233
234 for (auto &f : fNtplReader->GetDescriptor().GetTopLevelFields())
235 ids.emplace_back(f.GetId());
236
237 if (ids.size() == 0) return nullptr;
238 return std::make_unique<RFieldsIterator>(fNtplReader, std::move(ids));
239}
240
241
242// ==============================================================================================
243
244/** \class RNTupleBrowseProvider
245\ingroup rbrowser
246\brief Provider for browsing RNTuple classes
247\author Sergey Linev <S.Linev@gsi.de>
248\date 2021-03-08
249\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
250*/
251
253
254public:
255
257 {
258 RegisterNTupleFunc([](const std::string &tuple_name, const std::string &filename) -> std::shared_ptr<RElement> {
259 auto elem = std::make_shared<RNTupleElement>(tuple_name, filename);
260 return elem->IsNull() ? nullptr : elem;
261 });
262 }
263
265 {
266 RegisterNTupleFunc(nullptr);
267 }
268
270
RNTupleBrowseProvider newRNTupleBrowseProvider
#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.
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 filename
Browsing element representing of RField.
RFieldElement(std::shared_ptr< ROOT::Experimental::RNTupleReader > ntplReader, const std::string &parent_name, const ROOT::DescriptorId_t id)
std::shared_ptr< ROOT::Experimental::RNTupleReader > fNtplReader
EActionKind GetDefaultAction() const override
Get default action.
bool IsCapable(EActionKind kind) const override
Check if want to perform action.
virtual ~RFieldElement()=default
std::unique_ptr< RHolder > GetObject() override
Access object.
std::string GetTitle() const override
Title of RField.
std::unique_ptr< RLevelIter > GetChildsIter() override
Create iterator for childs elements if any.
const TClass * GetClass() const
Return class of field - for a moment using RNTuple class as dummy.
std::string GetName() const override
Name of RField.
ROOT::DescriptorId_t fFieldId
Iterator over RNTuple fields.
std::unique_ptr< RItem > CreateItem() override
Create element for the browser.
std::vector< ROOT::DescriptorId_t > fFieldIds
std::shared_ptr< RElement > GetElement() override
Create RElement for current entry - may take much time to load object or open file.
bool Next() override
Shift to next entry.
std::shared_ptr< ROOT::Experimental::RNTupleReader > fNtplReader
virtual ~RFieldsIterator()=default
RFieldsIterator(std::shared_ptr< ROOT::Experimental::RNTupleReader > ntplReader, std::vector< ROOT::DescriptorId_t > &&ids, const std::string &parent_name=""s)
std::string GetItemName() const override
Returns current entry name
bool CanItemHaveChilds() const override
Returns true if current item can have childs.
Provider for browsing RNTuple classes.
Browsing element representing of RNTuple.
std::string GetName() const override
Name of NTuple.
std::unique_ptr< RLevelIter > GetChildsIter() override
Create iterator for childs elements if any.
virtual ~RNTupleElement()=default
bool IsNull() const
Returns true if no ntuple found.
std::shared_ptr< ROOT::Experimental::RNTupleReader > fNtplReader
std::unique_ptr< RItem > CreateItem() const override
Returns item with element description.
RNTupleElement(const std::string &ntplName, const std::string &filename)
std::string GetTitle() const override
Title of NTuple.
const TClass * GetClass() const
Basic element of browsable hierarchy.
Definition RElement.hxx:34
EActionKind
Possible actions on double-click.
Definition RElement.hxx:50
@ kActDraw6
can be drawn inside ROOT6 canvas
Definition RElement.hxx:55
@ kActDraw7
can be drawn inside ROOT7 canvas
Definition RElement.hxx:56
Iterator over single level hierarchy like any array, keys list, ...
Provider of different browsing methods for supported classes.
Definition RProvider.hxx:37
void RegisterNTupleFunc(BrowseNTupleFunc_t func)
static std::unique_ptr< RNTupleReader > Open(std::string_view ntupleName, std::string_view storage, const ROOT::RNTupleReadOptions &options=ROOT::RNTupleReadOptions())
Open an RNTuple for reading.
const_iterator begin() const
const_iterator end() const
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:84
std::uint64_t DescriptorId_t
Distriniguishes elements of the same type within a descriptor, e.g. different fields.