JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwRootTreeBranchVector Class Reference

A helper class to manage a vector of branch entries for ROOT trees. More...

#include <QwRootFile.h>

+ Collaboration diagram for QwRootTreeBranchVector:

Data Structures

struct  Entry
 

Public Types

using size_type = std::size_t
 

Public Member Functions

 QwRootTreeBranchVector ()=default
 
void reserve (size_type count)
 
void shrink_to_fit ()
 
void clear ()
 
size_type size () const noexcept
 
bool empty () const noexcept
 
template<typename T = uint8_t>
const T & operator[] (size_type index) const
 
template<typename T = uint8_t>
T & operator[] (size_type index)
 
template<typename T>
T & value (size_type index)
 
template<typename T>
const T & value (size_type index) const
 
void SetValue (size_type index, Double_t val)
 
void SetValue (size_type index, Float_t val)
 
void SetValue (size_type index, Int_t val)
 
void SetValue (size_type index, Long64_t val)
 
void SetValue (size_type index, Short_t val)
 
void SetValue (size_type index, UShort_t val)
 
void SetValue (size_type index, UInt_t val)
 
void SetValue (size_type index, ULong64_t val)
 
void * data () noexcept
 
const void * data () const noexcept
 
size_type data_size () const noexcept
 
template<typename T>
T & back ()
 
template<typename T>
const T & back () const
 
void push_back (const std::string &name, const char type='D')
 
void push_back (const TString &name, const char type='D')
 
void push_back (const char *name, const char type='D')
 
std::string LeafList (size_type start_index=0) const
 
std::string Dump (size_type start_index=0, size_type end_index=0) const
 

Private Member Functions

std::string FormatValue (const Entry &entry, size_type index) const
 

Static Private Member Functions

static std::size_t GetTypeSize (char type)
 
static std::size_t AlignOffset (std::size_t offset)
 
template<typename T>
static std::string FormatNumeric (T input)
 

Private Attributes

std::vector< Entrym_entries
 
std::vector< std::uint8_t > m_buffer
 

Detailed Description

A helper class to manage a vector of branch entries for ROOT trees.

This class provides functionality to manage a collection of branch entries, including their names, types, offsets, and sizes. It supports adding new entries, accessing entries by index or name, and generating leaf lists for ROOT trees.

Definition at line 55 of file QwRootFile.h.

Member Typedef Documentation

◆ size_type

using QwRootTreeBranchVector::size_type = std::size_t

Definition at line 64 of file QwRootFile.h.

Constructor & Destructor Documentation

◆ QwRootTreeBranchVector()

QwRootTreeBranchVector::QwRootTreeBranchVector ( )
default

Member Function Documentation

◆ AlignOffset()

static std::size_t QwRootTreeBranchVector::AlignOffset ( std::size_t offset)
inlinestaticprivate

Definition at line 320 of file QwRootFile.h.

320 {
321 const std::size_t alignment = 4u;
322 return (offset + (alignment - 1u)) & ~(alignment - 1u);
323 }

Referenced by push_back().

+ Here is the caller graph for this function:

◆ back() [1/2]

template<typename T>
T & QwRootTreeBranchVector::back ( )
inline

Definition at line 180 of file QwRootFile.h.

180 {
181 if (m_entries.empty()) {
182 throw std::out_of_range("QwRootTreeBranchVector::back() called on empty container");
183 }
184 const auto& last_entry = m_entries.back();
185 return *reinterpret_cast<T*>(m_buffer.data() + last_entry.offset);
186 }
std::vector< Entry > m_entries
Definition QwRootFile.h:325
std::vector< std::uint8_t > m_buffer
Definition QwRootFile.h:326

References m_buffer, and m_entries.

Referenced by QwHelicityBase::ConstructBranchAndVector(), QwHelicityDecoder::ConstructBranchAndVector(), and QwSubsystemArrayParity::ConstructBranchAndVector().

+ Here is the caller graph for this function:

◆ back() [2/2]

