Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RNTupleDrawVisitor.hxx
Go to the documentation of this file.
1/// \file ROOT/RNTupleDrawVisitor.hxx
2/// \ingroup NTuple
3/// \author Sergey Linev <S.Linev@gsi.de>, Jakob Blomer <jblomer@cern.ch>
4/// \date 2025-07-24
5
6/*************************************************************************
7 * Copyright (C) 1995-2025, Rene Brun and Fons Rademakers. *
8 * All rights reserved. *
9 * *
10 * For the licensing terms see $ROOTSYS/LICENSE. *
11 * For the list of contributors see $ROOTSYS/README/CREDITS. *
12 *************************************************************************/
13
14#ifndef ROOT_RNTupleDrawVisitor
15#define ROOT_RNTupleDrawVisitor
16
17#include <ROOT/RField.hxx>
20#include <ROOT/RNTupleView.hxx>
21
22#include <TH1F.h>
23
24#include <map>
25#include <memory>
26#include <string>
27#include <utility>
28
29namespace ROOT {
30
31namespace Internal {
32
34private:
35 std::shared_ptr<ROOT::RNTupleReader> fNtplReader;
36 std::unique_ptr<TH1> fHist;
37 std::string fTitle;
38
39 /** Test collected entries if it looks like integer values and one can use better binning */
40 void TestHistBuffer();
41
42 template <typename ViewT>
44 {
45 fHist = std::make_unique<TH1F>("hdraw", fTitle.c_str(), 100, 0, 0);
46 fHist->SetDirectory(nullptr);
47
48 auto bufsize = (fHist->GetBufferSize() - 1) / 2;
49 int cnt = 0;
50 if (bufsize > 10) {
51 bufsize -= 3;
52 } else {
53 bufsize = -1;
54 }
55
56 for (auto i : view.GetFieldRange()) {
57 fHist->Fill(view(i));
58 if (++cnt == bufsize) {
60 ++cnt;
61 }
62 }
63 if (cnt < bufsize)
65
66 fHist->BufferEmpty();
67 }
68
69 template <typename T>
71 {
72 auto view = fNtplReader->GetDirectAccessView<T>(field.GetOnDiskId());
74 }
75
76 template <typename T>
78 {
79 auto view = fNtplReader->GetView<T>(field.GetOnDiskId());
81 }
82
84 {
85 std::map<std::string, std::uint64_t> values;
86
87 std::uint64_t nentries = 0;
88
89 auto view = fNtplReader->GetView<std::string>(field.GetOnDiskId());
90 for (auto i : view.GetFieldRange()) {
91 std::string v = view(i);
92 nentries++;
93 auto iter = values.find(v);
94 if (iter != values.end())
95 iter->second++;
96 else if (values.size() >= 50)
97 return;
98 else
99 values[v] = 0;
100 }
101
102 // now create histogram with labels
103 fHist = std::make_unique<TH1F>("h", fTitle.c_str(), 3, 0, 3);
104 fHist->SetDirectory(nullptr);
105 fHist->SetStats(0);
106 fHist->SetEntries(nentries);
107 fHist->SetCanExtend(TH1::kAllAxes);
108 for (auto &entry : values)
109 fHist->Fill(entry.first.c_str(), entry.second);
110 fHist->LabelsDeflate();
111 fHist->Sumw2(false);
112 }
113
114public:
115 RNTupleDrawVisitor(std::shared_ptr<ROOT::RNTupleReader> ntplReader, const std::string &title)
116 : fNtplReader(ntplReader), fTitle(title)
117 {
118 }
119
120 TH1 *MoveHist() { return fHist.release(); }
121
122 void VisitField(const ROOT::RFieldBase & /* field */) final {}
137 {
138 if (const auto f32 = field.As32Bit()) {
140 } else if (const auto f64 = field.As64Bit()) {
142 }
143 }
144}; // class RDrawVisitor
145
146} // namespace Internal
147} // namespace ROOT
148
149#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:319
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:292
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
@ kAllAxes
Definition TH1.h:126
Namespace for new ROOT classes and functions.