Project
Loading...
Searching...
No Matches
TimeFrame.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.
15
16#include <numeric>
17#include <sstream>
18
19#include "Framework/Logger.h"
30
31namespace
32{
33struct ClusterHelper {
34 float phi;
35 float r;
36 int bin;
37 int ind;
38};
39} // namespace
40
41namespace o2::its
42{
43
48
49template <int nLayers>
50void TimeFrame<nLayers>::addPrimaryVertices(const bounded_vector<Vertex>& vertices, const int iteration)
51{
52 for (const auto& vertex : vertices) {
53 mPrimaryVertices.emplace_back(vertex); // put a copy in the present
54 mTotVertPerIteration[iteration]++;
55 if (!isBeamPositionOverridden) { // beam position is updated only at first occurrence of the vertex. A bit sketchy if we have past/future vertices, it should not impact too much.
56 const float w = vertex.getNContributors();
57 mBeamPos[0] = (mBeamPos[0] * mBeamPosWeight + vertex.getX() * w) / (mBeamPosWeight + w);
58 mBeamPos[1] = (mBeamPos[1] * mBeamPosWeight + vertex.getY() * w) / (mBeamPosWeight + w);
59 mBeamPosWeight += w;
60 }
61 }
62 mROFramesPV.push_back(mPrimaryVertices.size()); // current rof must have number of vertices up to present
63}
64
65template <int nLayers>
67{
68 mVerticesMCRecInfo.insert(mVerticesMCRecInfo.end(), labels.begin(), labels.end());
69}
70
71template <int nLayers>
73{
74 mVerticesContributorLabels.insert(mVerticesContributorLabels.end(), labels.begin(), labels.end());
75}
76
77template <int nLayers>
78void TimeFrame<nLayers>::addPrimaryVerticesInROF(const bounded_vector<Vertex>& vertices, const int rofId, const int iteration)
79{
80 mPrimaryVertices.insert(mPrimaryVertices.begin() + mROFramesPV[rofId], vertices.begin(), vertices.end());
81 for (int i = rofId + 1; i < mROFramesPV.size(); ++i) {
82 mROFramesPV[i] += vertices.size();
83 }
84 mTotVertPerIteration[iteration] += vertices.size();
87template <int nLayers>
88void TimeFrame<nLayers>::addPrimaryVerticesLabelsInROF(const bounded_vector<std::pair<MCCompLabel, float>>& labels, const int rofId)
89{
90 mVerticesMCRecInfo.insert(mVerticesMCRecInfo.begin() + mROFramesPV[rofId], labels.begin(), labels.end());
91}
93template <int nLayers>
95{
96 // count the number of cont. in rofs before and including the target rof
97 unsigned int n{0};
98 const auto& pvs = getPrimaryVertices(0, rofId);
99 for (const auto& pv : pvs) {
100 n += pv.getNContributors();
101 }
102 mVerticesContributorLabels.insert(mVerticesContributorLabels.begin() + n, labels.begin(), labels.end());
103}
104
105template <int nLayers>
106int TimeFrame<nLayers>::loadROFrameData(gsl::span<const o2::itsmft::ROFRecord> rofs,
107 gsl::span<const itsmft::CompClusterExt> clusters,
108 gsl::span<const unsigned char>::iterator& pattIt,
109 const itsmft::TopologyDictionary* dict,
111{
114
115 resetROFrameData(rofs.size());
116 prepareROFrameData(rofs, clusters);
117
118 for (size_t iRof{0}; iRof < rofs.size(); ++iRof) {
119 const auto& rof = rofs[iRof];
120 for (int clusterId{rof.getFirstEntry()}; clusterId < rof.getFirstEntry() + rof.getNEntries(); ++clusterId) {
121 const auto& c = clusters[clusterId];
122
123 int layer = geom->getLayer(c.getSensorID());
124
125 auto pattID = c.getPatternID();
127 float sigmaY2 = DefClusError2Row, sigmaZ2 = DefClusError2Col, sigmaYZ = 0; // Dummy COG errors (about half pixel size)
128 unsigned int clusterSize{0};
130 sigmaY2 = dict->getErr2X(pattID);
131 sigmaZ2 = dict->getErr2Z(pattID);
132 if (!dict->isGroup(pattID)) {
133 locXYZ = dict->getClusterCoordinates(c);
134 clusterSize = dict->getNpixels(pattID);
135 } else {
136 o2::itsmft::ClusterPattern patt(pattIt);
137 locXYZ = dict->getClusterCoordinates(c, patt);
138 clusterSize = patt.getNPixels();
139 }
140 } else {
141 o2::itsmft::ClusterPattern patt(pattIt);
142 locXYZ = dict->getClusterCoordinates(c, patt, false);
143 clusterSize = patt.getNPixels();
144 }
145 mClusterSize[clusterId] = std::clamp(clusterSize, 0u, 255u);
146 auto sensorID = c.getSensorID();
147 // Inverse transformation to the local --> tracking
148 auto trkXYZ = geom->getMatrixT2L(sensorID) ^ locXYZ;
149 // Transformation to the local --> global
150 auto gloXYZ = geom->getMatrixL2G(sensorID) * locXYZ;
151
152 addTrackingFrameInfoToLayer(layer, gloXYZ.x(), gloXYZ.y(), gloXYZ.z(), trkXYZ.x(), geom->getSensorRefAlpha(sensorID),
153 std::array<float, 2>{trkXYZ.y(), trkXYZ.z()},
154 std::array<float, 3>{sigmaY2, sigmaYZ, sigmaZ2});
156 addClusterToLayer(layer, gloXYZ.x(), gloXYZ.y(), gloXYZ.z(), mUnsortedClusters[layer].size());
157 addClusterExternalIndexToLayer(layer, clusterId);
158 }
159 for (unsigned int iL{0}; iL < mUnsortedClusters.size(); ++iL) {
160 mROFramesClusters[iL][iRof + 1] = mUnsortedClusters[iL].size(); // effectively calculating and exclusive sum
161 }
162 }
163
164 for (auto i = 0; i < mNTrackletsPerCluster.size(); ++i) {
165 mNTrackletsPerCluster[i].resize(mUnsortedClusters[1].size());
166 mNTrackletsPerClusterSum[i].resize(mUnsortedClusters[1].size() + 1); // Exc sum "prepends" a 0
167 }
168
169 if (mcLabels != nullptr) {
170 mClusterLabels = mcLabels;
171 }
172
173 return mNrof;
174}
175
176template <int nLayers>
178{
179 for (int iLayer{0}; iLayer < nLayers; ++iLayer) {
180 deepVectorClear(mUnsortedClusters[iLayer], getMaybeFrameworkHostResource());
181 deepVectorClear(mTrackingFrameInfo[iLayer], getMaybeFrameworkHostResource());
182 clearResizeBoundedVector(mROFramesClusters[iLayer], nRofs + 1, getMaybeFrameworkHostResource());
183 deepVectorClear(mClusterExternalIndices[iLayer], mMemoryPool.get());
184
185 if (iLayer < 2) {
186 deepVectorClear(mTrackletsIndexROF[iLayer], mMemoryPool.get());
187 deepVectorClear(mNTrackletsPerCluster[iLayer], mMemoryPool.get());
188 deepVectorClear(mNTrackletsPerClusterSum[iLayer], mMemoryPool.get());
189 }
190 }
191}
193template <int nLayers>
194void TimeFrame<nLayers>::prepareROFrameData(gsl::span<const o2::itsmft::ROFRecord> rofs,
195 gsl::span<const itsmft::CompClusterExt> clusters)
198 mNrof = rofs.size();
199 clearResizeBoundedVector(mClusterSize, clusters.size(), mMemoryPool.get());
200 std::array<int, nLayers> clusterCountPerLayer{};
201 for (const auto& clus : clusters) {
202 ++clusterCountPerLayer[geom->getLayer(clus.getSensorID())];
203 }
204 for (int iLayer{0}; iLayer < nLayers; ++iLayer) {
205 mUnsortedClusters[iLayer].reserve(clusterCountPerLayer[iLayer]);
206 mTrackingFrameInfo[iLayer].reserve(clusterCountPerLayer[iLayer]);
207 mClusterExternalIndices[iLayer].reserve(clusterCountPerLayer[iLayer]);
208 }
209}
210
211template <int nLayers>
212void TimeFrame<nLayers>::prepareClusters(const TrackingParameters& trkParam, const int maxLayers)
213{
214 const int numBins{trkParam.PhiBins * trkParam.ZBins};
215 const int stride{numBins + 1};
216 bounded_vector<ClusterHelper> cHelper(mMemoryPool.get());
217 bounded_vector<int> clsPerBin(numBins, 0, mMemoryPool.get());
218 bounded_vector<int> lutPerBin(numBins, 0, mMemoryPool.get());
219 for (int rof{0}; rof < mNrof; ++rof) {
220 if ((int)mMultiplicityCutMask.size() == mNrof && !mMultiplicityCutMask[rof]) {
221 continue;
222 }
223 for (int iLayer{0}, stopLayer = std::min(trkParam.NLayers, maxLayers); iLayer < stopLayer; ++iLayer) {
224 const auto& unsortedClusters{getUnsortedClustersOnLayer(rof, iLayer)};
225 const int clustersNum{static_cast<int>(unsortedClusters.size())};
226 auto* tableBase = mIndexTables[iLayer].data() + rof * stride;
227
228 cHelper.resize(clustersNum);
229
230 for (int iCluster{0}; iCluster < clustersNum; ++iCluster) {
231 const Cluster& c = unsortedClusters[iCluster];
232 ClusterHelper& h = cHelper[iCluster];
233
234 const float x = c.xCoordinate - mBeamPos[0];
235 const float y = c.yCoordinate - mBeamPos[1];
236 const float z = c.zCoordinate;
237
238 float phi = math_utils::computePhi(x, y);
239 int zBin{mIndexTableUtils.getZBinIndex(iLayer, z)};
240 if (zBin < 0 || zBin >= trkParam.ZBins) {
241 zBin = std::clamp(zBin, 0, trkParam.ZBins - 1);
242 mBogusClusters[iLayer]++;
243 }
244 int bin = mIndexTableUtils.getBinIndex(zBin, mIndexTableUtils.getPhiBinIndex(phi));
245 h.phi = phi;
246 h.r = math_utils::hypot(x, y);
247 mMinR[iLayer] = o2::gpu::GPUCommonMath::Min(h.r, mMinR[iLayer]);
248 mMaxR[iLayer] = o2::gpu::GPUCommonMath::Max(h.r, mMaxR[iLayer]);
249 h.bin = bin;
250 h.ind = clsPerBin[bin]++;
251 }
252 std::exclusive_scan(clsPerBin.begin(), clsPerBin.end(), lutPerBin.begin(), 0);
253
254 auto clusters2beSorted{getClustersOnLayer(rof, iLayer)};
255 for (int iCluster{0}; iCluster < clustersNum; ++iCluster) {
256 const ClusterHelper& h = cHelper[iCluster];
257 Cluster& c = clusters2beSorted[lutPerBin[h.bin] + h.ind];
259 c = unsortedClusters[iCluster];
260 c.phi = h.phi;
261 c.radius = h.r;
262 c.indexTableBinIndex = h.bin;
264 std::copy_n(lutPerBin.data(), clsPerBin.size(), tableBase);
265 std::fill_n(tableBase + clsPerBin.size(), stride - clsPerBin.size(), clustersNum);
266
267 std::fill(clsPerBin.begin(), clsPerBin.end(), 0);
268 cHelper.clear();
269 }
270 }
271}
272
273template <int nLayers>
274void TimeFrame<nLayers>::initialise(const int iteration, const TrackingParameters& trkParam, const int maxLayers, bool resetVertices)
275{
276 if (iteration == 0) {
277 if (maxLayers < trkParam.NLayers && resetVertices) {
278 resetRofPV();
279 deepVectorClear(mTotVertPerIteration);
280 }
281 deepVectorClear(mTracks);
282 deepVectorClear(mTracksLabel);
283 deepVectorClear(mLines);
284 deepVectorClear(mLinesLabels);
285 if (resetVertices) {
286 deepVectorClear(mVerticesMCRecInfo);
287 deepVectorClear(mVerticesContributorLabels);
288 }
289 clearResizeBoundedVector(mTracks, mNrof, mMemoryPool.get());
290 clearResizeBoundedVector(mTracksLabel, mNrof, mMemoryPool.get());
291 clearResizeBoundedVector(mLinesLabels, mNrof, mMemoryPool.get());
292 clearResizeBoundedVector(mCells, trkParam.CellsPerRoad(), mMemoryPool.get());
293 clearResizeBoundedVector(mCellsLookupTable, trkParam.CellsPerRoad() - 1, mMemoryPool.get());
294 clearResizeBoundedVector(mCellsNeighbours, trkParam.CellsPerRoad() - 1, mMemoryPool.get());
295 clearResizeBoundedVector(mCellsNeighboursLUT, trkParam.CellsPerRoad() - 1, mMemoryPool.get());
296 clearResizeBoundedVector(mCellLabels, trkParam.CellsPerRoad(), mMemoryPool.get());
297 clearResizeBoundedVector(mTracklets, std::min(trkParam.TrackletsPerRoad(), maxLayers - 1), mMemoryPool.get());
298 clearResizeBoundedVector(mTrackletLabels, trkParam.TrackletsPerRoad(), mMemoryPool.get());
299 clearResizeBoundedVector(mTrackletsLookupTable, trkParam.TrackletsPerRoad(), mMemoryPool.get());
300 mIndexTableUtils.setTrackingParameters(trkParam);
301 clearResizeBoundedVector(mPositionResolution, trkParam.NLayers, mMemoryPool.get());
302 clearResizeBoundedVector(mBogusClusters, trkParam.NLayers, mMemoryPool.get());
303 deepVectorClear(mTrackletClusters);
304 for (unsigned int iLayer{0}; iLayer < std::min((int)mClusters.size(), maxLayers); ++iLayer) {
305 clearResizeBoundedVector(mClusters[iLayer], mUnsortedClusters[iLayer].size(), getMaybeFrameworkHostResource(maxLayers != nLayers));
306 clearResizeBoundedVector(mUsedClusters[iLayer], mUnsortedClusters[iLayer].size(), getMaybeFrameworkHostResource(maxLayers != nLayers));
307 mPositionResolution[iLayer] = o2::gpu::CAMath::Sqrt(0.5f * (trkParam.SystErrorZ2[iLayer] + trkParam.SystErrorY2[iLayer]) + trkParam.LayerResolution[iLayer] * trkParam.LayerResolution[iLayer]);
308 }
309 clearResizeBoundedArray(mIndexTables, mNrof * (trkParam.ZBins * trkParam.PhiBins + 1), getMaybeFrameworkHostResource(maxLayers != nLayers));
310 clearResizeBoundedVector(mLines, mNrof, mMemoryPool.get());
311 clearResizeBoundedVector(mTrackletClusters, mNrof, mMemoryPool.get());
312
313 for (int iLayer{0}; iLayer < trkParam.NLayers; ++iLayer) {
314 if (trkParam.SystErrorY2[iLayer] > 0.f || trkParam.SystErrorZ2[iLayer] > 0.f) {
315 for (auto& tfInfo : mTrackingFrameInfo[iLayer]) {
317 tfInfo.covarianceTrackingFrame[0] += trkParam.SystErrorY2[iLayer];
318 tfInfo.covarianceTrackingFrame[2] += trkParam.SystErrorZ2[iLayer];
319 }
320 }
321 }
322 mMinR.fill(10000.);
323 mMaxR.fill(-1.);
324 }
325 mNTrackletsPerROF.resize(2);
326 for (auto& v : mNTrackletsPerROF) {
327 v = bounded_vector<int>(mNrof + 1, 0, mMemoryPool.get());
328 }
329 if (iteration == 0 || iteration == 3) {
330 prepareClusters(trkParam, maxLayers);
331 }
332 mTotalTracklets = {0, 0};
333 if (maxLayers < trkParam.NLayers) { // Vertexer only, but in both iterations
334 for (size_t iLayer{0}; iLayer < maxLayers; ++iLayer) {
335 deepVectorClear(mUsedClusters[iLayer]);
336 clearResizeBoundedVector(mUsedClusters[iLayer], mUnsortedClusters[iLayer].size(), mMemoryPool.get());
337 }
338 }
339
340 mTotVertPerIteration.resize(1 + iteration);
341 mNoVertexROF = 0;
342 deepVectorClear(mRoads);
343 deepVectorClear(mRoadLabels);
344
345 mMSangles.resize(trkParam.NLayers);
346 mPhiCuts.resize(mClusters.size() - 1, 0.f);
347 float oneOverR{0.001f * 0.3f * std::abs(mBz) / trkParam.TrackletMinPt};
348 for (unsigned int iLayer{0}; iLayer < nLayers; ++iLayer) {
349 mMSangles[iLayer] = math_utils::MSangle(0.14f, trkParam.TrackletMinPt, trkParam.LayerxX0[iLayer]);
350 mPositionResolution[iLayer] = o2::gpu::CAMath::Sqrt(0.5f * (trkParam.SystErrorZ2[iLayer] + trkParam.SystErrorY2[iLayer]) + trkParam.LayerResolution[iLayer] * trkParam.LayerResolution[iLayer]);
351 if (iLayer < mClusters.size() - 1) {
352 const float& r1 = trkParam.LayerRadii[iLayer];
353 const float& r2 = trkParam.LayerRadii[iLayer + 1];
354 oneOverR = (0.5 * oneOverR >= 1.f / r2) ? 2.f / r2 - o2::constants::math::Almost0 : oneOverR;
355 const float res1 = o2::gpu::CAMath::Hypot(trkParam.PVres, mPositionResolution[iLayer]);
356 const float res2 = o2::gpu::CAMath::Hypot(trkParam.PVres, mPositionResolution[iLayer + 1]);
357 const float cosTheta1half = o2::gpu::CAMath::Sqrt(1.f - math_utils::Sq(0.5f * r1 * oneOverR));
358 const float cosTheta2half = o2::gpu::CAMath::Sqrt(1.f - math_utils::Sq(0.5f * r2 * oneOverR));
359 float x = r2 * cosTheta1half - r1 * cosTheta2half;
360 float delta = o2::gpu::CAMath::Sqrt(1.f / (1.f - 0.25f * math_utils::Sq(x * oneOverR)) * (math_utils::Sq(0.25f * r1 * r2 * math_utils::Sq(oneOverR) / cosTheta2half + cosTheta1half) * math_utils::Sq(res1) + math_utils::Sq(0.25f * r1 * r2 * math_utils::Sq(oneOverR) / cosTheta1half + cosTheta2half) * math_utils::Sq(res2)));
362 mPhiCuts[iLayer] = std::min(o2::gpu::CAMath::ASin(0.5f * x * oneOverR) + 2.f * mMSangles[iLayer] + delta, o2::constants::math::PI * 0.5f);
363 }
364 }
365
366 for (int iLayer{0}; iLayer < std::min((int)mTracklets.size(), maxLayers); ++iLayer) {
367 deepVectorClear(mTracklets[iLayer]);
368 deepVectorClear(mTrackletLabels[iLayer]);
369 if (iLayer < (int)mCells.size()) {
370 deepVectorClear(mCells[iLayer]);
371 deepVectorClear(mTrackletsLookupTable[iLayer]);
372 mTrackletsLookupTable[iLayer].resize(mClusters[iLayer + 1].size() + 1, 0);
373 deepVectorClear(mCellLabels[iLayer]);
374 }
375
376 if (iLayer < (int)mCells.size() - 1) {
377 deepVectorClear(mCellsLookupTable[iLayer]);
378 deepVectorClear(mCellsNeighbours[iLayer]);
379 deepVectorClear(mCellsNeighboursLUT[iLayer]);
380 }
381 }
382}
383
384template <int nLayers>
386{
387 unsigned long size{0};
388 for (const auto& trkl : mTracklets) {
389 size += sizeof(Tracklet) * trkl.size();
390 }
391 for (const auto& cells : mCells) {
392 size += sizeof(CellSeedN) * cells.size();
393 }
394 for (const auto& cellsN : mCellsNeighbours) {
395 size += sizeof(int) * cellsN.size();
396 }
397 return size + sizeof(Road<nLayers - 2>) * mRoads.size();
398}
399
400template <int nLayers>
402{
403 LOGP(info, "TimeFrame: Artefacts occupy {:.2f} MB", getArtefactsMemory() / constants::MB);
404}
405
406template <int nLayers>
408{
409 deepVectorClear(mPValphaX);
410 mPValphaX.reserve(mPrimaryVertices.size());
411 for (auto& pv : mPrimaryVertices) {
412 mPValphaX.emplace_back(std::array<float, 2>{o2::gpu::CAMath::Hypot(pv.getX(), pv.getY()), math_utils::computePhi(pv.getX(), pv.getY())});
413 }
414}
415
416template <int nLayers>
418{
419 for (ushort iLayer = 0; iLayer < 2; ++iLayer) {
420 for (unsigned int iRof{0}; iRof < mNrof; ++iRof) {
421 if (mMultiplicityCutMask[iRof]) {
422 mTotalTracklets[iLayer] += mNTrackletsPerROF[iLayer][iRof];
423 }
424 }
425 std::exclusive_scan(mNTrackletsPerROF[iLayer].begin(), mNTrackletsPerROF[iLayer].end(), mNTrackletsPerROF[iLayer].begin(), 0);
426 std::exclusive_scan(mNTrackletsPerCluster[iLayer].begin(), mNTrackletsPerCluster[iLayer].end(), mNTrackletsPerClusterSum[iLayer].begin(), 0);
427 }
428}
429
430template <int nLayers>
432{
433 for (uint32_t iLayer{0}; iLayer < getTracklets().size(); ++iLayer) {
434 int prev{-1};
435 int count{0};
436 for (uint32_t iTracklet{0}; iTracklet < getTracklets()[iLayer].size(); ++iTracklet) {
437 auto& trk = getTracklets()[iLayer][iTracklet];
438 int currentId{trk.firstClusterIndex};
439 if (currentId < prev) {
440 LOG(info) << "First Cluster Index not increasing monotonically on L:T:ID:Prev " << iLayer << "\t" << iTracklet << "\t" << currentId << "\t" << prev;
441 } else if (currentId == prev) {
442 count++;
443 } else {
444 if (iLayer > 0) {
445 auto& lut{getTrackletsLookupTable()[iLayer - 1]};
446 if (count != lut[prev + 1] - lut[prev]) {
447 LOG(info) << "LUT count broken " << iLayer - 1 << "\t" << prev << "\t" << count << "\t" << lut[prev + 1] << "\t" << lut[prev];
448 }
449 }
450 count = 1;
451 }
452 prev = currentId;
453 if (iLayer > 0) {
454 auto& lut{getTrackletsLookupTable()[iLayer - 1]};
455 if (iTracklet >= (uint32_t)(lut[currentId + 1]) || iTracklet < (uint32_t)(lut[currentId])) {
456 LOG(info) << "LUT broken: " << iLayer - 1 << "\t" << currentId << "\t" << iTracklet;
457 }
458 }
459 }
460 }
461}
462
463template <int nLayers>
465{
466 LOG(info) << "-------- Tracklet LUT " << i;
467 std::stringstream s;
468 for (int j : mTrackletsLookupTable[i]) {
469 s << j << "\t";
470 }
471 LOG(info) << s.str();
472 LOG(info) << "--------";
473}
474
475template <int nLayers>
477{
478 LOG(info) << "-------- Cell LUT " << i;
479 std::stringstream s;
480 for (int j : mCellsLookupTable[i]) {
481 s << j << "\t";
482 }
483 LOG(info) << s.str();
484 LOG(info) << "--------";
485}
486
487template <int nLayers>
489{
490 for (unsigned int i{0}; i < mTrackletsLookupTable.size(); ++i) {
491 printTrackletLUTonLayer(i);
492 }
493}
494
495template <int nLayers>
497{
498 for (unsigned int i{0}; i < mCellsLookupTable.size(); ++i) {
499 printCellLUTonLayer(i);
500 }
501}
502
503template <int nLayers>
505{
506 LOG(info) << "Vertices in ROF (nROF = " << mNrof << ", lut size = " << mROFramesPV.size() << ")";
507 for (unsigned int iR{0}; iR < mROFramesPV.size(); ++iR) {
508 LOG(info) << mROFramesPV[iR] << "\t";
509 }
510 LOG(info) << "\n\n Vertices:";
511 for (unsigned int iV{0}; iV < mPrimaryVertices.size(); ++iV) {
512 LOG(info) << mPrimaryVertices[iV].getX() << "\t" << mPrimaryVertices[iV].getY() << "\t" << mPrimaryVertices[iV].getZ();
513 }
514 LOG(info) << "--------";
515}
516
517template <int nLayers>
519{
520 LOG(info) << "--------";
521 for (unsigned int iLayer{0}; iLayer < mROFramesClusters.size(); ++iLayer) {
522 LOG(info) << "Layer " << iLayer;
523 std::stringstream s;
524 for (auto value : mROFramesClusters[iLayer]) {
525 s << value << "\t";
526 }
527 LOG(info) << s.str();
528 }
529}
530
531template <int nLayers>
533{
534 LOG(info) << "--------";
535 for (unsigned int iLayer{0}; iLayer < mNClustersPerROF.size(); ++iLayer) {
536 LOG(info) << "Layer " << iLayer;
537 std::stringstream s;
538 for (auto& value : mNClustersPerROF[iLayer]) {
539 s << value << "\t";
540 }
541 LOG(info) << s.str();
542 }
543}
544
545template <int nLayers>
546void TimeFrame<nLayers>::printSliceInfo(const int startROF, const int sliceSize)
547{
548 LOG(info) << "Dumping slice of " << sliceSize << " rofs:";
549 for (int iROF{startROF}; iROF < startROF + sliceSize; ++iROF) {
550 LOG(info) << "ROF " << iROF << " dump:";
551 for (unsigned int iLayer{0}; iLayer < mClusters.size(); ++iLayer) {
552 LOG(info) << "Layer " << iLayer << " has: " << getClustersOnLayer(iROF, iLayer).size() << " clusters.";
553 }
554 LOG(info) << "Number of seeding vertices: " << getPrimaryVertices(iROF).size();
555 int iVertex{0};
556 for (auto& v : getPrimaryVertices(iROF)) {
557 LOG(info) << "\t vertex " << iVertex++ << ": x=" << v.getX() << " "
558 << " y=" << v.getY() << " z=" << v.getZ() << " has " << v.getNContributors() << " contributors.";
559 }
560 }
561}
562
563template <int nLayers>
564void TimeFrame<nLayers>::setMemoryPool(std::shared_ptr<BoundedMemoryResource> pool)
565{
566 mMemoryPool = pool;
567
568 auto initVector = [&]<typename T>(bounded_vector<T>& vec, bool useExternal = false) {
569 std::pmr::memory_resource* mr = (useExternal) ? mExtMemoryPool.get() : mMemoryPool.get();
570 deepVectorClear(vec, mr);
571 };
572
573 auto initContainers = [&]<typename Container>(Container& container, bool useExternal = false) {
574 for (auto& v : container) {
575 initVector(v, useExternal);
576 }
577 };
578
579 // these will only reside on the host for the cpu part
580 initVector(mTotVertPerIteration);
581 initContainers(mClusterExternalIndices);
582 initContainers(mNTrackletsPerCluster);
583 initContainers(mNTrackletsPerClusterSum);
584 initContainers(mNClustersPerROF);
585 initVector(mROFramesPV);
586 initVector(mPrimaryVertices);
587 initVector(mRoads);
588 initVector(mMSangles);
589 initVector(mPhiCuts);
590 initVector(mPositionResolution);
591 initVector(mClusterSize);
592 initVector(mPValphaX);
593 initVector(mBogusClusters);
594 initContainers(mTrackletsIndexROF);
595 initContainers(mTracks);
596 initContainers(mTracklets);
597 initContainers(mCells);
598 initContainers(mCellsNeighbours);
599 initContainers(mCellsLookupTable);
600 // MC info (we don't know if we have MC)
601 initVector(mVerticesContributorLabels);
602 initContainers(mLinesLabels);
603 initContainers(mTrackletLabels);
604 initContainers(mCellLabels);
605 initVector(mRoadLabels);
606 initContainers(mTracksLabel);
607 // these will use possibly an externally provided allocator
608 initContainers(mClusters, hasFrameworkAllocator());
609 initContainers(mUsedClusters, hasFrameworkAllocator());
610 initContainers(mUnsortedClusters, hasFrameworkAllocator());
611 initContainers(mIndexTables, hasFrameworkAllocator());
612 initContainers(mTrackingFrameInfo, hasFrameworkAllocator());
613 initContainers(mROFramesClusters, hasFrameworkAllocator());
614}
615
616template <int nLayers>
618{
619 mExternalAllocator = ext;
620 mExtMemoryPool = std::make_shared<BoundedMemoryResource>(mExternalAllocator);
621}
622
623template <int nLayers>
625{
626 deepVectorClear(mTracks);
627 deepVectorClear(mTracklets);
628 deepVectorClear(mCells);
629 deepVectorClear(mRoads);
630 deepVectorClear(mCellsNeighbours);
631 deepVectorClear(mCellsLookupTable);
632 deepVectorClear(mTotVertPerIteration);
633 deepVectorClear(mPrimaryVertices);
634 deepVectorClear(mTrackletsLookupTable);
635 deepVectorClear(mClusterExternalIndices);
636 deepVectorClear(mNTrackletsPerCluster);
637 deepVectorClear(mNTrackletsPerClusterSum);
638 deepVectorClear(mNClustersPerROF);
639 deepVectorClear(mROFramesPV);
640 deepVectorClear(mMSangles);
641 deepVectorClear(mPhiCuts);
642 deepVectorClear(mPositionResolution);
643 deepVectorClear(mClusterSize);
644 deepVectorClear(mPValphaX);
645 deepVectorClear(mBogusClusters);
646 deepVectorClear(mTrackletsIndexROF);
647 deepVectorClear(mTrackletClusters);
648 deepVectorClear(mLines);
649 // if we use the external host allocator then the assumption is that we
650 // don't clear the memory ourself
651 if (!hasFrameworkAllocator()) {
652 deepVectorClear(mClusters);
653 deepVectorClear(mUsedClusters);
654 deepVectorClear(mUnsortedClusters);
655 deepVectorClear(mIndexTables);
656 deepVectorClear(mTrackingFrameInfo);
657 deepVectorClear(mROFramesClusters);
658 }
659 // only needed to clear if we have MC info
660 if (hasMCinformation()) {
661 deepVectorClear(mLinesLabels);
662 deepVectorClear(mVerticesContributorLabels);
663 deepVectorClear(mTrackletLabels);
664 deepVectorClear(mCellLabels);
665 deepVectorClear(mRoadLabels);
666 deepVectorClear(mTracksLabel);
667 }
668}
669
670template class TimeFrame<7>;
671// ALICE3 upgrade
672#ifdef ENABLE_UPGRADES
673template class TimeFrame<11>;
674#endif
675
676} // namespace o2::its
std::vector< std::string > labels
Definition of the ITSMFT compact cluster.
Definition of the ITSMFT cluster.
Definition of the ClusterTopology class.
uint64_t vertex
Definition RawEventData.h:9
int32_t i
Definition of the GeometryTGeo class.
Definition of the ITSMFT ROFrame (trigger) record.
uint32_t j
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
Definition of the SegmentationAlpide class.
int clusterSize
Class for time synchronization of RawReader instances.
A container to hold and manage MC truth information/labels.
const Mat3D & getMatrixL2G(int sensID) const
const Mat3D & getMatrixT2L(int lay, int hba, int sta, int det) const
float getSensorRefAlpha(int isn) const
static GeometryTGeo * Instance()
int getLayer(int index) const final
Get chip layer, from 0.
void fillMatrixCache(int mask) override
int getNPixels() const
Returns the number of fired pixels.
static constexpr unsigned short InvalidPatternID
Definition CompCluster.h:46
static constexpr float PitchCol
static constexpr float PitchRow
math_utils::Point3D< T > getClusterCoordinates(const CompCluster &cl) const
float getErr2X(int n) const
Returns the error^2 on the x position of the COG for the n_th element.
int getNpixels(int n) const
Returns the number of fired pixels of the n_th element.
bool isGroup(int n) const
Returns true if the element corresponds to a group of rare topologies.
float getErr2Z(int n) const
Returns the error^2 on the z position of the COG for the n_th element.
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLint GLsizei count
Definition glcorearb.h:399
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
const GLdouble * v
Definition glcorearb.h:832
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLint GLenum GLboolean GLsizei stride
Definition glcorearb.h:867
GLenum GLuint GLint GLint layer
Definition glcorearb.h:1310
GLboolean r
Definition glcorearb.h:1233
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
constexpr float Almost0
constexpr float PI
constexpr float MB
Definition Constants.h:29
void deepVectorClear(std::vector< T > &vec)
std::pmr::vector< T > bounded_vector
constexpr float DefClusError2Col
Definition TimeFrame.cxx:47
constexpr float DefClusError2Row
Definition TimeFrame.cxx:46
constexpr float DefClusErrorRow
Definition TimeFrame.cxx:44
constexpr float DefClusErrorCol
Definition TimeFrame.cxx:45
void clearResizeBoundedArray(std::array< bounded_vector< T >, S > &arr, size_t size, std::pmr::memory_resource *mr=nullptr, T def=T())
void clearResizeBoundedVector(bounded_vector< T > &vec, size_t sz, std::pmr::memory_resource *mr=nullptr, T def=T())
constexpr int nLayers
Definition Specs.h:45
void setFrameworkAllocator(ExternalAllocator *ext)
void checkTrackletLUTs()
Debug and printing.
void printTrackletLUTonLayer(int i)
void initialise(const int iteration, const TrackingParameters &trkParam, const int maxLayers=7, bool resetVertices=true)
void addPrimaryVertices(const bounded_vector< Vertex > &vertices, const int iteration)
Definition TimeFrame.cxx:50
virtual void wipe()
void addPrimaryVerticesInROF(const bounded_vector< Vertex > &vertices, const int rofId, const int iteration)
Definition TimeFrame.cxx:78
void setMemoryPool(std::shared_ptr< BoundedMemoryResource > pool)
memory management
void printSliceInfo(const int, const int)
void prepareROFrameData(gsl::span< const o2::itsmft::ROFRecord > rofs, gsl::span< const itsmft::CompClusterExt > clusters)
void computeTrackletsPerROFScans()
void addPrimaryVerticesLabelsInROF(const bounded_vector< std::pair< MCCompLabel, float > > &labels, const int rofId)
Definition TimeFrame.cxx:88
void addPrimaryVerticesContributorLabelsInROF(const bounded_vector< MCCompLabel > &labels, const int rofId)
Definition TimeFrame.cxx:94
void addPrimaryVerticesLabels(bounded_vector< std::pair< MCCompLabel, float > > &labels)
Definition TimeFrame.cxx:66
int loadROFrameData(const o2::itsmft::ROFRecord &rof, gsl::span< const itsmft::Cluster > clusters, const dataformats::MCTruthContainer< MCCompLabel > *mcLabels=nullptr)
void fillPrimaryVerticesXandAlpha()
unsigned long getArtefactsMemory() const
void printCellLUTonLayer(int i)
void printArtefactsMemory() const
void addPrimaryVerticesContributorLabels(bounded_vector< MCCompLabel > &labels)
Definition TimeFrame.cxx:72
void resetROFrameData(size_t nROFs)
void prepareClusters(const TrackingParameters &trkParam, const int maxLayers=nLayers)
std::vector< float > LayerRadii
std::vector< float > SystErrorY2
std::vector< float > SystErrorZ2
float TrackletMinPt
Trackleting cuts.
std::vector< float > LayerResolution
std::vector< float > LayerxX0
int CellsPerRoad() const noexcept
int TrackletsPerRoad() const noexcept
static constexpr int L2G
Definition Cartesian.h:54
static constexpr int T2L
Definition Cartesian.h:55
std::vector< o2::ctf::BufferType > vec
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< Cluster > clusters
std::vector< Cell > cells