Project
Loading...
Searching...
No Matches
TrackParametrizationWithError.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
16#include <GPUCommonLogger.h>
17
18#ifndef GPUCA_GPUCODE_DEVICE
19#include <iostream>
20#endif
21
22#ifndef GPUCA_ALIGPUCODE
23#include <fmt/printf.h>
24#endif
25
26using namespace o2::track;
27using namespace o2::gpu;
28
29//______________________________________________________________
30template <typename value_T>
32{
33 // Transform this track to the local coord. system rotated by 180 deg.
34 this->invertParam();
35 // since the fP1 and fP2 are not inverted, their covariances with others change sign
36 mC[kSigZY] = -mC[kSigZY];
37 mC[kSigSnpY] = -mC[kSigSnpY];
38 mC[kSigTglZ] = -mC[kSigTglZ];
39 mC[kSigTglSnp] = -mC[kSigTglSnp];
40 mC[kSigQ2PtZ] = -mC[kSigQ2PtZ];
41 mC[kSigQ2PtSnp] = -mC[kSigQ2PtSnp];
42}
43
44//______________________________________________________________
45template <typename value_T>
46GPUd() bool TrackParametrizationWithError<value_T>::propagateTo(value_t xk, value_t b)
47{
48 //----------------------------------------------------------------
49 // propagate this track to the plane X=xk (cm) in the field "b" (kG)
50 //----------------------------------------------------------------
51 value_t dx = xk - this->getX();
52 if (gpu::CAMath::Abs(dx) < constants::math::Almost0) {
53 return true;
54 }
55 value_t crv = this->getCurvature(b);
56 value_t x2r = crv * dx;
57 value_t f1 = this->getSnp(), f2 = f1 + x2r;
58 if ((gpu::CAMath::Abs(f1) > constants::math::Almost1) || (gpu::CAMath::Abs(f2) > constants::math::Almost1)) {
59 return false;
60 }
61 value_t r1 = gpu::CAMath::Sqrt((1.f - f1) * (1.f + f1));
62 if (gpu::CAMath::Abs(r1) < constants::math::Almost0) {
63 return false;
64 }
65 value_t r2 = gpu::CAMath::Sqrt((1.f - f2) * (1.f + f2));
66 if (gpu::CAMath::Abs(r2) < constants::math::Almost0) {
67 return false;
68 }
69 double dy2dx = (f1 + f2) / (r1 + r2);
70 bool arcz = gpu::CAMath::Abs(x2r) > 0.05f;
71 params_t dP{0.f};
72 if (arcz) {
73 // for small dx/R the linear apporximation of the arc by the segment is OK,
74 // but at large dx/R the error is very large and leads to incorrect Z propagation
75 // angle traversed delta = 2*asin(dist_start_end / R / 2), hence the arc is: R*deltaPhi
76 // The dist_start_end is obtained from sqrt(dx^2+dy^2) = x/(r1+r2)*sqrt(2+f1*f2+r1*r2)
77 // double chord = dx*TMath::Sqrt(1+dy2dx*dy2dx); // distance from old position to new one
78 // double rot = 2*TMath::ASin(0.5*chord*crv); // angular difference seen from the circle center
79 // track1 += rot/crv*track3;
80 //
81 auto arg = r1 * f2 - r2 * f1;
82 if (gpu::CAMath::Abs(arg) > constants::math::Almost1) {
83 return false;
84 }
85 value_t rot = gpu::CAMath::ASin(arg); // more economic version from Yura.
86 if (f1 * f1 + f2 * f2 > 1.f && f1 * f2 < 0.f) { // special cases of large rotations or large abs angles
87 if (f2 > 0.f) {
88 rot = constants::math::PI - rot; //
89 } else {
90 rot = -constants::math::PI - rot;
91 }
92 }
93 dP[kZ] = this->getTgl() / crv * rot;
94 } else {
95 dP[kZ] = dx * (r2 + f2 * dy2dx) * this->getTgl();
96 }
97 this->setX(xk);
98 dP[kY] = dx * dy2dx;
99 dP[kSnp] = x2r;
100
101 this->updateParams(dP); // apply corrections
102
103 value_t &c00 = mC[kSigY2], &c10 = mC[kSigZY], &c11 = mC[kSigZ2], &c20 = mC[kSigSnpY], &c21 = mC[kSigSnpZ],
104 &c22 = mC[kSigSnp2], &c30 = mC[kSigTglY], &c31 = mC[kSigTglZ], &c32 = mC[kSigTglSnp], &c33 = mC[kSigTgl2],
105 &c40 = mC[kSigQ2PtY], &c41 = mC[kSigQ2PtZ], &c42 = mC[kSigQ2PtSnp], &c43 = mC[kSigQ2PtTgl],
106 &c44 = mC[kSigQ2Pt2];
107
108 // evaluate matrix in double prec.
109 double rinv = 1. / r1;
110 double r3inv = rinv * rinv * rinv;
111 double f24 = dx * b * constants::math::B2C; // x2r/mP[kQ2Pt];
112 double f02 = dx * r3inv;
113 double f04 = 0.5 * f24 * f02;
114 double f12 = f02 * this->getTgl() * f1;
115 double f14 = 0.5 * f24 * f12; // 0.5*f24*f02*getTgl()*f1;
116 double f13 = dx * rinv;
117
118 // b = C*ft
119 double b00 = f02 * c20 + f04 * c40, b01 = f12 * c20 + f14 * c40 + f13 * c30;
120 double b02 = f24 * c40;
121 double b10 = f02 * c21 + f04 * c41, b11 = f12 * c21 + f14 * c41 + f13 * c31;
122 double b12 = f24 * c41;
123 double b20 = f02 * c22 + f04 * c42, b21 = f12 * c22 + f14 * c42 + f13 * c32;
124 double b22 = f24 * c42;
125 double b40 = f02 * c42 + f04 * c44, b41 = f12 * c42 + f14 * c44 + f13 * c43;
126 double b42 = f24 * c44;
127 double b30 = f02 * c32 + f04 * c43, b31 = f12 * c32 + f14 * c43 + f13 * c33;
128 double b32 = f24 * c43;
129
130 // a = f*b = f*C*ft
131 double a00 = f02 * b20 + f04 * b40, a01 = f02 * b21 + f04 * b41, a02 = f02 * b22 + f04 * b42;
132 double a11 = f12 * b21 + f14 * b41 + f13 * b31, a12 = f12 * b22 + f14 * b42 + f13 * b32;
133 double a22 = f24 * b42;
134
135 // F*C*Ft = C + (b + bt + a)
136 c00 += b00 + b00 + a00;
137 c10 += b10 + b01 + a01;
138 c20 += b20 + b02 + a02;
139 c30 += b30;
140 c40 += b40;
141 c11 += b11 + b11 + a11;
142 c21 += b21 + b12 + a12;
143 c31 += b31;
144 c41 += b41;
145 c22 += b22 + b22 + a22;
146 c32 += b32;
147 c42 += b42;
148
149 checkCovariance();
150
151 return true;
152}
153
154//______________________________________________________________
155template <typename value_T>
156GPUd() bool TrackParametrizationWithError<value_T>::testRotate(value_t) const
157{
158 // no ops
159 return true;
160}
161//______________________________________________________________
162template <typename value_T>
163GPUd() bool TrackParametrizationWithError<value_T>::rotate(value_t alpha)
164{
165 // rotate to alpha frame
166 if (gpu::CAMath::Abs(this->getSnp()) > constants::math::Almost1) {
167 LOGP(debug, "Precondition is not satisfied: |sin(phi)|>1 ! {:f}", this->getSnp());
168 return false;
169 }
170 //
171 math_utils::detail::bringToPMPi<value_t>(alpha);
172 //
173 value_t ca = 0, sa = 0;
174 math_utils::detail::sincos(alpha - this->getAlpha(), sa, ca);
175 value_t snp = this->getSnp(), csp = gpu::CAMath::Sqrt((1.f - snp) * (1.f + snp)); // Improve precision
176 // RS: check if rotation does no invalidate track model (cos(local_phi)>=0, i.e. particle
177 // direction in local frame is along the X axis
178 if ((csp * ca + snp * sa) < 0) {
179 // LOGP(warning,"Rotation failed: local cos(phi) would become {:.2f}", csp * ca + snp * sa);
180 return false;
181 }
182 //
183
184 value_t updSnp = snp * ca - csp * sa;
185 if (gpu::CAMath::Abs(updSnp) > constants::math::Almost1) {
186 LOGP(debug, "Rotation failed: new snp {:.2f}", updSnp);
187 return false;
188 }
189 value_t xold = this->getX(), yold = this->getY();
190 this->setAlpha(alpha);
191 this->setX(xold * ca + yold * sa);
192 this->setY(-xold * sa + yold * ca);
193 this->setSnp(updSnp);
194
195 if (gpu::CAMath::Abs(csp) < constants::math::Almost0) {
196 LOGP(debug, "Too small cosine value {:f}", csp);
197 csp = constants::math::Almost0;
198 }
199
200 value_t rr = (ca + snp / csp * sa);
201
202 mC[kSigY2] *= (ca * ca);
203 mC[kSigZY] *= ca;
204 mC[kSigSnpY] *= ca * rr;
205 mC[kSigSnpZ] *= rr;
206 mC[kSigSnp2] *= rr * rr;
207 mC[kSigTglY] *= ca;
208 mC[kSigTglSnp] *= rr;
209 mC[kSigQ2PtY] *= ca;
210 mC[kSigQ2PtSnp] *= rr;
211
212 checkCovariance();
213 return true;
214}
215
216//_______________________________________________________________________
217template <typename value_T>
218GPUd() bool TrackParametrizationWithError<value_T>::propagateToDCA(const o2::dataformats::VertexBase& vtx, value_t b, o2::dataformats::DCA* dca, value_t maxD)
219{
220 // propagate track to DCA to the vertex
221 value_t sn, cs, alp = this->getAlpha();
222 o2::math_utils::detail::sincos(alp, sn, cs);
223 value_t x = this->getX(), y = this->getY(), snp = this->getSnp(), csp = gpu::CAMath::Sqrt((1.f - snp) * (1.f + snp));
224 value_t xv = vtx.getX() * cs + vtx.getY() * sn, yv = -vtx.getX() * sn + vtx.getY() * cs, zv = vtx.getZ();
225 x -= xv;
226 y -= yv;
227 // Estimate the impact parameter neglecting the track curvature
228 value_t d = gpu::CAMath::Abs(x * snp - y * csp);
229 if (d > maxD) {
230 return false;
231 }
232 value_t crv = this->getCurvature(b);
233 value_t tgfv = -(crv * x - snp) / (crv * y + csp);
234 sn = tgfv / gpu::CAMath::Sqrt(1.f + tgfv * tgfv);
235 cs = gpu::CAMath::Sqrt((1.f - sn) * (1.f + sn));
236 cs = (gpu::CAMath::Abs(tgfv) > constants::math::Almost0) ? sn / tgfv : constants::math::Almost1;
237
238 x = xv * cs + yv * sn;
239 yv = -xv * sn + yv * cs;
240 xv = x;
241
242 auto tmpT(*this); // operate on the copy to recover after the failure
243 alp += gpu::CAMath::ASin(sn);
244 if (!tmpT.rotate(alp) || !tmpT.propagateTo(xv, b)) {
245#if !defined(GPUCA_ALIGPUCODE)
246 LOG(debug) << "failed to propagate to alpha=" << alp << " X=" << xv << vtx << " | Track is: " << tmpT.asString();
247#endif
248 return false;
249 }
250 *this = tmpT;
251 if (dca) {
252 o2::math_utils::detail::sincos(alp, sn, cs);
253 auto s2ylocvtx = vtx.getSigmaX2() * sn * sn + vtx.getSigmaY2() * cs * cs - 2. * vtx.getSigmaXY() * cs * sn;
254 dca->set(this->getY() - yv, this->getZ() - zv, getSigmaY2() + s2ylocvtx, getSigmaZY(), getSigmaZ2() + vtx.getSigmaZ2());
255 }
256 return true;
257}
258
259//______________________________________________________________
260template <typename value_T>
261GPUd() TrackParametrizationWithError<value_T>::TrackParametrizationWithError(const dim3_t& xyz, const dim3_t& pxpypz,
262 const gpu::gpustd::array<value_t, kLabCovMatSize>& cv, int charge, bool sectorAlpha, const PID pid)
263{
264 // construct track param and covariance from kinematics and lab errors
265 set(xyz, pxpypz, cv, charge, sectorAlpha, pid);
266}
267
268//______________________________________________________________
269template <typename value_T>
270GPUd() void TrackParametrizationWithError<value_T>::set(const dim3_t& xyz, const dim3_t& pxpypz,
271 const gpu::gpustd::array<value_t, kLabCovMatSize>& cv, int charge, bool sectorAlpha, const PID pid)
272{
273 // set track param and covariance from kinematics and lab errors
274
275 // Alpha of the frame is defined as:
276 // sectorAlpha == false : -> angle of pt direction
277 // sectorAlpha == true : -> angle of the sector from X,Y coordinate for r>1
278 // angle of pt direction for r==0
279 //
280 //
281 constexpr value_t kSafe = 1e-5f;
282 value_t radPos2 = xyz[0] * xyz[0] + xyz[1] * xyz[1];
283 value_t alp = 0;
284 if (sectorAlpha || radPos2 < 1) {
285 alp = gpu::CAMath::ATan2(pxpypz[1], pxpypz[0]);
286 } else {
287 alp = gpu::CAMath::ATan2(xyz[1], xyz[0]);
288 }
289 if (sectorAlpha) {
290 alp = math_utils::detail::angle2Alpha<value_t>(alp);
291 }
292 //
293 value_t sn, cs;
294 math_utils::detail::sincos(alp, sn, cs);
295 // protection against cosp<0
296 if (cs * pxpypz[0] + sn * pxpypz[1] < 0) {
297 LOG(debug) << "alpha from phiPos() will invalidate this track parameters, overriding to alpha from phi()";
298 alp = gpu::CAMath::ATan2(pxpypz[1], pxpypz[0]);
299 if (sectorAlpha) {
300 alp = math_utils::detail::angle2Alpha<value_t>(alp);
301 }
302 math_utils::detail::sincos(alp, sn, cs);
303 }
304 // protection: avoid alpha being too close to 0 or +-pi/2
305 if (gpu::CAMath::Abs(sn) < 2.f * kSafe) {
306 if (alp > 0) {
307 alp += alp < constants::math::PIHalf ? 2.f * kSafe : -2.f * kSafe;
308 } else {
309 alp += alp > -constants::math::PIHalf ? -2.f * kSafe : 2.f * kSafe;
310 }
311 math_utils::detail::sincos(alp, sn, cs);
312 } else if (gpu::CAMath::Abs(cs) < 2.f * kSafe) {
313 if (alp > 0) {
314 alp += alp > constants::math::PIHalf ? 2.f * kSafe : -2.f * kSafe;
315 } else {
316 alp += alp > -constants::math::PIHalf ? 2.f * kSafe : -2.f * kSafe;
317 }
318 math_utils::detail::sincos(alp, sn, cs);
319 }
320 // get the vertex of origin and the momentum
321 dim3_t ver{xyz[0], xyz[1], xyz[2]};
322 dim3_t mom{pxpypz[0], pxpypz[1], pxpypz[2]};
323 //
324 // Rotate to the local coordinate system
325 math_utils::detail::rotateZ<value_t>(ver, -alp);
326 math_utils::detail::rotateZ<value_t>(mom, -alp);
327 //
328 value_t pt = gpu::CAMath::Sqrt(mom[0] * mom[0] + mom[1] * mom[1]);
329 value_t ptI = 1.f / pt;
330 this->setX(ver[0]);
331 this->setAlpha(alp);
332 this->setY(ver[1]);
333 this->setZ(ver[2]);
334 this->setSnp(mom[1] * ptI); // cos(phi)
335 this->setTgl(mom[2] * ptI); // tg(lambda)
336 this->setAbsCharge(gpu::CAMath::Abs(charge));
337 this->setQ2Pt(charge ? ptI * charge : ptI);
338 this->setPID(pid);
339 //
340 if (gpu::CAMath::Abs(1.f - this->getSnp()) < kSafe) {
341 this->setSnp(1.f - kSafe); // Protection
342 } else if (gpu::CAMath::Abs(-1.f - this->getSnp()) < kSafe) {
343 this->setSnp(-1.f + kSafe); // Protection
344 }
345 //
346 // Covariance matrix (formulas to be simplified)
347 value_t r = mom[0] * ptI; // cos(phi)
348 value_t cv34 = gpu::CAMath::Sqrt(cv[3] * cv[3] + cv[4] * cv[4]);
349 //
350 int special = 0;
351 value_t sgcheck = r * sn + this->getSnp() * cs;
352 if (gpu::CAMath::Abs(sgcheck) > 1 - kSafe) { // special case: lab phi is +-pi/2
353 special = 1;
354 sgcheck = sgcheck < 0 ? -1.f : 1.f;
355 } else if (gpu::CAMath::Abs(sgcheck) < kSafe) {
356 sgcheck = cs < 0 ? -1.0f : 1.0f;
357 special = 2; // special case: lab phi is 0
358 }
359 //
360 mC[kSigY2] = cv[0] + cv[2];
361 mC[kSigZY] = (-cv[3] * sn) < 0 ? -cv34 : cv34;
362 mC[kSigZ2] = cv[5];
363 //
364 value_t ptI2 = ptI * ptI;
365 value_t tgl2 = this->getTgl() * this->getTgl();
366 if (special == 1) {
367 mC[kSigSnpY] = cv[6] * ptI;
368 mC[kSigSnpZ] = -sgcheck * cv[8] * r * ptI;
369 mC[kSigSnp2] = gpu::CAMath::Abs(cv[9] * r * r * ptI2);
370 mC[kSigTglY] = (cv[10] * this->getTgl() - sgcheck * cv[15]) * ptI / r;
371 mC[kSigTglZ] = (cv[17] - sgcheck * cv[12] * this->getTgl()) * ptI;
372 mC[kSigTglSnp] = (-sgcheck * cv[18] + cv[13] * this->getTgl()) * r * ptI2;
373 mC[kSigTgl2] = gpu::CAMath::Abs(cv[20] - 2 * sgcheck * cv[19] * mC[4] + cv[14] * tgl2) * ptI2;
374 mC[kSigQ2PtY] = cv[10] * ptI2 / r * charge;
375 mC[kSigQ2PtZ] = -sgcheck * cv[12] * ptI2 * charge;
376 mC[kSigQ2PtSnp] = cv[13] * r * ptI * ptI2 * charge;
377 mC[kSigQ2PtTgl] = (-sgcheck * cv[19] + cv[14] * this->getTgl()) * r * ptI2 * ptI;
378 mC[kSigQ2Pt2] = gpu::CAMath::Abs(cv[14] * ptI2 * ptI2);
379 } else if (special == 2) {
380 mC[kSigSnpY] = -cv[10] * ptI * cs / sn;
381 mC[kSigSnpZ] = cv[12] * cs * ptI;
382 mC[kSigSnp2] = gpu::CAMath::Abs(cv[14] * cs * cs * ptI2);
383 mC[kSigTglY] = (sgcheck * cv[6] * this->getTgl() - cv[15]) * ptI / sn;
384 mC[kSigTglZ] = (cv[17] - sgcheck * cv[8] * this->getTgl()) * ptI;
385 mC[kSigTglSnp] = (cv[19] - sgcheck * cv[13] * this->getTgl()) * cs * ptI2;
386 mC[kSigTgl2] = gpu::CAMath::Abs(cv[20] - 2 * sgcheck * cv[18] * this->getTgl() + cv[9] * tgl2) * ptI2;
387 mC[kSigQ2PtY] = sgcheck * cv[6] * ptI2 / sn * charge;
388 mC[kSigQ2PtZ] = -sgcheck * cv[8] * ptI2 * charge;
389 mC[kSigQ2PtSnp] = -sgcheck * cv[13] * cs * ptI * ptI2 * charge;
390 mC[kSigQ2PtTgl] = (-sgcheck * cv[18] + cv[9] * this->getTgl()) * ptI2 * ptI * charge;
391 mC[kSigQ2Pt2] = gpu::CAMath::Abs(cv[9] * ptI2 * ptI2);
392 } else {
393 double m00 = -sn; // m10=cs;
394 double m23 = -pt * (sn + this->getSnp() * cs / r), m43 = -pt * pt * (r * cs - this->getSnp() * sn);
395 double m24 = pt * (cs - this->getSnp() * sn / r), m44 = -pt * pt * (r * sn + this->getSnp() * cs);
396 double m35 = pt, m45 = -pt * pt * this->getTgl();
397 //
398 if (charge) { // RS: this is a hack, proper treatment to be implemented
399 m43 *= charge;
400 m44 *= charge;
401 m45 *= charge;
402 }
403 //
404 double a1 = cv[13] - cv[9] * (m23 * m44 + m43 * m24) / m23 / m43;
405 double a2 = m23 * m24 - m23 * (m23 * m44 + m43 * m24) / m43;
406 double a3 = m43 * m44 - m43 * (m23 * m44 + m43 * m24) / m23;
407 double a4 = cv[14] + 2. * cv[9];
408 double a5 = m24 * m24 - 2. * m24 * m44 * m23 / m43;
409 double a6 = m44 * m44 - 2. * m24 * m44 * m43 / m23;
410 //
411 mC[kSigSnpY] = (cv[10] * m43 - cv[6] * m44) / (m24 * m43 - m23 * m44) / m00;
412 mC[kSigQ2PtY] = (cv[6] / m00 - mC[kSigSnpY] * m23) / m43;
413 mC[kSigTglY] = (cv[15] / m00 - mC[kSigQ2PtY] * m45) / m35;
414 mC[kSigSnpZ] = (cv[12] * m43 - cv[8] * m44) / (m24 * m43 - m23 * m44);
415 mC[kSigQ2PtZ] = (cv[8] - mC[kSigSnpZ] * m23) / m43;
416 mC[kSigTglZ] = cv[17] / m35 - mC[kSigQ2PtZ] * m45 / m35;
417 mC[kSigSnp2] = gpu::CAMath::Abs((a4 * a3 - a6 * a1) / (a5 * a3 - a6 * a2));
418 mC[kSigQ2Pt2] = gpu::CAMath::Abs((a1 - a2 * mC[kSigSnp2]) / a3);
419 mC[kSigQ2PtSnp] = (cv[9] - mC[kSigSnp2] * m23 * m23 - mC[kSigQ2Pt2] * m43 * m43) / m23 / m43;
420 double b1 = cv[18] - mC[kSigQ2PtSnp] * m23 * m45 - mC[kSigQ2Pt2] * m43 * m45;
421 double b2 = m23 * m35;
422 double b3 = m43 * m35;
423 double b4 = cv[19] - mC[kSigQ2PtSnp] * m24 * m45 - mC[kSigQ2Pt2] * m44 * m45;
424 double b5 = m24 * m35;
425 double b6 = m44 * m35;
426 mC[kSigTglSnp] = (b4 - b6 * b1 / b3) / (b5 - b6 * b2 / b3);
427 mC[kSigQ2PtTgl] = b1 / b3 - b2 * mC[kSigTglSnp] / b3;
428 mC[kSigTgl2] = gpu::CAMath::Abs((cv[20] - mC[kSigQ2Pt2] * (m45 * m45) - mC[kSigQ2PtTgl] * 2.f * m35 * m45) / (m35 * m35));
429 }
430 checkCovariance();
431}
432
433//____________________________________________________________
434template <typename value_T>
435GPUd() bool TrackParametrizationWithError<value_T>::propagateTo(value_t xk, const dim3_t& b)
436{
437 //----------------------------------------------------------------
438 // Extrapolate this track to the plane X=xk in the field b[].
439 //
440 // X [cm] is in the "tracking coordinate system" of this track.
441 // b[]={Bx,By,Bz} [kG] is in the Global coordidate system.
442 //----------------------------------------------------------------
443
444 value_t dx = xk - this->getX();
445 if (gpu::CAMath::Abs(dx) < constants::math::Almost0) {
446 return true;
447 }
448 // Do not propagate tracks outside the ALICE detector
449 if (gpu::CAMath::Abs(dx) > 1e5 || gpu::CAMath::Abs(this->getY()) > 1e5 || gpu::CAMath::Abs(this->getZ()) > 1e5) {
450 LOG(warning) << "Anomalous track, target X:" << xk;
451 // print();
452 return false;
453 }
454 value_t crv = (gpu::CAMath::Abs(b[2]) < constants::math::Almost0) ? 0.f : this->getCurvature(b[2]);
455 if (gpu::CAMath::Abs(crv) < constants::math::Almost0) {
456 return propagateTo(xk, 0.);
457 }
458 value_t x2r = crv * dx;
459 value_t f1 = this->getSnp(), f2 = f1 + x2r;
460 if ((gpu::CAMath::Abs(f1) > constants::math::Almost1) || (gpu::CAMath::Abs(f2) > constants::math::Almost1)) {
461 return false;
462 }
463 value_t r1 = gpu::CAMath::Sqrt((1.f - f1) * (1.f + f1));
464 if (gpu::CAMath::Abs(r1) < constants::math::Almost0) {
465 return false;
466 }
467 value_t r2 = gpu::CAMath::Sqrt((1.f - f2) * (1.f + f2));
468 if (gpu::CAMath::Abs(r2) < constants::math::Almost0) {
469 return false;
470 }
471
472 value_t dy2dx = (f1 + f2) / (r1 + r2);
473 value_t step = (gpu::CAMath::Abs(x2r) < 0.05f) ? dx * gpu::CAMath::Abs(r2 + f2 * dy2dx) // chord
474 : 2.f * gpu::CAMath::ASin(0.5f * dx * gpu::CAMath::Sqrt(1.f + dy2dx * dy2dx) * crv) / crv; // arc
475 step *= gpu::CAMath::Sqrt(1.f + this->getTgl() * this->getTgl());
476 //
477 // get the track x,y,z,px/p,py/p,pz/p,p,sinAlpha,cosAlpha in the Global System
479 if (!this->getPosDirGlo(vecLab)) {
480 return false;
481 }
482 //
483 // matrix transformed with Bz component only
484 value_t &c00 = mC[kSigY2], &c10 = mC[kSigZY], &c11 = mC[kSigZ2], &c20 = mC[kSigSnpY], &c21 = mC[kSigSnpZ],
485 &c22 = mC[kSigSnp2], &c30 = mC[kSigTglY], &c31 = mC[kSigTglZ], &c32 = mC[kSigTglSnp], &c33 = mC[kSigTgl2],
486 &c40 = mC[kSigQ2PtY], &c41 = mC[kSigQ2PtZ], &c42 = mC[kSigQ2PtSnp], &c43 = mC[kSigQ2PtTgl],
487 &c44 = mC[kSigQ2Pt2];
488 // evaluate matrix in double prec.
489 double rinv = 1. / r1;
490 double r3inv = rinv * rinv * rinv;
491 double f24 = dx * b[2] * constants::math::B2C; // x2r/track[kQ2Pt];
492 double f02 = dx * r3inv;
493 double f04 = 0.5 * f24 * f02;
494 double f12 = f02 * this->getTgl() * f1;
495 double f14 = 0.5 * f24 * f12; // 0.5*f24*f02*getTgl()*f1;
496 double f13 = dx * rinv;
497
498 // b = C*ft
499 double b00 = f02 * c20 + f04 * c40, b01 = f12 * c20 + f14 * c40 + f13 * c30;
500 double b02 = f24 * c40;
501 double b10 = f02 * c21 + f04 * c41, b11 = f12 * c21 + f14 * c41 + f13 * c31;
502 double b12 = f24 * c41;
503 double b20 = f02 * c22 + f04 * c42, b21 = f12 * c22 + f14 * c42 + f13 * c32;
504 double b22 = f24 * c42;
505 double b40 = f02 * c42 + f04 * c44, b41 = f12 * c42 + f14 * c44 + f13 * c43;
506 double b42 = f24 * c44;
507 double b30 = f02 * c32 + f04 * c43, b31 = f12 * c32 + f14 * c43 + f13 * c33;
508 double b32 = f24 * c43;
509
510 // a = f*b = f*C*ft
511 double a00 = f02 * b20 + f04 * b40, a01 = f02 * b21 + f04 * b41, a02 = f02 * b22 + f04 * b42;
512 double a11 = f12 * b21 + f14 * b41 + f13 * b31, a12 = f12 * b22 + f14 * b42 + f13 * b32;
513 double a22 = f24 * b42;
514
515 // F*C*Ft = C + (b + bt + a)
516 c00 += b00 + b00 + a00;
517 c10 += b10 + b01 + a01;
518 c20 += b20 + b02 + a02;
519 c30 += b30;
520 c40 += b40;
521 c11 += b11 + b11 + a11;
522 c21 += b21 + b12 + a12;
523 c31 += b31;
524 c41 += b41;
525 c22 += b22 + b22 + a22;
526 c32 += b32;
527 c42 += b42;
528
529 checkCovariance();
530
531 // Rotate to the system where Bx=By=0.
532 value_t bxy2 = b[0] * b[0] + b[1] * b[1];
533 value_t bt = gpu::CAMath::Sqrt(bxy2);
534 value_t cosphi = 1.f, sinphi = 0.f;
535 if (bt > constants::math::Almost0) {
536 cosphi = b[0] / bt;
537 sinphi = b[1] / bt;
538 }
539 value_t bb = gpu::CAMath::Sqrt(bxy2 + b[2] * b[2]);
540 value_t costet = 1., sintet = 0.;
541 if (bb > constants::math::Almost0) {
542 costet = b[2] / bb;
543 sintet = bt / bb;
544 }
545 gpu::gpustd::array<value_t, 7> vect{costet * cosphi * vecLab[0] + costet * sinphi * vecLab[1] - sintet * vecLab[2],
546 -sinphi * vecLab[0] + cosphi * vecLab[1],
547 sintet * cosphi * vecLab[0] + sintet * sinphi * vecLab[1] + costet * vecLab[2],
548 costet * cosphi * vecLab[3] + costet * sinphi * vecLab[4] - sintet * vecLab[5],
549 -sinphi * vecLab[3] + cosphi * vecLab[4],
550 sintet * cosphi * vecLab[3] + sintet * sinphi * vecLab[4] + costet * vecLab[5],
551 vecLab[6]};
552
553 // Do the helix step
554 value_t q = this->getCharge();
555 g3helx3(q * bb, step, vect);
556
557 // Rotate back to the Global System
558 vecLab[0] = cosphi * costet * vect[0] - sinphi * vect[1] + cosphi * sintet * vect[2];
559 vecLab[1] = sinphi * costet * vect[0] + cosphi * vect[1] + sinphi * sintet * vect[2];
560 vecLab[2] = -sintet * vect[0] + costet * vect[2];
561
562 vecLab[3] = cosphi * costet * vect[3] - sinphi * vect[4] + cosphi * sintet * vect[5];
563 vecLab[4] = sinphi * costet * vect[3] + cosphi * vect[4] + sinphi * sintet * vect[5];
564 vecLab[5] = -sintet * vect[3] + costet * vect[5];
565
566 // Rotate back to the Tracking System
567 value_t sinalp = -vecLab[7], cosalp = vecLab[8];
568 value_t t = cosalp * vecLab[0] - sinalp * vecLab[1];
569 vecLab[1] = sinalp * vecLab[0] + cosalp * vecLab[1];
570 vecLab[0] = t;
571 t = cosalp * vecLab[3] - sinalp * vecLab[4];
572 vecLab[4] = sinalp * vecLab[3] + cosalp * vecLab[4];
573 vecLab[3] = t;
574
575 // Do the final correcting step to the target plane (linear approximation)
576 value_t x = vecLab[0], y = vecLab[1], z = vecLab[2];
577 if (gpu::CAMath::Abs(dx) > constants::math::Almost0) {
578 if (gpu::CAMath::Abs(vecLab[3]) < constants::math::Almost0) {
579 return false;
580 }
581 dx = xk - vecLab[0];
582 x += dx;
583 y += vecLab[4] / vecLab[3] * dx;
584 z += vecLab[5] / vecLab[3] * dx;
585 }
586
587 // Calculate the track parameters
588 t = 1.f / gpu::CAMath::Sqrt(vecLab[3] * vecLab[3] + vecLab[4] * vecLab[4]);
589 this->setX(xk);
590 this->setY(y);
591 this->setZ(z);
592 this->setSnp(vecLab[4] * t);
593 this->setTgl(vecLab[5] * t);
594 this->setQ2Pt(q * t / vecLab[6]);
595
596 return true;
597}
598
599//______________________________________________
600template <typename value_T>
601GPUd() void TrackParametrizationWithError<value_T>::checkCorrelations()
602{
603 // This function forces the abs of correlation coefficients to be <1.
604 constexpr float MaxCorr = 0.99;
605 for (int i = kNParams; i--;) {
606 for (int j = i; j--;) {
607 auto sig2 = mC[DiagMap[i]] * mC[DiagMap[j]];
608 auto& cov = mC[CovarMap[i][j]];
609 if (cov * cov >= MaxCorr * sig2) { // constrain correlation
610 cov = gpu::CAMath::Sqrt(sig2) * (cov > 0. ? MaxCorr : -MaxCorr);
611 }
612 }
613 }
614}
615
616//______________________________________________
617template <typename value_T>
618GPUd() void TrackParametrizationWithError<value_T>::checkCovariance()
619{
620 // This function forces the diagonal elements of the covariance matrix to be positive and abs of correlation coefficients to be <1.
621 // In case the diagonal element is bigger than the maximal allowed value, it is set to
622 // the limit and the off-diagonal elements that correspond to it are set to zero.
623
624 mC[kSigY2] = gpu::CAMath::Abs(mC[kSigY2]);
625 if (mC[kSigY2] > kCY2max) {
626 value_t scl = gpu::CAMath::Sqrt(kCY2max / mC[kSigY2]);
627 mC[kSigY2] = kCY2max;
628 mC[kSigZY] *= scl;
629 mC[kSigSnpY] *= scl;
630 mC[kSigTglY] *= scl;
631 mC[kSigQ2PtY] *= scl;
632 }
633 mC[kSigZ2] = gpu::CAMath::Abs(mC[kSigZ2]);
634 if (mC[kSigZ2] > kCZ2max) {
635 value_t scl = gpu::CAMath::Sqrt(kCZ2max / mC[kSigZ2]);
636 mC[kSigZ2] = kCZ2max;
637 mC[kSigZY] *= scl;
638 mC[kSigSnpZ] *= scl;
639 mC[kSigTglZ] *= scl;
640 mC[kSigQ2PtZ] *= scl;
641 }
642 mC[kSigSnp2] = gpu::CAMath::Abs(mC[kSigSnp2]);
643 if (mC[kSigSnp2] > kCSnp2max) {
644 value_t scl = gpu::CAMath::Sqrt(kCSnp2max / mC[kSigSnp2]);
645 mC[kSigSnp2] = kCSnp2max;
646 mC[kSigSnpY] *= scl;
647 mC[kSigSnpZ] *= scl;
648 mC[kSigTglSnp] *= scl;
649 mC[kSigQ2PtSnp] *= scl;
650 }
651 mC[kSigTgl2] = gpu::CAMath::Abs(mC[kSigTgl2]);
652 if (mC[kSigTgl2] > kCTgl2max) {
653 value_t scl = gpu::CAMath::Sqrt(kCTgl2max / mC[kSigTgl2]);
654 mC[kSigTgl2] = kCTgl2max;
655 mC[kSigTglY] *= scl;
656 mC[kSigTglZ] *= scl;
657 mC[kSigTglSnp] *= scl;
658 mC[kSigQ2PtTgl] *= scl;
659 }
660 mC[kSigQ2Pt2] = gpu::CAMath::Abs(mC[kSigQ2Pt2]);
661 if (mC[kSigQ2Pt2] > kC1Pt2max) {
662 value_t scl = gpu::CAMath::Sqrt(kC1Pt2max / mC[kSigQ2Pt2]);
663 mC[kSigQ2Pt2] = kC1Pt2max;
664 mC[kSigQ2PtY] *= scl;
665 mC[kSigQ2PtZ] *= scl;
666 mC[kSigQ2PtSnp] *= scl;
667 mC[kSigQ2PtTgl] *= scl;
668 }
669}
670
671//______________________________________________
672template <typename value_T>
673GPUd() void TrackParametrizationWithError<value_T>::resetCovariance(value_t s2)
674{
675 // Reset the covarince matrix to "something big"
676 double d0(kCY2max), d1(kCZ2max), d2(kCSnp2max), d3(kCTgl2max), d4(kC1Pt2max);
677 if (s2 > constants::math::Almost0) {
678 d0 = getSigmaY2() * s2;
679 d1 = getSigmaZ2() * s2;
680 d2 = getSigmaSnp2() * s2;
681 d3 = getSigmaTgl2() * s2;
682 d4 = getSigma1Pt2() * s2;
683 if (d0 > kCY2max) {
684 d0 = kCY2max;
685 }
686 if (d1 > kCZ2max) {
687 d1 = kCZ2max;
688 }
689 if (d2 > kCSnp2max) {
690 d2 = kCSnp2max;
691 }
692 if (d3 > kCTgl2max) {
693 d3 = kCTgl2max;
694 }
695 if (d4 > kC1Pt2max) {
696 d4 = kC1Pt2max;
697 }
698 }
699 for (int i = 0; i < kCovMatSize; i++) {
700 mC[i] = 0;
701 }
702 mC[kSigY2] = d0;
703 mC[kSigZ2] = d1;
704 mC[kSigSnp2] = d2;
705 mC[kSigTgl2] = d3;
706 mC[kSigQ2Pt2] = d4;
707}
708
709//______________________________________________
710template <typename value_T>
711GPUd() auto TrackParametrizationWithError<value_T>::getPredictedChi2(const value_t* p, const value_t* cov) const -> value_t
712{
713 // Estimate the chi2 of the space point "p" with the cov. matrix "cov"
714 auto sdd = static_cast<double>(getSigmaY2()) + static_cast<double>(cov[0]);
715 auto sdz = static_cast<double>(getSigmaZY()) + static_cast<double>(cov[1]);
716 auto szz = static_cast<double>(getSigmaZ2()) + static_cast<double>(cov[2]);
717 auto det = sdd * szz - sdz * sdz;
718
719 if (gpu::CAMath::Abs(det) < constants::math::Almost0) {
720 return constants::math::VeryBig;
721 }
722
723 value_t d = this->getY() - p[0];
724 value_t z = this->getZ() - p[1];
725 auto chi2 = (d * (szz * d - sdz * z) + z * (sdd * z - d * sdz)) / det;
726 if (chi2 < 0.) {
727#ifndef GPUCA_ALIGPUCODE
728 LOGP(warning, "Negative chi2={}, Cluster: {} {} {} Dy:{} Dz:{} | sdd:{} sdz:{} szz:{} det:{}", chi2, cov[0], cov[1], cov[2], d, z, sdd, sdz, szz, det);
729 LOGP(warning, "Track: {}", asString());
730#endif
731 }
732 return chi2;
733}
734
735//______________________________________________
736template <typename value_T>
737GPUd() auto TrackParametrizationWithError<value_T>::getPredictedChi2Quiet(const value_t* p, const value_t* cov) const -> value_t
738{
739 // Estimate the chi2 of the space point "p" with the cov. matrix "cov"
740 auto sdd = static_cast<double>(getSigmaY2()) + static_cast<double>(cov[0]);
741 auto sdz = static_cast<double>(getSigmaZY()) + static_cast<double>(cov[1]);
742 auto szz = static_cast<double>(getSigmaZ2()) + static_cast<double>(cov[2]);
743 auto det = sdd * szz - sdz * sdz;
744
745 if (gpu::CAMath::Abs(det) < constants::math::Almost0) {
746 return constants::math::VeryBig;
747 }
748
749 value_t d = this->getY() - p[0];
750 value_t z = this->getZ() - p[1];
751
752 return (d * (szz * d - sdz * z) + z * (sdd * z - d * sdz)) / det;
753}
754
755//______________________________________________
756template <typename value_T>
757GPUd() auto TrackParametrizationWithError<value_T>::getPredictedChi2(const TrackParametrizationWithError<value_T>& rhs) const -> value_t
758{
759 MatrixDSym5 cov; // perform matrix operations in double!
760 return getPredictedChi2(rhs, cov);
761}
762
763//______________________________________________
764template <typename value_T>
765GPUd() void TrackParametrizationWithError<value_T>::buildCombinedCovMatrix(const TrackParametrizationWithError<value_T>& rhs, MatrixDSym5& cov) const
766{
767 // fill combined cov.matrix (NOT inverted)
768 cov(kY, kY) = static_cast<double>(getSigmaY2()) + static_cast<double>(rhs.getSigmaY2());
769 cov(kZ, kY) = static_cast<double>(getSigmaZY()) + static_cast<double>(rhs.getSigmaZY());
770 cov(kZ, kZ) = static_cast<double>(getSigmaZ2()) + static_cast<double>(rhs.getSigmaZ2());
771 cov(kSnp, kY) = static_cast<double>(getSigmaSnpY()) + static_cast<double>(rhs.getSigmaSnpY());
772 cov(kSnp, kZ) = static_cast<double>(getSigmaSnpZ()) + static_cast<double>(rhs.getSigmaSnpZ());
773 cov(kSnp, kSnp) = static_cast<double>(getSigmaSnp2()) + static_cast<double>(rhs.getSigmaSnp2());
774 cov(kTgl, kY) = static_cast<double>(getSigmaTglY()) + static_cast<double>(rhs.getSigmaTglY());
775 cov(kTgl, kZ) = static_cast<double>(getSigmaTglZ()) + static_cast<double>(rhs.getSigmaTglZ());
776 cov(kTgl, kSnp) = static_cast<double>(getSigmaTglSnp()) + static_cast<double>(rhs.getSigmaTglSnp());
777 cov(kTgl, kTgl) = static_cast<double>(getSigmaTgl2()) + static_cast<double>(rhs.getSigmaTgl2());
778 cov(kQ2Pt, kY) = static_cast<double>(getSigma1PtY()) + static_cast<double>(rhs.getSigma1PtY());
779 cov(kQ2Pt, kZ) = static_cast<double>(getSigma1PtZ()) + static_cast<double>(rhs.getSigma1PtZ());
780 cov(kQ2Pt, kSnp) = static_cast<double>(getSigma1PtSnp()) + static_cast<double>(rhs.getSigma1PtSnp());
781 cov(kQ2Pt, kTgl) = static_cast<double>(getSigma1PtTgl()) + static_cast<double>(rhs.getSigma1PtTgl());
782 cov(kQ2Pt, kQ2Pt) = static_cast<double>(getSigma1Pt2()) + static_cast<double>(rhs.getSigma1Pt2());
783}
784
785//______________________________________________
786template <typename value_T>
787GPUd() auto TrackParametrizationWithError<value_T>::getPredictedChi2(const TrackParametrizationWithError<value_T>& rhs, MatrixDSym5& covToSet) const -> value_t
788{
789 // get chi2 wrt other track, which must be defined at the same parameters X,alpha
790 // Supplied non-initialized covToSet matrix is filled by inverse combined matrix for further use
791
792 if (gpu::CAMath::Abs(this->getAlpha() - rhs.getAlpha()) > o2::constants::math::Epsilon) {
793 LOG(error) << "The reference Alpha of the tracks differ: " << this->getAlpha() << " : " << rhs.getAlpha();
794 return 2.f * HugeF;
795 }
796 if (gpu::CAMath::Abs(this->getX() - rhs.getX()) > o2::constants::math::Epsilon) {
797 LOG(error) << "The reference X of the tracks differ: " << this->getX() << " : " << rhs.getX();
798 return 2.f * HugeF;
799 }
800 buildCombinedCovMatrix(rhs, covToSet);
801 if (!covToSet.Invert()) {
802 LOG(warning) << "Cov.matrix inversion failed: " << covToSet;
803 return 2.f * HugeF;
804 }
805 double chi2diag = 0., chi2ndiag = 0., diff[kNParams];
806 for (int i = kNParams; i--;) {
807 diff[i] = this->getParam(i) - rhs.getParam(i);
808 chi2diag += diff[i] * diff[i] * covToSet(i, i);
809 }
810 for (int i = kNParams; i--;) {
811 for (int j = i; j--;) {
812 chi2ndiag += diff[i] * diff[j] * covToSet(i, j);
813 }
814 }
815 return chi2diag + 2. * chi2ndiag;
816}
817
818//______________________________________________
819template <typename value_T>
820GPUd() bool TrackParametrizationWithError<value_T>::update(const TrackParametrizationWithError<value_T>& rhs, const MatrixDSym5& covInv)
821{
822 // update track with other track, the inverted combined cov matrix should be supplied
823
824 // consider skipping this check, since it is usually already done upstream
825 if (gpu::CAMath::Abs(this->getAlpha() - rhs.getAlpha()) > o2::constants::math::Epsilon) {
826 LOG(error) << "The reference Alpha of the tracks differ: " << this->getAlpha() << " : " << rhs.getAlpha();
827 return false;
828 }
829 if (gpu::CAMath::Abs(this->getX() - rhs.getX()) > o2::constants::math::Epsilon) {
830 LOG(error) << "The reference X of the tracks differ: " << this->getX() << " : " << rhs.getX();
831 return false;
832 }
833
834 // gain matrix K = Cov0*H*(Cov0+Cov0)^-1 (for measurement matrix H=I)
835 MatrixDSym5 matC0;
836 matC0(kY, kY) = getSigmaY2();
837 matC0(kZ, kY) = getSigmaZY();
838 matC0(kZ, kZ) = getSigmaZ2();
839 matC0(kSnp, kY) = getSigmaSnpY();
840 matC0(kSnp, kZ) = getSigmaSnpZ();
841 matC0(kSnp, kSnp) = getSigmaSnp2();
842 matC0(kTgl, kY) = getSigmaTglY();
843 matC0(kTgl, kZ) = getSigmaTglZ();
844 matC0(kTgl, kSnp) = getSigmaTglSnp();
845 matC0(kTgl, kTgl) = getSigmaTgl2();
846 matC0(kQ2Pt, kY) = getSigma1PtY();
847 matC0(kQ2Pt, kZ) = getSigma1PtZ();
848 matC0(kQ2Pt, kSnp) = getSigma1PtSnp();
849 matC0(kQ2Pt, kTgl) = getSigma1PtTgl();
850 matC0(kQ2Pt, kQ2Pt) = getSigma1Pt2();
851 MatrixD5 matK = matC0 * covInv;
852
853 // updated state vector: x = K*(x1-x0)
854 // RS: why SMatix, SVector does not provide multiplication operators ???
855 double diff[kNParams];
856 for (int i = kNParams; i--;) {
857 diff[i] = rhs.getParam(i) - this->getParam(i);
858 }
859 for (int i = kNParams; i--;) {
860 for (int j = kNParams; j--;) {
861 this->updateParam(matK(i, j) * diff[j], i);
862 }
863 }
864
865 // updated covariance: Cov0 = Cov0 - K*Cov0
867 mC[kSigY2] -= matK(kY, kY);
868 mC[kSigZY] -= matK(kZ, kY);
869 mC[kSigZ2] -= matK(kZ, kZ);
870 mC[kSigSnpY] -= matK(kSnp, kY);
871 mC[kSigSnpZ] -= matK(kSnp, kZ);
872 mC[kSigSnp2] -= matK(kSnp, kSnp);
873 mC[kSigTglY] -= matK(kTgl, kY);
874 mC[kSigTglZ] -= matK(kTgl, kZ);
875 mC[kSigTglSnp] -= matK(kTgl, kSnp);
876 mC[kSigTgl2] -= matK(kTgl, kTgl);
877 mC[kSigQ2PtY] -= matK(kQ2Pt, kY);
878 mC[kSigQ2PtZ] -= matK(kQ2Pt, kZ);
879 mC[kSigQ2PtSnp] -= matK(kQ2Pt, kSnp);
880 mC[kSigQ2PtTgl] -= matK(kQ2Pt, kTgl);
881 mC[kSigQ2Pt2] -= matK(kQ2Pt, kQ2Pt);
882
883 return true;
884}
885
886//______________________________________________
887template <typename value_T>
888GPUd() bool TrackParametrizationWithError<value_T>::update(const TrackParametrizationWithError<value_T>& rhs)
889{
890 // update track with other track
891 MatrixDSym5 covI; // perform matrix operations in double!
892 buildCombinedCovMatrix(rhs, covI);
893 if (!covI.Invert()) {
894 LOG(warning) << "Cov.matrix inversion failed: " << covI;
895 return false;
896 }
897 return update(rhs, covI);
898}
899
900//______________________________________________
901template <typename value_T>
902GPUd() bool TrackParametrizationWithError<value_T>::update(const value_t* p, const value_t* cov)
903{
904 // Update the track parameters with the space point "p" having
905 // the covariance matrix "cov"
906
907 value_t &cm00 = mC[kSigY2], &cm10 = mC[kSigZY], &cm11 = mC[kSigZ2], &cm20 = mC[kSigSnpY], &cm21 = mC[kSigSnpZ],
908 &cm22 = mC[kSigSnp2], &cm30 = mC[kSigTglY], &cm31 = mC[kSigTglZ], &cm32 = mC[kSigTglSnp], &cm33 = mC[kSigTgl2],
909 &cm40 = mC[kSigQ2PtY], &cm41 = mC[kSigQ2PtZ], &cm42 = mC[kSigQ2PtSnp], &cm43 = mC[kSigQ2PtTgl],
910 &cm44 = mC[kSigQ2Pt2];
911
912 // use double precision?
913 double r00 = static_cast<double>(cov[0]) + static_cast<double>(cm00);
914 double r01 = static_cast<double>(cov[1]) + static_cast<double>(cm10);
915 double r11 = static_cast<double>(cov[2]) + static_cast<double>(cm11);
916 double det = r00 * r11 - r01 * r01;
917
918 if (gpu::CAMath::Abs(det) < constants::math::Almost0) {
919 return false;
920 }
921 double detI = 1. / det;
922 double tmp = r00;
923 r00 = r11 * detI;
924 r11 = tmp * detI;
925 r01 = -r01 * detI;
926
927 double k00 = cm00 * r00 + cm10 * r01, k01 = cm00 * r01 + cm10 * r11;
928 double k10 = cm10 * r00 + cm11 * r01, k11 = cm10 * r01 + cm11 * r11;
929 double k20 = cm20 * r00 + cm21 * r01, k21 = cm20 * r01 + cm21 * r11;
930 double k30 = cm30 * r00 + cm31 * r01, k31 = cm30 * r01 + cm31 * r11;
931 double k40 = cm40 * r00 + cm41 * r01, k41 = cm40 * r01 + cm41 * r11;
932
933 value_t dy = p[kY] - this->getY(), dz = p[kZ] - this->getZ();
934 value_t dsnp = k20 * dy + k21 * dz;
935 if (gpu::CAMath::Abs(this->getSnp() + dsnp) > constants::math::Almost1) {
936 return false;
937 }
938
939 const params_t dP{value_t(k00 * dy + k01 * dz), value_t(k10 * dy + k11 * dz), dsnp, value_t(k30 * dy + k31 * dz),
940 value_t(k40 * dy + k41 * dz)};
941 this->updateParams(dP);
942
943 double c01 = cm10, c02 = cm20, c03 = cm30, c04 = cm40;
944 double c12 = cm21, c13 = cm31, c14 = cm41;
945
946 cm00 -= k00 * cm00 + k01 * cm10;
947 cm10 -= k00 * c01 + k01 * cm11;
948 cm20 -= k00 * c02 + k01 * c12;
949 cm30 -= k00 * c03 + k01 * c13;
950 cm40 -= k00 * c04 + k01 * c14;
951
952 cm11 -= k10 * c01 + k11 * cm11;
953 cm21 -= k10 * c02 + k11 * c12;
954 cm31 -= k10 * c03 + k11 * c13;
955 cm41 -= k10 * c04 + k11 * c14;
956
957 cm22 -= k20 * c02 + k21 * c12;
958 cm32 -= k20 * c03 + k21 * c13;
959 cm42 -= k20 * c04 + k21 * c14;
960
961 cm33 -= k30 * c03 + k31 * c13;
962 cm43 -= k30 * c04 + k31 * c14;
963
964 cm44 -= k40 * c04 + k41 * c14;
965
966 checkCovariance();
967
968 return true;
969}
970
971//______________________________________________
972template <typename value_T>
973GPUd() value_T TrackParametrizationWithError<value_T>::update(const o2::dataformats::VertexBase& vtx, value_T maxChi2)
974{
975 // Update track with vertex if the track-vertex chi2 does not exceed maxChi2. Track must be already propagated to the DCA to vertex
976 // return update chi2 or -chi2 if chi2 if chi2 exceeds maxChi2
977 auto vtLoc = this->getVertexInTrackFrame(vtx);
978 value_T chi2 = getPredictedChi2(vtLoc.yz, vtLoc.yzerr);
979 return chi2 < maxChi2 && update(vtLoc.yz, vtLoc.yzerr) ? chi2 : -chi2;
980}
981
982//______________________________________________
983template <typename value_T>
984GPUd() bool TrackParametrizationWithError<value_T>::correctForMaterial(value_t x2x0, value_t xrho, bool anglecorr)
985{
986 //------------------------------------------------------------------
987 // This function corrects the track parameters for the crossed material.
988 // "x2x0" - X/X0, the thickness in units of the radiation length.
989 // "xrho" - is the product length*density (g/cm^2).
990 // It should be passed as negative when propagating tracks
991 // from the intreaction point to the outside of the central barrel.
992 // "dedx" - mean enery loss (GeV/(g/cm^2), if <=kCalcdEdxAuto : calculate on the fly
993 // "anglecorr" - switch for the angular correction
994 //------------------------------------------------------------------
995 constexpr value_t kMSConst2 = 0.0136f * 0.0136f;
996 constexpr value_t kMinP = 0.01f; // kill below this momentum
997
998 value_t csp2 = (1.f - this->getSnp()) * (1.f + this->getSnp()); // cos(phi)^2
999 value_t cst2I = (1.f + this->getTgl() * this->getTgl()); // 1/cos(lambda)^2
1000 if (anglecorr) { // Apply angle correction, if requested
1001 value_t angle = gpu::CAMath::Sqrt(cst2I / csp2);
1002 x2x0 *= angle;
1003 xrho *= angle;
1004 }
1005 auto m = this->getPID().getMass();
1006 int charge2 = this->getAbsCharge() * this->getAbsCharge();
1007 value_t p = this->getP(), p0 = p, p02 = p * p, e2 = p02 + this->getPID().getMass2(), massInv = 1. / m, bg = p * massInv, dETot = 0.;
1008 value_t e = gpu::CAMath::Sqrt(e2), e0 = e;
1009 if (m > 0 && xrho != 0.f) {
1010 value_t ekin = e - m, dedx = this->getdEdxBBOpt(bg);
1011#ifdef _BB_NONCONST_CORR_
1012 value_t dedxDer = 0., dedx1 = dedx;
1013#endif
1014 if (charge2 != 1) {
1015 dedx *= charge2;
1016 }
1017 value_t dE = dedx * xrho;
1018 int na = 1 + int(gpu::CAMath::Abs(dE) / ekin * ELoss2EKinThreshInv);
1019 if (na > MaxELossIter) {
1020 na = MaxELossIter;
1021 }
1022 if (na > 1) {
1023 dE /= na;
1024 xrho /= na;
1025#ifdef _BB_NONCONST_CORR_
1026 dedxDer = this->getBetheBlochSolidDerivativeApprox(dedx1, bg); // require correction for non-constantness of dedx vs betagamma
1027 if (charge2 != 1) {
1028 dedxDer *= charge2;
1029 }
1030#endif
1031 }
1032 while (na--) {
1033#ifdef _BB_NONCONST_CORR_
1034 if (dedxDer != 0.) { // correction for non-constantness of dedx vs beta*gamma (in linear approximation): for a single step dE -> dE * [(exp(dedxDer) - 1)/dedxDer]
1035 if (xrho < 0) {
1036 dedxDer = -dedxDer; // E.loss ( -> positive derivative)
1037 }
1038 auto corrC = (gpu::CAMath::Exp(dedxDer) - 1.) / dedxDer;
1039 dE *= corrC;
1040 }
1041#endif
1042 e += dE;
1043 if (e > m) { // stopped
1044 p = gpu::CAMath::Sqrt(e * e - this->getPID().getMass2());
1045 } else {
1046 return false;
1047 }
1048 if (na) {
1049 bg = p * massInv;
1050 dedx = this->getdEdxBBOpt(bg);
1051#ifdef _BB_NONCONST_CORR_
1052 dedxDer = this->getBetheBlochSolidDerivativeApprox(dedx, bg);
1053#endif
1054 if (charge2 != 1) {
1055 dedx *= charge2;
1056#ifdef _BB_NONCONST_CORR_
1057 dedxDer *= charge2;
1058#endif
1059 }
1060 dE = dedx * xrho;
1061 }
1062 }
1063
1064 if (p < kMinP) {
1065 return false;
1066 }
1067 dETot = e - e0;
1068 } // end of e.loss correction
1069
1070 // Calculating the multiple scattering corrections******************
1071 value_t& fC22 = mC[kSigSnp2];
1072 value_t& fC33 = mC[kSigTgl2];
1073 value_t& fC43 = mC[kSigQ2PtTgl];
1074 value_t& fC44 = mC[kSigQ2Pt2];
1075 //
1076 value_t cC22(0.f), cC33(0.f), cC43(0.f), cC44(0.f);
1077 if (x2x0 != 0.f) {
1078 value_t beta2 = p02 / e2, theta2 = kMSConst2 / (beta2 * p02) * gpu::CAMath::Abs(x2x0);
1079 value_t fp34 = this->getTgl();
1080 if (charge2 != 1) {
1081 theta2 *= charge2;
1082 fp34 *= this->getCharge2Pt();
1083 }
1084 if (theta2 > constants::math::PI * constants::math::PI) {
1085 return false;
1086 }
1087 value_t t2c2I = theta2 * cst2I;
1088 cC22 = t2c2I * csp2;
1089 cC33 = t2c2I * cst2I;
1090 cC43 = t2c2I * fp34;
1091 cC44 = theta2 * fp34 * fp34;
1092 // optimize this
1093 // cC22 = theta2*((1.-getSnp())*(1.+getSnp()))*(1. + this->getTgl()*getTgl());
1094 // cC33 = theta2*(1. + this->getTgl()*getTgl())*(1. + this->getTgl()*getTgl());
1095 // cC43 = theta2*getTgl()*this->getQ2Pt()*(1. + this->getTgl()*getTgl());
1096 // cC44 = theta2*getTgl()*this->getQ2Pt()*getTgl()*this->getQ2Pt();
1097 }
1098
1099 // the energy loss correction contribution to cov.matrix: approximate energy loss fluctuation (M.Ivanov)
1100 constexpr value_t knst = 0.0007f; // To be tuned.
1101 value_t sigmadE = knst * gpu::CAMath::Sqrt(gpu::CAMath::Abs(dETot)) * e0 / p02 * this->getCharge2Pt();
1102 cC44 += sigmadE * sigmadE;
1103
1104 // Applying the corrections*****************************
1105 fC22 += cC22;
1106 fC33 += cC33;
1107 fC43 += cC43;
1108 fC44 += cC44;
1109 this->setQ2Pt(this->getQ2Pt() * p0 / p);
1110
1111 checkCovariance();
1112
1113 return true;
1114}
1115
1116//______________________________________________________________
1117template <typename value_T>
1118GPUd() bool TrackParametrizationWithError<value_T>::getCovXYZPxPyPzGlo(gpu::gpustd::array<value_t, kLabCovMatSize>& cv) const
1119{
1120 //---------------------------------------------------------------------
1121 // This function returns the global covariance matrix of the track params
1122 //
1123 // Cov(x,x) ... : cv[0]
1124 // Cov(y,x) ... : cv[1] cv[2]
1125 // Cov(z,x) ... : cv[3] cv[4] cv[5]
1126 // Cov(px,x)... : cv[6] cv[7] cv[8] cv[9]
1127 // Cov(py,x)... : cv[10] cv[11] cv[12] cv[13] cv[14]
1128 // Cov(pz,x)... : cv[15] cv[16] cv[17] cv[18] cv[19] cv[20]
1129 //
1130 // Results for (nearly) straight tracks are meaningless !
1131 //---------------------------------------------------------------------
1132 if (gpu::CAMath::Abs(this->getQ2Pt()) <= constants::math::Almost0 || gpu::CAMath::Abs(this->getSnp()) > constants::math::Almost1) {
1133 for (int i = 0; i < 21; i++) {
1134 cv[i] = 0.;
1135 }
1136 return false;
1137 }
1138
1139 auto pt = this->getPt();
1140 value_t sn, cs;
1141 o2::math_utils::detail::sincos(this->getAlpha(), sn, cs);
1142 auto r = gpu::CAMath::Sqrt((1. - this->getSnp()) * (1. + this->getSnp()));
1143 auto m00 = -sn, m10 = cs;
1144 auto m23 = -pt * (sn + this->getSnp() * cs / r), m43 = -pt * pt * (r * cs - this->getSnp() * sn);
1145 auto m24 = pt * (cs - this->getSnp() * sn / r), m44 = -pt * pt * (r * sn + this->getSnp() * cs);
1146 auto m35 = pt, m45 = -pt * pt * this->getTgl();
1147
1148 if (this->getSign() < 0) {
1149 m43 = -m43;
1150 m44 = -m44;
1151 m45 = -m45;
1152 }
1153
1154 cv[0] = mC[0] * m00 * m00;
1155 cv[1] = mC[0] * m00 * m10;
1156 cv[2] = mC[0] * m10 * m10;
1157 cv[3] = mC[1] * m00;
1158 cv[4] = mC[1] * m10;
1159 cv[5] = mC[2];
1160 cv[6] = m00 * (mC[3] * m23 + mC[10] * m43);
1161 cv[7] = m10 * (mC[3] * m23 + mC[10] * m43);
1162 cv[8] = mC[4] * m23 + mC[11] * m43;
1163 cv[9] = m23 * (mC[5] * m23 + mC[12] * m43) + m43 * (mC[12] * m23 + mC[14] * m43);
1164 cv[10] = m00 * (mC[3] * m24 + mC[10] * m44);
1165 cv[11] = m10 * (mC[3] * m24 + mC[10] * m44);
1166 cv[12] = mC[4] * m24 + mC[11] * m44;
1167 cv[13] = m23 * (mC[5] * m24 + mC[12] * m44) + m43 * (mC[12] * m24 + mC[14] * m44);
1168 cv[14] = m24 * (mC[5] * m24 + mC[12] * m44) + m44 * (mC[12] * m24 + mC[14] * m44);
1169 cv[15] = m00 * (mC[6] * m35 + mC[10] * m45);
1170 cv[16] = m10 * (mC[6] * m35 + mC[10] * m45);
1171 cv[17] = mC[7] * m35 + mC[11] * m45;
1172 cv[18] = m23 * (mC[8] * m35 + mC[12] * m45) + m43 * (mC[13] * m35 + mC[14] * m45);
1173 cv[19] = m24 * (mC[8] * m35 + mC[12] * m45) + m44 * (mC[13] * m35 + mC[14] * m45);
1174 cv[20] = m35 * (mC[9] * m35 + mC[13] * m45) + m45 * (mC[13] * m35 + mC[14] * m45);
1175
1176 return true;
1177}
1178
1179#ifndef GPUCA_ALIGPUCODE
1180//______________________________________________________________
1181template <typename value_T>
1183{
1185 fmt::format(" Cov: [{:+.3e}] [{:+.3e} {:+.3e}] [{:+.3e} {:+.3e} {:+.3e}] [{:+.3e} {:+.3e} {:+.3e} {:+.3e}] [{:+.3e} {:+.3e} {:+.3e} {:+.3e} {:+.3e}]",
1186 mC[kSigY2], mC[kSigZY], mC[kSigZ2], mC[kSigSnpY], mC[kSigSnpZ], mC[kSigSnp2], mC[kSigTglY],
1187 mC[kSigTglZ], mC[kSigTglSnp], mC[kSigTgl2], mC[kSigQ2PtY], mC[kSigQ2PtZ], mC[kSigQ2PtSnp], mC[kSigQ2PtTgl],
1188 mC[kSigQ2Pt2]);
1189}
1190
1191template <typename value_T>
1193{
1195 fmt::format(" Cov: [{:x}] [{:x} {:x}] [{:x} {:x} {:x}] [{:x} {:x} {:x} {:x}] [{:x} {:x} {:x} {:x} {:x}]",
1196 reinterpret_cast<const unsigned int&>(mC[kSigY2]), reinterpret_cast<const unsigned int&>(mC[kSigZY]), reinterpret_cast<const unsigned int&>(mC[kSigZ2]),
1197 reinterpret_cast<const unsigned int&>(mC[kSigSnpY]), reinterpret_cast<const unsigned int&>(mC[kSigSnpZ]), reinterpret_cast<const unsigned int&>(mC[kSigSnp2]),
1198 reinterpret_cast<const unsigned int&>(mC[kSigTglY]), reinterpret_cast<const unsigned int&>(mC[kSigTglZ]), reinterpret_cast<const unsigned int&>(mC[kSigTglSnp]),
1199 reinterpret_cast<const unsigned int&>(mC[kSigTgl2]), reinterpret_cast<const unsigned int&>(mC[kSigQ2PtY]), reinterpret_cast<const unsigned int&>(mC[kSigQ2PtZ]),
1200 reinterpret_cast<const unsigned int&>(mC[kSigQ2PtSnp]), reinterpret_cast<const unsigned int&>(mC[kSigQ2PtTgl]), reinterpret_cast<const unsigned int&>(mC[kSigQ2Pt2]));
1201}
1202#endif
1203
1204//______________________________________________________________
1205template <typename value_T>
1206GPUd() void TrackParametrizationWithError<value_T>::print() const
1207{
1208 // print parameters
1209#ifndef GPUCA_ALIGPUCODE
1210 printf("%s\n", asString().c_str());
1211#elif !defined(GPUCA_GPUCODE_DEVICE) || (!defined(__OPENCL__) && defined(GPUCA_GPU_DEBUG_PRINT))
1213 printf(
1214 " Cov: [%+.3e] [%+.3e %+.3e] [%+.3e %+.3e %+.3e] [%+.3e %+.3e %+.3e %+.3e] [%+.3e %+.3e %+.3e %+.3e %+.3e]\n",
1215 mC[kSigY2], mC[kSigZY], mC[kSigZ2], mC[kSigSnpY], mC[kSigSnpZ], mC[kSigSnp2], mC[kSigTglY],
1216 mC[kSigTglZ], mC[kSigTglSnp], mC[kSigTgl2], mC[kSigQ2PtY], mC[kSigQ2PtZ], mC[kSigQ2PtSnp], mC[kSigQ2PtTgl],
1217 mC[kSigQ2Pt2]);
1218#endif
1219}
1220
1221//______________________________________________________________
1222template <typename value_T>
1223GPUd() void TrackParametrizationWithError<value_T>::printHexadecimal()
1224{
1225 // print parameters
1226#ifndef GPUCA_ALIGPUCODE
1227 printf("%s\n", asStringHexadecimal().c_str());
1228#elif !defined(GPUCA_GPUCODE_DEVICE) || (!defined(__OPENCL__) && defined(GPUCA_GPU_DEBUG_PRINT))
1230 printf(
1231 " Cov: [%x] [%x %x] [%x %x %x] [%x %x %x %x] [%x %x %x %x %x]\n",
1232 gpu::CAMath::Float2UIntReint(mC[kSigY2]),
1233 gpu::CAMath::Float2UIntReint(mC[kSigZY]), gpu::CAMath::Float2UIntReint(mC[kSigZ2]),
1234 gpu::CAMath::Float2UIntReint(mC[kSigSnpY]), gpu::CAMath::Float2UIntReint(mC[kSigSnpZ]), gpu::CAMath::Float2UIntReint(mC[kSigSnp2]),
1235 gpu::CAMath::Float2UIntReint(mC[kSigTglY]), gpu::CAMath::Float2UIntReint(mC[kSigTglZ]), gpu::CAMath::Float2UIntReint(mC[kSigTglSnp]), gpu::CAMath::Float2UIntReint(mC[kSigTgl2]),
1236 gpu::CAMath::Float2UIntReint(mC[kSigQ2PtY]), gpu::CAMath::Float2UIntReint(mC[kSigQ2PtZ]), gpu::CAMath::Float2UIntReint(mC[kSigQ2PtSnp]), gpu::CAMath::Float2UIntReint(mC[kSigQ2PtTgl]), gpu::CAMath::Float2UIntReint(mC[kSigQ2Pt2]));
1237#endif
1238}
1239
1240namespace o2::track
1241{
1242#if !defined(GPUCA_GPUCODE) || defined(GPUCA_GPUCODE_DEVICE) // FIXME: DR: WORKAROUND to avoid CUDA bug creating host symbols for device code.
1244#endif
1245#ifndef GPUCA_GPUCODE
1247#endif
1248} // namespace o2::track
std::string asString(TDataMember const &dm, char *pointer)
int16_t charge
Definition RawEventData.h:5
void print() const
int32_t i
const GPUTPCGMMerger::trackCluster & b1
const GPUTPCGMMerger::trackCluster & a1
const int16_t bb
useful math constants
uint32_t j
Definition RawData.h:0
uint16_t pid
Definition RawData.h:2
std::ostringstream debug
std::string asString() const
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
const GLfloat * m
Definition glcorearb.h:4066
GLenum array
Definition glcorearb.h:4274
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLfloat angle
Definition glcorearb.h:4071
GLboolean r
Definition glcorearb.h:1233
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
GLboolean invert
Definition glcorearb.h:543
void set(T *h, size_t v)
Definition PageParser.h:107
typename trackParam_t::params_t params_t
Definition utils.h:33
typename trackParam_t::dim3_t dim3_t
Definition utils.h:32
constexpr float Epsilon
constexpr float Almost1
std::array< T, N > array
D const SVectorGPU< T, D > & rhs
Definition SMatrixGPU.h:191
std::array< int, 24 > p0
double * getX(double *xyDxy, int N)
double * getY(double *xyDxy, int N)
return * this
value_T bg
Definition TrackUtils.h:194
constexpr float kCTgl2max
constexpr float HugeF
constexpr int kCovMatSize
constexpr float kCSnp2max
constexpr int kNParams
constexpr int kLabCovMatSize
value_T step
Definition TrackUtils.h:42
value_T f1
Definition TrackUtils.h:91
const value_T x
Definition TrackUtils.h:136
constexpr float kCY2max
value_T d2
Definition TrackUtils.h:135
value_T gpu::gpustd::array< value_T, 7 > & vect
Definition TrackUtils.h:42
kp1 *kp2 *value_T beta2
Definition TrackUtils.h:131
GPUd() value_T BetheBlochSolid(value_T bg
Definition TrackUtils.h:150
constexpr float kC1Pt2max
value_T f2
Definition TrackUtils.h:92
constexpr float kCZ2max
constexpr int MaxELossIter
constexpr float ELoss2EKinThreshInv
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"