template<typename T>
const T & QwRootTreeBranchVector::back ( ) const
inline

Definition at line 189 of file QwRootFile.h.

189 {
190 if (m_entries.empty()) {
191 throw std::out_of_range("QwRootTreeBranchVector::back() called on empty container");
192 }
193 const auto& last_entry = m_entries.back();
194 return *reinterpret_cast<const T*>(m_buffer.data() + last_entry.offset);
195 }

References m_buffer, and m_entries.

◆ clear()

void QwRootTreeBranchVector::clear ( )
inline

Definition at line 78 of file QwRootFile.h.

78 {
79 m_entries.clear();
80 m_buffer.clear();
81 }

References m_buffer, and m_entries.

◆ data() [1/2]

const void * QwRootTreeBranchVector::data ( ) const
inlinenoexcept

Definition at line 176 of file QwRootFile.h.

176{ return m_buffer.data(); }

References m_buffer.

◆ data() [2/2]

void * QwRootTreeBranchVector::data ( )
inlinenoexcept

Definition at line 175 of file QwRootFile.h.

175{ return m_buffer.data(); }

References m_buffer.

◆ data_size()

size_type QwRootTreeBranchVector::data_size ( ) const
inlinenoexcept

Definition at line 177 of file QwRootFile.h.

177{ return m_buffer.size(); }

References m_buffer.

◆ Dump()

std::string QwRootTreeBranchVector::Dump ( size_type start_index = 0,
size_type end_index = 0 ) const
inline

Definition at line 245 of file QwRootFile.h.

245 {
246 std::ostringstream stream;
247 stream << "QwRootTreeBranchVector: " << m_entries.size() << " entries, "
248 << m_buffer.size() << " bytes\n";
249 size_t end_offset = (end_index == 0 || end_index > m_entries.size()) ?
250 m_buffer.size()
251 : m_entries[end_index - 1].offset + m_entries[end_index - 1].size;
252 stream << "QwRootTreeBranchVector: buffer at 0x" << std::hex << (void*) &m_buffer[0] << '\n';
253 stream << "QwRootTreeBranchVector: entries at 0x" << std::hex << (void*) &m_entries[0] << '\n';
254 for (size_t offset = m_entries[start_index].offset; offset < end_offset; offset += 4) {
255 stream << std::dec
256 << " [" << offset << "] "
257 << std::hex
258 << " offset=0x" << offset
259 << " (0x" << std::setw(4) << std::setfill('0')
260 << offset - m_entries[start_index].offset << ")"
261 << " buff=";
262 // Little-endian
263 for (std::size_t byte = 0; byte < 4; ++byte) {
264 stream << std::hex << std::setw(2) << std::setfill('0')
265 << static_cast<unsigned int>(m_buffer[offset + byte])
266 << " ";
267 }
268 stream << '\n';
269 }
270 end_index = (end_index == 0 || end_index > m_entries.size()) ?
271 m_entries.size()
272 : end_index;
273 for (size_type index = start_index; index < end_index; ++index) {
274 const auto& entry = m_entries[index];
275 stream << std::dec
276 << " [" << index << "] "
277 << std::hex
278 << " offset=0x" << entry.offset
279 << " (0x" << std::setw(4) << std::setfill('0')
280 << entry.offset - m_entries[start_index].offset << ")"
281 << " size=0x" << entry.size
282 << " buff=0x";
283 // Little-endian
284 for (std::size_t byte = GetTypeSize(entry.type); byte > 0; --byte) {
285 stream << std::hex << std::setw(2) << std::setfill('0')
286 << static_cast<unsigned int>((m_buffer.data() + entry.offset)[byte - 1]);
287 }
288 stream << std::dec
289 << " name=" << entry.name << "/" << entry.type
290 << " value=" << FormatValue(entry, index);
291 stream << '\n';
292 }
293 return stream.str();
294 }
std::size_t size_type
Definition QwRootFile.h:64
std::string FormatValue(const Entry &entry, size_type index) const
Definition QwRootFile.h:328
static std::size_t GetTypeSize(char type)
Definition QwRootFile.h:297

