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
32// Find the populated interior layer closest to the radial midpoint.
33// If no layer can be found, return constants::UnusedIndex.
34// Should minimize the sagitta bias.
35template <int NLayers>
36GPUdi() int selectReseedMidLayer(int minLayer, int maxLayer, const float* layerRadii, const TrackSeed<NLayers>& seed)
37{
39 float distanceToMidR = layerRadii[NLayers - 1]; // midpoint cannot be last layer
40 const float midR = 0.5f * (layerRadii[maxLayer] + layerRadii[minLayer]);
41 for (int iLayer = minLayer + 1; iLayer < maxLayer; ++iLayer) {
42 if (seed.getCluster(iLayer) != constants::UnusedIndex) {
43 const float distance = o2::gpu::CAMath::Abs(midR - layerRadii[iLayer]);
44 if (distance < distanceToMidR) { // keep the smaller-radius layer on ties
45 midLayer = iLayer;
47 }
48 }
49 }
50 return midLayer;
51}
52
54{
55 track.resetCovariance();
56 track.setCov(track.getQ2Pt() * track.getQ2Pt() * track.getCov()[o2::track::CovLabels::kSigQ2Pt2], o2::track::CovLabels::kSigQ2Pt2);
57}
58
59GPUdi() o2::track::TrackParCov buildTrackSeed(const Cluster& cluster1,
60 const Cluster& cluster2,
61 const TrackingFrameInfo& tf3,
62 const float bz,
63 const bool reverse = false)
64{
66 o2::gpu::CAMath::SinCos(tf3.alphaTrackingFrame, sa, ca);
67 const float sign = reverse ? -1.f : 1.f;
68 const float x1 = (cluster1.xCoordinate * ca) + (cluster1.yCoordinate * sa);
69 const float y1 = (-cluster1.xCoordinate * sa) + (cluster1.yCoordinate * ca);
70 const float x2 = (cluster2.xCoordinate * ca) + (cluster2.yCoordinate * sa);
71 const float y2 = (-cluster2.xCoordinate * sa) + (cluster2.yCoordinate * ca);
72 const float x3 = tf3.xTrackingFrame;
73 const float y3 = tf3.positionTrackingFrame[0];
74 if (o2::gpu::CAMath::Abs(bz) < 0.01f) { // zero field
75 const float dx = x3 - x1;
76 const float dy = y3 - y1;
77 snp = sign * dy / o2::gpu::CAMath::Hypot(dx, dy);
78 q2pt = 1.f / o2::track::kMostProbablePt;
79 q2pt2 = 1.f;
80 } else {
81 const float crv = math_utils::computeCurvature(x3, y3, x2, y2, x1, y1);
82 snp = sign * crv * (x3 - math_utils::computeCurvatureCentreX(x3, y3, x2, y2, x1, y1));
83 q2pt = sign * crv / (bz * o2::constants::math::B2C);
84 q2pt2 = crv * crv;
85 }
86 const float tgl = -0.5f * sign * (math_utils::computeTanDipAngle(x1, y1, x2, y2, cluster1.zCoordinate, cluster2.zCoordinate) + math_utils::computeTanDipAngle(x2, y2, x3, y3, cluster2.zCoordinate, tf3.positionTrackingFrame[1]));
87 const float sg2q2pt = o2::track::kC1Pt2max * o2::gpu::CAMath::Clamp(q2pt2, 0.0005f, 1.0f);
88 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}};
89}
90
91template <int NLayers>
92GPUdi() TrackITSExt seedTrackForRefit(const TrackSeed<NLayers>& seed,
95 const float* layerRadii,
96 const float bz,
97 const int reseedIfShorter)
98{
100 int lrMin = NLayers;
101 int lrMax = 0;
102 for (int iL{0}; iL < NLayers; ++iL) {
103 const int idx = seed.getCluster(iL);
104 temporaryTrack.setExternalClusterIndex(iL, idx, idx != constants::UnusedIndex);
105 if (idx != constants::UnusedIndex) {
106 lrMin = o2::gpu::CAMath::Min(lrMin, iL);
107 lrMax = o2::gpu::CAMath::Max(lrMax, iL);
108 }
109 }
110
111 const int ncl = temporaryTrack.getNClusters();
112 if (ncl < reseedIfShorter && ncl > 2) {
113 const int lrMid = selectReseedMidLayer<NLayers>(lrMin, lrMax, layerRadii, seed);
114 if (lrMid != constants::UnusedIndex) {
115 const auto& cluster0TF = foundTrackingFrameInfo[lrMin][seed.getCluster(lrMin)];
116 const auto& cluster1GL = unsortedClusters[lrMid][seed.getCluster(lrMid)];
117 const auto& cluster2GL = unsortedClusters[lrMax][seed.getCluster(lrMax)];
118 temporaryTrack.getParamIn() = buildTrackSeed(cluster2GL, cluster1GL, cluster0TF, bz, true);
119 }
120 }
121
124}
125
126GPUdi() bool fitTrack(TrackITSExt& trk,
127 int start,
128 int end,
129 int step,
133 int nCl,
134 const float bz,
136 const float* layerxX0,
137 const o2::base::Propagator* propagator,
138 const o2::base::PropagatorF::MatCorrType matCorrType,
139 o2::track::TrackPar* linRef = nullptr,
140 const bool shiftRefToCluster = false)
141{
142 for (int iLayer{start}; iLayer != end; iLayer += step) {
143 if (trk.getClusterIndex(iLayer) == constants::UnusedIndex) {
144 continue;
145 }
146
147 const TrackingFrameInfo& trackingHit = tfInfos[iLayer][trk.getClusterIndex(iLayer)];
148 if (linRef) {
149 if (!trk.o2::track::TrackParCovF::rotate(trackingHit.alphaTrackingFrame, *linRef, bz)) {
150 return false;
151 }
152 if (!propagator->propagateToX(trk, *linRef, trackingHit.xTrackingFrame, bz,
155 matCorrType)) {
156 return false;
157 }
158 if (matCorrType == o2::base::PropagatorF::MatCorrType::USEMatCorrNONE) {
159 if (!trk.correctForMaterial(*linRef, layerxX0[iLayer], layerxX0[iLayer] * constants::Radl * constants::Rho, true)) {
160 continue;
161 }
162 }
163 } else {
164 if (!trk.o2::track::TrackParCovF::rotate(trackingHit.alphaTrackingFrame)) {
165 return false;
166 }
167 if (!propagator->propagateToX(trk, trackingHit.xTrackingFrame, bz,
170 matCorrType)) {
171 return false;
172 }
173 if (matCorrType == o2::base::PropagatorF::MatCorrType::USEMatCorrNONE) {
174 if (!trk.correctForMaterial(layerxX0[iLayer], layerxX0[iLayer] * constants::Radl * constants::Rho, true)) {
175 continue;
176 }
177 }
178 }
179
180 const auto predChi2{trk.getPredictedChi2Quiet(trackingHit.positionTrackingFrame, trackingHit.covarianceTrackingFrame)};
181 if ((nCl >= 3 && predChi2 > chi2clcut) || predChi2 < 0.f) {
182 return false;
183 }
184 trk.setChi2(trk.getChi2() + predChi2);
185 if (!trk.o2::track::TrackParCov::update(trackingHit.positionTrackingFrame, trackingHit.covarianceTrackingFrame)) {
186 return false;
187 }
188 if (linRef && shiftRefToCluster) {
189 linRef->setY(trackingHit.positionTrackingFrame[0]);
190 linRef->setZ(trackingHit.positionTrackingFrame[1]);
191 }
192 nCl++;
193 }
194
195 return o2::gpu::CAMath::Abs(trk.getQ2Pt()) < maxQoverPt && trk.getChi2() < chi2ndfcut * (float)((nCl * 2) - 5);
196}
197
198template <int NLayers>
199GPUdi() bool refitTrack(const TrackSeed<NLayers>& trackSeed,
201 float chi2clcut,
202 float chi2ndfcut,
203 const float bz,
204 const TrackingFrameInfo* const* tfInfos,
205 const Cluster* const* clusters,
206 const float* layerxX0,
207 const float* layerRadii,
208 const float* minPt,
209 const o2::base::Propagator* propagator,
210 const o2::base::PropagatorF::MatCorrType matCorrType,
211 const int reseedIfShorter,
212 const bool shiftRefToCluster,
213 const bool repeatRefitOut)
214{
215 temporaryTrack = seedTrackForRefit(trackSeed,
216 tfInfos,
217 clusters,
219 bz,
222 bool fitSuccess = fitTrack(temporaryTrack,
223 0,
224 NLayers,
225 1,
226 chi2clcut,
229 0,
230 bz,
231 tfInfos,
232 layerxX0,
235 &linRef,
236 shiftRefToCluster);
237 if (!fitSuccess) {
238 return false;
239 }
240 temporaryTrack.getParamOut() = temporaryTrack.getParamIn();
241 linRef = temporaryTrack.getParamOut(); // use refitted track as lin.reference
243 temporaryTrack.setChi2(0);
244 fitSuccess = fitTrack(temporaryTrack,
245 NLayers - 1,
246 -1,
247 -1,
248 chi2clcut,
250 50.f,
251 0,
252 bz,
253 tfInfos,
254 layerxX0,
257 &linRef,
258 shiftRefToCluster);
259 if (!fitSuccess || temporaryTrack.getPt() < minPt[NLayers - temporaryTrack.getNClusters()]) {
260 return false;
261 }
262 if (repeatRefitOut) { // repeat outward refit seeding and linearizing with the stable inward fit result
264 linRef = saveInw; // use refitted track as lin.reference
265 float saveChi2 = temporaryTrack.getChi2();
267 temporaryTrack.setChi2(0);
268 fitSuccess = o2::its::track::fitTrack(temporaryTrack,
269 0,
270 NLayers,
271 1,
272 chi2clcut,
275 0,
276 bz,
277 tfInfos,
278 layerxX0,
281 &linRef,
282 shiftRefToCluster);
283 if (!fitSuccess) {
284 return false;
285 }
286 temporaryTrack.getParamOut() = temporaryTrack.getParamIn();
287 temporaryTrack.getParamIn() = saveInw;
288 temporaryTrack.setChi2(saveChi2);
289 }
290 return true;
291}
292
293} // namespace o2::its::track
294
295#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:32
constexpr float Radl
Definition Constants.h:34
constexpr float Rho
Definition Constants.h:35
constexpr float UnsetValue
Definition Constants.h:33
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 const float * layerRadii
return temporaryTrack
int int int float float float int const float const TrackingFrameInfo *const const float * layerxX0
int const float const TrackSeed< NLayers > & seed
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
GPUdi() int selectReseedMidLayer(int minLayer
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