File indexing completed on 2025-02-22 10:53:09
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #ifndef ROOT_RDF_RTREECOLUMNREADER
0012 #define ROOT_RDF_RTREECOLUMNREADER
0013
0014 #include "RColumnReaderBase.hxx"
0015 #include <ROOT/RVec.hxx>
0016 #include <Rtypes.h> // Long64_t, R__CLING_PTRCHECK
0017 #include <TTreeReader.h>
0018 #include <TTreeReaderValue.h>
0019 #include <TTreeReaderArray.h>
0020
0021 #include <memory>
0022 #include <string>
0023
0024 namespace ROOT {
0025 namespace Internal {
0026 namespace RDF {
0027
0028
0029 template <typename T>
0030 class R__CLING_PTRCHECK(off) RTreeColumnReader final : public ROOT::Detail::RDF::RColumnReaderBase {
0031 std::unique_ptr<TTreeReaderValue<T>> fTreeValue;
0032
0033 void *GetImpl(Long64_t) final { return fTreeValue->Get(); }
0034 public:
0035
0036 RTreeColumnReader(TTreeReader &r, const std::string &colName)
0037 : fTreeValue(std::make_unique<TTreeReaderValue<T>>(r, colName.c_str()))
0038 {
0039 }
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049 ~RTreeColumnReader() override { fTreeValue.reset(); }
0050 };
0051
0052
0053
0054
0055 template <typename T>
0056 class R__CLING_PTRCHECK(off) RTreeColumnReader<RVec<T>> final : public ROOT::Detail::RDF::RColumnReaderBase {
0057 std::unique_ptr<TTreeReaderArray<T>> fTreeArray;
0058
0059
0060 enum class EStorageType : char { kContiguous, kUnknown, kSparse };
0061
0062
0063 RVec<T> fRVec;
0064
0065
0066
0067 EStorageType fStorageType = EStorageType::kUnknown;
0068 Long64_t fLastEntry = -1;
0069
0070
0071 bool fCopyWarningPrinted = false;
0072
0073 void *GetImpl(Long64_t entry) final
0074 {
0075 if (entry == fLastEntry)
0076 return &fRVec;
0077
0078 auto &readerArray = *fTreeArray;
0079
0080
0081
0082
0083 const auto readerArraySize = readerArray.GetSize();
0084 if (EStorageType::kUnknown == fStorageType && readerArraySize > 1) {
0085
0086 fStorageType = EStorageType::kContiguous;
0087 for (auto i = 0u; i < readerArraySize - 1; ++i) {
0088 if ((char *)&readerArray[i + 1] - (char *)&readerArray[i] != sizeof(T)) {
0089 fStorageType = EStorageType::kSparse;
0090 break;
0091 }
0092 }
0093 }
0094
0095 if (EStorageType::kContiguous == fStorageType ||
0096 (EStorageType::kUnknown == fStorageType && readerArray.GetSize() < 2)) {
0097 if (readerArraySize > 0) {
0098
0099
0100
0101 auto readerArrayAddr = &readerArray.At(0);
0102 RVec<T> rvec(readerArrayAddr, readerArraySize);
0103 swap(fRVec, rvec);
0104 } else {
0105 RVec<T> emptyVec{};
0106 swap(fRVec, emptyVec);
0107 }
0108 } else {
0109
0110 #ifndef NDEBUG
0111 if (!fCopyWarningPrinted) {
0112 Warning("RTreeColumnReader::Get",
0113 "Branch %s hangs from a non-split branch. A copy is being performed in order "
0114 "to properly read the content.",
0115 readerArray.GetBranchName());
0116 fCopyWarningPrinted = true;
0117 }
0118 #else
0119 (void)fCopyWarningPrinted;
0120 #endif
0121 if (readerArraySize > 0) {
0122 RVec<T> rvec(readerArray.begin(), readerArray.end());
0123 swap(fRVec, rvec);
0124 } else {
0125 RVec<T> emptyVec{};
0126 swap(fRVec, emptyVec);
0127 }
0128 }
0129 fLastEntry = entry;
0130 return &fRVec;
0131 }
0132
0133 public:
0134 RTreeColumnReader(TTreeReader &r, const std::string &colName)
0135 : fTreeArray(std::make_unique<TTreeReaderArray<T>>(r, colName.c_str()))
0136 {
0137 }
0138
0139
0140 ~RTreeColumnReader() override { fTreeArray.reset(); }
0141 };
0142
0143
0144
0145
0146 template <>
0147 class R__CLING_PTRCHECK(off) RTreeColumnReader<RVec<bool>> final : public ROOT::Detail::RDF::RColumnReaderBase {
0148
0149 std::unique_ptr<TTreeReaderArray<bool>> fTreeArray;
0150
0151
0152 RVec<bool> fRVec;
0153
0154
0155
0156
0157
0158
0159 void *GetImpl(Long64_t) final
0160 {
0161 auto &readerArray = *fTreeArray;
0162 const auto readerArraySize = readerArray.GetSize();
0163 if (readerArraySize > 0) {
0164
0165 RVec<bool> rvec(readerArray.begin(), readerArray.end());
0166 swap(fRVec, rvec);
0167 } else {
0168 RVec<bool> emptyVec{};
0169 swap(fRVec, emptyVec);
0170 }
0171 return &fRVec;
0172 }
0173
0174 public:
0175 RTreeColumnReader(TTreeReader &r, const std::string &colName)
0176 : fTreeArray(std::make_unique<TTreeReaderArray<bool>>(r, colName.c_str()))
0177 {
0178 }
0179
0180
0181 ~RTreeColumnReader() override { fTreeArray.reset(); }
0182 };
0183
0184 }
0185 }
0186 }
0187
0188 #endif