References FormatValue(), GetTypeSize(), m_buffer, and m_entries.

+ Here is the call graph for this function:

◆ empty()

bool QwRootTreeBranchVector::empty ( ) const
inlinenoexcept

Definition at line 84 of file QwRootFile.h.

84{ return m_entries.empty(); }

References m_entries.

◆ FormatNumeric()

template<typename T>
static std::string QwRootTreeBranchVector::FormatNumeric ( T input)
inlinestaticprivate

Definition at line 352 of file QwRootFile.h.

352 {
353 std::ostringstream stream;
354 stream << input;
355 return stream.str();
356 }

Referenced by FormatValue().

+ Here is the caller graph for this function:

◆ FormatValue()

std::string QwRootTreeBranchVector::FormatValue ( const Entry & entry,
size_type index ) const
inlineprivate

Definition at line 328 of file QwRootFile.h.

328 {
329 switch (entry.type) {
330 case 'D':
331 return FormatNumeric(value<double>(index));
332 case 'F':
333 return FormatNumeric(value<float>(index));
334 case 'L':
335 return FormatNumeric(value<long long>(index));
336 case 'l':
338 case 'I':
339 return FormatNumeric(value<int>(index));
340 case 'i':
341 return FormatNumeric(value<unsigned int>(index));
342 case 'S':
343 return FormatNumeric(value<short>(index));
344 case 's':
346 default:
347 return "<unknown>";
348 }
349 }
T & value(size_type index)
Definition QwRootFile.h:97
static std::string FormatNumeric(T input)
Definition QwRootFile.h:352

References FormatNumeric(), QwRootTreeBranchVector::Entry::type, and value().

Referenced by Dump().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetTypeSize()

static std::size_t QwRootTreeBranchVector::GetTypeSize ( char type)
inlinestaticprivate

Definition at line 297 of file QwRootFile.h.

297 {
298 switch (type) {
299 case 'D':
300 return sizeof(double);
301 case 'F':
302 return sizeof(float);
303 case 'L':
304 return sizeof(long long);
305 case 'l':
306 return sizeof(unsigned long long);
307 case 'I':
308 return sizeof(int);
309 case 'i':
310 return sizeof(unsigned int);
311 case 'S':
312 return sizeof(short);
313 case 's':
314 return sizeof(unsigned short);
315 default:
316 throw std::invalid_argument("Unsupported branch type code: " + std::string(1, type));
317 }
318 }

Referenced by Dump(), and push_back().

+ Here is the caller graph for this function:

◆ LeafList()

std::string QwRootTreeBranchVector::LeafList ( size_type start_index = 0) const
inline

Definition at line 230 of file QwRootFile.h.

230 {
231 static const std::string separator = ":";
232 std::ostringstream stream;
233 bool first = true;
234 for (size_type index = start_index; index < m_entries.size(); ++index) {
235 const auto& entry = m_entries[index];
236 if (!first) {
237 stream << separator;
238 }
239 stream << entry.name << "/" << entry.type;
240 first = false;
241 }
242 return stream.str();
243 }

References m_entries.

Referenced by QwADC18_Channel::ConstructBranchAndVector(), QwBeamMod::ConstructBranchAndVector(), QwEPICSEvent::ConstructBranchAndVector(), QwMollerADC_Channel::ConstructBranchAndVector(), QwPMT_Channel::ConstructBranchAndVector(), and QwVQWK_Channel::ConstructBranchAndVector().

+ Here is the caller graph for this function:

◆ operator[]() [1/2]

template<typename T = uint8_t>
T & QwRootTreeBranchVector::operator[] ( size_type index)
inline

Definition at line 92 of file QwRootFile.h.

92 {
93 return value<T>(index);
94 }

References value().

+ Here is the call graph for this function:

◆ operator[]() [2/2]

