Logo ROOT  
Reference Guide
 
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
17
18#include "TClass.h"
19#include "RFieldHolder.hxx"
20
21
22using namespace std::string_literals;
23
24using namespace ROOT::Browsable;
25
26
27// ==============================================================================================
28
29/** \class RFieldElement
30\ingroup rbrowser
31\brief Browsing element representing of RField
32\author Sergey Linev <S.Linev@gsi.de>
33\date 2021-03-08
34\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
35*/
36
37class RFieldElement : public RElement {
38protected:
39 std::shared_ptr<ROOT::RNTupleReader> fNtplReader;
40 std::string fParentName;
42 std::string fDisplayName;
43
44public:
45 RFieldElement(std::shared_ptr<ROOT::RNTupleReader> ntplReader, const std::string &parent_name,
46 const ROOT::DescriptorId_t id, const std::string &displayName)
48 {
49 }
50
51 ~RFieldElement() override = 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, fDisplayName);
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::RNTupleReader> fNtplReader;
106
107public:
108 RNTupleElement(const std::string &ntplName, const std::string &filename)
109 {
111 }
112
113 ~RNTupleElement() override = 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::RNTupleReader> fNtplReader;
156 std::vector<ROOT::DescriptorId_t> fProvidedFieldIds;
157 std::vector<ROOT::DescriptorId_t> fActualFieldIds;
158 std::string fParentName;
159 int fCounter{-1};
160
161public:
162 RFieldsIterator(std::shared_ptr<ROOT::RNTupleReader> ntplReader, std::vector<ROOT::DescriptorId_t> &&ids,
163 const std::string &parent_name = "")
165 {
166 const auto &desc = fNtplReader->GetDescriptor();
167 fActualFieldIds.reserve(fProvidedFieldIds.size());
168 for (auto fid : fProvidedFieldIds) {
170 }
171 }
172
173 ~RFieldsIterator() override = default;
174
175 bool Next() override { return ++fCounter < (int)fProvidedFieldIds.size(); }
176
177 std::string GetItemName() const override
178 {
179 return fNtplReader->GetDescriptor().GetFieldDescriptor(fProvidedFieldIds[fCounter]).GetFieldName();
180 }
181
182 bool CanItemHaveChilds() const override
183 {
184 auto subrange = fNtplReader->GetDescriptor().GetFieldIterable(fActualFieldIds[fCounter]);
185 return subrange.begin() != subrange.end();
186 }
187
188 /** Create element for the browser */
189 std::unique_ptr<RItem> CreateItem() override
190 {
191 int nchilds = 0;
192 for (auto &sub : fNtplReader->GetDescriptor().GetFieldIterable(fActualFieldIds[fCounter])) {
193 (void)sub;
194 nchilds++;
195 }
196
197 const auto &field = fNtplReader->GetDescriptor().GetFieldDescriptor(fProvidedFieldIds[fCounter]);
198
199 auto item =
200 std::make_unique<RItem>(field.GetFieldName(), nchilds, nchilds > 0 ? "sap-icon://split" : "sap-icon://e-care");
201
202 item->SetTitle("RField name "s + field.GetFieldName() + " type "s + field.GetTypeName());
203
204 return item;
205 }
206
207 std::shared_ptr<RElement> GetElement() override
208 {
209 const auto name = fNtplReader->GetDescriptor().GetFieldDescriptor(fProvidedFieldIds[fCounter]).GetFieldName();
210 return std::make_shared<RFieldElement>(fNtplReader, fParentName, fActualFieldIds[fCounter], name);
211 }
212};
213
214
215std::unique_ptr<RLevelIter> RFieldElement::GetChildsIter()
216{
217 std::vector<ROOT::DescriptorId_t> ids;
218 std::string prefix;
219
220 const auto &desc = fNtplReader->GetDescriptor();
221 for (auto &f : fNtplReader->GetDescriptor().GetFieldIterable(ROOT::Internal::GetNextBrowsableField(fFieldId, desc)))
222 ids.emplace_back(f.GetId());
223
224 if (ids.size() == 0)
225 return nullptr;
226
227 prefix = fParentName;
228 const auto &fld = desc.GetFieldDescriptor(fFieldId);
229 prefix.append(fld.GetFieldName());
230 prefix.append(".");
231
232 return std::make_unique<RFieldsIterator>(fNtplReader, std::move(ids), prefix);
233}
234
235std::unique_ptr<RLevelIter> RNTupleElement::GetChildsIter()
236{
237 std::vector<ROOT::DescriptorId_t> ids;
238
239 for (auto &f : fNtplReader->GetDescriptor().GetTopLevelFields())
240 ids.emplace_back(f.GetId());
241
242 if (ids.size() == 0) return nullptr;
243 return std::make_unique<RFieldsIterator>(fNtplReader, std::move(ids));
244}
245
246
247// ==============================================================================================
248
249/** \class RNTupleBrowseProvider
250\ingroup rbrowser
251\brief Provider for browsing RNTuple classes
252\author Sergey Linev <S.Linev@gsi.de>
253\date 2021-03-08
254\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
255*/
256
258
259public:
260
262 {
263 RegisterNTupleFunc([](const std::string &tuple_name, const std::string &filename) -> std::shared_ptr<RElement> {
264 auto elem = std::make_shared<RNTupleElement>(tuple_name, filename);
265 return elem->IsNull() ? nullptr : elem;
266 });
267 }
268
270
272
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
char name[80]
Definition TGX11.cxx:110
Browsing element representing of RField.
EActionKind GetDefaultAction() const override
Get default action.
bool IsCapable(EActionKind kind) const override
Check if want to perform action.
std::shared_ptr< ROOT::RNTupleReader > fNtplReader
RFieldElement(std::shared_ptr< ROOT::RNTupleReader > ntplReader, const std::string &parent_name, const ROOT::DescriptorId_t id, const std::string &displayName)
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
~RFieldElement() override=default
Iterator over RNTuple fields.
std::unique_ptr< RItem > CreateItem() override
Create element for the browser.
std::shared_ptr< RElement > GetElement() override
Create RElement for current entry - may take much time to load object or open file.
~RFieldsIterator() override=default
std::vector< ROOT::DescriptorId_t > fActualFieldIds
bool Next() override
Shift to next entry.
std::vector< ROOT::DescriptorId_t > fProvidedFieldIds
RFieldsIterator(std::shared_ptr< ROOT::RNTupleReader > ntplReader, std::vector< ROOT::DescriptorId_t > &&ids, const std::string &parent_name="")
std::string GetItemName() const override
Returns current entry name
bool CanItemHaveChilds() const override
Returns true if current item can have childs.
std::shared_ptr< ROOT::RNTupleReader > fNtplReader
Provider for browsing RNTuple classes.
Browsing element representing of RNTuple.
std::string GetName() const override
Name of NTuple.
~RNTupleElement() override=default
std::unique_ptr< RLevelIter > GetChildsIter() override
Create iterator for childs elements if any.
std::shared_ptr< ROOT::RNTupleReader > fNtplReader
bool IsNull() const
Returns true if no ntuple found.
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
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.