File indexing completed on 2025-12-14 10:31:17
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #ifndef ROOT_RNTupleDrawVisitor
0015 #define ROOT_RNTupleDrawVisitor
0016
0017 #include <ROOT/RField.hxx>
0018 #include <ROOT/RFieldVisitor.hxx>
0019 #include <ROOT/RNTupleReader.hxx>
0020 #include <ROOT/RNTupleView.hxx>
0021
0022 #include <TH1F.h>
0023
0024 #include <map>
0025 #include <memory>
0026 #include <string>
0027 #include <utility>
0028
0029 namespace ROOT {
0030
0031 namespace Internal {
0032
0033 class RNTupleDrawVisitor : public ROOT::Detail::RFieldVisitor {
0034 private:
0035 std::shared_ptr<ROOT::RNTupleReader> fNtplReader;
0036 std::unique_ptr<TH1> fHist;
0037 std::string fTitle;
0038
0039
0040 void TestHistBuffer();
0041
0042 template <typename ViewT>
0043 void FillHistogramImpl(ViewT &view)
0044 {
0045 fHist = std::make_unique<TH1F>("hdraw", fTitle.c_str(), 100, 0, 0);
0046 fHist->SetDirectory(nullptr);
0047
0048 auto bufsize = (fHist->GetBufferSize() - 1) / 2;
0049 int cnt = 0;
0050 if (bufsize > 10) {
0051 bufsize -= 3;
0052 } else {
0053 bufsize = -1;
0054 }
0055
0056 for (auto i : view.GetFieldRange()) {
0057 fHist->Fill(view(i));
0058 if (++cnt == bufsize) {
0059 TestHistBuffer();
0060 ++cnt;
0061 }
0062 }
0063 if (cnt < bufsize)
0064 TestHistBuffer();
0065
0066 fHist->BufferEmpty();
0067 }
0068
0069 template <typename T>
0070 void FillHistogram(const ROOT::RIntegralField<T> &field)
0071 {
0072 auto view = fNtplReader->GetDirectAccessView<T>(field.GetOnDiskId());
0073 FillHistogramImpl(view);
0074 }
0075
0076 template <typename T>
0077 void FillHistogram(const ROOT::RField<T> &field)
0078 {
0079 auto view = fNtplReader->GetView<T>(field.GetOnDiskId());
0080 FillHistogramImpl(view);
0081 }
0082
0083 void FillStringHistogram(const ROOT::RField<std::string> &field)
0084 {
0085 std::map<std::string, std::uint64_t> values;
0086
0087 std::uint64_t nentries = 0;
0088
0089 auto view = fNtplReader->GetView<std::string>(field.GetOnDiskId());
0090 for (auto i : view.GetFieldRange()) {
0091 std::string v = view(i);
0092 nentries++;
0093 auto iter = values.find(v);
0094 if (iter != values.end())
0095 iter->second++;
0096 else if (values.size() >= 50)
0097 return;
0098 else
0099 values[v] = 0;
0100 }
0101
0102
0103 fHist = std::make_unique<TH1F>("h", fTitle.c_str(), 3, 0, 3);
0104 fHist->SetDirectory(nullptr);
0105 fHist->SetStats(0);
0106 fHist->SetEntries(nentries);
0107 fHist->SetCanExtend(TH1::kAllAxes);
0108 for (auto &entry : values)
0109 fHist->Fill(entry.first.c_str(), entry.second);
0110 fHist->LabelsDeflate();
0111 fHist->Sumw2(false);
0112 }
0113
0114 public:
0115 RNTupleDrawVisitor(std::shared_ptr<ROOT::RNTupleReader> ntplReader, const std::string &title)
0116 : fNtplReader(ntplReader), fTitle(title)
0117 {
0118 }
0119
0120 TH1 *MoveHist() { return fHist.release(); }
0121
0122 void VisitField(const ROOT::RFieldBase & ) final {}
0123 void VisitBoolField(const ROOT::RField<bool> &field) final { FillHistogram(field); }
0124 void VisitFloatField(const ROOT::RField<float> &field) final { FillHistogram(field); }
0125 void VisitDoubleField(const ROOT::RField<double> &field) final { FillHistogram(field); }
0126 void VisitCharField(const ROOT::RField<char> &field) final { FillHistogram(field); }
0127 void VisitInt8Field(const ROOT::RIntegralField<std::int8_t> &field) final { FillHistogram(field); }
0128 void VisitInt16Field(const ROOT::RIntegralField<std::int16_t> &field) final { FillHistogram(field); }
0129 void VisitInt32Field(const ROOT::RIntegralField<std::int32_t> &field) final { FillHistogram(field); }
0130 void VisitInt64Field(const ROOT::RIntegralField<std::int64_t> &field) final { FillHistogram(field); }
0131 void VisitStringField(const ROOT::RField<std::string> &field) final { FillStringHistogram(field); }
0132 void VisitUInt16Field(const ROOT::RIntegralField<std::uint16_t> &field) final { FillHistogram(field); }
0133 void VisitUInt32Field(const ROOT::RIntegralField<std::uint32_t> &field) final { FillHistogram(field); }
0134 void VisitUInt64Field(const ROOT::RIntegralField<std::uint64_t> &field) final { FillHistogram(field); }
0135 void VisitUInt8Field(const ROOT::RIntegralField<std::uint8_t> &field) final { FillHistogram(field); }
0136 void VisitCardinalityField(const ROOT::RCardinalityField &field) final
0137 {
0138 if (const auto f32 = field.As32Bit()) {
0139 FillHistogram(*f32);
0140 } else if (const auto f64 = field.As64Bit()) {
0141 FillHistogram(*f64);
0142 }
0143 }
0144 };
0145
0146 }
0147 }
0148
0149 #endif