template<typename T = uint8_t>
const T & QwRootTreeBranchVector::operator[] ( size_type index) const
inline

Definition at line 87 of file QwRootFile.h.

87 {
88 return value<T>(index);
89 }

References value().

+ Here is the call graph for this function:

◆ push_back() [1/3]

void QwRootTreeBranchVector::push_back ( const char * name,
const char type = 'D' )
inline

Definition at line 226 of file QwRootFile.h.

226 {
227 push_back(std::string(name), type);
228 }
void push_back(const std::string &name, const char type='D')
Definition QwRootFile.h:197

References push_back().

+ Here is the call graph for this function:

◆ push_back() [2/3]

void QwRootTreeBranchVector::push_back ( const std::string & name,
const char type = 'D' )
inline

Definition at line 197 of file QwRootFile.h.

197 {
198 const std::size_t entry_size = GetTypeSize(type);
199 const std::size_t offset = AlignOffset(m_buffer.size());
200
201 if (offset > m_buffer.capacity()) {
202 throw std::out_of_range("QwRootTreeBranchVector::push_back() requires buffer resize beyond reserved capacity");
203 }
204 if (offset > m_buffer.size()) {
205 m_buffer.resize(offset, 0u);
206 }
207
208 Entry entry{name, offset, entry_size, type};
209 m_entries.push_back(entry);
210
211 const std::size_t required = offset + entry_size;
212 if (required > m_buffer.capacity()) {
213 throw std::out_of_range("QwRootTreeBranchVector::push_back() requires buffer resize beyond reserved capacity");
214 }
215 if (required > m_buffer.size()) {
216 m_buffer.resize(required, 0u);
217 }
218 }
static std::size_t AlignOffset(std::size_t offset)
Definition QwRootFile.h:320

References AlignOffset(), GetTypeSize(), m_buffer, and m_entries.

Referenced by QwADC18_Channel::ConstructBranchAndVector(), QwBeamMod::ConstructBranchAndVector(), QwEPICSEvent::ConstructBranchAndVector(), QwHelicityBase::ConstructBranchAndVector(), QwHelicityDecoder::ConstructBranchAndVector(), QwMollerADC_Channel::ConstructBranchAndVector(), QwPMT_Channel::ConstructBranchAndVector(), QwScaler_Channel< data_mask, data_shift >::ConstructBranchAndVector(), QwSubsystemArray::ConstructBranchAndVector(), QwSubsystemArrayParity::ConstructBranchAndVector(), QwVQWK_Channel::ConstructBranchAndVector(), push_back(), and push_back().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ push_back() [3/3]

void QwRootTreeBranchVector::push_back ( const TString & name,
const char type = 'D' )
inline

Definition at line 221 of file QwRootFile.h.

221 {
222 push_back(std::string(name.Data()), type);
223 }

References push_back().

+ Here is the call graph for this function:

◆ reserve()

void QwRootTreeBranchVector::reserve ( size_type count)
inline

Definition at line 68 of file QwRootFile.h.

68 {
69 m_entries.reserve(count);
70 m_buffer.reserve(sizeof(double)*count);
71 }

References m_buffer, and m_entries.

◆ SetValue() [1/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
Double_t val )
inline

Definition at line 110 of file QwRootFile.h.

110 {
111 const auto& entry = m_entries.at(index);
112 if (entry.type != 'D') {
113 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store double value '" + entry.name + "'");
114 }
115 this->value<Double_t>(index) = val;
116 }

References m_entries, and value().

Referenced by QwADC18_Channel::FillTreeVector(), QwBeamMod::FillTreeVector(), QwEPICSEvent::FillTreeVector(), QwHelicityBase::FillTreeVector(), QwHelicityDecoder::FillTreeVector(), QwMollerADC_Channel::FillTreeVector(), QwPMT_Channel::FillTreeVector(), QwScaler_Channel< data_mask, data_shift >::FillTreeVector(), QwSubsystemArray::FillTreeVector(), QwSubsystemArrayParity::FillTreeVector(), and QwVQWK_Channel::FillTreeVector().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetValue() [2/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
Float_t val )
inline

