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
18#include "Framework/Logger.h"
27
28namespace
29{
30struct ClusterHelper {
31 float phi;
32 float r;
33 int bin;
34 int ind;
35};
36} // namespace
37
38namespace o2::its
39{
40
45
46template <int NLayers>
48{
49 mPrimaryVertices.emplace_back(vert);
50 if (!isBeamPositionOverridden) {
51 const float w = vert.getNContributors();
52 mBeamPos[0] = (mBeamPos[0] * mBeamPosWeight + vert.getX() * w) / (mBeamPosWeight + w);
53 mBeamPos[1] = (mBeamPos[1] * mBeamPosWeight + vert.getY() * w) / (mBeamPosWeight + w);
54 mBeamPosWeight += w;
55 }
56}
57
58template <int NLayers>
59void TimeFrame<NLayers>::loadROFrameData(gsl::span<const o2::itsmft::ROFRecord> rofs,
60 gsl::span<const itsmft::CompClusterExt> clusters,
61 gsl::span<const unsigned char>::iterator& pattIt,
63 int layer,
65{
68 resetROFrameData(layer);
69 prepareROFrameData(clusters, layer);
70
71 // check for missing/empty/unset rofs
72 // the code requires consistent monotonically increasing input without gaps
73 const auto& timing = mROFOverlapTableView.getLayer(layer >= 0 ? layer : 0);
74 if (timing.mNROFsTF != rofs.size()) {
75 LOGP(fatal, "Received inconsistent number of rofs on layer:{} expected:{} received:{}", layer, timing.mNROFsTF, rofs.size());
76 }
77
78 for (int32_t iRof{0}; iRof < rofs.size(); ++iRof) {
79 const auto& rof = rofs[iRof];
80 for (int clusterId{rof.getFirstEntry()}; clusterId < rof.getFirstEntry() + rof.getNEntries(); ++clusterId) {
81 const auto& c = clusters[clusterId];
82 int lay = geom->getLayer(c.getSensorID());
83 auto pattID = c.getPatternID();
85 float sigmaY2 = DefClusError2Row, sigmaZ2 = DefClusError2Col, sigmaYZ = 0; // Dummy COG errors (about half pixel size)
86 unsigned int clusterSize{0};
88 sigmaY2 = dict->getErr2X(pattID);
89 sigmaZ2 = dict->getErr2Z(pattID);
90 if (!dict->isGroup(pattID)) {
91 locXYZ = dict->getClusterCoordinates(c);
92 clusterSize = dict->getNpixels(pattID);
93 } else {
95 locXYZ = dict->getClusterCoordinates(c, patt);
96 clusterSize = patt.getNPixels();
97 }
98 } else {
99 o2::itsmft::ClusterPattern patt(pattIt);
100 locXYZ = dict->getClusterCoordinates(c, patt, false);
101 clusterSize = patt.getNPixels();
102 }
103 mClusterSize[layer >= 0 ? layer : 0][clusterId] = std::clamp(clusterSize, 0u, 255u);
104 auto sensorID = c.getSensorID();
105 // Inverse transformation to the local --> tracking
106 auto trkXYZ = geom->getMatrixT2L(sensorID) ^ locXYZ;
107 // Transformation to the local --> global
108 auto gloXYZ = geom->getMatrixL2G(sensorID) * locXYZ;
109 addTrackingFrameInfoToLayer(lay, gloXYZ.x(), gloXYZ.y(), gloXYZ.z(), trkXYZ.x(), geom->getSensorRefAlpha(sensorID),
110 std::array<float, 2>{trkXYZ.y(), trkXYZ.z()},
111 std::array<float, 3>{sigmaY2, sigmaYZ, sigmaZ2});
113 addClusterToLayer(lay, gloXYZ.x(), gloXYZ.y(), gloXYZ.z(), mUnsortedClusters[lay].size());
114 addClusterExternalIndexToLayer(lay, clusterId);
115 }
116 // effectively calculating an exclusive sum
117 if (layer >= 0) {
118 mROFramesClusters[layer][iRof + 1] = mUnsortedClusters[layer].size();
119 } else {
120 for (unsigned int iL{0}; iL < mUnsortedClusters.size(); ++iL) {
121 mROFramesClusters[iL][iRof + 1] = mUnsortedClusters[iL].size();
122 }
123 }
124 }
125
126 if (layer == 1 || layer == -1) {
127 for (auto i = 0; i < mNTrackletsPerCluster.size(); ++i) {
128 mNTrackletsPerCluster[i].resize(mUnsortedClusters[1].size());
129 mNTrackletsPerClusterSum[i].resize(mUnsortedClusters[1].size() + 1);
130 }
131 }
132
133 if (mcLabels != nullptr) {
134 mClusterLabels[layer >= 0 ? layer : 0] = mcLabels;
135 } else {
136 mClusterLabels[layer >= 0 ? layer : 0] = nullptr;
137 }
138}
139
140template <int NLayers>
142{
143 if (layer >= 0) {
144 deepVectorClear(mUnsortedClusters[layer], getMaybeFrameworkHostResource());
145 deepVectorClear(mTrackingFrameInfo[layer], getMaybeFrameworkHostResource());
146 deepVectorClear(mClusterExternalIndices[layer], mMemoryPool.get());
147 clearResizeBoundedVector(mROFramesClusters[layer], mROFOverlapTableView.getLayer(layer).mNROFsTF + 1, getMaybeFrameworkHostResource());
148 } else {
149 for (int iLayer{0}; iLayer < NLayers; ++iLayer) {
150 deepVectorClear(mUnsortedClusters[iLayer], getMaybeFrameworkHostResource());
151 deepVectorClear(mTrackingFrameInfo[iLayer], getMaybeFrameworkHostResource());
152 deepVectorClear(mClusterExternalIndices[iLayer], mMemoryPool.get());
153 clearResizeBoundedVector(mROFramesClusters[iLayer], mROFOverlapTableView.getLayer(iLayer).mNROFsTF + 1, getMaybeFrameworkHostResource());
154 }
155 }
156}
157
158template <int NLayers>
159void TimeFrame<NLayers>::prepareROFrameData(gsl::span<const itsmft::CompClusterExt> clusters, int layer)
160{
161 if (layer >= 0) {
162 mUnsortedClusters[layer].reserve(clusters.size());
163 mTrackingFrameInfo[layer].reserve(clusters.size());
164 mClusterExternalIndices[layer].reserve(clusters.size());
165 clearResizeBoundedVector(mClusterSize[layer], clusters.size(), mMemoryPool.get());
166 } else {
167 auto* geom = GeometryTGeo::Instance();
168 clearResizeBoundedVector(mClusterSize[0], clusters.size(), mMemoryPool.get());
169 std::array<size_t, NLayers> clusterCountPerLayer{0};
170 for (const auto& cls : clusters) {
171 ++clusterCountPerLayer[geom->getLayer(cls.getChipID())];
172 }
173 for (int iLayer{0}; iLayer < NLayers; ++iLayer) {
174 mUnsortedClusters[iLayer].reserve(clusterCountPerLayer[iLayer]);
175 mTrackingFrameInfo[iLayer].reserve(clusterCountPerLayer[iLayer]);
176 mClusterExternalIndices[iLayer].reserve(clusterCountPerLayer[iLayer]);
177 }
178 }
179}
180
181template <int NLayers>
182void TimeFrame<NLayers>::prepareClusters(const TrackingParameters& trkParam, const int maxLayers)
183{
184 const int numBins{trkParam.PhiBins * trkParam.ZBins};
185 const int stride{numBins + 1};
186 bounded_vector<ClusterHelper> cHelper(mMemoryPool.get());
187 bounded_vector<int> clsPerBin(numBins, 0, mMemoryPool.get());
188 bounded_vector<int> lutPerBin(numBins, 0, mMemoryPool.get());
189 for (int iLayer{0}, stopLayer = std::min(trkParam.NLayers, maxLayers); iLayer < stopLayer; ++iLayer) {
190 for (int rof{0}; rof < getNrof(iLayer); ++rof) {
191 if (!mROFMaskView.isROFEnabled(iLayer, rof)) {
192 continue;
193 }
194 const auto& unsortedClusters{getUnsortedClustersOnLayer(rof, iLayer)};
195 const int clustersNum{static_cast<int>(unsortedClusters.size())};
196 auto* tableBase = mIndexTables[iLayer].data() + rof * stride;
197
198 cHelper.resize(clustersNum);
199
200 for (int iCluster{0}; iCluster < clustersNum; ++iCluster) {
201 const Cluster& c = unsortedClusters[iCluster];
202 ClusterHelper& h = cHelper[iCluster];
203
204 const float x = c.xCoordinate - mBeamPos[0];
205 const float y = c.yCoordinate - mBeamPos[1];
206 const float z = c.zCoordinate;
207
208 float phi = math_utils::computePhi(x, y);
209 int zBin{mIndexTableUtils.getZBinIndex(iLayer, z)};
210 if (zBin < 0 || zBin >= trkParam.ZBins) {
211 zBin = std::clamp(zBin, 0, trkParam.ZBins - 1);
212 mBogusClusters[iLayer]++;
213 }
214 int bin = mIndexTableUtils.getBinIndex(zBin, mIndexTableUtils.getPhiBinIndex(phi));
215 h.phi = phi;
216 h.r = math_utils::hypot(x, y);
217 mMinR[iLayer] = o2::gpu::GPUCommonMath::Min(h.r, mMinR[iLayer]);
218 mMaxR[iLayer] = o2::gpu::GPUCommonMath::Max(h.r, mMaxR[iLayer]);
219 h.bin = bin;
220 h.ind = clsPerBin[bin]++;
222 std::exclusive_scan(clsPerBin.begin(), clsPerBin.end(), lutPerBin.begin(), 0);
223
224 auto clusters2beSorted{getClustersOnLayer(rof, iLayer)};
225 for (int iCluster{0}; iCluster < clustersNum; ++iCluster) {
226 const ClusterHelper& h = cHelper[iCluster];
227 Cluster& c = clusters2beSorted[lutPerBin[h.bin] + h.ind];
228
229 c = unsortedClusters[iCluster];
230 c.phi = h.phi;
231 c.radius = h.r;
232 c.indexTableBinIndex = h.bin;
233 }
234 std::copy_n(lutPerBin.data(), clsPerBin.size(), tableBase);
235 std::fill_n(tableBase + clsPerBin.size(), stride - clsPerBin.size(), clustersNum);
236
237 std::fill(clsPerBin.begin(), clsPerBin.end(), 0);
238 cHelper.clear();
239 }
240 }
241}
242
243template <int NLayers>
245{
246 mVertexingTopology.init(3, trkParam.MaxHoles, trkParam.HoleLayerMask);
247}
248
249template <int NLayers>
251{
252 mDefaultTrackingTopology.init(maxLayers, trkParam.MaxHoles, trkParam.HoleLayerMask);
253}
254
255template <int NLayers>
256void TimeFrame<NLayers>::initTrackerTopologies(gsl::span<const TrackingParameters> trkParams, const int maxLayers)
257{
258 mTrackerTopologies.resize(trkParams.size());
259 for (size_t iteration = 0; iteration < trkParams.size(); ++iteration) {
260 const int iterationMaxLayers = std::min(maxLayers, trkParams[iteration].NLayers);
261 mTrackerTopologies[iteration].init(iterationMaxLayers, trkParams[iteration].MaxHoles, trkParams[iteration].HoleLayerMask);
262 }
263}
264
265template <int NLayers>
266void TimeFrame<NLayers>::initialise(const TrackingParameters& trkParam, const int maxLayers, const int iteration)
267{
268 mTrackingTopologyView = iteration != constants::UnusedIndex ? mTrackerTopologies[iteration].getView() : (maxLayers == 3 ? mVertexingTopology.getView() : mDefaultTrackingTopology.getView());
269
270 if (trkParam.PassFlags[IterationStep::FirstPass]) {
271 deepVectorClear(mTracks);
272 deepVectorClear(mTracksLabel);
273 deepVectorClear(mLines);
274 deepVectorClear(mLinesLabels);
276 deepVectorClear(mPrimaryVertices);
277 deepVectorClear(mPrimaryVerticesLabels);
278 }
279 clearResizeBoundedVector(mLinesLabels, getNrof(1), mMemoryPool.get());
280 mIndexTableUtils.setTrackingParameters(trkParam);
281 clearResizeBoundedVector(mPositionResolution, trkParam.NLayers, mMemoryPool.get());
282 clearResizeBoundedVector(mBogusClusters, trkParam.NLayers, mMemoryPool.get());
283 deepVectorClear(mTrackletClusters);
284 for (unsigned int iLayer{0}; iLayer < std::min((int)mClusters.size(), maxLayers); ++iLayer) {
285 clearResizeBoundedVector(mClusters[iLayer], mUnsortedClusters[iLayer].size(), getMaybeFrameworkHostResource(maxLayers != NLayers));
286 clearResizeBoundedVector(mUsedClusters[iLayer], mUnsortedClusters[iLayer].size(), getMaybeFrameworkHostResource(maxLayers != NLayers));
287 mPositionResolution[iLayer] = o2::gpu::CAMath::Sqrt((0.5f * (trkParam.SystErrorZ2[iLayer] + trkParam.SystErrorY2[iLayer])) + (trkParam.LayerResolution[iLayer] * trkParam.LayerResolution[iLayer]));
288 }
289 clearResizeBoundedVector(mLines, getNrof(1), mMemoryPool.get());
290 clearResizeBoundedVector(mTrackletClusters, getNrof(1), mMemoryPool.get());
291
292 for (int iLayer{0}; iLayer < NLayers; ++iLayer) {
293 clearResizeBoundedVector(mIndexTables[iLayer], getNrof(iLayer) * ((trkParam.ZBins * trkParam.PhiBins) + 1), getMaybeFrameworkHostResource());
294 }
295 for (int iLayer{0}; iLayer < trkParam.NLayers; ++iLayer) {
296 if (trkParam.SystErrorY2[iLayer] > 0.f || trkParam.SystErrorZ2[iLayer] > 0.f) {
297 for (auto& tfInfo : mTrackingFrameInfo[iLayer]) {
299 tfInfo.covarianceTrackingFrame[0] += trkParam.SystErrorY2[iLayer];
300 tfInfo.covarianceTrackingFrame[2] += trkParam.SystErrorZ2[iLayer];
301 }
302 }
303 }
304
305 mMinR.fill(std::numeric_limits<float>::max());
306 mMaxR.fill(std::numeric_limits<float>::min());
307 }
308 clearResizeBoundedVector(mCells, mTrackingTopologyView.nCells, mMemoryPool.get());
309 clearResizeBoundedVector(mCellsLookupTable, mTrackingTopologyView.nCells, mMemoryPool.get());
310 clearResizeBoundedVector(mCellsNeighbours, mTrackingTopologyView.nCells, mMemoryPool.get());
311 clearResizeBoundedVector(mCellsNeighboursTopology, mTrackingTopologyView.nCells, mMemoryPool.get());
312 clearResizeBoundedVector(mCellsNeighboursLUT, mTrackingTopologyView.nCells, mMemoryPool.get());
313 clearResizeBoundedVector(mCellLabels, mTrackingTopologyView.nCells, mMemoryPool.get());
314 clearResizeBoundedVector(mTracklets, mTrackingTopologyView.nTransitions, mMemoryPool.get());
315 clearResizeBoundedVector(mTrackletLabels, mTrackingTopologyView.nTransitions, mMemoryPool.get());
316 clearResizeBoundedVector(mTrackletsLookupTable, mTrackingTopologyView.nTransitions, mMemoryPool.get());
317 clearResizeBoundedVector(mTransitionPhiCuts, mTrackingTopologyView.nTransitions, mMemoryPool.get());
318 clearResizeBoundedVector(mTransitionMSAngles, mTrackingTopologyView.nTransitions, mMemoryPool.get());
319 mNTrackletsPerROF.resize(2);
320 for (auto& v : mNTrackletsPerROF) {
321 v = bounded_vector<int>(getNrof(1) + 1, 0, mMemoryPool.get());
322 }
324 prepareClusters(trkParam, maxLayers);
325 }
326 mTotalTracklets = {0, 0};
327 if (maxLayers < trkParam.NLayers) { // Vertexer only, but in both iterations
328 for (size_t iLayer{0}; iLayer < maxLayers; ++iLayer) {
329 deepVectorClear(mUsedClusters[iLayer]);
330 clearResizeBoundedVector(mUsedClusters[iLayer], mUnsortedClusters[iLayer].size(), mMemoryPool.get());
331 }
332 }
333
334 // estimate MS per layer
335 std::array<float, NLayers> msAngles{};
336 for (unsigned int iLayer{0}; iLayer < NLayers; ++iLayer) {
337 msAngles[iLayer] = math_utils::MSangle(0.14f, trkParam.TrackletMinPt, trkParam.LayerxX0[iLayer]);
338 mPositionResolution[iLayer] = o2::gpu::CAMath::Sqrt((0.5f * (trkParam.SystErrorZ2[iLayer] + trkParam.SystErrorY2[iLayer])) + (trkParam.LayerResolution[iLayer] * trkParam.LayerResolution[iLayer]));
339 }
340
341 // for each transition calculate the phi-cuts + integrated MS
342 float oneOverR{0.001f * 0.3f * std::abs(mBz) / trkParam.TrackletMinPt};
343 for (int transitionId{0}; transitionId < (int)mTracklets.size(); ++transitionId) {
344 const auto& transition = mTrackingTopologyView.getTransition(transitionId);
345 float ms2 = 0.;
346 for (int layer = transition.fromLayer; layer < transition.toLayer; ++layer) {
347 ms2 += math_utils::Sq(msAngles[layer]);
348 }
349 mTransitionMSAngles[transitionId] = o2::gpu::CAMath::Sqrt(ms2);
350 const float& r1 = trkParam.LayerRadii[transition.fromLayer];
351 const float& r2 = trkParam.LayerRadii[transition.toLayer];
352 oneOverR = (0.5 * oneOverR >= 1.f / r2) ? (2.f / r2) - o2::constants::math::Almost0 : oneOverR;
353 const float res1 = o2::gpu::CAMath::Hypot(trkParam.PVres, mPositionResolution[transition.fromLayer]);
354 const float res2 = o2::gpu::CAMath::Hypot(trkParam.PVres, mPositionResolution[transition.toLayer]);
355 const float cosTheta1half = o2::gpu::CAMath::Sqrt(1.f - math_utils::Sq(0.5f * r1 * oneOverR));
356 const float cosTheta2half = o2::gpu::CAMath::Sqrt(1.f - math_utils::Sq(0.5f * r2 * oneOverR));
357 float x = (r2 * cosTheta1half) - (r1 * cosTheta2half);
358 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)));
360 mTransitionPhiCuts[transitionId] = o2::gpu::CAMath::Min(o2::gpu::CAMath::ASin(0.5f * x * oneOverR) + 2.f * mTransitionMSAngles[transitionId] + delta, o2::constants::math::PI * 0.5f);
361
362 // some cleanup
363 deepVectorClear(mTracklets[transitionId]);
364 deepVectorClear(mTrackletLabels[transitionId]);
365 deepVectorClear(mTrackletsLookupTable[transitionId]);
366 mTrackletsLookupTable[transitionId].resize(mClusters[transition.fromLayer].size() + 1, 0);
367 }
368
369 for (int cellId{0}; cellId < (int)mCells.size(); ++cellId) {
370 deepVectorClear(mCells[cellId]);
371 deepVectorClear(mCellsLookupTable[cellId]);
372 deepVectorClear(mCellsNeighbours[cellId]);
373 deepVectorClear(mCellsNeighboursTopology[cellId]);
374 deepVectorClear(mCellsNeighboursLUT[cellId]);
375 deepVectorClear(mCellLabels[cellId]);
376 }
377}
378
379template <int NLayers>
381{
382 unsigned long size{0};
383 for (const auto& trkl : mTracklets) {
384 size += sizeof(Tracklet) * trkl.size();
385 }
386 for (const auto& cells : mCells) {
387 size += sizeof(CellSeed) * cells.size();
388 }
389 for (const auto& cellsN : mCellsNeighbours) {
390 size += sizeof(int) * cellsN.size();
391 }
392 for (const auto& cellsN : mCellsNeighboursTopology) {
393 size += sizeof(int) * cellsN.size();
394 }
395 return size;
396}
397
398template <int NLayers>
400{
401 LOGP(info, "TimeFrame: Artefacts occupy {:.2f} MB", getArtefactsMemory() / constants::MB);
402}
403
404template <int NLayers>
406{
407 for (ushort iLayer = 0; iLayer < 2; ++iLayer) {
408 for (unsigned int iRof{0}; iRof < getNrof(1); ++iRof) {
409 if (mROFMaskView.isROFEnabled(1, iRof)) {
410 mTotalTracklets[iLayer] += mNTrackletsPerROF[iLayer][iRof];
411 }
412 }
413 std::exclusive_scan(mNTrackletsPerROF[iLayer].begin(), mNTrackletsPerROF[iLayer].end(), mNTrackletsPerROF[iLayer].begin(), 0);
414 std::exclusive_scan(mNTrackletsPerCluster[iLayer].begin(), mNTrackletsPerCluster[iLayer].end(), mNTrackletsPerClusterSum[iLayer].begin(), 0);
415 }
416}
417
418template <int NLayers>
419void TimeFrame<NLayers>::setMemoryPool(std::shared_ptr<BoundedMemoryResource> pool)
420{
421 mMemoryPool = pool;
422
423 auto initVector = [&]<typename T>(bounded_vector<T>& vec, bool useExternal = false) {
424 std::pmr::memory_resource* mr = (useExternal) ? mExtMemoryPool.get() : mMemoryPool.get();
425 deepVectorClear(vec, mr);
426 };
427
428 auto initContainers = [&]<typename Container>(Container& container, bool useExternal = false) {
429 for (auto& v : container) {
430 initVector(v, useExternal);
431 }
432 };
433
434 // these will only reside on the host for the cpu part
435 initContainers(mClusterExternalIndices);
436 initContainers(mNTrackletsPerCluster);
437 initContainers(mNTrackletsPerClusterSum);
438 initContainers(mNClustersPerROF);
439 initVector(mPrimaryVertices);
440 initVector(mTransitionPhiCuts);
441 initVector(mTransitionMSAngles);
442 initVector(mPositionResolution);
443 initContainers(mClusterSize);
444 initVector(mPValphaX);
445 initVector(mBogusClusters);
446 initContainers(mTrackletsIndexROF);
447 initVector(mTracks);
448 initContainers(mTracklets);
449 initContainers(mCells);
450 initContainers(mCellsNeighbours);
451 initContainers(mCellsLookupTable);
452 // MC info (we don't know if we have MC)
453 initVector(mPrimaryVerticesLabels);
454 initContainers(mLinesLabels);
455 initContainers(mTrackletLabels);
456 initContainers(mCellLabels);
457 initVector(mTracksLabel);
458 // these will use possibly an externally provided allocator
459 initContainers(mClusters, hasFrameworkAllocator());
460 initContainers(mUsedClusters, hasFrameworkAllocator());
461 initContainers(mUnsortedClusters, hasFrameworkAllocator());
462 initContainers(mIndexTables, hasFrameworkAllocator());
463 initContainers(mTrackingFrameInfo, hasFrameworkAllocator());
464 initContainers(mROFramesClusters, hasFrameworkAllocator());
465}
466
467template <int NLayers>
469{
470 mExternalAllocator = ext;
471 mExtMemoryPool = std::make_shared<BoundedMemoryResource>(mExternalAllocator);
472}
473
474template <int NLayers>
476{
477 deepVectorClear(mTracks);
478 deepVectorClear(mTracklets);
479 deepVectorClear(mCells);
480 deepVectorClear(mCellsNeighbours);
481 deepVectorClear(mCellsNeighboursTopology);
482 deepVectorClear(mCellsLookupTable);
483 deepVectorClear(mPrimaryVertices);
484 deepVectorClear(mTrackletsLookupTable);
485 deepVectorClear(mClusterExternalIndices);
486 deepVectorClear(mNTrackletsPerCluster);
487 deepVectorClear(mNTrackletsPerClusterSum);
488 deepVectorClear(mNClustersPerROF);
489 deepVectorClear(mTransitionPhiCuts);
490 deepVectorClear(mTransitionMSAngles);
491 deepVectorClear(mPositionResolution);
492 deepVectorClear(mClusterSize);
493 deepVectorClear(mPValphaX);
494 deepVectorClear(mBogusClusters);
495 deepVectorClear(mTrackletsIndexROF);
496 deepVectorClear(mTrackletClusters);
497 deepVectorClear(mLines);
498 // if we use the external host allocator then the assumption is that we
499 // don't clear the memory ourself
500 if (!hasFrameworkAllocator()) {
501 deepVectorClear(mClusters);
502 deepVectorClear(mUsedClusters);
503 deepVectorClear(mUnsortedClusters);
504 deepVectorClear(mIndexTables);
505 deepVectorClear(mTrackingFrameInfo);
506 deepVectorClear(mROFramesClusters);
507 }
508 // only needed to clear if we have MC info
509 if (hasMCinformation()) {
510 deepVectorClear(mLinesLabels);
511 deepVectorClear(mPrimaryVerticesLabels);
512 deepVectorClear(mTrackletLabels);
513 deepVectorClear(mCellLabels);
514 deepVectorClear(mTracksLabel);
515 }
516}
517
518template class TimeFrame<7>;
519// ALICE3 upgrade
520#ifdef ENABLE_UPGRADES
521template class TimeFrame<11>;
522#endif
523
524} // namespace o2::its
Definition of the ITSMFT compact cluster.
Definition of the ClusterTopology class.
int32_t i
Definition of the GeometryTGeo class.
Definition of the ITSMFT ROFrame (trigger) record.
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
HMPID cluster implementation.
Definition Cluster.h:27
CellSeed: connections of three clusters.
Definition Cell.h:81
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
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.
GLint GLenum GLint x
Definition glcorearb.h:403
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
const GLdouble * v
Definition glcorearb.h:832
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 int UnusedIndex
Definition Constants.h:32
constexpr float MB
Definition Constants.h:26
const TrackingFrameInfo *const const Cluster *const * unsortedClusters
void deepVectorClear(std::vector< T > &vec)
std::pmr::vector< T > bounded_vector
constexpr float DefClusError2Col
Definition TimeFrame.cxx:44
constexpr float DefClusError2Row
Definition TimeFrame.cxx:43
constexpr float DefClusErrorRow
Definition TimeFrame.cxx:41
constexpr float DefClusErrorCol
Definition TimeFrame.cxx:42
void clearResizeBoundedVector(bounded_vector< T > &vec, size_t sz, std::pmr::memory_resource *mr=nullptr, T def=T())
void prepareClusters(const TrackingParameters &trkParam, const int maxLayers=NLayers)
void initDefaultTrackingTopology(const TrackingParameters &trkParam, const int maxLayers=NLayers)
void setMemoryPool(std::shared_ptr< BoundedMemoryResource > pool)
memory management
void initVertexingTopology(const TrackingParameters &trkParam)
void printArtefactsMemory() const
void setFrameworkAllocator(ExternalAllocator *ext)
void initialise(const TrackingParameters &trkParam, const int maxLayers=NLayers, const int iteration=constants::UnusedIndex)
void prepareROFrameData(gsl::span< const itsmft::CompClusterExt > clusters, int layer)
void computeTrackletsPerROFScans()
void addPrimaryVertex(const Vertex &vertex)
Definition TimeFrame.cxx:47
void loadROFrameData(gsl::span< const o2::itsmft::ROFRecord > rofs, gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, const itsmft::TopologyDictionary *dict, int layer, const dataformats::MCTruthContainer< MCCompLabel > *mcLabels=nullptr)
Definition TimeFrame.cxx:59
void resetROFrameData(int iLayer)
void initTrackerTopologies(gsl::span< const TrackingParameters > trkParams, const int maxLayers=NLayers)
virtual void wipe()
unsigned long getArtefactsMemory() const
std::vector< float > LayerRadii
std::vector< float > SystErrorY2
std::vector< float > SystErrorZ2
float TrackletMinPt
Trackleting cuts.
std::vector< float > LayerResolution
std::vector< float > LayerxX0
static constexpr int L2G
Definition Cartesian.h:54
static constexpr int T2L
Definition Cartesian.h:55
std::vector< o2::ctf::BufferType > vec
std::vector< Cluster > clusters
std::vector< Cell > cells