Project
Loading...
Searching...
No Matches
MCTruthContainer.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15
16#ifndef ALICEO2_DATAFORMATS_MCTRUTH_H_
17#define ALICEO2_DATAFORMATS_MCTRUTH_H_
18
19#include "GPUCommonRtypes.h" // to have the ClassDef macros
20#include <cstdint> // uint8_t etc
21#include <cassert>
22#include <stdexcept>
23#include <gsl/span> // for guideline support library span
24#include <type_traits>
25#include <cstring> // memmove, memcpy
26#include <memory>
27#include <vector>
28
29// type traits are needed for the compile time consistency check
30// maybe to be moved out of Framework first
31//#include "Framework/TypeTraits.h"
32
33namespace o2
34{
35class MCCompLabel;
36namespace dataformats
37{
38
43 MCTruthHeaderElement() = default; // for ROOT IO
44
45 MCTruthHeaderElement(uint32_t i) : index(i) {}
46 uint32_t index = (uint32_t)-1; // the index into the actual MC track storage (-1 if invalid)
48};
49
81template <typename TruthElement>
83{
84 private:
85 // for the moment we require the truth element to be messageable in order to simply flatten the object
86 // if it turnes out that other types are required this needs to be extended and method flatten needs to
87 // be conditionally added
88 // TODO: activate this check
89 // static_assert(o2::framework::is_messageable<TruthElement>::value, "truth element type must be messageable");
90
91 std::vector<MCTruthHeaderElement> mHeaderArray; // the header structure array serves as an index into the actual storage
92 std::vector<TruthElement> mTruthArray; // the buffer containing the actual truth information
96 std::vector<char> mStreamerData; // buffer used for streaming a flat raw buffer
97
98 size_t getSize(uint32_t dataindex) const
99 {
100 // calculate size / number of labels from a difference in pointed indices
101 const auto size = (dataindex < getIndexedSize() - 1)
102 ? getMCTruthHeader(dataindex + 1).index - getMCTruthHeader(dataindex).index
103 : getNElements() - getMCTruthHeader(dataindex).index;
104 return size;
105 }
106
107 public:
108 // constructor
109 MCTruthContainer() = default;
110 // destructor
111 ~MCTruthContainer() = default;
112 // copy constructor
114 // move constructor
116 // construct from raw data
117 MCTruthContainer(std::vector<MCTruthHeaderElement>& header, std::vector<TruthElement>& truthArray)
118 {
119 setFrom(header, truthArray);
120 }
121 // assignment operator
123 // move assignment operator
125
127 struct FlatHeader {
128 uint8_t version = 1;
130 uint8_t sizeofTruthElement = sizeof(TruthElement);
131 uint8_t reserved = 0;
134 };
135
136 // access
137 MCTruthHeaderElement const& getMCTruthHeader(uint32_t dataindex) const { return mHeaderArray[dataindex]; }
138 // access the element directly (can be encapsulated better away)... needs proper element index
139 // which can be obtained from the MCTruthHeader startposition and size
140 TruthElement const& getElement(uint32_t elementindex) const { return mTruthArray[elementindex]; }
141 // return the number of original data indexed here
142 size_t getIndexedSize() const { return mHeaderArray.size(); }
143 // return the number of elements (labels) pointed to in this container
144 size_t getNElements() const { return mTruthArray.size(); }
145 // return unterlaying vector of elements
146 const std::vector<TruthElement>& getTruthArray() const
147 {
148 return mTruthArray;
149 }
150
151 // get individual "view" container for a given data index
152 // the caller can do modifications on this view (such as sorting)
153 gsl::span<TruthElement> getLabels(uint32_t dataindex)
154 {
155 if (dataindex >= getIndexedSize()) {
156 return gsl::span<TruthElement>();
157 }
158 return gsl::span<TruthElement>(&mTruthArray[getMCTruthHeader(dataindex).index], getSize(dataindex));
159 }
160
161 // get individual const "view" container for a given data index
162 // the caller can't do modifications on this view
163 gsl::span<const TruthElement> getLabels(uint32_t dataindex) const
164 {
165 if (dataindex >= getIndexedSize()) {
166 return gsl::span<const TruthElement>();
167 }
168 return gsl::span<const TruthElement>(&mTruthArray[getMCTruthHeader(dataindex).index], getSize(dataindex));
169 }
170
171 void clear()
172 {
173 mHeaderArray.clear();
174 mTruthArray.clear();
175 }
176
177 // clear and force freeing the memory
179 {
180 clear();
181 // this forces the desctructor being called on existing buffers
182 mHeaderArray = std::vector<MCTruthHeaderElement>();
183 mTruthArray = std::vector<TruthElement>();
184 mStreamerData = std::vector<char>();
185 }
186
187 // add element for a particular dataindex
188 // at the moment only strictly consecutive modes are supported
189 // noElement indicates that actually no label/element will be added for this dataindex
190 void addElement(uint32_t dataindex, TruthElement const& element, bool noElement = false)
191 {
192 if (dataindex < mHeaderArray.size()) {
193 // look if we have something for this dataindex already
194 // must currently be the last one
195 if (dataindex != (mHeaderArray.size() - 1)) {
196 throw std::runtime_error("MCTruthContainer: unsupported code path");
197 }
198 } else {
199 // assert(dataindex == mHeaderArray.size());
200
201 // add empty holes
202 int holes = dataindex - mHeaderArray.size();
203 assert(holes >= 0);
204 for (int i = 0; i < holes; ++i) {
205 mHeaderArray.emplace_back(mTruthArray.size());
206 }
207 // add a new one
208 mHeaderArray.emplace_back(mTruthArray.size());
209 }
210 if (!noElement) {
211 mTruthArray.emplace_back(element);
212 }
213 }
214
216 void addNoLabelIndex(uint32_t dataindex)
217 {
218 addElement(dataindex, TruthElement(), true);
219 }
220
221 // convenience interface to add multiple labels at once
222 // can use elements of any assignable type or sub-type
223 template <typename CompatibleLabel>
224 void addElements(uint32_t dataindex, gsl::span<CompatibleLabel> elements)
225 {
226 static_assert(std::is_same<TruthElement, CompatibleLabel>::value ||
227 std::is_assignable<TruthElement, CompatibleLabel>::value ||
228 std::is_base_of<TruthElement, CompatibleLabel>::value,
229 "Need to add compatible labels");
230 addNoLabelIndex(dataindex);
231 for (auto& e : elements) {
232 addElement(dataindex, e);
233 }
234 }
235
236 template <typename CompatibleLabel>
237 void addElements(uint32_t dataindex, const std::vector<CompatibleLabel>& v)
238 {
239 using B = typename std::remove_const<CompatibleLabel>::type;
240 auto s = gsl::span<CompatibleLabel>(const_cast<B*>(&v[0]), v.size());
241 addElements(dataindex, s);
242 }
243
244 // Add element at last position or for a previous index
245 // (at random access position).
246 // This might be a slow process since data has to be moved internally
247 // so this function should be used with care.
248 void addElementRandomAccess(uint32_t dataindex, TruthElement const& element)
249 {
250 if (dataindex >= mHeaderArray.size()) {
251 // a new dataindex -> push element at back
252
253 // we still forbid to leave holes
254 assert(dataindex == mHeaderArray.size());
255
256 mHeaderArray.resize(dataindex + 1);
257 mHeaderArray[dataindex] = mTruthArray.size();
258 mTruthArray.emplace_back(element);
259 } else {
260 // if appending at end use fast function
261 if (dataindex == mHeaderArray.size() - 1) {
262 addElement(dataindex, element);
263 return;
264 }
265
266 // existing dataindex
267 // have to:
268 // a) move data;
269 // b) insert new element;
270 // c) adjust indices of all headers right to this
271 auto currentindex = mHeaderArray[dataindex].index;
272 auto lastindex = currentindex + getSize(dataindex);
273 assert(currentindex >= 0);
274
275 // resize truth array
276 mTruthArray.resize(mTruthArray.size() + 1);
277 // move data (have to do from right to left)
278 for (int i = mTruthArray.size() - 1; i > lastindex; --i) {
279 mTruthArray[i] = mTruthArray[i - 1];
280 }
281 // insert new element
282 mTruthArray[lastindex] = element;
283
284 // fix headers
285 for (uint32_t i = dataindex + 1; i < mHeaderArray.size(); ++i) {
286 auto oldindex = mHeaderArray[i].index;
287 mHeaderArray[i].index = (oldindex != (uint32_t)-1) ? oldindex + 1 : oldindex;
288 }
289 }
290 }
291
292 // Set container directly from header + truthArray.
293 void setFrom(std::vector<MCTruthHeaderElement>& header, std::vector<TruthElement>& truthArray)
294 {
295 mHeaderArray = std::move(header);
296 mTruthArray = std::move(truthArray);
297 }
298
299 // merge another container to the back of this one
301 {
302 const auto oldtruthsize = mTruthArray.size();
303 const auto oldheadersize = mHeaderArray.size();
304
305 // copy from other
306 std::copy(other.mHeaderArray.begin(), other.mHeaderArray.end(), std::back_inserter(mHeaderArray));
307 std::copy(other.mTruthArray.begin(), other.mTruthArray.end(), std::back_inserter(mTruthArray));
308
309 // adjust information of newly attached part
310 for (uint32_t i = oldheadersize; i < mHeaderArray.size(); ++i) {
311 mHeaderArray[i].index += oldtruthsize;
312 }
313 }
314
315 // merge part of another container ("n" entries starting from "from") to the back of this one
316 void mergeAtBack(MCTruthContainer<TruthElement> const& other, size_t from, size_t n)
317 {
318 const auto oldtruthsize = mTruthArray.size();
319 const auto oldheadersize = mHeaderArray.size();
320 auto endIdx = from + n;
321 assert(endIdx <= other.mHeaderArray.size());
322 const auto* headBeg = &other.mHeaderArray[from];
323 const auto* headEnd = headBeg + n;
324 const auto* trtArrBeg = &other.mTruthArray[other.getMCTruthHeader(from).index];
325 const auto* trtArrEnd = (endIdx == other.mHeaderArray.size()) ? (&other.mTruthArray.back()) + 1 : &other.mTruthArray[other.getMCTruthHeader(endIdx).index];
326
327 // copy from other
328 std::copy(headBeg, headEnd, std::back_inserter(mHeaderArray));
329 std::copy(trtArrBeg, trtArrEnd, std::back_inserter(mTruthArray));
330 long offset = long(oldtruthsize) - other.getMCTruthHeader(from).index;
331 // adjust information of newly attached part
332 for (uint32_t i = oldheadersize; i < mHeaderArray.size(); ++i) {
333 mHeaderArray[i].index += offset;
334 }
335 }
336
341 template <typename ContainerType>
342 size_t flatten_to(ContainerType& container) const
343 {
344 size_t bufferSize = sizeof(FlatHeader) + sizeof(MCTruthHeaderElement) * mHeaderArray.size() + sizeof(TruthElement) * mTruthArray.size();
345 container.resize((bufferSize / sizeof(typename ContainerType::value_type)) + ((bufferSize % sizeof(typename ContainerType::value_type)) > 0 ? 1 : 0));
346 char* target = reinterpret_cast<char*>(container.data());
347 auto& flatheader = *reinterpret_cast<FlatHeader*>(target);
348 target += sizeof(FlatHeader);
349 flatheader.version = 1;
350 flatheader.sizeofHeaderElement = sizeof(MCTruthHeaderElement);
351 flatheader.sizeofTruthElement = sizeof(TruthElement);
352 flatheader.reserved = 0;
353 flatheader.nofHeaderElements = mHeaderArray.size();
354 flatheader.nofTruthElements = mTruthArray.size();
355 size_t copySize = flatheader.sizeofHeaderElement * flatheader.nofHeaderElements;
356 memcpy(target, mHeaderArray.data(), copySize);
357 target += copySize;
358 copySize = flatheader.sizeofTruthElement * flatheader.nofTruthElements;
359 memcpy(target, mTruthArray.data(), copySize);
360 return bufferSize;
361 }
362
366 void restore_from(const char* buffer, size_t bufferSize)
367 {
368 if (buffer == nullptr || bufferSize < sizeof(FlatHeader)) {
369 return;
370 }
371 auto* source = buffer;
372 auto& flatheader = *reinterpret_cast<FlatHeader const*>(source);
373 source += sizeof(FlatHeader);
374 if (bufferSize < sizeof(FlatHeader) + flatheader.sizeofHeaderElement * flatheader.nofHeaderElements + flatheader.sizeofTruthElement * flatheader.nofTruthElements) {
375 throw std::runtime_error("inconsistent buffer size: too small");
376 return;
377 }
378 if (flatheader.sizeofHeaderElement != sizeof(MCTruthHeaderElement) || flatheader.sizeofTruthElement != sizeof(TruthElement)) {
379 // not yet handled
380 throw std::runtime_error("member element sizes don't match");
381 }
382 // TODO: with a spectator memory ressource the vectors can be built directly
383 // over the original buffer, there is the implementation for a memory ressource
384 // working on a FairMQ message, here we would need two memory resources over
385 // the two ranges of the input buffer
386 // for now doing a copy
387 mHeaderArray.resize(flatheader.nofHeaderElements);
388 mTruthArray.resize(flatheader.nofTruthElements);
389 size_t copySize = flatheader.sizeofHeaderElement * flatheader.nofHeaderElements;
390 memcpy(mHeaderArray.data(), source, copySize);
391 source += copySize;
392 copySize = flatheader.sizeofTruthElement * flatheader.nofTruthElements;
393 memcpy(mTruthArray.data(), source, copySize);
394 }
395
397 template <typename Stream>
399 {
400 stream << "MCTruthContainer index = " << getIndexedSize() << " for " << getNElements() << " elements(s), flat buffer size " << mStreamerData.size() << "\n";
401 }
402
408 void inflate()
409 {
410 if (mHeaderArray.size() > 0) {
411 mStreamerData.clear();
412 return;
413 }
414 restore_from(mStreamerData.data(), mStreamerData.size());
415 mStreamerData.clear();
416 }
417
423 void deflate()
424 {
425 if (mStreamerData.size() > 0) {
426 clear();
427 return;
428 }
429 mStreamerData.clear();
430 flatten_to(mStreamerData);
431 clear();
432 }
433
435}; // end class
436
438
439} // namespace dataformats
440} // namespace o2
441
442#endif
int32_t i
uint32_t version
Definition RawData.h:8
Definition B.h:16
A container to hold and manage MC truth information/labels.
gsl::span< TruthElement > getLabels(uint32_t dataindex)
MCTruthContainer & operator=(const MCTruthContainer &other)=default
void addNoLabelIndex(uint32_t dataindex)
adds a data index that has no label
TruthElement const & getElement(uint32_t elementindex) const
ClassDefNV(MCTruthContainer, 2)
void addElements(uint32_t dataindex, const std::vector< CompatibleLabel > &v)
void mergeAtBack(MCTruthContainer< TruthElement > const &other)
gsl::span< const TruthElement > getLabels(uint32_t dataindex) const
void addElements(uint32_t dataindex, gsl::span< CompatibleLabel > elements)
MCTruthContainer & operator=(MCTruthContainer &&other)=default
void addElement(uint32_t dataindex, TruthElement const &element, bool noElement=false)
const std::vector< TruthElement > & getTruthArray() const
void restore_from(const char *buffer, size_t bufferSize)
MCTruthContainer(const MCTruthContainer &other)=default
void setFrom(std::vector< MCTruthHeaderElement > &header, std::vector< TruthElement > &truthArray)
MCTruthContainer(std::vector< MCTruthHeaderElement > &header, std::vector< TruthElement > &truthArray)
MCTruthHeaderElement const & getMCTruthHeader(uint32_t dataindex) const
void mergeAtBack(MCTruthContainer< TruthElement > const &other, size_t from, size_t n)
void addElementRandomAccess(uint32_t dataindex, TruthElement const &element)
MCTruthContainer(MCTruthContainer &&other)=default
size_t flatten_to(ContainerType &container) const
void print(Stream &stream)
Print some info.
GLdouble n
Definition glcorearb.h:1982
GLuint buffer
Definition glcorearb.h:655
GLsizeiptr size
Definition glcorearb.h:659
const GLdouble * v
Definition glcorearb.h:832
GLuint index
Definition glcorearb.h:781
GLsizei GLsizei GLchar * source
Definition glcorearb.h:798
GLenum target
Definition glcorearb.h:1641
GLintptr offset
Definition glcorearb.h:660
GLuint GLuint stream
Definition glcorearb.h:1806
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Simple struct having information about truth elements for particular indices: how many associations w...
ClassDefNV(MCTruthHeaderElement, 1)
VectorOfTObjectPtrs other