Project
Loading...
Searching...
No Matches
TrackHelpers.h
Go to the documentation of this file.
1// Copyright 2019-2026 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#ifndef O2_ITS_TRACKING_TRACKHELPERS_H_
17#define O2_ITS_TRACKING_TRACKHELPERS_H_
18
19#include <cmath>
20
22#include "ITStracking/Cell.h"
23#include "ITStracking/Cluster.h"
28
30{
31
32GPUdi() int selectReseedMidLayer(int minLayer, int maxLayer, int nLayers, const float* layerRadii)
33{
34 if (maxLayer - minLayer == nLayers - 1) {
35 return (minLayer + maxLayer) / 2;
36 }
37 int midLayer = minLayer + 1;
38 const float midR = 0.5f * (layerRadii[maxLayer] + layerRadii[minLayer]);
39 float distanceToMidR = o2::gpu::CAMath::Abs(midR - layerRadii[midLayer]);
40 for (int iLayer = midLayer + 1; iLayer < maxLayer; ++iLayer) { // find the midpoint as closest to the midR
41 const float distance = o2::gpu::CAMath::Abs(midR - layerRadii[iLayer]);
43 midLayer = iLayer;
45 }
46 }
47 return midLayer;
48}
49
51{
52 track.resetCovariance();
53 track.setCov(track.getQ2Pt() * track.getQ2Pt() * track.getCov()[o2::track::CovLabels::kSigQ2Pt2], o2::track::CovLabels::kSigQ2Pt2);
54}
55
56GPUdi() o2::track::TrackParCov buildTrackSeed(const Cluster& cluster1,
57 const Cluster& cluster2,
58 const TrackingFrameInfo& tf3,
59 const float bz,
60 const bool reverse = false)
61{
62 float ca = NAN, sa = NAN, snp = NAN, q2pt = NAN, q2pt2 = NAN;
63 o2::gpu::CAMath::SinCos(tf3.alphaTrackingFrame, sa, ca);
64 const float sign = reverse ? -1.f : 1.f;
65 const float x1 = (cluster1.xCoordinate * ca) + (cluster1.yCoordinate * sa);
66 const float y1 = (-cluster1.xCoordinate * sa) + (cluster1.yCoordinate * ca);
67 const float x2 = (cluster2.xCoordinate * ca) + (cluster2.yCoordinate * sa);
68 const float y2 = (-cluster2.xCoordinate * sa) + (cluster2.yCoordinate * ca);
69 const float x3 = tf3.xTrackingFrame;
70 const float y3 = tf3.positionTrackingFrame[0];
71 if (o2::gpu::CAMath::Abs(bz) < 0.01f) { // zero field
72 const float dx = x3 - x1;
73 const float dy = y3 - y1;
74 snp = sign * dy / o2::gpu::CAMath::Hypot(dx, dy);
75 q2pt = 1.f / o2::track::kMostProbablePt;
76 q2pt2 = 1.f;
77 } else {
78 const float crv = math_utils::computeCurvature(x3, y3, x2, y2, x1, y1);
79 snp = sign * crv * (x3 - math_utils::computeCurvatureCentreX(x3, y3, x2, y2, x1, y1));
80 q2pt = sign * crv / (bz * o2::constants::math::B2C);
81 q2pt2 = crv * crv;
82 }
83 const float tgl = 0.5f * (math_utils::computeTanDipAngle(x1, y1, x2, y2, cluster1.zCoordinate, cluster2.zCoordinate) +
84 math_utils::computeTanDipAngle(x2, y2, x3, y3, cluster2.zCoordinate, tf3.positionTrackingFrame[1]));
85 const float sg2q2pt = o2::track::kC1Pt2max * o2::gpu::CAMath::Clamp(q2pt2, 0.0005f, 1.0f);
86 return {x3, tf3.alphaTrackingFrame, {y3, tf3.positionTrackingFrame[1], snp, tgl, q2pt}, {tf3.covarianceTrackingFrame[0], tf3.covarianceTrackingFrame[1], tf3.covarianceTrackingFrame[2], 0.f, 0.f, o2::track::kCSnp2max, 0.f, 0.f, 0.f, o2::track::kCTgl2max, 0.f, 0.f, 0.f, 0.f, sg2q2pt}};
87}
88
89template <int NLayers>
90GPUdi() TrackITSExt seedTrackForRefit(const TrackSeed<NLayers>& seed,
93 const float* layerRadii,
94 const float bz,
95 const int reseedIfShorter)
96{
98 int lrMin = NLayers;
99 int lrMax = 0;
100 for (int iL{0}; iL < NLayers; ++iL) {
101 const int idx = seed.getCluster(iL);
102 temporaryTrack.setExternalClusterIndex(iL, idx, idx != constants::UnusedIndex);
103 if (idx != constants::UnusedIndex) {
104 lrMin = o2::gpu::CAMath::Min(lrMin, iL);
105 lrMax = o2::gpu::CAMath::Max(lrMax, iL);
106 }
107 }
108
109 const int ncl = temporaryTrack.getNClusters();
110 if (ncl < reseedIfShorter && ncl > 1) {
111 const int lrMid = selectReseedMidLayer(lrMin, lrMax, NLayers, layerRadii);
112 const auto& cluster0TF = foundTrackingFrameInfo[lrMin][seed.getCluster(lrMin)];
113 const auto& cluster1GL = unsortedClusters[lrMid][seed.getCluster(lrMid)];
114 const auto& cluster2GL = unsortedClusters[lrMax][seed.getCluster(lrMax)];
115 temporaryTrack.getParamIn() = buildTrackSeed(cluster2GL, cluster1GL, cluster0TF, bz, true);
116 }
117
120}
121
122GPUdi() bool fitTrack(TrackITSExt& trk,
123 int start,
124 int end,
125 int step,
129 int nCl,
130 const float bz,
132 const float* layerxX0,
133 const o2::base::Propagator* propagator,
134 const o2::base::PropagatorF::MatCorrType matCorrType,
135 o2::track::TrackPar* linRef = nullptr,
136 const bool shiftRefToCluster = false)
137{
138 for (int iLayer{start}; iLayer != end; iLayer += step) {
139 if (trk.getClusterIndex(iLayer) == constants::UnusedIndex) {
140 continue;
141 }
142
143 const TrackingFrameInfo& trackingHit = tfInfos[iLayer][trk.getClusterIndex(iLayer)];
144 if (linRef) {
145 if (!trk.o2::track::TrackParCovF::rotate(trackingHit.alphaTrackingFrame, *linRef, bz)) {
146 return false;
147 }
148 if (!propagator->propagateToX(trk, *linRef, trackingHit.xTrackingFrame, bz,
151 matCorrType)) {
152 return false;
153 }
154 if (matCorrType == o2::base::PropagatorF::MatCorrType::USEMatCorrNONE) {
155 if (!trk.correctForMaterial(*linRef, layerxX0[iLayer], layerxX0[iLayer] * constants::Radl * constants::Rho, true)) {
156 continue;
157 }
158 }
159 } else {
160 if (!trk.o2::track::TrackParCovF::rotate(trackingHit.alphaTrackingFrame)) {
161 return false;
162 }
163 if (!propagator->propagateToX(trk, trackingHit.xTrackingFrame, bz,
166 matCorrType)) {
167 return false;
168 }
169 if (matCorrType == o2::base::PropagatorF::MatCorrType::USEMatCorrNONE) {
170 if (!trk.correctForMaterial(layerxX0[iLayer], layerxX0[iLayer] * constants::Radl * constants::Rho, true)) {
171 continue;
172 }
173 }
174 }
175
176 const auto predChi2{trk.getPredictedChi2Quiet(trackingHit.positionTrackingFrame, trackingHit.covarianceTrackingFrame)};
177 if ((nCl >= 3 && predChi2 > chi2clcut) || predChi2 < 0.f) {
178 return false;
179 }
180 trk.setChi2(trk.getChi2() + predChi2);
181 if (!trk.o2::track::TrackParCov::update(trackingHit.positionTrackingFrame, trackingHit.covarianceTrackingFrame)) {
182 return false;
183 }
184 if (linRef && shiftRefToCluster) {
185 linRef->setY(trackingHit.positionTrackingFrame[0]);
186 linRef->setZ(trackingHit.positionTrackingFrame[1]);
187 }
188 nCl++;
189 }
190
191 return o2::gpu::CAMath::Abs(trk.getQ2Pt()) < maxQoverPt && trk.getChi2() < chi2ndfcut * (float)((nCl * 2) - 5);
192}
193
194template <int NLayers>
195GPUdi() bool refitTrack(const TrackSeed<NLayers>& trackSeed,
197 float chi2clcut,
198 float chi2ndfcut,
199 const float bz,
200 const TrackingFrameInfo* const* tfInfos,
201 const Cluster* const* clusters,
202 const float* layerxX0,
203 const float* layerRadii,
204 const float* minPt,
205 const o2::base::Propagator* propagator,
206 const o2::base::PropagatorF::MatCorrType matCorrType,
207 const int reseedIfShorter,
208 const bool shiftRefToCluster,
209 const bool repeatRefitOut)
210{
211 temporaryTrack = seedTrackForRefit(trackSeed,
212 tfInfos,
213 clusters,
215 bz,
218 bool fitSuccess = fitTrack(temporaryTrack,
219 0,
220 NLayers,
221 1,
222 chi2clcut,
225 0,
226 bz,
227 tfInfos,
228 layerxX0,
231 &linRef,
232 shiftRefToCluster);
233 if (!fitSuccess) {
234 return false;
235 }
236 temporaryTrack.getParamOut() = temporaryTrack.getParamIn();
237 linRef = temporaryTrack.getParamOut(); // use refitted track as lin.reference
239 temporaryTrack.setChi2(0);
240 fitSuccess = fitTrack(temporaryTrack,
241 NLayers - 1,
242 -1,
243 -1,
244 chi2clcut,
246 50.f,
247 0,
248 bz,
249 tfInfos,
250 layerxX0,
253 &linRef,
254 shiftRefToCluster);
255 if (!fitSuccess || temporaryTrack.getPt() < minPt[NLayers - temporaryTrack.getNClusters()]) {
256 return false;
257 }
258 if (repeatRefitOut) { // repeat outward refit seeding and linearizing with the stable inward fit result
260 linRef = saveInw; // use refitted track as lin.reference
261 float saveChi2 = temporaryTrack.getChi2();
263 temporaryTrack.setChi2(0);
264 fitSuccess = o2::its::track::fitTrack(temporaryTrack,
265 0,
266 NLayers,
267 1,
268 chi2clcut,
271 0,
272 bz,
273 tfInfos,
274 layerxX0,
277 &linRef,
278 shiftRefToCluster);
279 if (!fitSuccess) {
280 return false;
281 }
282 temporaryTrack.getParamOut() = temporaryTrack.getParamIn();
283 temporaryTrack.getParamIn() = saveInw;
284 temporaryTrack.setChi2(saveChi2);
285 }
286 return true;
287}
288
289} // namespace o2::its::track
290
291#endif
Base track model for the Barrel, params only, w/o covariance.
o2::track::TrackParCov TrackParCov
Definition Recon.h:39
std::array< int, 64 > reverse(std::array< int, 64 > a)
Definition of the ITS track.
GLuint GLfloat GLfloat GLfloat GLfloat y1
Definition glcorearb.h:5034
GLuint GLuint end
Definition glcorearb.h:469
GLuint GLfloat GLfloat GLfloat x1
Definition glcorearb.h:5034
GLsizei GLsizei GLfloat distance
Definition glcorearb.h:5506
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint start
Definition glcorearb.h:469
constexpr float B2C
constexpr float VeryBig
constexpr int UnusedIndex
Definition Constants.h:35
constexpr float Radl
Definition Constants.h:36
constexpr float Rho
Definition Constants.h:37
int int int float chi2clcut
const TrackingFrameInfo *const * foundTrackingFrameInfo
int int int float float float int const float const TrackingFrameInfo *const const float const o2::base::Propagator const o2::base::PropagatorF::MatCorrType matCorrType
int int int float float float int nCl
int int int float float float maxQoverPt
const TrackingFrameInfo *const const Cluster *const const float const float bz
resetTrackCovariance(temporaryTrack)
int int int step
const TrackingFrameInfo *const const Cluster *const * unsortedClusters
int int const float * layerRadii
GPUdi() int selectReseedMidLayer(int minLayer
return temporaryTrack
int int int float float float int const float const TrackingFrameInfo *const const float * layerxX0
int int int float float float int const float const TrackingFrameInfo *const const float const o2::base::Propagator const o2::base::PropagatorF::MatCorrType o2::track::TrackPar * linRef
int int int float float float int const float const TrackingFrameInfo *const * tfInfos
int int int float float chi2ndfcut
const TrackingFrameInfo *const const Cluster *const const float const float const int reseedIfShorter
int int int float float float int const float const TrackingFrameInfo *const const float const o2::base::Propagator * propagator
const float midR
constexpr float kCTgl2max
TrackParCovF TrackParCov
Definition Track.h:33
constexpr float kCSnp2max
constexpr float kMostProbablePt
constexpr float kC1Pt2max
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::array< float, 2 > positionTrackingFrame
Definition Cluster.h:75
std::array< float, 3 > covarianceTrackingFrame
Definition Cluster.h:76
std::vector< Cluster > clusters