Project
Loading...
Searching...
No Matches
testMCTruthContainer.cxx
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
12#define BOOST_TEST_MODULE Test MCTruthContainer class
13#define BOOST_TEST_MAIN
14#define BOOST_TEST_DYN_LINK
15#include <boost/test/unit_test.hpp>
20#include <algorithm>
21#include <iostream>
22#include <TFile.h>
23#include <TTree.h>
24
25namespace o2
26{
28{
29 using TruthElement = long;
31 container.addElement(0, TruthElement(1));
32 container.addElement(0, TruthElement(2));
33 container.addElement(1, TruthElement(1));
34 container.addElement(2, TruthElement(10));
35 container.addNoLabelIndex(3);
36 container.addElement(4, TruthElement(4));
37
38 // not supported, must throw
39 BOOST_CHECK_THROW(container.addElement(0, TruthElement(0)), std::runtime_error);
40
41 // check header/index information
42 BOOST_CHECK(container.getMCTruthHeader(0).index == 0);
43 BOOST_CHECK(container.getMCTruthHeader(1).index == 2);
44 BOOST_CHECK(container.getMCTruthHeader(2).index == 3);
45
46 // check MC truth information
47 BOOST_CHECK(container.getElement(0) == 1);
48 BOOST_CHECK(container.getElement(1) == 2);
49 BOOST_CHECK(container.getElement(2) == 1);
50 BOOST_CHECK(container.getElement(3) == 10);
51
52 // get iterable container view on labels for index 0
53 auto view = container.getLabels(0);
54 BOOST_CHECK(view.size() == 2);
55 BOOST_CHECK(view[0] == 1);
56 BOOST_CHECK(view[1] == 2);
57 // try to sort the view
58 std::sort(view.begin(), view.end(), [](TruthElement a, TruthElement b) { return a > b; });
59 BOOST_CHECK(view[0] == 2);
60 BOOST_CHECK(view[1] == 1);
61 // verify sorting took effect inside container
62 auto view2 = container.getLabels(0);
63 BOOST_CHECK(view2[0] == 2);
64 BOOST_CHECK(view2[1] == 1);
65
66 // same for another data index
67 view = container.getLabels(2);
68 BOOST_CHECK(view.size() == 1);
69 BOOST_CHECK(view[0] == 10);
70
71 // try to get something invalid
72 view = container.getLabels(10);
73 BOOST_CHECK(view.size() == 0);
74
75 // test assignment/copy
76 auto copy = container;
77 view = copy.getLabels(2);
78 BOOST_CHECK(view.size() == 1);
79 BOOST_CHECK(view[0] == 10);
80
81 // add multiple labels (to last index which already had a label)
82 std::vector<TruthElement> newlabels = {101, 102, 103};
83 container.addElements(4, newlabels);
84 view = container.getLabels(4);
85 BOOST_CHECK(view.size() == 4);
86 BOOST_CHECK(view[0] == 4);
87 BOOST_CHECK(view[1] == 101);
88 BOOST_CHECK(view[2] == 102);
89 BOOST_CHECK(view[3] == 103);
90
91 // check empty labels case
92 view = container.getLabels(3);
93 BOOST_CHECK(view.size() == 0);
94
95 // add empty label vector
96 std::vector<TruthElement> newlabels2 = {};
97 container.addElements(5, newlabels2);
98 view = container.getLabels(5);
99 BOOST_CHECK(view.size() == 0);
100
101 // test merging
102 {
104 container1.addElement(0, TruthElement(1));
105 container1.addElement(0, TruthElement(2));
106 container1.addElement(1, TruthElement(1));
107 container1.addElement(2, TruthElement(10));
108
110 container2.addElement(0, TruthElement(11));
111 container2.addElement(0, TruthElement(12));
112 container2.addElement(1, TruthElement(1));
113 container2.addElement(2, TruthElement(10));
114
116
117 container1.mergeAtBack(container2);
118
119 containerA.mergeAtBack(container1, 0, 2);
120 containerA.mergeAtBack(container1, 2, 2);
121
122 auto lview = container1.getLabels(3); //
123 auto lviewA = containerA.getLabels(3);
124 BOOST_CHECK(lview.size() == 2);
125 BOOST_CHECK(lview[0] == 11);
126 BOOST_CHECK(lview[1] == 12);
127 BOOST_CHECK(container1.getIndexedSize() == 6);
128 BOOST_CHECK(container1.getNElements() == 8);
129 BOOST_CHECK(lview.size() == lviewA.size());
130 BOOST_CHECK(lview[0] == lviewA[0] && lview[1] == lviewA[1]);
131 }
132}
133
134BOOST_AUTO_TEST_CASE(MCTruth_RandomAccess)
135{
136 using TruthElement = long;
138 container.addElementRandomAccess(0, TruthElement(1));
139 container.addElementRandomAccess(0, TruthElement(2));
140 container.addElementRandomAccess(1, TruthElement(1));
141 container.addElementRandomAccess(2, TruthElement(10));
142 container.addElementRandomAccess(1, TruthElement(5));
143 container.addElementRandomAccess(0, TruthElement(5));
144 // add element at end
145 container.addElement(3, TruthElement(20));
146 container.addElement(3, TruthElement(21));
147
148 // check header/index information
149 BOOST_CHECK(container.getMCTruthHeader(0).index == 0);
150 BOOST_CHECK(container.getMCTruthHeader(1).index == 3);
151 BOOST_CHECK(container.getMCTruthHeader(2).index == 5);
152 BOOST_CHECK(container.getMCTruthHeader(3).index == 6);
153
154 // get iterable container view on labels
155 {
156 auto view = container.getLabels(1);
157 BOOST_CHECK(view.size() == 2);
158 BOOST_CHECK(view[0] == 1);
159 BOOST_CHECK(view[1] == 5);
160 }
161
162 {
163 auto view = container.getLabels(3);
164 BOOST_CHECK(view.size() == 2);
165 BOOST_CHECK(view[0] == 20);
166 BOOST_CHECK(view[1] == 21);
167 }
168}
169
170BOOST_AUTO_TEST_CASE(MCTruthContainer_flatten)
171{
172 using TruthElement = long;
173 using TruthContainer = dataformats::MCTruthContainer<TruthElement>;
174 TruthContainer container;
175 container.addElement(0, TruthElement(1));
176 container.addElement(0, TruthElement(2));
177 container.addElement(1, TruthElement(1));
178 container.addElement(2, TruthElement(10));
179
180 std::vector<char> buffer;
181 container.flatten_to(buffer);
182 BOOST_REQUIRE(buffer.size() > sizeof(TruthContainer::FlatHeader));
183 auto& header = *reinterpret_cast<TruthContainer::FlatHeader*>(buffer.data());
184 BOOST_CHECK(header.nofHeaderElements == container.getIndexedSize());
185 BOOST_CHECK(header.nofTruthElements == container.getNElements());
186
187 TruthContainer restoredContainer;
188 restoredContainer.restore_from(buffer.data(), buffer.size());
189
190 // check header/index information
191 BOOST_CHECK(restoredContainer.getMCTruthHeader(0).index == 0);
192 BOOST_CHECK(restoredContainer.getMCTruthHeader(1).index == 2);
193 BOOST_CHECK(restoredContainer.getMCTruthHeader(2).index == 3);
194
195 // check MC truth information
196 BOOST_CHECK(restoredContainer.getElement(0) == 1);
197 BOOST_CHECK(restoredContainer.getElement(1) == 2);
198 BOOST_CHECK(restoredContainer.getElement(2) == 1);
199 BOOST_CHECK(restoredContainer.getElement(3) == 10);
200
201 // check the special version ConstMCTruthContainer
202 using ConstMCTruthContainer = dataformats::ConstMCTruthContainer<TruthElement>;
203 ConstMCTruthContainer cc;
204 container.flatten_to(cc);
205
206 BOOST_CHECK(cc.getIndexedSize() == container.getIndexedSize());
207 BOOST_CHECK(cc.getNElements() == container.getNElements());
208 BOOST_CHECK(cc.getLabels(0).size() == container.getLabels(0).size());
209 BOOST_CHECK(cc.getLabels(1).size() == container.getLabels(1).size());
210 BOOST_CHECK(cc.getLabels(2).size() == container.getLabels(2).size());
211 BOOST_CHECK(cc.getLabels(2)[0] == container.getLabels(2)[0]);
212 BOOST_CHECK(cc.getLabels(2)[0] == 10);
213}
214
215BOOST_AUTO_TEST_CASE(LabelContainer_noncont)
216{
217 using TruthElement = long;
218 // creates a container where labels might not be contiguous
220 container.addLabel(0, TruthElement(1));
221 container.addLabel(1, TruthElement(1));
222 container.addLabel(0, TruthElement(10));
223 container.addLabel(2, TruthElement(20));
224
225 auto view = container.getLabels(0);
226 BOOST_CHECK(view.size() == 2);
227 for (auto& e : view) {
228 std::cerr << e << "\n";
229 }
230
231 {
232 auto view2 = container.getLabels(10);
233 BOOST_CHECK(view2.size() == 0);
234 for (auto& e : view2) {
235 // should not come here
236 BOOST_CHECK(false);
237 }
238 }
239
240 {
241 auto view2 = container.getLabels(2);
242 BOOST_CHECK(view2.size() == 1);
243 std::cout << "ELEMENTS OF LABEL 2\n";
244 for (auto& e : view2) {
245 // should not come here
246 std::cout << e << "\n";
247 }
248 std::cout << "------\n";
249 }
250
251 std::vector<TruthElement> v;
252 container.fillVectorOfLabels(0, v);
253 BOOST_CHECK(v.size() == 2);
254 std::sort(v.begin(), v.end(), [](TruthElement a, TruthElement b) { return a > b; });
255 for (auto& e : v) {
256 std::cerr << e << "\n";
257 }
258
259 const int R = 3;
260 const int L = 5;
261 // test with more elements
263 for (int run = 0; run < R; ++run) {
264 for (int i = 0; i < L; ++i) {
265 container2.addLabel(i, TruthElement(run));
266 }
267 }
268 // testing stage
269 for (int i = 0; i < L; ++i) {
270 auto labelview = container2.getLabels(i);
271 BOOST_CHECK(labelview.size() == R);
272 // count elements when iterating over view
273 int counter = 0;
274 std::cout << "CHECK CONTENT FOR INDEX " << i << "\n";
275 std::cout << "----- \n";
276 for (auto& l : labelview) {
277 counter++;
278 std::cout << l << "\n";
279 }
280 std::cout << "#### " << i << "\n";
281 std::cout << counter << "\n";
282 BOOST_CHECK(labelview.size() == counter);
283 }
284
285 // in this case we have to add the elements contiguously per dataindex:
287 cont2.addLabel(0, TruthElement(1));
288 cont2.addLabel(0, TruthElement(10));
289 cont2.addLabel(1, TruthElement(1));
290 cont2.addLabel(2, TruthElement(20));
291 {
292 auto view2 = cont2.getLabels(0);
293 BOOST_CHECK(view2.size() == 2);
294 for (auto& e : view2) {
295 std::cerr << e << "\n";
296 }
297 }
298 {
299 auto view2 = cont2.getLabels(1);
300 BOOST_CHECK(view2.size() == 1);
301 for (auto& e : view2) {
302 std::cerr << e << "\n";
303 }
304 }
305 {
306 auto view2 = cont2.getLabels(2);
307 BOOST_CHECK(view2.size() == 1);
308 for (auto& e : view2) {
309 std::cerr << e << "\n";
310 }
311 }
312
313 // get labels for nonexisting dataelement
314 BOOST_CHECK(cont2.getLabels(100).size() == 0);
315}
316
317BOOST_AUTO_TEST_CASE(MCTruthContainer_move)
318{
319 using TruthElement = long;
321 Container container;
322 container.addElement(0, TruthElement(1));
323 container.addElement(0, TruthElement(2));
324 container.addElement(1, TruthElement(1));
325 container.addElement(2, TruthElement(10));
326
327 Container container2 = std::move(container);
328 BOOST_CHECK(container.getIndexedSize() == 0);
329 BOOST_CHECK(container.getNElements() == 0);
330
331 std::swap(container, container2);
332 BOOST_CHECK(container2.getIndexedSize() == 0);
333 BOOST_CHECK(container2.getNElements() == 0);
334 BOOST_CHECK(container.getIndexedSize() == 3);
335 BOOST_CHECK(container.getNElements() == 4);
336}
337
338BOOST_AUTO_TEST_CASE(MCTruthContainer_ROOTIO)
339{
340 using TruthElement = o2::MCCompLabel;
342 Container container;
343 const size_t BIGSIZE{1000000};
344 container.addNoLabelIndex(0); // the first index does not have a label
345 for (int i = 1; i < BIGSIZE; ++i) {
346 container.addElement(i, TruthElement(i, i, i));
347 container.addElement(i, TruthElement(i + 1, i, i));
348 }
349 std::vector<char> buffer;
350 container.flatten_to(buffer);
351
352 // We use the special IO split container to stream to a file and back
354 {
355 TFile f("tmp2.root", "RECREATE");
356 TTree tree("o2sim", "o2sim");
357 auto br = tree.Branch("Labels", &io, 32000, 2);
358 tree.Branch("LabelsOriginal", &container, 32000, 2);
359 tree.Fill();
360 tree.Write();
361 f.Close();
362 }
363
364 // read back
365 TFile f2("tmp2.root", "OPEN");
366 auto tree2 = (TTree*)f2.Get("o2sim");
368 auto br2 = tree2->GetBranch("Labels");
369 BOOST_CHECK(br2 != nullptr);
370 br2->SetAddress(&io2);
371 br2->GetEntry(0);
372
373 // make a const MC label container out of it
374 using ConstMCTruthContainer = dataformats::ConstMCTruthContainer<TruthElement>;
375 ConstMCTruthContainer cc;
376 io2->copyandflatten(cc);
377
378 BOOST_CHECK(cc.getNElements() == (BIGSIZE - 1) * 2);
379 BOOST_CHECK(cc.getIndexedSize() == BIGSIZE);
380 BOOST_CHECK(cc.getLabels(0).size() == 0);
381 BOOST_CHECK(cc.getLabels(1).size() == 2);
382 BOOST_CHECK(cc.getLabels(1)[0] == TruthElement(1, 1, 1));
383 BOOST_CHECK(cc.getLabels(1)[1] == TruthElement(2, 1, 1));
384 BOOST_CHECK(cc.getLabels(BIGSIZE - 1).size() == 2);
385 BOOST_CHECK(cc.getLabels(BIGSIZE - 1)[0] == TruthElement(BIGSIZE - 1, BIGSIZE - 1, BIGSIZE - 1));
386 BOOST_CHECK(cc.getLabels(BIGSIZE - 1)[1] == TruthElement(BIGSIZE, BIGSIZE - 1, BIGSIZE - 1));
387
388 // testing convenience API to retrieve a constant label container from a ROOT file, entry 0
389 auto cont = o2::dataformats::MCLabelIOHelper::loadFromTTree(tree2, "Labels", 0);
390 auto cont2 = o2::dataformats::MCLabelIOHelper::loadFromTTree(tree2, "LabelsOriginal", 0);
391
392 BOOST_CHECK(cont);
393 BOOST_CHECK(cont2);
394 BOOST_CHECK(cont->getNElements() == (BIGSIZE - 1) * 2);
395 BOOST_CHECK(cont2->getNElements() == (BIGSIZE - 1) * 2);
396 BOOST_CHECK(cont->getLabels(0).size() == 0);
397 BOOST_CHECK(cont2->getLabels(0).size() == 0);
398 BOOST_CHECK(cont->getLabels(BIGSIZE - 1)[0] == TruthElement(BIGSIZE - 1, BIGSIZE - 1, BIGSIZE - 1));
399 BOOST_CHECK(cont->getLabels(BIGSIZE - 1)[1] == TruthElement(BIGSIZE, BIGSIZE - 1, BIGSIZE - 1));
400 BOOST_CHECK(cont2->getLabels(BIGSIZE - 1)[0] == TruthElement(BIGSIZE - 1, BIGSIZE - 1, BIGSIZE - 1));
401 BOOST_CHECK(cont2->getLabels(BIGSIZE - 1)[1] == TruthElement(BIGSIZE, BIGSIZE - 1, BIGSIZE - 1));
402}
403
404} // namespace o2
A const (ready only) version of MCTruthContainer.
int32_t i
A special IO container - splitting a given vector to enable ROOT IO.
A read-only version of MCTruthContainer allowing for storage optimisation.
void copyandflatten(std::vector< char, Alloc > &output) const
void addLabel(unsigned int dataindex, LabelType const &label)
add a label for a dataindex
void fillVectorOfLabels(int dataindex, std::vector< LabelType > &v)
LabelView getLabels(int dataindex)
get a container view on labels allowing use standard forward iteration in user code
static ConstMCTruthContainer< o2::MCCompLabel > * loadFromTTree(TTree *tree, std::string const &brname, int entry)
A container to hold and manage MC truth information/labels.
gsl::span< TruthElement > getLabels(uint32_t dataindex)
void addNoLabelIndex(uint32_t dataindex)
adds a data index that has no label
TruthElement const & getElement(uint32_t elementindex) const
void mergeAtBack(MCTruthContainer< TruthElement > const &other)
void addElements(uint32_t dataindex, gsl::span< CompatibleLabel > elements)
void addElement(uint32_t dataindex, TruthElement const &element, bool noElement=false)
MCTruthHeaderElement const & getMCTruthHeader(uint32_t dataindex) const
void addElementRandomAccess(uint32_t dataindex, TruthElement const &element)
GLuint buffer
Definition glcorearb.h:655
const GLdouble * v
Definition glcorearb.h:832
GLdouble f
Definition glcorearb.h:310
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLuint counter
Definition glcorearb.h:3987
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
BOOST_AUTO_TEST_CASE(FlatHisto)
std::vector< o2::mch::ChannelCode > cc
BOOST_CHECK(tree)
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))