Definition at line 118 of file QwRootFile.h.

118 {
119 const auto& entry = m_entries.at(index);
120 if (entry.type != 'F') {
121 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store float value '" + entry.name + "'");
122 }
123 this->value<Float_t>(index) = val;
124 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ SetValue() [3/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
Int_t val )
inline

Definition at line 126 of file QwRootFile.h.

126 {
127 const auto& entry = m_entries.at(index);
128 if (entry.type != 'I') {
129 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store int value '" + entry.name + "'");
130 }
131 this->value<Int_t>(index) = val;
132 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ SetValue() [4/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
Long64_t val )
inline

Definition at line 134 of file QwRootFile.h.

134 {
135 const auto& entry = m_entries.at(index);
136 if (entry.type != 'L') {
137 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store long long value '" + entry.name + "'");
138 }
139 this->value<Long64_t>(index) = val;
140 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ SetValue() [5/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
Short_t val )
inline

Definition at line 142 of file QwRootFile.h.

142 {
143 const auto& entry = m_entries.at(index);
144 if (entry.type != 'S') {
145 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store short value '" + entry.name + "'");
146 }
147 this->value<Short_t>(index) = val;
148 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ SetValue() [6/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
UInt_t val )
inline

Definition at line 159 of file QwRootFile.h.

159 {
160 const auto& entry = m_entries.at(index);
161 if (entry.type != 'i') {
162 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store unsigned int value '" + entry.name + "'");
163 }
164 this->value<UInt_t>(index) = val;
165 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ SetValue() [7/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
ULong64_t val )
inline

Definition at line 167 of file QwRootFile.h.

167 {
168 const auto& entry = m_entries.at(index);
169 if (entry.type != 'l') {
170 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store long long value '" + entry.name + "'");
171 }
172 this->value<ULong64_t>(index) = val;
173 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ SetValue() [8/8]

void QwRootTreeBranchVector::SetValue ( size_type index,
UShort_t val )
inline

Definition at line 151 of file QwRootFile.h.

151 {
152 const auto& entry = m_entries.at(index);
153 if (entry.type != 's') {
154 throw std::invalid_argument("Type mismatch: entry type '" + std::string(1, entry.type) + "' cannot store short value '" + entry.name + "'");
155 }
156 this->value<UShort_t>(index) = val;
157 }

References m_entries, and value().

+ Here is the call graph for this function:

◆ shrink_to_fit()

void QwRootTreeBranchVector::shrink_to_fit ( )
inline

Definition at line 73 of file QwRootFile.h.

73 {
74 m_entries.shrink_to_fit();
75 m_buffer.shrink_to_fit();
76 }

References m_buffer, and m_entries.

◆ size()

◆ value() [1/2]

template<typename T>
T & QwRootTreeBranchVector::value ( size_type index)
inline

Definition at line 97 of file QwRootFile.h.

97 {
98 auto& entry = m_entries.at(index);
99 return *reinterpret_cast<T*>(m_buffer.data() + entry.offset);
100 }

References m_buffer, and m_entries.

Referenced by FormatValue(), operator[](), operator[](), SetValue(), SetValue(), SetValue(), SetValue(), SetValue(), SetValue(), SetValue(), and SetValue().

+ Here is the caller graph for this function:

◆ value() [2/2]

template<typename T>
const T & QwRootTreeBranchVector::value ( size_type index) const
inline

Definition at line 103 of file QwRootFile.h.

103 {
104 const auto& entry = m_entries.at(index);
105 return *reinterpret_cast<const T*>(m_buffer.data() + entry.offset);
106 }

References m_buffer, and m_entries.

Field Documentation

◆ m_buffer

std::vector<std::uint8_t> QwRootTreeBranchVector::m_buffer
private

◆ m_entries

std::vector<Entry> QwRootTreeBranchVector::m_entries
private

The documentation for this class was generated from the following file: