Project
Loading...
Searching...
No Matches
ConstMCTruthContainer.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 O2_CONSTMCTRUTHCONTAINER_H
17#define O2_CONSTMCTRUTHCONTAINER_H
18
20#ifndef GPUCA_STANDALONE
21#include <Framework/Traits.h>
22#endif
23
24class TTree;
25
26namespace o2
27{
28namespace dataformats
29{
30
38template <typename TruthElement>
39class ConstMCTruthContainer : public std::vector<char>
40{
41 public:
42 // (unfortunately we need these constructors for DPL)
43 using std::vector<char>::vector;
45
46 // const data access
47 // get individual const "view" container for a given data index
48 // the caller can't do modifications on this view
49 MCTruthHeaderElement const& getMCTruthHeader(uint32_t dataindex) const
50 {
51 return getHeaderStart()[dataindex];
52 }
53
54 gsl::span<const TruthElement> getLabels(uint32_t dataindex) const
55 {
56 if (dataindex >= getIndexedSize()) {
57 return gsl::span<const TruthElement>();
58 }
59 const auto start = getMCTruthHeader(dataindex).index;
60 const auto labelsptr = getLabelStart();
61 return gsl::span<const TruthElement>(&labelsptr[start], getSize(dataindex));
62 }
63
64 // return the number of original data indexed here
65 size_t getIndexedSize() const { return size() >= sizeof(FlatHeader) ? getHeader().nofHeaderElements : 0; }
66
67 // return the number of labels managed in this container
68 size_t getNElements() const { return size() >= sizeof(FlatHeader) ? getHeader().nofTruthElements : 0; }
69
70 private:
71 using FlatHeader = typename MCTruthContainer<TruthElement>::FlatHeader;
72
73 size_t getSize(uint32_t dataindex) const
74 {
75 // calculate size / number of labels from a difference in pointed indices
76 const auto size = (dataindex < getIndexedSize() - 1)
77 ? getMCTruthHeader(dataindex + 1).index - getMCTruthHeader(dataindex).index
78 : getNElements() - getMCTruthHeader(dataindex).index;
79 return size;
80 }
81
85 TruthElement const* getLabelStart() const
86 {
87 auto* source = &(*this)[0];
88 auto flatheader = getHeader();
89 source += sizeof(FlatHeader);
90 const size_t headerSize = flatheader.sizeofHeaderElement * flatheader.nofHeaderElements;
91 source += headerSize;
92 return (TruthElement const*)source;
93 }
94
95 FlatHeader const& getHeader() const
96 {
97 const auto* source = &(*this)[0];
98 const auto& flatheader = *reinterpret_cast<FlatHeader const*>(source);
99 return flatheader;
100 }
101
102 MCTruthHeaderElement const* getHeaderStart() const
103 {
104 auto* source = &(*this)[0];
105 source += sizeof(FlatHeader);
106 return (MCTruthHeaderElement const*)source;
107 }
108};
109} // namespace dataformats
110} // namespace o2
111
112// This is done so that DPL treats this container as a vector.
113// In particular in enables
114// a) --> snapshot without ROOT dictionary (as a flat buffer)
115// b) --> requesting the resource in shared mem using make<T>
116#ifndef GPUCA_STANDALONE
117namespace o2::framework
118{
119template <typename T>
120struct is_specialization<o2::dataformats::ConstMCTruthContainer<T>, std::vector> : std::true_type {
121};
122} // namespace o2::framework
123#endif
124
125namespace o2
126{
127namespace dataformats
128{
129
130// A "view" label container without owning the storage (similar to gsl::span)
131template <typename TruthElement>
133{
134 public:
135 ConstMCTruthContainerView(gsl::span<const char> const bufferview) : mStorage(bufferview){};
136 ConstMCTruthContainerView(ConstMCTruthContainer<TruthElement> const& cont) : mStorage(gsl::span<const char>(cont)){};
137 // be explicit that we want nullptr / 0 for an uninitialized container (needs (void)0 to avoid false codechecker warning)
138 ConstMCTruthContainerView() : mStorage{nullptr, static_cast<gsl::span<const char>::size_type>(0)}
139 {
140 (void)0;
141 }
142 // ConstMCTruthContainerView(const ConstMCTruthContainerView&) = default;
143
144 // const data access
145 // get individual const "view" container for a given data index
146 // the caller can't do modifications on this view
147 MCTruthHeaderElement const& getMCTruthHeader(uint32_t dataindex) const
148 {
149 return getHeaderStart()[dataindex];
150 }
151
152 gsl::span<const TruthElement> getLabels(uint32_t dataindex) const
153 {
154 if (dataindex >= getIndexedSize()) {
155 return gsl::span<const TruthElement>();
156 }
157 const auto start = getMCTruthHeader(dataindex).index;
158 const auto labelsptr = getLabelStart();
159 return gsl::span<const TruthElement>(&labelsptr[start], getSize(dataindex));
160 }
161
162 // return the number of original data indexed here
163 size_t getIndexedSize() const { return (size_t)mStorage.size() >= sizeof(FlatHeader) ? getHeader().nofHeaderElements : 0; }
164
165 // return the number of labels managed in this container
166 size_t getNElements() const { return (size_t)mStorage.size() >= sizeof(FlatHeader) ? getHeader().nofTruthElements : 0; }
167
168 // return underlying buffer
169 const gsl::span<const char>& getBuffer() const { return mStorage; }
170
171 private:
172 gsl::span<const char> mStorage;
173
174 using FlatHeader = typename MCTruthContainer<TruthElement>::FlatHeader;
175
176 size_t getSize(uint32_t dataindex) const
177 {
178 // calculate size / number of labels from a difference in pointed indices
179 const auto size = (dataindex < getIndexedSize() - 1)
180 ? getMCTruthHeader(dataindex + 1).index - getMCTruthHeader(dataindex).index
181 : getNElements() - getMCTruthHeader(dataindex).index;
182 return size;
183 }
184
188 TruthElement const* getLabelStart() const
189 {
190 auto* source = &(mStorage)[0];
191 auto flatheader = getHeader();
192 source += sizeof(FlatHeader);
193 const size_t headerSize = flatheader.sizeofHeaderElement * flatheader.nofHeaderElements;
194 source += headerSize;
195 return (TruthElement const*)source;
196 }
197
198 FlatHeader const& getHeader() const
199 {
200 const auto* source = &(mStorage)[0];
201 const auto& flatheader = *reinterpret_cast<FlatHeader const*>(source);
202 return flatheader;
203 }
204
205 MCTruthHeaderElement const* getHeaderStart() const
206 {
207 auto* source = &(mStorage)[0];
208 source += sizeof(FlatHeader);
209 return (MCTruthHeaderElement const*)source;
210 }
211};
212
215
217{
218 public:
221 static ConstMCTruthContainer<o2::MCCompLabel>* loadFromTTree(TTree* tree, std::string const& brname, int entry);
222};
223
224} // namespace dataformats
225} // namespace o2
226
227#endif //O2_CONSTMCTRUTHCONTAINER_H
Definition of a container to keep Monte Carlo truth external to simulation objects.
gsl::span< const TruthElement > getLabels(uint32_t dataindex) const
ConstMCTruthContainerView(ConstMCTruthContainer< TruthElement > const &cont)
MCTruthHeaderElement const & getMCTruthHeader(uint32_t dataindex) const
ConstMCTruthContainerView(gsl::span< const char > const bufferview)
const gsl::span< const char > & getBuffer() const
A read-only version of MCTruthContainer allowing for storage optimisation.
gsl::span< const TruthElement > getLabels(uint32_t dataindex) const
MCTruthHeaderElement const & getMCTruthHeader(uint32_t dataindex) const
static ConstMCTruthContainer< o2::MCCompLabel > * loadFromTTree(TTree *tree, std::string const &brname, int entry)
A container to hold and manage MC truth information/labels.
GLuint entry
Definition glcorearb.h:5735
GLsizeiptr size
Definition glcorearb.h:659
GLsizei GLsizei GLchar * source
Definition glcorearb.h:798
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint start
Definition glcorearb.h:469
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
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...
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))