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