Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
31// ==============================================================================================
32
33/** \class RFieldProvider
34\ingroup rbrowser
35\brief Base class for provider of RNTuple drawing
36\author Sergey Linev <S.Linev@gsi.de>
37\date 2021-03-09
38\warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
39*/
40
41class RFieldProvider : public RProvider {
43 private:
44 std::shared_ptr<ROOT::RNTupleReader> fNtplReader;
45 std::unique_ptr<TH1> fHist;
46
47 /** Test collected entries if it looks like integer values and one can use better binning */
49 {
50 auto len = fHist->GetBufferLength();
51 auto buf = fHist->GetBuffer();
52
53 if (!buf || (len < 5))
54 return;
55
56 Double_t min = buf[1], max = buf[1];
58
59 for (Int_t n = 0; n < len; ++n) {
60 Double_t v = buf[2 + 2*n];
61 if (v > max) max = v;
62 if (v < min) min = v;
63 if (TMath::Abs(v - TMath::Nint(v)) > 1e-5) { is_integer = kFALSE; break; }
64 }
65
66 // special case when only integer values in short range - better binning
67 if (is_integer && (max-min < 100)) {
68 max += 2;
69 if (min > 1) min -= 2;
70 int npoints = TMath::Nint(max - min);
71 std::unique_ptr<TH1> h1 = std::make_unique<TH1F>(fHist->GetName(), fHist->GetTitle(), npoints, min, max);
72 h1->SetDirectory(nullptr);
73 for (Int_t n = 0; n < len; ++n)
74 h1->Fill(buf[2 + 2*n], buf[1 + 2*n]);
75 std::swap(fHist, h1);
76 }
77 }
78
79 template <typename T>
81 {
82 std::string title = "Drawing of RField "s + field.GetFieldName();
83
84 fHist = std::make_unique<TH1F>("hdraw", title.c_str(), 100, 0, 0);
85 fHist->SetDirectory(nullptr);
86
87 auto bufsize = (fHist->GetBufferSize() - 1) / 2;
88 int cnt = 0;
89 if (bufsize > 10) bufsize-=3; else bufsize = -1;
90
91 for (auto i : view.GetFieldRange()) {
92 fHist->Fill(view(i));
93 if (++cnt == bufsize) {
95 ++cnt;
96 }
97 }
98 if (cnt < bufsize)
100
101 fHist->BufferEmpty();
102 }
103
104 template <typename T>
106 {
107 auto view = fNtplReader->GetView<T>(field.GetOnDiskId());
109 }
110
111 template <typename T>
113 {
114 auto view = fNtplReader->GetView<T>(field.GetOnDiskId());
116 }
117
119 {
120 std::map<std::string, int> values;
121
122 int nentries = 0;
123
124 auto view = fNtplReader->GetView<std::string>(field.GetOnDiskId());
125 for (auto i : view.GetFieldRange()) {
126 std::string v = view(i);
127 nentries++;
128 auto iter = values.find(v);
129 if (iter != values.end())
130 iter->second++;
131 else if (values.size() >= 50)
132 return;
133 else
134 values[v] = 0;
135 }
136
137 // now create histogram with labels
138
139 std::string title = "Drawing of RField "s + field.GetFieldName();
140 fHist = std::make_unique<TH1F>("h",title.c_str(),3,0,3);
141 fHist->SetDirectory(nullptr);
142 fHist->SetStats(0);
143 fHist->SetEntries(nentries);
144 fHist->SetCanExtend(TH1::kAllAxes);
145 for (auto &entry : values)
146 fHist->Fill(entry.first.c_str(), entry.second);
147 fHist->LabelsDeflate();
148 fHist->Sumw2(kFALSE);
149 }
150
151 public:
152 explicit RDrawVisitor(std::shared_ptr<ROOT::RNTupleReader> ntplReader) : fNtplReader(ntplReader) {}
153
155 return fHist.release();
156 }
157
158 void VisitField(const ROOT::RFieldBase & /* field */) final {}
173 {
174 if (const auto f32 = field.As32Bit()) {
176 } else if (const auto f64 = field.As64Bit()) {
178 }
179 }
180 }; // class RDrawVisitor
181
182public:
183 // virtual ~RFieldProvider() = default;
184
186 {
187 if (!holder) return nullptr;
188
189 auto ntplReader = holder->GetNtplReader();
190 std::string name = holder->GetParentName();
191
192 const auto fieldName = ntplReader->GetDescriptor().GetFieldDescriptor(holder->GetId()).GetFieldName();
193 const auto qualifiedFieldName = ntplReader->GetDescriptor().GetQualifiedFieldName(holder->GetId());
194 auto view = ntplReader->GetView<void>(qualifiedFieldName);
195 name.append(fieldName);
196
198 view.GetField().AcceptVisitor(drawVisitor);
199 return drawVisitor.MoveHist();
200 }
201};
202
203#endif
#define e(i)
Definition RSha256.hxx:103
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
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 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
void FillHistogram(const ROOT::RField< T > &field)
std::unique_ptr< TH1 > fHist
void TestHistBuffer()
Test collected entries if it looks like integer values and one can use better binning.
void VisitUInt64Field(const ROOT::RIntegralField< std::uint64_t > &field) final
void VisitStringField(const ROOT::RField< std::string > &field) final
void VisitField(const ROOT::RFieldBase &) final
void VisitBoolField(const ROOT::RField< bool > &field) final
void FillStringHistogram(const ROOT::RField< std::string > &field)
void VisitCharField(const ROOT::RField< char > &field) final
void VisitCardinalityField(const ROOT::RCardinalityField &field) final
void VisitUInt32Field(const ROOT::RIntegralField< std::uint32_t > &field) final
void VisitUInt8Field(const ROOT::RIntegralField< std::uint8_t > &field) final
std::shared_ptr< ROOT::RNTupleReader > fNtplReader
void VisitFloatField(const ROOT::RField< float > &field) final
void VisitInt32Field(const ROOT::RIntegralField< std::int32_t > &field) final
void VisitInt8Field(const ROOT::RIntegralField< std::int8_t > &field) final
void VisitInt16Field(const ROOT::RIntegralField< std::int16_t > &field) final
RDrawVisitor(std::shared_ptr< ROOT::RNTupleReader > ntplReader)
void VisitInt64Field(const ROOT::RIntegralField< std::int64_t > &field) final
void VisitDoubleField(const ROOT::RField< double > &field) final
void VisitUInt16Field(const ROOT::RIntegralField< std::uint16_t > &field) final
void FillHistogramImpl(const ROOT::RFieldBase &field, ROOT::RNTupleView< T > &view)
void FillHistogram(const ROOT::RIntegralField< T > &field)
Base class for provider of RNTuple drawing.
TH1 * DrawField(RFieldHolder *holder)
Provider of different browsing methods for supported classes.
Definition RProvider.hxx:37
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:310
A field translates read and write calls from/to underlying columns to/from tree values.
Classes with dictionaries that can be inspected by TClass.
Definition RField.hxx:283
ROOT::RNTupleGlobalRange GetFieldRange() const
Returns the global field range of this view.
An RNTupleView for a known type.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
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:8933
@ kAllAxes
Definition TH1.h:76
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3316
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:697
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123