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 53 of file QwRootFile.h.

Member Typedef Documentation

◆ size_type

using QwRootTreeBranchVector::size_type = std::size_t

Definition at line 62 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 318 of file QwRootFile.h.

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

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 178 of file QwRootFile.h.

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

References m_buffer, and m_entries.

Referenced by QwHelicity::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 187 of file QwRootFile.h.

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

References m_buffer, and m_entries.

◆ clear()

void QwRootTreeBranchVector::clear ( )
inline

Definition at line 76 of file QwRootFile.h.

76 {
77 m_entries.clear();
78 m_buffer.clear();
79 }

References m_buffer, and m_entries.

◆ data() [1/2]

const void * QwRootTreeBranchVector::data ( ) const
inlinenoexcept

Definition at line 174 of file QwRootFile.h.

174{ return m_buffer.data(); }

References m_buffer.

◆ data() [2/2]

void * QwRootTreeBranchVector::data ( )
inlinenoexcept

Definition at line 173 of file QwRootFile.h.

173{ return m_buffer.data(); }

References m_buffer.

◆ data_size()

size_type QwRootTreeBranchVector::data_size ( ) const
inlinenoexcept

Definition at line 175 of file QwRootFile.h.

175{ 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 243 of file QwRootFile.h.

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

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 82 of file QwRootFile.h.

82{ return m_entries.empty(); }

References m_entries.

◆ FormatNumeric()

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

Definition at line 350 of file QwRootFile.h.

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

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 326 of file QwRootFile.h.

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

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 295 of file QwRootFile.h.

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

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 228 of file QwRootFile.h.

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

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 90 of file QwRootFile.h.

90 {
91 return value<T>(index);
92 }

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 85 of file QwRootFile.h.

85 {
86 return value<T>(index);
87 }

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 224 of file QwRootFile.h.

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

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 195 of file QwRootFile.h.

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

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

Referenced by QwADC18_Channel::ConstructBranchAndVector(), QwBeamMod::ConstructBranchAndVector(), QwEPICSEvent::ConstructBranchAndVector(), QwHelicity::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 219 of file QwRootFile.h.

219 {
220 push_back(std::string(name.Data()), type);
221 }

References push_back().

+ Here is the call graph for this function:

◆ reserve()

void QwRootTreeBranchVector::reserve ( size_type count)
inline

Definition at line 66 of file QwRootFile.h.

66 {
67 m_entries.reserve(count);
68 m_buffer.reserve(sizeof(double)*count);
69 }

References m_buffer, and m_entries.

◆ SetValue() [1/8]

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

Definition at line 108 of file QwRootFile.h.

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

References m_entries, and value().

Referenced by QwADC18_Channel::FillTreeVector(), QwBeamMod::FillTreeVector(), QwEPICSEvent::FillTreeVector(), QwHelicity::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 116 of file QwRootFile.h.

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

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 124 of file QwRootFile.h.

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

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 132 of file QwRootFile.h.

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

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 140 of file QwRootFile.h.

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

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 157 of file QwRootFile.h.

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

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 165 of file QwRootFile.h.

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

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 149 of file QwRootFile.h.

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

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 71 of file QwRootFile.h.

71 {
72 m_entries.shrink_to_fit();
73 m_buffer.shrink_to_fit();
74 }

References m_buffer, and m_entries.

◆ size()

◆ value() [1/2]

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

Definition at line 95 of file QwRootFile.h.

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

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 101 of file QwRootFile.h.

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

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: