Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RNTupleDrawVisitor.hxx
Go to the documentation of this file.
1/// \file ROOT/RNTupleDrawVisitor.hxx
2/// \author Sergey Linev <S.Linev@gsi.de>, Jakob Blomer <jblomer@cern.ch>
3/// \date 2025-07-24
4
5/*************************************************************************
6 * Copyright (C) 1995-2025, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13#ifndef ROOT_RNTupleDrawVisitor
14#define ROOT_RNTupleDrawVisitor
15
16#include <ROOT/RField.hxx>
19#include <ROOT/RNTupleView.hxx>
20
21#include <TH1F.h>
22
23#include <map>
24#include <memory>
25#include <string>
26#include <utility>
27
28namespace ROOT {
29
30namespace Internal {
31
33private:
34 std::shared_ptr<ROOT::RNTupleReader> fNtplReader;
35 std::unique_ptr<TH1> fHist;
36 std::string fTitle;
37
38 /** Test collected entries if it looks like integer values and one can use better binning */
39 void TestHistBuffer();
40
41 template <typename ViewT>
42 void FillHistogramImpl(ViewT &view)
43 {
44 fHist = std::make_unique<TH1F>("hdraw", fTitle.c_str(), 100, 0, 0);
45 fHist->SetDirectory(nullptr);
46
47 auto bufsize = (fHist->GetBufferSize() - 1) / 2;
48 int cnt = 0;
49 if (bufsize > 10) {
50 bufsize -= 3;
51 } else {
52 bufsize = -1;
53 }
54
55 for (auto i : view.GetFieldRange()) {
56 fHist->Fill(view(i));
57 if (++cnt == bufsize) {
59 ++cnt;
60 }
61 }
62 if (cnt < bufsize)
64
65 fHist->BufferEmpty();
66 }
67
68 template <typename T>
70 {
71 auto view = fNtplReader->GetDirectAccessView<T>(field.GetOnDiskId());
73 }
74
75 template <typename T>
76 void FillHistogram(const ROOT::RField<T> &field)
77 {
78 auto view = fNtplReader->GetView<T>(field.GetOnDiskId());
80 }
81
83 {
84 std::map<std::string, std::uint64_t> values;
85
86 std::uint64_t nentries = 0;
87
88 auto view = fNtplReader->GetView<std::string>(field.GetOnDiskId());
89 for (auto i : view.GetFieldRange()) {
90 std::string v = view(i);
91 nentries++;
92 auto iter = values.find(v);
93 if (iter != values.end())
94 iter->second++;
95 else if (values.size() >= 50)
96 return;
97 else
98 values[v] = 0;
99 }
100
101 // now create histogram with labels
102 fHist = std::make_unique<TH1F>("h", fTitle.c_str(), 3, 0, 3);
103 fHist->SetDirectory(nullptr);
104 fHist->SetStats(0);
105 fHist->SetEntries(nentries);
106 fHist->SetCanExtend(TH1::kAllAxes);
107 for (auto &entry : values)
108 fHist->Fill(entry.first.c_str(), entry.second);
109 fHist->LabelsDeflate();
110 fHist->Sumw2(false);
111 }
112
113public:
114 RNTupleDrawVisitor(std::shared_ptr<ROOT::RNTupleReader> ntplReader, const std::string &title)
115 : fNtplReader(ntplReader), fTitle(title)
116 {
117 }
118
119 TH1 *MoveHist() { return fHist.release(); }
120
121 void VisitField(const ROOT::RFieldBase & /* field */) final {}
122 void VisitBoolField(const ROOT::RField<bool> &field) final { FillHistogram(field); }
123 void VisitFloatField(const ROOT::RField<float> &field) final { FillHistogram(field); }
124 void VisitDoubleField(const ROOT::RField<double> &field) final { FillHistogram(field); }
125 void VisitCharField(const ROOT::RField<char> &field) final { FillHistogram(field); }
136 {
137 if (const auto f32 = field.As32Bit()) {
138 FillHistogram(*f32);
139 } else if (const auto f64 = field.As64Bit()) {
140 FillHistogram(*f64);
141 }
142 }
143}; // class RDrawVisitor
144
145} // namespace Internal
146} // namespace ROOT
147
148#endif
int nentries
Abstract base class for classes implementing the visitor design pattern.
RNTupleDrawVisitor(std::shared_ptr< ROOT::RNTupleReader > ntplReader, const std::string &title)
void VisitDoubleField(const ROOT::RField< double > &field) final
void VisitStringField(const ROOT::RField< std::string > &field) final
void VisitUInt32Field(const ROOT::RIntegralField< std::uint32_t > &field) final
void VisitField(const ROOT::RFieldBase &) final
void VisitInt64Field(const ROOT::RIntegralField< std::int64_t > &field) final
void VisitInt8Field(const ROOT::RIntegralField< std::int8_t > &field) final
void FillStringHistogram(const ROOT::RField< std::string > &field)
void VisitUInt8Field(const ROOT::RIntegralField< std::uint8_t > &field) final
void VisitInt32Field(const ROOT::RIntegralField< std::int32_t > &field) final
void VisitInt16Field(const ROOT::RIntegralField< std::int16_t > &field) final
void VisitCharField(const ROOT::RField< char > &field) final
void VisitUInt64Field(const ROOT::RIntegralField< std::uint64_t > &field) final
void VisitUInt16Field(const ROOT::RIntegralField< std::uint16_t > &field) final
void VisitBoolField(const ROOT::RField< bool > &field) final
void VisitCardinalityField(const ROOT::RCardinalityField &field) final
void VisitFloatField(const ROOT::RField< float > &field) final
std::shared_ptr< ROOT::RNTupleReader > fNtplReader
void FillHistogram(const ROOT::RIntegralField< T > &field)
void TestHistBuffer()
Test collected entries if it looks like integer values and one can use better binning.
void FillHistogram(const ROOT::RField< T > &field)
An artificial field that transforms an RNTuple column that contains the offset of collections into co...
Definition RField.hxx:347
A field translates read and write calls from/to underlying columns to/from tree values.
ROOT::DescriptorId_t GetOnDiskId() const
Classes with dictionaries that can be inspected by TClass.
Definition RField.hxx:320
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
@ kAllAxes
Definition TH1.h:126