Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RFieldProvider.hxx
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
9#ifndef ROOT_Browsable_RFieldProvider
10#define ROOT_Browsable_RFieldProvider
11
12#include "TH1.h"
13#include "TMath.h"
14#include <map>
15#include <memory>
16#include <string>
17#include <utility>
18
20
22#include <ROOT/RPageStorage.hxx>
23#include <ROOT/RNTupleView.hxx>
24
25#include "RFieldHolder.hxx"
26
27using namespace ROOT::Browsable;
28
29using namespace std::string_literals;
30
31template<typename T>
33
34// ==============================================================================================
35
36/** \class RFieldProvider
37\ingroup rbrowser
38\brief Base class for provider of RNTuple drawing
39\author Sergey Linev <S.Linev@gsi.de>
40\date 2021-03-09
41\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
42*/
43
44class RFieldProvider : public RProvider {
46 private:
47 std::shared_ptr<ROOT::Experimental::Detail::RPageSource> fNtplSource;
48 std::unique_ptr<TH1> fHist;
49
50 /** Test collected entries if it looks like integer values and one can use better binning */
52 {
53 auto len = fHist->GetBufferLength();
54 auto buf = fHist->GetBuffer();
55
56 if (!buf || (len < 5))
57 return;
58
59 Double_t min = buf[1], max = buf[1];
60 Bool_t is_integer = kTRUE;
61
62 for (Int_t n = 0; n < len; ++n) {
63 Double_t v = buf[2 + 2*n];
64 if (v > max) max = v;
65 if (v < min) min = v;
66 if (TMath::Abs(v - TMath::Nint(v)) > 1e-5) { is_integer = kFALSE; break; }
67 }
68
69 // special case when only integer values in short range - better binning
70 if (is_integer && (max-min < 100)) {
71 max += 2;
72 if (min > 1) min -= 2;
73 int npoints = TMath::Nint(max - min);
74 std::unique_ptr<TH1> h1 = std::make_unique<TH1F>(fHist->GetName(), fHist->GetTitle(), npoints, min, max);
75 h1->SetDirectory(nullptr);
76 for (Int_t n = 0; n < len; ++n)
77 h1->Fill(buf[2 + 2*n], buf[1 + 2*n]);
78 std::swap(fHist, h1);
79 }
80 }
81
82 template<typename T>
83 void FillHistogram(const RField<T> &field)
84 {
85 std::string title = "Drawing of RField "s + field.GetName();
86
87 fHist = std::make_unique<TH1F>("hdraw", title.c_str(), 100, 0, 0);
88 fHist->SetDirectory(nullptr);
89
90 auto bufsize = (fHist->GetBufferSize() - 1) / 2;
91 int cnt = 0;
92 if (bufsize > 10) bufsize-=3; else bufsize = -1;
93
95 for (auto i : view.GetFieldRange()) {
96 fHist->Fill(view(i));
97 if (++cnt == bufsize) {
99 ++cnt;
100 }
101 }
102 if (cnt < bufsize)
104
105 fHist->BufferEmpty();
106 }
107
109 {
110 std::map<std::string, int> values;
111
112 int nentries = 0;
113
115 for (auto i : view.GetFieldRange()) {
116 std::string v = view(i);
117 nentries++;
118 auto iter = values.find(v);
119 if (iter != values.end())
120 iter->second++;
121 else if (values.size() >= 50)
122 return;
123 else
124 values[v] = 0;
125 }
126
127 // now create histogram with labels
128
129 std::string title = "Drawing of RField "s + field.GetName();
130 fHist = std::make_unique<TH1F>("h",title.c_str(),3,0,3);
131 fHist->SetDirectory(nullptr);
132 fHist->SetStats(0);
133 fHist->SetEntries(nentries);
134 fHist->SetCanExtend(TH1::kAllAxes);
135 for (auto &entry : values)
136 fHist->Fill(entry.first.c_str(), entry.second);
137 fHist->LabelsDeflate();
138 fHist->Sumw2(kFALSE);
139 }
140
141 public:
142 explicit RDrawVisitor(std::shared_ptr<ROOT::Experimental::Detail::RPageSource> ntplSource)
143 : fNtplSource(ntplSource)
144 {
145 }
146
148 return fHist.release();
149 }
150
151 void VisitField(const ROOT::Experimental::Detail::RFieldBase & /* field */) final {}
152 void VisitBoolField(const RField<bool> &field) final { FillHistogram(field); }
153 void VisitFloatField(const RField<float> &field) final { FillHistogram(field); }
154 void VisitDoubleField(const RField<double> &field) final { FillHistogram(field); }
155 void VisitCharField(const RField<char> &field) final { FillHistogram(field); }
156 void VisitInt8Field(const RField<std::int8_t> &field) final { FillHistogram(field); }
157 void VisitInt16Field(const RField<std::int16_t> &field) final { FillHistogram(field); }
158 void VisitIntField(const RField<int> &field) final { FillHistogram(field); }
159 void VisitInt64Field(const RField<std::int64_t> &field) final { FillHistogram(field); }
160 void VisitStringField(const RField<std::string> &field) final { FillStringHistogram(field); }
161 void VisitUInt16Field(const RField<std::uint16_t> &field) final { FillHistogram(field); }
162 void VisitUInt32Field(const RField<std::uint32_t> &field) final { FillHistogram(field); }
163 void VisitUInt64Field(const RField<std::uint64_t> &field) final { FillHistogram(field); }
164 void VisitUInt8Field(const RField<std::uint8_t> &field) final { FillHistogram(field); }
166 {
167 if (const auto f32 = field.As32Bit()) {
168 FillHistogram(*f32);
169 } else if (const auto f64 = field.As64Bit()) {
170 FillHistogram(*f64);
171 }
172 }
173 }; // class RDrawVisitor
174
175public:
176 // virtual ~RFieldProvider() = default;
177
179 {
180 if (!holder) return nullptr;
181
182 auto ntplSource = holder->GetNtplSource();
183 std::string name = holder->GetParentName();
184
185 std::unique_ptr<ROOT::Experimental::Detail::RFieldBase> field;
186 {
187 auto descriptorGuard = ntplSource->GetSharedDescriptorGuard();
188 field = descriptorGuard->GetFieldDescriptor(holder->GetId()).CreateField(descriptorGuard.GetRef());
189 }
190 name.append(field->GetName());
191
192 RDrawVisitor drawVisitor(ntplSource);
193 field->AcceptVisitor(drawVisitor);
194 return drawVisitor.MoveHist();
195 }
196};
197
198#endif
#define e(i)
Definition RSha256.hxx:103
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
char name[80]
Definition TGX11.cxx:110
int nentries
auto GetParentName() const
auto GetNtplSource() const
auto GetId() const
void VisitCardinalityField(const ROOT::Experimental::RCardinalityField &field) final
void VisitUInt16Field(const RField< std::uint16_t > &field) final
void FillStringHistogram(const RField< std::string > &field)
std::unique_ptr< TH1 > fHist
void VisitUInt32Field(const RField< std::uint32_t > &field) final
void TestHistBuffer()
Test collected entries if it looks like integer values and one can use better binning.
void FillHistogram(const RField< T > &field)
void VisitDoubleField(const RField< double > &field) final
std::shared_ptr< ROOT::Experimental::Detail::RPageSource > fNtplSource
void VisitInt8Field(const RField< std::int8_t > &field) final
void VisitCharField(const RField< char > &field) final
RDrawVisitor(std::shared_ptr< ROOT::Experimental::Detail::RPageSource > ntplSource)
void VisitUInt64Field(const RField< std::uint64_t > &field) final
void VisitStringField(const RField< std::string > &field) final
void VisitInt64Field(const RField< std::int64_t > &field) final
void VisitBoolField(const RField< bool > &field) final
void VisitUInt8Field(const RField< std::uint8_t > &field) final
void VisitInt16Field(const RField< std::int16_t > &field) final
void VisitField(const ROOT::Experimental::Detail::RFieldBase &) final
void VisitIntField(const RField< int > &field) final
void VisitFloatField(const RField< float > &field) final
Base class for provider of RNTuple drawing.
TH1 * DrawField(RFieldHolder *holder)
Provider of different browsing methods for supported classes.
Definition RProvider.hxx:37
A field translates read and write calls from/to underlying columns to/from tree values.
Definition RField.hxx:83
DescriptorId_t GetOnDiskId() const
Definition RField.hxx:610
Abstract base class for classes implementing the visitor design pattern.
An artificial field that transforms an RNTuple column that contains the offset of collections into co...
Definition RField.hxx:1422
Classes with dictionaries that can be inspected by TClass.
Definition RField.hxx:1224
An RNTupleView provides read-only access to a single field of the ntuple.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual void SetDirectory(TDirectory *dir)
By default, when a histogram is created, it is added to the list of histogram objects in the current ...
Definition TH1.cxx:8854
@ kAllAxes
Definition TH1.h:75
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3345
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:693
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123